Changeset 1441 for trunk/MagicSoft/Mars/mhist/MHHillas.cc
- Timestamp:
- 07/25/02 13:58:19 (22 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/MagicSoft/Mars/mhist/MHHillas.cc
r1283 r1441 17 17 ! 18 18 ! Author(s): Thomas Bretz 2001 <mailto:tbretz@uni-sw.gwdg.de> 19 ! Wolfgang Wittek 2002 <mailto:wittek@mppmu.mpg.de> 19 20 ! 20 21 ! Copyright: MAGIC Software Development, 2000-2002 … … 27 28 // MHHillas 28 29 // 29 // This class contains histograms for every Hillas parameter30 // This class contains histograms for the source independent image parameters 30 31 // 31 32 ///////////////////////////////////////////////////////////////////////////// … … 36 37 37 38 #include <TH1.h> 39 #include <TH2.h> 38 40 #include <TPad.h> 41 #include <TStyle.h> 39 42 #include <TCanvas.h> 40 43 … … 42 45 #include "MLogManip.h" 43 46 47 #include "MParList.h" 48 49 #include "MHillas.h" 44 50 #include "MGeomCam.h" 45 #include "M Hillas.h"46 #include "MParList.h" 51 #include "MBinning.h" 52 47 53 48 54 ClassImp(MHHillas); … … 53 59 // 54 60 MHHillas::MHHillas(const char *name, const char *title) 55 : fMm2Deg( -1), fUseMmScale(kTRUE)61 : fMm2Deg(1), fUseMmScale(kFALSE) 56 62 { 57 63 // … … 59 65 // 60 66 fName = name ? name : "MHHillas"; 61 fTitle = title ? title : " Container for Hillas histograms";62 63 //64 // loop over all Pixels and create two histograms65 // one for the Low and one for the High gain66 // connect all the histogram with the container fHist67 // 67 fTitle = title ? title : "Source independent image parameters"; 68 69 fLength = new TH1F("Length", "Length of Ellipse", 100, 0, 300); 70 fLength->SetDirectory(NULL); 71 fLength->SetXTitle("Length [mm]"); 72 fLength->SetYTitle("Counts"); 73 68 74 fWidth = new TH1F("Width", "Width of Ellipse", 100, 0, 300); 69 fLength = new TH1F("Length", "Length of Ellipse", 100, 0, 300);70 71 fLength->SetDirectory(NULL);72 75 fWidth->SetDirectory(NULL); 73 74 fLength->SetXTitle("Length [mm]");75 76 fWidth->SetXTitle("Width [mm]"); 76 77 fLength->SetYTitle("Counts");78 77 fWidth->SetYTitle("Counts"); 79 } 80 81 // -------------------------------------------------------------------------- 82 // 83 // Delete the four histograms 78 79 fDistC = new TH1F("DistC", "Distance from center of camera", 100, 0, 600); 80 fDistC->SetDirectory(NULL); 81 fDistC->SetXTitle("Distance [mm]"); 82 fDistC->SetYTitle("Counts"); 83 84 fDelta = new TH1F("Delta", "Angle (Main axis - x-axis)", 100, -90, 90); 85 fDelta->SetDirectory(NULL); 86 fDelta->SetXTitle("Delta [\\circ]"); 87 fDelta->SetYTitle("Counts"); 88 89 MBinning bins; 90 bins.SetEdgesLog(50, 1, 1e7); 91 92 fSize = new TH1F; 93 fSize->SetName("Size"); 94 fSize->SetTitle("Number of Photons"); 95 fSize->SetDirectory(NULL); 96 fSize->SetXTitle("Size"); 97 fSize->SetYTitle("Counts"); 98 fSize->GetXaxis()->SetTitleOffset(1.2); 99 fSize->GetXaxis()->SetLabelOffset(-0.015); 100 101 MH::SetBinning(fSize, &bins); 102 103 fCenter = new TH2F("Center", "Center of Ellipse", 60, -600, 600, 60, -600, 600); 104 fCenter->SetDirectory(NULL); 105 fCenter->SetXTitle("x [mm]"); 106 fCenter->SetYTitle("y [mm]"); 107 fCenter->SetZTitle("Counts"); 108 } 109 110 // -------------------------------------------------------------------------- 111 // 112 // Delete the histograms 84 113 // 85 114 MHHillas::~MHHillas() 86 115 { 116 delete fLength; 87 117 delete fWidth; 88 delete fLength; 118 119 delete fDistC; 120 delete fDelta; 121 122 delete fSize; 123 delete fCenter; 89 124 } 90 125 … … 101 136 Bool_t MHHillas::SetupFill(const MParList *plist) 102 137 { 103 const MBinning* binsw = (MBinning*)plist->FindObject("BinningWidth"); 104 const MBinning* binsl = (MBinning*)plist->FindObject("BinningLength"); 105 if (!binsw || !binsl) 138 const MBinning *binsw = (MBinning*)plist->FindObject("BinningWidth"); 139 const MBinning *binsl = (MBinning*)plist->FindObject("BinningLength"); 140 const MBinning *binsd = (MBinning*)plist->FindObject("BinningDist"); 141 const MBinning *binsc = (MBinning*)plist->FindObject("BinningCamera"); 142 143 if (!binsw || !binsl || !binsd || !binsc) 106 144 { 107 145 *fLog << err << dbginf << "At least one MBinning not found... aborting." << endl; … … 111 149 SetBinning(fWidth, binsw); 112 150 SetBinning(fLength, binsl); 151 SetBinning(fDistC, binsd); 152 SetBinning(fCenter, binsc, binsc); 113 153 114 154 const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam"); … … 119 159 } 120 160 161 fMm2Deg = geom->GetConvMm2Deg(); 162 fUseMmScale = kFALSE; 163 121 164 fLength->SetXTitle("Length [\\circ]"); 122 165 fWidth->SetXTitle("Width [\\circ]"); 123 124 fMm2Deg = geom->GetConvMm2Deg(); 125 fUseMmScale = kFALSE; 126 127 return kTRUE; 128 } 129 130 // -------------------------------------------------------------------------- 131 // 132 // Fill the four histograms with data from a MHillas-Container. 133 // Be careful: Only call this with an object of type MHillas 134 // 135 Bool_t MHHillas::Fill(const MParContainer *par) 136 { 137 const MHillas &h = *(MHillas*)par; 138 139 if (fUseMmScale) 140 { 141 fWidth ->Fill(h.GetWidth()); 142 fLength->Fill(h.GetLength()); 143 } 144 else 145 { 146 fWidth ->Fill(fMm2Deg*h.GetWidth()); 147 fLength->Fill(fMm2Deg*h.GetLength()); 148 } 166 fDistC->SetXTitle("Distance [\\circ]"); 167 fCenter->SetXTitle("x [\\circ]"); 168 fCenter->SetYTitle("y [\\circ]"); 149 169 150 170 return kTRUE; … … 190 210 if (fUseMmScale) 191 211 { 212 MH::ScaleAxis(fLength, 1./fMm2Deg); 213 MH::ScaleAxis(fWidth, 1./fMm2Deg); 214 MH::ScaleAxis(fDistC, 1./fMm2Deg); 215 MH::ScaleAxis(fCenter, 1./fMm2Deg, 1./fMm2Deg); 216 192 217 fLength->SetXTitle("Length [mm]"); 193 218 fWidth->SetXTitle("Width [mm]"); 194 195 fLength->Scale(1./fMm2Deg); 196 fWidth->Scale(1./fMm2Deg); 219 fCenter->SetXTitle("x [mm]"); 220 fCenter->SetYTitle("y [mm]"); 197 221 } 198 222 else 199 223 { 224 MH::ScaleAxis(fLength, fMm2Deg); 225 MH::ScaleAxis(fWidth, fMm2Deg); 226 MH::ScaleAxis(fDistC, fMm2Deg); 227 MH::ScaleAxis(fCenter, fMm2Deg, fMm2Deg); 228 200 229 fLength->SetXTitle("Length [\\circ]"); 201 230 fWidth->SetXTitle("Width [\\circ]"); 202 203 f Length->Scale(fMm2Deg);204 f Width->Scale(fMm2Deg);231 fDistC->SetXTitle("Distance [\\circ]"); 232 fCenter->SetXTitle("x [\\circ]"); 233 fCenter->SetYTitle("y [\\circ]"); 205 234 } 206 235 … … 210 239 // -------------------------------------------------------------------------- 211 240 // 212 // Draw clones of all four histograms. So that the object can be deleted 241 // Fill the histograms with data from a MHillas-Container. 242 // Be careful: Only call this with an object of type MHillas 243 // 244 Bool_t MHHillas::Fill(const MParContainer *par) 245 { 246 const MHillas &h = *(MHillas*)par; 247 248 const Double_t d = sqrt(h.GetMeanX()*h.GetMeanX() + h.GetMeanY()*h.GetMeanY()); 249 250 if (fUseMmScale) 251 { 252 fLength->Fill(h.GetLength()); 253 fWidth ->Fill(h.GetWidth()); 254 fDistC ->Fill(d); 255 fCenter->Fill(h.GetMeanX(), h.GetMeanY()); 256 } 257 else 258 { 259 fLength->Fill(fMm2Deg*h.GetLength()); 260 fWidth ->Fill(fMm2Deg*h.GetWidth()); 261 fDistC ->Fill(fMm2Deg*d); 262 fCenter->Fill(fMm2Deg*h.GetMeanX(), fMm2Deg*h.GetMeanY()); 263 } 264 265 fDelta->Fill(kRad2Deg*h.GetDelta()); 266 fSize->Fill(h.GetSize()); 267 268 return kTRUE; 269 } 270 271 // -------------------------------------------------------------------------- 272 // 273 // Setup a inversed deep blue sea palette for the fCenter histogram. 274 // 275 void MHHillas::SetColors() const 276 { 277 gStyle->SetPalette(51, NULL); 278 Int_t c[50]; 279 for (int i=0; i<50; i++) 280 c[49-i] = gStyle->GetColorPalette(i); 281 gStyle->SetPalette(50, c); 282 } 283 284 // -------------------------------------------------------------------------- 285 // 286 // Draw clones of four histograms. So that the object can be deleted 213 287 // and the histograms are still visible in the canvas. 214 288 // The cloned object are deleted together with the canvas if the canvas is … … 218 292 TObject *MHHillas::DrawClone(Option_t *opt) const 219 293 { 220 TCanvas *c = MakeDefCanvas("Hillas", "Histograms of Hillas Parameters", 221 350, 500); 222 c->Divide(1, 2); 294 TCanvas *c = MakeDefCanvas("Hillas", fTitle, 700, 750); 295 c->Divide(2,3); 223 296 224 297 gROOT->SetSelectedPad(NULL); 225 298 226 //227 // This is necessary to get the expected bahviour of DrawClone228 //229 299 c->cd(1); 230 300 fLength->DrawCopy(); … … 232 302 c->cd(2); 233 303 fWidth->DrawCopy(); 304 305 c->cd(3); 306 gPad->SetLogx(); 307 fSize->DrawCopy(); 308 309 c->cd(4); 310 fDelta->DrawCopy(); 311 312 c->cd(5); 313 fDistC->DrawCopy(); 314 315 c->cd(6); 316 SetColors(); 317 fCenter->DrawCopy("colz"); 234 318 235 319 c->Modified(); … … 248 332 { 249 333 if (!gPad) 250 MakeDefCanvas("Hillas", "Histograms of Hillas Parameters", 350, 500);251 252 gPad->Divide( 1, 2);334 MakeDefCanvas("Hillas", fTitle, 700, 750); 335 336 gPad->Divide(2,3); 253 337 254 338 gPad->cd(1); … … 258 342 fWidth->Draw(); 259 343 344 gPad->cd(3); 345 gPad->SetLogx(); 346 fSize->Draw(); 347 348 gPad->cd(4); 349 fDelta->Draw(); 350 351 gPad->cd(5); 352 fDistC->Draw(); 353 354 gPad->cd(6); 355 SetColors(); 356 fCenter->Draw("colz"); 357 260 358 gPad->Modified(); 261 359 gPad->Update();
Note:
See TracChangeset
for help on using the changeset viewer.