| 1 | /* ======================================================================== *\ | 
|---|
| 2 | ! | 
|---|
| 3 | ! * | 
|---|
| 4 | ! * This file is part of MARS, the MAGIC Analysis and Reconstruction | 
|---|
| 5 | ! * Software. It is distributed to you in the hope that it can be a useful | 
|---|
| 6 | ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes. | 
|---|
| 7 | ! * It is distributed WITHOUT ANY WARRANTY. | 
|---|
| 8 | ! * | 
|---|
| 9 | ! * Permission to use, copy, modify and distribute this software and its | 
|---|
| 10 | ! * documentation for any purpose is hereby granted without fee, | 
|---|
| 11 | ! * provided that the above copyright notice appear in all copies and | 
|---|
| 12 | ! * that both that copyright notice and this permission notice appear | 
|---|
| 13 | ! * in supporting documentation. It is provided "as is" without express | 
|---|
| 14 | ! * or implied warranty. | 
|---|
| 15 | ! * | 
|---|
| 16 | ! | 
|---|
| 17 | ! | 
|---|
| 18 | !   Author(s): Thomas Bretz  2001 <mailto:tbretz@uni-sw.gwdg.de> | 
|---|
| 19 | !              Wolfgang Wittek  2002 <mailto:wittek@mppmu.mpg.de> | 
|---|
| 20 | ! | 
|---|
| 21 | !   Copyright: MAGIC Software Development, 2000-2002 | 
|---|
| 22 | ! | 
|---|
| 23 | ! | 
|---|
| 24 | \* ======================================================================== */ | 
|---|
| 25 |  | 
|---|
| 26 | ///////////////////////////////////////////////////////////////////////////// | 
|---|
| 27 | // | 
|---|
| 28 | // MHHillas | 
|---|
| 29 | // | 
|---|
| 30 | // This class contains histograms for the source independent image parameters | 
|---|
| 31 | // | 
|---|
| 32 | ///////////////////////////////////////////////////////////////////////////// | 
|---|
| 33 | #include "MHHillas.h" | 
|---|
| 34 |  | 
|---|
| 35 | #include <math.h> | 
|---|
| 36 |  | 
|---|
| 37 | #include <TH1.h> | 
|---|
| 38 | #include <TH2.h> | 
|---|
| 39 | #include <TPad.h> | 
|---|
| 40 | #include <TStyle.h> | 
|---|
| 41 | #include <TCanvas.h> | 
|---|
| 42 |  | 
|---|
| 43 | #include "MLog.h" | 
|---|
| 44 | #include "MLogManip.h" | 
|---|
| 45 |  | 
|---|
| 46 | #include "MParList.h" | 
|---|
| 47 |  | 
|---|
| 48 | #include "MHillas.h" | 
|---|
| 49 | #include "MGeomCam.h" | 
|---|
| 50 | #include "MBinning.h" | 
|---|
| 51 |  | 
|---|
| 52 |  | 
|---|
| 53 | ClassImp(MHHillas); | 
|---|
| 54 |  | 
|---|
| 55 | using namespace std; | 
|---|
| 56 |  | 
|---|
| 57 | // -------------------------------------------------------------------------- | 
|---|
| 58 | // | 
|---|
| 59 | // Setup four histograms for Width, Length | 
|---|
| 60 | // | 
|---|
| 61 | MHHillas::MHHillas(const char *name, const char *title) | 
|---|
| 62 | : fMm2Deg(1), fUseMmScale(kTRUE) | 
|---|
| 63 | { | 
|---|
| 64 | // | 
|---|
| 65 | //   set the name and title of this object | 
|---|
| 66 | // | 
|---|
| 67 | fName  = name  ? name  : "MHHillas"; | 
|---|
| 68 | fTitle = title ? title : "Source independent image parameters"; | 
|---|
| 69 |  | 
|---|
| 70 | fLength  = new TH1F("Length",  "Length of Ellipse",               100,   0, 296.7); | 
|---|
| 71 | fWidth   = new TH1F("Width",   "Width of Ellipse",                100,   0, 296.7); | 
|---|
| 72 | fDistC   = new TH1F("DistC",   "Distance from center of camera",  100,   0, 445); | 
|---|
| 73 | fDelta   = new TH1F("Delta",   "Angle (Main axis - x-axis)",      101, -90,  90); | 
|---|
| 74 |  | 
|---|
| 75 | fLength->SetLineColor(kBlue); | 
|---|
| 76 |  | 
|---|
| 77 | fLength->SetDirectory(NULL); | 
|---|
| 78 | fWidth->SetDirectory(NULL); | 
|---|
| 79 | fDistC->SetDirectory(NULL); | 
|---|
| 80 | fDelta->SetDirectory(NULL); | 
|---|
| 81 |  | 
|---|
| 82 | fLength->SetXTitle("Length [mm]"); | 
|---|
| 83 | fWidth->SetXTitle("Width [mm]"); | 
|---|
| 84 | fDistC->SetXTitle("Distance [mm]"); | 
|---|
| 85 | fDelta->SetXTitle("Delta [\\circ]"); | 
|---|
| 86 |  | 
|---|
| 87 | fLength->SetYTitle("Counts"); | 
|---|
| 88 | fWidth->SetYTitle("Counts"); | 
|---|
| 89 | fDistC->SetYTitle("Counts"); | 
|---|
| 90 | fDelta->SetYTitle("Counts"); | 
|---|
| 91 |  | 
|---|
| 92 | MBinning bins; | 
|---|
| 93 | bins.SetEdgesLog(50, 1, 1e7); | 
|---|
| 94 |  | 
|---|
| 95 | fSize  = new TH1F; | 
|---|
| 96 | fSize->SetName("Size"); | 
|---|
| 97 | fSize->SetTitle("Number of Photons"); | 
|---|
| 98 | fSize->SetDirectory(NULL); | 
|---|
| 99 | fSize->SetXTitle("Size"); | 
|---|
| 100 | fSize->SetYTitle("Counts"); | 
|---|
| 101 | fSize->GetXaxis()->SetTitleOffset(1.2); | 
|---|
| 102 | fSize->GetXaxis()->SetLabelOffset(-0.015); | 
|---|
| 103 | fSize->SetFillStyle(4000); | 
|---|
| 104 | fSize->UseCurrentStyle(); | 
|---|
| 105 |  | 
|---|
| 106 | bins.Apply(*fSize); | 
|---|
| 107 |  | 
|---|
| 108 | fCenter = new TH2F("Center", "Center of Ellipse", 51, -445, 445, 51, -445, 445); | 
|---|
| 109 | fCenter->SetDirectory(NULL); | 
|---|
| 110 | fCenter->SetXTitle("x [mm]"); | 
|---|
| 111 | fCenter->SetYTitle("y [mm]"); | 
|---|
| 112 | fCenter->SetZTitle("Counts"); | 
|---|
| 113 | } | 
|---|
| 114 |  | 
|---|
| 115 | // -------------------------------------------------------------------------- | 
|---|
| 116 | // | 
|---|
| 117 | // Delete the histograms | 
|---|
| 118 | // | 
|---|
| 119 | MHHillas::~MHHillas() | 
|---|
| 120 | { | 
|---|
| 121 | delete fLength; | 
|---|
| 122 | delete fWidth; | 
|---|
| 123 |  | 
|---|
| 124 | delete fDistC; | 
|---|
| 125 | delete fDelta; | 
|---|
| 126 |  | 
|---|
| 127 | delete fSize; | 
|---|
| 128 | delete fCenter; | 
|---|
| 129 | } | 
|---|
| 130 |  | 
|---|
| 131 | // -------------------------------------------------------------------------- | 
|---|
| 132 | // | 
|---|
| 133 | // Setup the Binning for the histograms automatically if the correct | 
|---|
| 134 | // instances of MBinning (with the names 'BinningWidth' and 'BinningLength') | 
|---|
| 135 | // are found in the parameter list | 
|---|
| 136 | // Use this function if you want to set the conversion factor which | 
|---|
| 137 | // is used to convert the mm-scale in the camera plain into the deg-scale | 
|---|
| 138 | // used for histogram presentations. The conversion factor is part of | 
|---|
| 139 | // the camera geometry. Please create a corresponding MGeomCam container. | 
|---|
| 140 | // | 
|---|
| 141 | Bool_t MHHillas::SetupFill(const MParList *plist) | 
|---|
| 142 | { | 
|---|
| 143 | const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam"); | 
|---|
| 144 | if (!geom) | 
|---|
| 145 | *fLog << warn << GetDescriptor() << ": No Camera Geometry available. Using mm-scale for histograms." << endl; | 
|---|
| 146 | else | 
|---|
| 147 | { | 
|---|
| 148 | fMm2Deg = geom->GetConvMm2Deg(); | 
|---|
| 149 | SetMmScale(kFALSE); | 
|---|
| 150 | } | 
|---|
| 151 |  | 
|---|
| 152 | ApplyBinning(*plist, "Width",  fWidth); | 
|---|
| 153 | ApplyBinning(*plist, "Length", fLength); | 
|---|
| 154 | ApplyBinning(*plist, "Dist",   fDistC); | 
|---|
| 155 | ApplyBinning(*plist, "Delta",  fDelta); | 
|---|
| 156 | ApplyBinning(*plist, "Size",   fSize); | 
|---|
| 157 |  | 
|---|
| 158 | const MBinning *bins = (MBinning*)plist->FindObject("BinningCamera"); | 
|---|
| 159 | if (!bins) | 
|---|
| 160 | { | 
|---|
| 161 | float r = geom ? geom->GetMaxRadius() : 600; | 
|---|
| 162 | r *= 0.9; | 
|---|
| 163 | if (!fUseMmScale) | 
|---|
| 164 | r *= fMm2Deg; | 
|---|
| 165 |  | 
|---|
| 166 | MBinning b; | 
|---|
| 167 | b.SetEdges(61, -r, r); | 
|---|
| 168 | SetBinning(fCenter, &b, &b); | 
|---|
| 169 | } | 
|---|
| 170 | else | 
|---|
| 171 | SetBinning(fCenter, bins, bins); | 
|---|
| 172 |  | 
|---|
| 173 |  | 
|---|
| 174 | return kTRUE; | 
|---|
| 175 | } | 
|---|
| 176 |  | 
|---|
| 177 | // -------------------------------------------------------------------------- | 
|---|
| 178 | // | 
|---|
| 179 | // Use this function to setup your own conversion factor between degrees | 
|---|
| 180 | // and millimeters. The conversion factor should be the one calculated in | 
|---|
| 181 | // MGeomCam. Use this function with Caution: You could create wrong values | 
|---|
| 182 | // by setting up your own scale factor. | 
|---|
| 183 | // | 
|---|
| 184 | void MHHillas::SetMm2Deg(Float_t mmdeg) | 
|---|
| 185 | { | 
|---|
| 186 | if (mmdeg<0) | 
|---|
| 187 | { | 
|---|
| 188 | *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl; | 
|---|
| 189 | return; | 
|---|
| 190 | } | 
|---|
| 191 |  | 
|---|
| 192 | if (fMm2Deg>=0) | 
|---|
| 193 | *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl; | 
|---|
| 194 |  | 
|---|
| 195 | fMm2Deg = mmdeg; | 
|---|
| 196 | } | 
|---|
| 197 |  | 
|---|
| 198 | // -------------------------------------------------------------------------- | 
|---|
| 199 | // | 
|---|
| 200 | // With this function you can convert the histogram ('on the fly') between | 
|---|
| 201 | // degrees and millimeters. | 
|---|
| 202 | // | 
|---|
| 203 | void MHHillas::SetMmScale(Bool_t mmscale) | 
|---|
| 204 | { | 
|---|
| 205 | if (fUseMmScale == mmscale) | 
|---|
| 206 | return; | 
|---|
| 207 |  | 
|---|
| 208 | if (fMm2Deg<0) | 
|---|
| 209 | { | 
|---|
| 210 | *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl; | 
|---|
| 211 | return; | 
|---|
| 212 | } | 
|---|
| 213 |  | 
|---|
| 214 | const Double_t scale = mmscale ? 1./fMm2Deg : fMm2Deg; | 
|---|
| 215 | MH::ScaleAxis(fLength, scale); | 
|---|
| 216 | MH::ScaleAxis(fWidth,  scale); | 
|---|
| 217 | MH::ScaleAxis(fDistC,  scale); | 
|---|
| 218 | MH::ScaleAxis(fCenter, scale, scale); | 
|---|
| 219 |  | 
|---|
| 220 | if (mmscale) | 
|---|
| 221 | { | 
|---|
| 222 | fLength->SetXTitle("Length [mm]"); | 
|---|
| 223 | fWidth->SetXTitle("Width [mm]"); | 
|---|
| 224 | fDistC->SetXTitle("Distance [mm]"); | 
|---|
| 225 | fCenter->SetXTitle("x [mm]"); | 
|---|
| 226 | fCenter->SetYTitle("y [mm]"); | 
|---|
| 227 | } | 
|---|
| 228 | else | 
|---|
| 229 | { | 
|---|
| 230 | fLength->SetXTitle("Length [\\circ]"); | 
|---|
| 231 | fWidth->SetXTitle("Width [\\circ]"); | 
|---|
| 232 | fDistC->SetXTitle("Distance [\\circ]"); | 
|---|
| 233 | fCenter->SetXTitle("x [\\circ]"); | 
|---|
| 234 | fCenter->SetYTitle("y [\\circ]"); | 
|---|
| 235 | } | 
|---|
| 236 |  | 
|---|
| 237 | fUseMmScale = mmscale; | 
|---|
| 238 | } | 
|---|
| 239 |  | 
|---|
| 240 | // -------------------------------------------------------------------------- | 
|---|
| 241 | // | 
|---|
| 242 | // Fill the histograms with data from a MHillas-Container. | 
|---|
| 243 | // Be careful: Only call this with an object of type MHillas | 
|---|
| 244 | // | 
|---|
| 245 | Bool_t MHHillas::Fill(const MParContainer *par, const Stat_t w) | 
|---|
| 246 | { | 
|---|
| 247 | if (!par) | 
|---|
| 248 | { | 
|---|
| 249 | *fLog << err << "MHHillas::Fill: Pointer (!=NULL) expected." << endl; | 
|---|
| 250 | return kFALSE; | 
|---|
| 251 | } | 
|---|
| 252 |  | 
|---|
| 253 | const MHillas &h = *(MHillas*)par; | 
|---|
| 254 |  | 
|---|
| 255 | const Double_t d = sqrt(h.GetMeanX()*h.GetMeanX() + h.GetMeanY()*h.GetMeanY()); | 
|---|
| 256 | const Double_t scale = fUseMmScale ? 1 : fMm2Deg; | 
|---|
| 257 |  | 
|---|
| 258 | fLength->Fill(scale*h.GetLength(), w); | 
|---|
| 259 | fWidth ->Fill(scale*h.GetWidth(), w); | 
|---|
| 260 | fDistC ->Fill(scale*d, w); | 
|---|
| 261 | fCenter->Fill(scale*h.GetMeanX(), scale*h.GetMeanY(), w); | 
|---|
| 262 | fDelta ->Fill(kRad2Deg*h.GetDelta(), w); | 
|---|
| 263 | fSize  ->Fill(h.GetSize(), w); | 
|---|
| 264 |  | 
|---|
| 265 | return kTRUE; | 
|---|
| 266 | } | 
|---|
| 267 |  | 
|---|
| 268 | // -------------------------------------------------------------------------- | 
|---|
| 269 | // | 
|---|
| 270 | // Setup a inversed deep blue sea palette for the fCenter histogram. | 
|---|
| 271 | // | 
|---|
| 272 | void MHHillas::SetColors() const | 
|---|
| 273 | { | 
|---|
| 274 | gStyle->SetPalette(51, NULL); | 
|---|
| 275 | Int_t c[50]; | 
|---|
| 276 | for (int i=0; i<50; i++) | 
|---|
| 277 | c[49-i] = gStyle->GetColorPalette(i); | 
|---|
| 278 | gStyle->SetPalette(50, c); | 
|---|
| 279 | } | 
|---|
| 280 |  | 
|---|
| 281 | // -------------------------------------------------------------------------- | 
|---|
| 282 | // | 
|---|
| 283 | // Creates a new canvas and draws the four histograms into it. | 
|---|
| 284 | // Be careful: The histograms belongs to this object and won't get deleted | 
|---|
| 285 | // together with the canvas. | 
|---|
| 286 | // | 
|---|
| 287 | void MHHillas::Draw(Option_t *) | 
|---|
| 288 | { | 
|---|
| 289 | TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this); | 
|---|
| 290 | pad->SetBorderMode(0); | 
|---|
| 291 |  | 
|---|
| 292 | AppendPad(""); | 
|---|
| 293 |  | 
|---|
| 294 | pad->Divide(2,3); | 
|---|
| 295 |  | 
|---|
| 296 | pad->cd(1); | 
|---|
| 297 | gPad->SetBorderMode(0); | 
|---|
| 298 | MH::DrawSame(*fWidth, *fLength, "Width'n'Length"); | 
|---|
| 299 |  | 
|---|
| 300 | pad->cd(2); | 
|---|
| 301 | gPad->SetBorderMode(0); | 
|---|
| 302 | fDistC->Draw(); | 
|---|
| 303 |  | 
|---|
| 304 | pad->cd(3); | 
|---|
| 305 | gPad->SetBorderMode(0); | 
|---|
| 306 | gPad->SetLogx(); | 
|---|
| 307 | fSize->Draw(); | 
|---|
| 308 |  | 
|---|
| 309 | pad->cd(4); | 
|---|
| 310 | gPad->SetBorderMode(0); | 
|---|
| 311 | gPad->SetPad(0.51, 0.01, 0.99, 0.65); | 
|---|
| 312 | SetColors(); | 
|---|
| 313 | fCenter->Draw("colz"); | 
|---|
| 314 |  | 
|---|
| 315 | pad->cd(5); | 
|---|
| 316 | gPad->SetBorderMode(0); | 
|---|
| 317 | fDelta->Draw(); | 
|---|
| 318 |  | 
|---|
| 319 | pad->cd(6); | 
|---|
| 320 | delete gPad; | 
|---|
| 321 |  | 
|---|
| 322 | pad->Modified(); | 
|---|
| 323 | pad->Update(); | 
|---|
| 324 | } | 
|---|
| 325 |  | 
|---|
| 326 | TH1 *MHHillas::GetHistByName(const TString name) | 
|---|
| 327 | { | 
|---|
| 328 | if (name.Contains("Width", TString::kIgnoreCase)) | 
|---|
| 329 | return fWidth; | 
|---|
| 330 | if (name.Contains("Length", TString::kIgnoreCase)) | 
|---|
| 331 | return fLength; | 
|---|
| 332 | if (name.Contains("Size", TString::kIgnoreCase)) | 
|---|
| 333 | return fSize; | 
|---|
| 334 | if (name.Contains("Delta", TString::kIgnoreCase)) | 
|---|
| 335 | return fDelta; | 
|---|
| 336 | if (name.Contains("DistC", TString::kIgnoreCase)) | 
|---|
| 337 | return fDistC; | 
|---|
| 338 | if (name.Contains("Center", TString::kIgnoreCase)) | 
|---|
| 339 | return fCenter; | 
|---|
| 340 |  | 
|---|
| 341 | return NULL; | 
|---|
| 342 | } | 
|---|
| 343 |  | 
|---|
| 344 | void MHHillas::Paint(Option_t *opt) | 
|---|
| 345 | { | 
|---|
| 346 | SetColors(); | 
|---|
| 347 | MH::Paint(); | 
|---|
| 348 | } | 
|---|