| 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, 05/2002 <mailto:tbretz@astro.uni-wuerzburg.de> | 
|---|
| 19 | !   Author(s): Harald Kornmayer, 1/2001 | 
|---|
| 20 | !   Author(s): Markus Gaug, 03/2004 <mailto:markus@ifae.es> | 
|---|
| 21 | ! | 
|---|
| 22 | !   Copyright: MAGIC Software Development, 2000-2004 | 
|---|
| 23 | ! | 
|---|
| 24 | ! | 
|---|
| 25 | \* ======================================================================== */ | 
|---|
| 26 |  | 
|---|
| 27 | ///////////////////////////////////////////////////////////////////////////// | 
|---|
| 28 | // | 
|---|
| 29 | // MHCamera | 
|---|
| 30 | // | 
|---|
| 31 | // Camera Display, based on a TH1D. Pleas be carefull using the | 
|---|
| 32 | // underlaying TH1D. | 
|---|
| 33 | // | 
|---|
| 34 | // To change the scale to a logarithmic scale SetLogy() of the Pad. | 
|---|
| 35 | // | 
|---|
| 36 | // You can correct for the abberation. Assuming that the distance | 
|---|
| 37 | // between the mean position of the light distribution and the position | 
|---|
| 38 | // of a perfect reflection on a perfect mirror in the distance r on | 
|---|
| 39 | // the camera plane is dr it is  d = a*dr  while a is the abberation | 
|---|
| 40 | // constant (for the MAGIC mirror it is roughly 0.0713). You can | 
|---|
| 41 | // set this constant by calling SetAbberation(a) which results in a | 
|---|
| 42 | // 'corrected' display (all outer pixels are shifted towards the center | 
|---|
| 43 | // of the camera to correct for this abberation) | 
|---|
| 44 | // | 
|---|
| 45 | // Be carefull: Entries in this context means Entries/bin or Events | 
|---|
| 46 | // | 
|---|
| 47 | // FIXME? Maybe MHCamera can take the fLog object from MGeomCam? | 
|---|
| 48 | // | 
|---|
| 49 | //////////////////////////////////////////////////////////////////////////// | 
|---|
| 50 | #include "MHCamera.h" | 
|---|
| 51 |  | 
|---|
| 52 | #include <fstream> | 
|---|
| 53 | #include <iostream> | 
|---|
| 54 |  | 
|---|
| 55 | #include <TBox.h> | 
|---|
| 56 | #include <TArrow.h> | 
|---|
| 57 | #include <TLatex.h> | 
|---|
| 58 | #include <TStyle.h> | 
|---|
| 59 | #include <TCanvas.h> | 
|---|
| 60 | #include <TArrayF.h> | 
|---|
| 61 | #include <TRandom.h> | 
|---|
| 62 | #include <TPaveText.h> | 
|---|
| 63 | #include <TPaveStats.h> | 
|---|
| 64 | #include <TClonesArray.h> | 
|---|
| 65 | #include <THistPainter.h> | 
|---|
| 66 | #include <THLimitsFinder.h> | 
|---|
| 67 | #include <TProfile.h> | 
|---|
| 68 | #include <TH1.h> | 
|---|
| 69 | #include <TF1.h> | 
|---|
| 70 | #include <TCanvas.h> | 
|---|
| 71 | #include <TLegend.h> | 
|---|
| 72 |  | 
|---|
| 73 | #include "MLog.h" | 
|---|
| 74 | #include "MLogManip.h" | 
|---|
| 75 |  | 
|---|
| 76 | #include "MH.h" | 
|---|
| 77 | #include "MBinning.h" | 
|---|
| 78 | #include "MHexagon.h" | 
|---|
| 79 |  | 
|---|
| 80 | #include "MGeomPix.h" | 
|---|
| 81 | #include "MGeomCam.h" | 
|---|
| 82 |  | 
|---|
| 83 | #include "MCamEvent.h" | 
|---|
| 84 |  | 
|---|
| 85 | #include "MArrayD.h" | 
|---|
| 86 | #include "MMath.h"    // MMath::GaussProb | 
|---|
| 87 |  | 
|---|
| 88 | #define kItemsLegend 48 // see SetPalette(1,0) | 
|---|
| 89 |  | 
|---|
| 90 | ClassImp(MHCamera); | 
|---|
| 91 |  | 
|---|
| 92 | using namespace std; | 
|---|
| 93 |  | 
|---|
| 94 | void MHCamera::Init() | 
|---|
| 95 | { | 
|---|
| 96 | Sumw2(); | 
|---|
| 97 |  | 
|---|
| 98 | UseCurrentStyle(); | 
|---|
| 99 |  | 
|---|
| 100 | SetDirectory(NULL); | 
|---|
| 101 |  | 
|---|
| 102 | SetLineColor(kGreen); | 
|---|
| 103 | SetMarkerStyle(kFullDotMedium); | 
|---|
| 104 | SetXTitle("Pixel Index"); | 
|---|
| 105 |  | 
|---|
| 106 | fNotify  = new TList; | 
|---|
| 107 | fNotify->SetBit(kMustCleanup); | 
|---|
| 108 | gROOT->GetListOfCleanups()->Add(fNotify); | 
|---|
| 109 |  | 
|---|
| 110 | TVirtualPad *save = gPad; | 
|---|
| 111 | gPad = 0; | 
|---|
| 112 | /* | 
|---|
| 113 | #if ROOT_VERSION_CODE < ROOT_VERSION(3,01,06) | 
|---|
| 114 | SetPalette(1, 0); | 
|---|
| 115 | #endif | 
|---|
| 116 | */ | 
|---|
| 117 | /* | 
|---|
| 118 | #if ROOT_VERSION_CODE < ROOT_VERSION(4,04,00) | 
|---|
| 119 | SetPrettyPalette(); | 
|---|
| 120 | #elese | 
|---|
| 121 | // WORAROUND - FIXME: Calling it many times becomes slower and slower | 
|---|
| 122 | SetInvDeepBlueSeaPalette(); | 
|---|
| 123 | #endif | 
|---|
| 124 | */ | 
|---|
| 125 | gPad = save; | 
|---|
| 126 | } | 
|---|
| 127 |  | 
|---|
| 128 | // ------------------------------------------------------------------------ | 
|---|
| 129 | // | 
|---|
| 130 | //  Default Constructor. To be used by the root system ONLY. | 
|---|
| 131 | // | 
|---|
| 132 | MHCamera::MHCamera() : TH1D(), fGeomCam(NULL), fAbberation(0) | 
|---|
| 133 | { | 
|---|
| 134 | Init(); | 
|---|
| 135 | } | 
|---|
| 136 |  | 
|---|
| 137 | // ------------------------------------------------------------------------ | 
|---|
| 138 | // | 
|---|
| 139 | //  Constructor. Makes a clone of MGeomCam. Removed the TH1D from the | 
|---|
| 140 | // current directory. Calls Sumw2(). Set the histogram line color | 
|---|
| 141 | // (for error bars) to Green and the marker style to kFullDotMedium. | 
|---|
| 142 | // | 
|---|
| 143 | MHCamera::MHCamera(const MGeomCam &geom, const char *name, const char *title) | 
|---|
| 144 | : fGeomCam(NULL), fAbberation(0) | 
|---|
| 145 | { | 
|---|
| 146 | //fGeomCam = (MGeomCam*)geom.Clone(); | 
|---|
| 147 | SetGeometry(geom, name, title); | 
|---|
| 148 | Init(); | 
|---|
| 149 |  | 
|---|
| 150 | // | 
|---|
| 151 | // root 3.02 | 
|---|
| 152 | //  * base/inc/TObject.h: | 
|---|
| 153 | //    register BIT(8) as kNoContextMenu. If an object has this bit set it will | 
|---|
| 154 | //    not get an automatic context menu when clicked with the right mouse button. | 
|---|
| 155 | } | 
|---|
| 156 |  | 
|---|
| 157 | void MHCamera::SetGeometry(const MGeomCam &geom, const char *name, const char *title) | 
|---|
| 158 | { | 
|---|
| 159 | SetNameTitle(name, title); | 
|---|
| 160 |  | 
|---|
| 161 | TAxis &x = *GetXaxis(); | 
|---|
| 162 |  | 
|---|
| 163 | SetBins(geom.GetNumPixels(), 0, 1); | 
|---|
| 164 | x.Set(geom.GetNumPixels(), -0.5, geom.GetNumPixels()-0.5); | 
|---|
| 165 |  | 
|---|
| 166 | //SetBins(geom.GetNumPixels(), -0.5, geom.GetNumPixels()-0.5); | 
|---|
| 167 | //Rebuild(); | 
|---|
| 168 |  | 
|---|
| 169 | if (fGeomCam) | 
|---|
| 170 | delete fGeomCam; | 
|---|
| 171 | fGeomCam = (MGeomCam*)geom.Clone(); | 
|---|
| 172 |  | 
|---|
| 173 | fUsed.Set(geom.GetNumPixels()); | 
|---|
| 174 | for (Int_t i=0; i<fNcells-2; i++) | 
|---|
| 175 | ResetUsed(i); | 
|---|
| 176 |  | 
|---|
| 177 | fBinEntries.Set(geom.GetNumPixels()+2); | 
|---|
| 178 | fBinEntries.Reset(); | 
|---|
| 179 | } | 
|---|
| 180 |  | 
|---|
| 181 | // ------------------------------------------------------------------------ | 
|---|
| 182 | // | 
|---|
| 183 | // Destructor. Deletes the cloned fGeomCam and the notification list. | 
|---|
| 184 | // | 
|---|
| 185 | MHCamera::~MHCamera() | 
|---|
| 186 | { | 
|---|
| 187 | if (fGeomCam) | 
|---|
| 188 | delete fGeomCam; | 
|---|
| 189 | if (fNotify) | 
|---|
| 190 | delete fNotify; | 
|---|
| 191 | } | 
|---|
| 192 |  | 
|---|
| 193 | // ------------------------------------------------------------------------ | 
|---|
| 194 | // | 
|---|
| 195 | // Return kTRUE for sector<0. Otherwise return kTRUE only if the specified | 
|---|
| 196 | // sector idx matches the sector of the pixel with index idx. | 
|---|
| 197 | // | 
|---|
| 198 | Bool_t MHCamera::MatchSector(Int_t idx, const TArrayI §or, const TArrayI &aidx) const | 
|---|
| 199 | { | 
|---|
| 200 | const MGeomPix &pix = (*fGeomCam)[idx]; | 
|---|
| 201 | return FindVal(sector, pix.GetSector()) && FindVal(aidx, pix.GetAidx()); | 
|---|
| 202 | } | 
|---|
| 203 |  | 
|---|
| 204 | // ------------------------------------------------------------------------ | 
|---|
| 205 | // | 
|---|
| 206 | // Taken from TH1D::Fill(). Uses the argument directly as bin index. | 
|---|
| 207 | // Doesn't increment the number of entries. | 
|---|
| 208 | // | 
|---|
| 209 | //   -*-*-*-*-*-*-*-*Increment bin with abscissa X by 1*-*-*-*-*-*-*-*-*-*-* | 
|---|
| 210 | //                   ================================== | 
|---|
| 211 | // | 
|---|
| 212 | //    if x is less than the low-edge of the first bin, the Underflow bin is incremented | 
|---|
| 213 | //    if x is greater than the upper edge of last bin, the Overflow bin is incremented | 
|---|
| 214 | // | 
|---|
| 215 | //    If the storage of the sum of squares of weights has been triggered, | 
|---|
| 216 | //    via the function Sumw2, then the sum of the squares of weights is incremented | 
|---|
| 217 | //    by 1 in the bin corresponding to x. | 
|---|
| 218 | // | 
|---|
| 219 | //   -*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* | 
|---|
| 220 | Int_t MHCamera::Fill(Axis_t x) | 
|---|
| 221 | { | 
|---|
| 222 | #if ROOT_VERSION_CODE > ROOT_VERSION(3,05,00) | 
|---|
| 223 | if (fBuffer) return BufferFill(x,1); | 
|---|
| 224 | #endif | 
|---|
| 225 | const Int_t bin = (Int_t)x+1; | 
|---|
| 226 | AddBinContent(bin); | 
|---|
| 227 | fBinEntries[bin]++; | 
|---|
| 228 | if (fSumw2.fN) | 
|---|
| 229 | fSumw2.fArray[bin]++; | 
|---|
| 230 |  | 
|---|
| 231 | if (bin<=0 || bin>fNcells-2) | 
|---|
| 232 | return -1; | 
|---|
| 233 |  | 
|---|
| 234 | fTsumw++; | 
|---|
| 235 | fTsumw2++; | 
|---|
| 236 | fTsumwx  += x; | 
|---|
| 237 | fTsumwx2 += x*x; | 
|---|
| 238 | return bin; | 
|---|
| 239 | } | 
|---|
| 240 |  | 
|---|
| 241 | // ------------------------------------------------------------------------ | 
|---|
| 242 | // | 
|---|
| 243 | // Taken from TH1D::Fill(). Uses the argument directly as bin index. | 
|---|
| 244 | // Doesn't increment the number of entries. | 
|---|
| 245 | // | 
|---|
| 246 | //   -*-*-*-*-*-*Increment bin with abscissa X with a weight w*-*-*-*-*-*-*-* | 
|---|
| 247 | //               ============================================= | 
|---|
| 248 | // | 
|---|
| 249 | //    if x is less than the low-edge of the first bin, the Underflow bin is incremented | 
|---|
| 250 | //    if x is greater than the upper edge of last bin, the Overflow bin is incremented | 
|---|
| 251 | // | 
|---|
| 252 | //    If the storage of the sum of squares of weights has been triggered, | 
|---|
| 253 | //    via the function Sumw2, then the sum of the squares of weights is incremented | 
|---|
| 254 | //    by w^2 in the bin corresponding to x. | 
|---|
| 255 | // | 
|---|
| 256 | //   -*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* | 
|---|
| 257 | Int_t MHCamera::Fill(Axis_t x, Stat_t w) | 
|---|
| 258 | { | 
|---|
| 259 | #if ROOT_VERSION_CODE > ROOT_VERSION(3,05,00) | 
|---|
| 260 | if (fBuffer) return BufferFill(x,w); | 
|---|
| 261 | #endif | 
|---|
| 262 | const Int_t bin = (Int_t)x+1; | 
|---|
| 263 | AddBinContent(bin, w); | 
|---|
| 264 | fBinEntries[bin]++; | 
|---|
| 265 | if (fSumw2.fN) | 
|---|
| 266 | fSumw2.fArray[bin] += w*w; | 
|---|
| 267 |  | 
|---|
| 268 | if (bin<=0 || bin>fNcells-2) | 
|---|
| 269 | return -1; | 
|---|
| 270 |  | 
|---|
| 271 | const Stat_t z = (w > 0 ? w : -w); | 
|---|
| 272 | fTsumw   += z; | 
|---|
| 273 | fTsumw2  += z*z; | 
|---|
| 274 | fTsumwx  += z*x; | 
|---|
| 275 | fTsumwx2 += z*x*x; | 
|---|
| 276 | return bin; | 
|---|
| 277 | } | 
|---|
| 278 |  | 
|---|
| 279 | // ------------------------------------------------------------------------ | 
|---|
| 280 | // | 
|---|
| 281 | // Use x and y in millimeters | 
|---|
| 282 | // | 
|---|
| 283 | Int_t MHCamera::Fill(Axis_t x, Axis_t y, Stat_t w) | 
|---|
| 284 | { | 
|---|
| 285 | if (fNcells<=1 || IsFreezed()) | 
|---|
| 286 | return -1; | 
|---|
| 287 |  | 
|---|
| 288 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 289 | { | 
|---|
| 290 | MHexagon hex((*fGeomCam)[idx]); | 
|---|
| 291 | if (hex.DistanceToPrimitive(x, y)>0) | 
|---|
| 292 | continue; | 
|---|
| 293 |  | 
|---|
| 294 | SetUsed(idx); | 
|---|
| 295 | return Fill(idx, w); | 
|---|
| 296 | } | 
|---|
| 297 | return -1; | 
|---|
| 298 | } | 
|---|
| 299 |  | 
|---|
| 300 | // ------------------------------------------------------------------------ | 
|---|
| 301 | // | 
|---|
| 302 | // Call this if you want to change the display status (displayed or not) | 
|---|
| 303 | // for all pixels. val==0 means that the pixel is not displayed. | 
|---|
| 304 | // | 
|---|
| 305 | void MHCamera::SetUsed(const TArrayC &arr) | 
|---|
| 306 | { | 
|---|
| 307 | if (fNcells-2 != arr.GetSize()) | 
|---|
| 308 | { | 
|---|
| 309 | gLog << warn << "WARNING - MHCamera::SetUsed: array size mismatch... ignored." << endl; | 
|---|
| 310 | return; | 
|---|
| 311 | } | 
|---|
| 312 |  | 
|---|
| 313 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 314 | arr[idx] ? SetUsed(idx) : ResetUsed(idx); | 
|---|
| 315 | } | 
|---|
| 316 |  | 
|---|
| 317 | // ------------------------------------------------------------------------ | 
|---|
| 318 | // | 
|---|
| 319 | // Return the median value on the y-axis (profile option is correctly taken | 
|---|
| 320 | // into account) | 
|---|
| 321 | // | 
|---|
| 322 | Stat_t MHCamera::GetMedian() const | 
|---|
| 323 | { | 
|---|
| 324 | // Just for speed reasons | 
|---|
| 325 | if (!TestBit(kProfile)) | 
|---|
| 326 | return TMath::Median(GetSize()-2, GetArray()+1); | 
|---|
| 327 |  | 
|---|
| 328 | // Copy profiled data into new array (FIXME: Should we take errors into account?) | 
|---|
| 329 | TArrayD arr(fNcells-2); | 
|---|
| 330 | for (int i=1; i<fNcells-1; i++) | 
|---|
| 331 | arr[i-1] = GetBinContent(i); | 
|---|
| 332 |  | 
|---|
| 333 | // return Median of the profile data | 
|---|
| 334 | return TMath::Median(arr.GetSize(), arr.GetArray()); | 
|---|
| 335 | } | 
|---|
| 336 |  | 
|---|
| 337 | // ------------------------------------------------------------------------ | 
|---|
| 338 | // | 
|---|
| 339 | // Return the median value (divided by MMath::GausProb(1.0)) of the | 
|---|
| 340 | // distribution of abs(y[i]-Median). This is my Median equivalent of the RMS | 
|---|
| 341 | // | 
|---|
| 342 | Stat_t MHCamera::GetMedianDev() const | 
|---|
| 343 | { | 
|---|
| 344 | // Just for speed reasons | 
|---|
| 345 | if (!TestBit(kProfile)) | 
|---|
| 346 | return MMath::MedianDev(GetSize()-2, GetArray()+1); | 
|---|
| 347 |  | 
|---|
| 348 | // Copy profiled data into new array (FIXME: Should we take errors into account?) | 
|---|
| 349 | TArrayD arr(fNcells-2); | 
|---|
| 350 | for (int i=1; i<fNcells-1; i++) | 
|---|
| 351 | arr[i-1] = GetBinContent(i); | 
|---|
| 352 |  | 
|---|
| 353 | // return MedianDev of the profile data | 
|---|
| 354 | return MMath::MedianDev(arr.GetSize(), arr.GetArray()); | 
|---|
| 355 | } | 
|---|
| 356 |  | 
|---|
| 357 | // ------------------------------------------------------------------------ | 
|---|
| 358 | // | 
|---|
| 359 | // Return the mean value of all entries which are used if all=kFALSE and | 
|---|
| 360 | // of all entries if all=kTRUE if sector<0. If sector>=0 only | 
|---|
| 361 | // entries with match the given sector are taken into account. | 
|---|
| 362 | // | 
|---|
| 363 | Stat_t MHCamera::GetMeanSectors(const TArrayI §or, const TArrayI &aidx, Bool_t all) const | 
|---|
| 364 | { | 
|---|
| 365 | if (fNcells<=1) | 
|---|
| 366 | return 0; | 
|---|
| 367 |  | 
|---|
| 368 | Int_t n=0; | 
|---|
| 369 |  | 
|---|
| 370 | Stat_t mean = 0; | 
|---|
| 371 | for (int i=0; i<fNcells-2; i++) | 
|---|
| 372 | { | 
|---|
| 373 | if ((all || IsUsed(i)) && MatchSector(i, sector, aidx)) | 
|---|
| 374 | { | 
|---|
| 375 | if (TestBit(kProfile) && fBinEntries[i+1]==0) | 
|---|
| 376 | continue; | 
|---|
| 377 |  | 
|---|
| 378 | mean += TestBit(kProfile) ? fArray[i+1]/fBinEntries[i+1] : fArray[i+1]; | 
|---|
| 379 | n++; | 
|---|
| 380 | } | 
|---|
| 381 | } | 
|---|
| 382 |  | 
|---|
| 383 | return n==0 ? 0 : mean/n; | 
|---|
| 384 | } | 
|---|
| 385 |  | 
|---|
| 386 | // ------------------------------------------------------------------------ | 
|---|
| 387 | // | 
|---|
| 388 | // Return the sqrt variance of all entries which are used if all=kFALSE and | 
|---|
| 389 | // of all entries if all=kTRUE if sector<0. If sector>=0 only | 
|---|
| 390 | // entries with match the given sector are taken into account. | 
|---|
| 391 | // | 
|---|
| 392 | Stat_t MHCamera::GetRmsSectors(const TArrayI §or, const TArrayI &aidx, Bool_t all) const | 
|---|
| 393 | { | 
|---|
| 394 | if (fNcells<=1) | 
|---|
| 395 | return -1; | 
|---|
| 396 |  | 
|---|
| 397 | Int_t n=0; | 
|---|
| 398 |  | 
|---|
| 399 | Stat_t sum = 0; | 
|---|
| 400 | Stat_t sq  = 0; | 
|---|
| 401 | for (int i=0; i<fNcells-2; i++) | 
|---|
| 402 | { | 
|---|
| 403 | if ((all || IsUsed(i)) && MatchSector(i, sector, aidx)) | 
|---|
| 404 | { | 
|---|
| 405 | if (TestBit(kProfile) && fBinEntries[i+1]==0) | 
|---|
| 406 | continue; | 
|---|
| 407 |  | 
|---|
| 408 | const Double_t val = TestBit(kProfile) ? fArray[i+1]/fBinEntries[i+1] : fArray[i+1]; | 
|---|
| 409 |  | 
|---|
| 410 | sum += val; | 
|---|
| 411 | sq  += val*val; | 
|---|
| 412 | n++; | 
|---|
| 413 | } | 
|---|
| 414 | } | 
|---|
| 415 |  | 
|---|
| 416 | if (n==0) | 
|---|
| 417 | return 0; | 
|---|
| 418 |  | 
|---|
| 419 | sum /= n; | 
|---|
| 420 | sq  /= n; | 
|---|
| 421 |  | 
|---|
| 422 | return TMath::Sqrt(sq-sum*sum); | 
|---|
| 423 | } | 
|---|
| 424 |  | 
|---|
| 425 | // ------------------------------------------------------------------------ | 
|---|
| 426 | // | 
|---|
| 427 | // Return the minimum contents of all pixels (if all is set, otherwise | 
|---|
| 428 | // only of all 'used' pixels), fMinimum if fMinimum set. If sector>=0 | 
|---|
| 429 | // only pixels with matching sector number are taken into account. | 
|---|
| 430 | // | 
|---|
| 431 | Double_t MHCamera::GetMinimumSectors(const TArrayI §or, const TArrayI &aidx, Bool_t all) const | 
|---|
| 432 | { | 
|---|
| 433 | if (fMinimum != -1111) | 
|---|
| 434 | return fMinimum; | 
|---|
| 435 |  | 
|---|
| 436 | if (fNcells<=1) | 
|---|
| 437 | return 0; | 
|---|
| 438 |  | 
|---|
| 439 | Double_t minimum=FLT_MAX; | 
|---|
| 440 |  | 
|---|
| 441 | for (Int_t i=0; i<fNcells-2; i++) | 
|---|
| 442 | { | 
|---|
| 443 | if (TestBit(kProfile) && fBinEntries[i+1]==0) | 
|---|
| 444 | continue; | 
|---|
| 445 |  | 
|---|
| 446 | const Double_t val = TestBit(kProfile) ? fArray[i+1]/fBinEntries[i+1] : fArray[i+1]; | 
|---|
| 447 | if (MatchSector(i, sector, aidx) && (all || IsUsed(i)) && val<minimum) | 
|---|
| 448 | minimum = val; | 
|---|
| 449 | } | 
|---|
| 450 |  | 
|---|
| 451 | return minimum; | 
|---|
| 452 | } | 
|---|
| 453 |  | 
|---|
| 454 | // ------------------------------------------------------------------------ | 
|---|
| 455 | // | 
|---|
| 456 | // Return the maximum contents of all pixels (if all is set, otherwise | 
|---|
| 457 | // only of all 'used' pixels), fMaximum if fMaximum set. If sector>=0 | 
|---|
| 458 | // only pixels with matching sector number are taken into account. | 
|---|
| 459 | // | 
|---|
| 460 | Double_t MHCamera::GetMaximumSectors(const TArrayI §or, const TArrayI &aidx, Bool_t all) const | 
|---|
| 461 | { | 
|---|
| 462 | if (fMaximum!=-1111) | 
|---|
| 463 | return fMaximum; | 
|---|
| 464 |  | 
|---|
| 465 | if (fNcells<=1) | 
|---|
| 466 | return 1; | 
|---|
| 467 |  | 
|---|
| 468 | Double_t maximum=-FLT_MAX; | 
|---|
| 469 | for (Int_t i=0; i<fNcells-2; i++) | 
|---|
| 470 | { | 
|---|
| 471 | if (TestBit(kProfile) && fBinEntries[i+1]==0) | 
|---|
| 472 | continue; | 
|---|
| 473 |  | 
|---|
| 474 | const Double_t val = TestBit(kProfile) ? fArray[i+1]/fBinEntries[i+1] : fArray[i+1]; | 
|---|
| 475 | if (MatchSector(i, sector, aidx) && (all || IsUsed(i)) && val>maximum) | 
|---|
| 476 | maximum = val; | 
|---|
| 477 | } | 
|---|
| 478 |  | 
|---|
| 479 | return maximum; | 
|---|
| 480 | } | 
|---|
| 481 |  | 
|---|
| 482 | // ------------------------------------------------------------------------ | 
|---|
| 483 | // | 
|---|
| 484 | // Call this function to draw the camera layout into your canvas. | 
|---|
| 485 | // Setup a drawing canvas. Add this object and all child objects | 
|---|
| 486 | // (hexagons, etc) to the current pad. If no pad exists a new one is | 
|---|
| 487 | // created. (To access the 'real' pad containing the camera you have | 
|---|
| 488 | // to do a cd(1) in the current layer. | 
|---|
| 489 | // | 
|---|
| 490 | // To draw a camera into its own pad do something like: | 
|---|
| 491 | // | 
|---|
| 492 | //   MGeomCamMagic m; | 
|---|
| 493 | //   MHCamera *d=new MHCamera(m); | 
|---|
| 494 | // | 
|---|
| 495 | //   TCanvas *c = new TCanvas; | 
|---|
| 496 | //   c->Divide(2,1); | 
|---|
| 497 | //   c->cd(1); | 
|---|
| 498 | // | 
|---|
| 499 | //   d->FillRandom(); | 
|---|
| 500 | //   d->Draw(); | 
|---|
| 501 | //   d->SetBit(kCanDelete); | 
|---|
| 502 | // | 
|---|
| 503 | // There are several drawing options: | 
|---|
| 504 | //   'hist'        Draw as a standard TH1 histogram (value vs. pixel index) | 
|---|
| 505 | //   'box'         Draw hexagons which size is in respect to its contents | 
|---|
| 506 | //   'nocol'       Leave the 'boxed' hexagons empty | 
|---|
| 507 | //   'pixelindex'  Display the pixel index in each pixel | 
|---|
| 508 | //   'sectorindex' Display the sector index in each pixel | 
|---|
| 509 | //   'content'     Display the relative content aligned to GetMaximum() and | 
|---|
| 510 | //                 GeMinimum() ((val-min)/(max-min)) | 
|---|
| 511 | //   'proj'        Display the y-projection of the histogram | 
|---|
| 512 | //   'pal0'        Use Pretty palette | 
|---|
| 513 | //   'pal1'        Use Deep Blue Sea palette | 
|---|
| 514 | //   'pal2'        Use Inverse Depp Blue Sea palette | 
|---|
| 515 | //   'same'        Draw trandparent pixels on top of an existing pad. This | 
|---|
| 516 | //                 makes it possible to draw the camera image on top of an | 
|---|
| 517 | //                 existing TH2, but also allows for distorted camera images | 
|---|
| 518 | // | 
|---|
| 519 | void MHCamera::Draw(Option_t *option) | 
|---|
| 520 | { | 
|---|
| 521 | const Bool_t hassame = TString(option).Contains("same", TString::kIgnoreCase) && gPad; | 
|---|
| 522 |  | 
|---|
| 523 | // root 3.02: | 
|---|
| 524 | // gPad->SetFixedAspectRatio() | 
|---|
| 525 | const Color_t col = gPad ? gPad->GetFillColor() : 16; | 
|---|
| 526 | TVirtualPad  *pad = gPad ? gPad : MH::MakeDefCanvas("CamDisplay", "Mars Camera Display", 656, 600); | 
|---|
| 527 |  | 
|---|
| 528 | if (!hassame) | 
|---|
| 529 | { | 
|---|
| 530 | pad->SetBorderMode(0); | 
|---|
| 531 | pad->SetFillColor(col); | 
|---|
| 532 |  | 
|---|
| 533 | // | 
|---|
| 534 | // Create an own pad for the MHCamera-Object which can be | 
|---|
| 535 | // resized in paint to keep the correct aspect ratio | 
|---|
| 536 | // | 
|---|
| 537 | // The margin != 0 is a workaround for a problem in root 4.02/00 | 
|---|
| 538 | pad->Divide(1, 1, 1e-10, 1e-10, col); | 
|---|
| 539 | pad->cd(1); | 
|---|
| 540 | gPad->SetBorderMode(0); | 
|---|
| 541 | } | 
|---|
| 542 |  | 
|---|
| 543 | AppendPad(option); | 
|---|
| 544 | //fGeomCam->AppendPad(); | 
|---|
| 545 |  | 
|---|
| 546 | // | 
|---|
| 547 | // Do not change gPad. The user should not see, that Draw | 
|---|
| 548 | // changes gPad... | 
|---|
| 549 | // | 
|---|
| 550 | if (!hassame) | 
|---|
| 551 | pad->cd(); | 
|---|
| 552 | } | 
|---|
| 553 |  | 
|---|
| 554 | // ------------------------------------------------------------------------ | 
|---|
| 555 | // | 
|---|
| 556 | // This is TObject::DrawClone but completely ignores | 
|---|
| 557 | // gROOT->GetSelectedPad(). tbretz had trouble with this in the past. | 
|---|
| 558 | // If this makes trouble please write a bug report. | 
|---|
| 559 | // | 
|---|
| 560 | TObject *MHCamera::DrawClone(Option_t *option) const | 
|---|
| 561 | { | 
|---|
| 562 | // Draw a clone of this object in the current pad | 
|---|
| 563 |  | 
|---|
| 564 | //TVirtualPad *pad = gROOT->GetSelectedPad(); | 
|---|
| 565 | TVirtualPad *padsav = gPad; | 
|---|
| 566 | //if (pad) pad->cd(); | 
|---|
| 567 |  | 
|---|
| 568 | TObject *newobj = Clone(); | 
|---|
| 569 |  | 
|---|
| 570 | if (!newobj) | 
|---|
| 571 | return 0; | 
|---|
| 572 |  | 
|---|
| 573 | /* | 
|---|
| 574 | if (pad) { | 
|---|
| 575 | if (strlen(option)) pad->GetListOfPrimitives()->Add(newobj,option); | 
|---|
| 576 | else                pad->GetListOfPrimitives()->Add(newobj,GetDrawOption()); | 
|---|
| 577 | pad->Modified(kTRUE); | 
|---|
| 578 | pad->Update(); | 
|---|
| 579 | if (padsav) padsav->cd(); | 
|---|
| 580 | return newobj; | 
|---|
| 581 | } | 
|---|
| 582 | */ | 
|---|
| 583 |  | 
|---|
| 584 | TString opt(option); | 
|---|
| 585 | opt.ToLower(); | 
|---|
| 586 |  | 
|---|
| 587 | newobj->Draw(opt.IsNull() ? GetDrawOption() : option); | 
|---|
| 588 |  | 
|---|
| 589 | if (padsav) | 
|---|
| 590 | padsav->cd(); | 
|---|
| 591 |  | 
|---|
| 592 | return newobj; | 
|---|
| 593 | } | 
|---|
| 594 |  | 
|---|
| 595 | // ------------------------------------------------------------------------ | 
|---|
| 596 | // | 
|---|
| 597 | // Creates a TH1D which contains the projection of the contents of the | 
|---|
| 598 | // MHCamera onto the y-axis. The maximum and minimum are calculated | 
|---|
| 599 | // such that a slighly wider range than (GetMinimum(), GetMaximum()) is | 
|---|
| 600 | // displayed using THLimitsFinder::OptimizeLimits. | 
|---|
| 601 | // | 
|---|
| 602 | // If no name is given the newly allocated histogram is removed from | 
|---|
| 603 | // the current directory calling SetDirectory(0) in any other case | 
|---|
| 604 | // the newly created histogram is removed from the current directory | 
|---|
| 605 | // and added to gROOT such the gROOT->FindObject can find the histogram. | 
|---|
| 606 | // | 
|---|
| 607 | // If the standard name "_py" is given "_py" is appended to the name | 
|---|
| 608 | // of the MHCamera and the corresponding histogram is searched using | 
|---|
| 609 | // gROOT->FindObject and updated with the present projection. | 
|---|
| 610 | // | 
|---|
| 611 | // It is the responsibility of the user to make sure, that the newly | 
|---|
| 612 | // created histogram is freed correctly. | 
|---|
| 613 | // | 
|---|
| 614 | // Currently the new histogram is restrictred to 50 bins. | 
|---|
| 615 | // Maybe a optimal number can be calulated from the number of | 
|---|
| 616 | // bins on the x-axis of the MHCamera? | 
|---|
| 617 | // | 
|---|
| 618 | // The code was taken mainly from TH2::ProjectX such the interface | 
|---|
| 619 | // is more or less the same than to TH2-projections. | 
|---|
| 620 | // | 
|---|
| 621 | // If sector>=0 only entries with matching sector index are taken | 
|---|
| 622 | // into account. | 
|---|
| 623 | // | 
|---|
| 624 | TH1D *MHCamera::ProjectionS(const TArrayI §or, const TArrayI &aidx, const char *name, const Int_t nbins) const | 
|---|
| 625 | { | 
|---|
| 626 |  | 
|---|
| 627 | // Create the projection histogram | 
|---|
| 628 | TString pname(name); | 
|---|
| 629 | if (name=="_py") | 
|---|
| 630 | { | 
|---|
| 631 | pname.Prepend(GetName()); | 
|---|
| 632 | if (sector.GetSize()>0) | 
|---|
| 633 | { | 
|---|
| 634 | pname += ";"; | 
|---|
| 635 | for (int i=0; i<sector.GetSize(); i++) | 
|---|
| 636 | pname += sector[i]; | 
|---|
| 637 | } | 
|---|
| 638 | if (aidx.GetSize()>0) | 
|---|
| 639 | { | 
|---|
| 640 | pname += ";"; | 
|---|
| 641 | for (int i=0; i<aidx.GetSize(); i++) | 
|---|
| 642 | pname += aidx[i]; | 
|---|
| 643 | } | 
|---|
| 644 | } | 
|---|
| 645 |  | 
|---|
| 646 | TH1D *h1=0; | 
|---|
| 647 |  | 
|---|
| 648 | //check if histogram with identical name exist | 
|---|
| 649 | TObject *h1obj = gROOT->FindObject(pname); | 
|---|
| 650 | if (h1obj && h1obj->InheritsFrom("TH1D")) { | 
|---|
| 651 | h1 = (TH1D*)h1obj; | 
|---|
| 652 | h1->Reset(); | 
|---|
| 653 | } | 
|---|
| 654 |  | 
|---|
| 655 | if (!h1) | 
|---|
| 656 | { | 
|---|
| 657 | h1 = new TH1D; | 
|---|
| 658 | h1->UseCurrentStyle(); | 
|---|
| 659 | h1->SetName(pname); | 
|---|
| 660 | h1->SetTitle(GetTitle()); | 
|---|
| 661 | h1->SetDirectory(0); | 
|---|
| 662 | h1->SetXTitle(GetYaxis()->GetTitle()); | 
|---|
| 663 | h1->SetYTitle("Counts"); | 
|---|
| 664 | //h1->Sumw2(); | 
|---|
| 665 | } | 
|---|
| 666 |  | 
|---|
| 667 | Double_t min = GetMinimumSectors(sector, aidx); | 
|---|
| 668 | Double_t max = GetMaximumSectors(sector, aidx); | 
|---|
| 669 |  | 
|---|
| 670 | if (min==max && max>0) | 
|---|
| 671 | min=0; | 
|---|
| 672 | if (min==max && min<0) | 
|---|
| 673 | max=0; | 
|---|
| 674 |  | 
|---|
| 675 | Int_t newbins=0; | 
|---|
| 676 | THLimitsFinder::OptimizeLimits(nbins, newbins, min, max, kFALSE); | 
|---|
| 677 |  | 
|---|
| 678 | MBinning bins(nbins, min, max); | 
|---|
| 679 | bins.Apply(*h1); | 
|---|
| 680 |  | 
|---|
| 681 | // Fill the projected histogram | 
|---|
| 682 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 683 | if (IsUsed(idx) && MatchSector(idx, sector, aidx)) | 
|---|
| 684 | h1->Fill(GetBinContent(idx+1)); | 
|---|
| 685 |  | 
|---|
| 686 | return h1; | 
|---|
| 687 | } | 
|---|
| 688 |  | 
|---|
| 689 | // ------------------------------------------------------------------------ | 
|---|
| 690 | // | 
|---|
| 691 | // Creates a TH1D which contains the projection of the contents of the | 
|---|
| 692 | // MHCamera onto the radius from the camera center. | 
|---|
| 693 | // The maximum and minimum are calculated | 
|---|
| 694 | // such that a slighly wider range than (GetMinimum(), GetMaximum()) is | 
|---|
| 695 | // displayed using THLimitsFinder::OptimizeLimits. | 
|---|
| 696 | // | 
|---|
| 697 | // If no name is given the newly allocated histogram is removed from | 
|---|
| 698 | // the current directory calling SetDirectory(0) in any other case | 
|---|
| 699 | // the newly created histogram is removed from the current directory | 
|---|
| 700 | // and added to gROOT such the gROOT->FindObject can find the histogram. | 
|---|
| 701 | // | 
|---|
| 702 | // If the standard name "_rad" is given "_rad" is appended to the name | 
|---|
| 703 | // of the MHCamera and the corresponding histogram is searched using | 
|---|
| 704 | // gROOT->FindObject and updated with the present projection. | 
|---|
| 705 | // | 
|---|
| 706 | // It is the responsibility of the user to make sure, that the newly | 
|---|
| 707 | // created histogram is freed correctly. | 
|---|
| 708 | // | 
|---|
| 709 | // Currently the new histogram is restrictred to 50 bins. | 
|---|
| 710 | // Maybe a optimal number can be calulated from the number of | 
|---|
| 711 | // bins on the x-axis of the MHCamera? | 
|---|
| 712 | // | 
|---|
| 713 | // The code was taken mainly from TH2::ProjectX such the interface | 
|---|
| 714 | // is more or less the same than to TH2-projections. | 
|---|
| 715 | // | 
|---|
| 716 | // If sector>=0 only entries with matching sector index are taken | 
|---|
| 717 | // into account. | 
|---|
| 718 | // | 
|---|
| 719 | TProfile *MHCamera::RadialProfileS(const TArrayI §or, const TArrayI &aidx, const char *name, const Int_t nbins) const | 
|---|
| 720 | { | 
|---|
| 721 | // Create the projection histogram | 
|---|
| 722 | TString pname(name); | 
|---|
| 723 | if (name=="_rad") | 
|---|
| 724 | { | 
|---|
| 725 | pname.Prepend(GetName()); | 
|---|
| 726 | if (sector.GetSize()>0) | 
|---|
| 727 | { | 
|---|
| 728 | pname += ";"; | 
|---|
| 729 | for (int i=0; i<sector.GetSize(); i++) | 
|---|
| 730 | pname += sector[i]; | 
|---|
| 731 | } | 
|---|
| 732 | if (aidx.GetSize()>0) | 
|---|
| 733 | { | 
|---|
| 734 | pname += ";"; | 
|---|
| 735 | for (int i=0; i<aidx.GetSize(); i++) | 
|---|
| 736 | pname += aidx[i]; | 
|---|
| 737 | } | 
|---|
| 738 | } | 
|---|
| 739 |  | 
|---|
| 740 | TProfile *h1=0; | 
|---|
| 741 |  | 
|---|
| 742 | //check if histogram with identical name exist | 
|---|
| 743 | TObject *h1obj = gROOT->FindObject(pname); | 
|---|
| 744 | if (h1obj && h1obj->InheritsFrom("TProfile")) { | 
|---|
| 745 | h1 = (TProfile*)h1obj; | 
|---|
| 746 | h1->Reset(); | 
|---|
| 747 | } | 
|---|
| 748 |  | 
|---|
| 749 | if (!h1) | 
|---|
| 750 | { | 
|---|
| 751 | h1 = new TProfile; | 
|---|
| 752 | h1->UseCurrentStyle(); | 
|---|
| 753 | h1->SetName(pname); | 
|---|
| 754 | h1->SetTitle(GetTitle()); | 
|---|
| 755 | h1->SetDirectory(0); | 
|---|
| 756 | h1->SetXTitle("Radius from camera center [mm]"); | 
|---|
| 757 | h1->SetYTitle(GetYaxis()->GetTitle()); | 
|---|
| 758 | } | 
|---|
| 759 |  | 
|---|
| 760 | const Double_t m2d = fGeomCam->GetConvMm2Deg(); | 
|---|
| 761 |  | 
|---|
| 762 | Double_t min = 0.; | 
|---|
| 763 | Double_t max = fGeomCam->GetMaxRadius()*m2d; | 
|---|
| 764 |  | 
|---|
| 765 | Int_t newbins=0; | 
|---|
| 766 |  | 
|---|
| 767 | THLimitsFinder::OptimizeLimits(nbins, newbins, min, max, kFALSE); | 
|---|
| 768 |  | 
|---|
| 769 | MBinning bins(nbins, min, max); | 
|---|
| 770 | bins.Apply(*h1); | 
|---|
| 771 |  | 
|---|
| 772 | // Fill the projected histogram | 
|---|
| 773 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 774 | if (IsUsed(idx) && MatchSector(idx, sector, aidx)) | 
|---|
| 775 | h1->Fill(TMath::Hypot((*fGeomCam)[idx].GetX(),(*fGeomCam)[idx].GetY())*m2d, | 
|---|
| 776 | GetBinContent(idx+1)); | 
|---|
| 777 | return h1; | 
|---|
| 778 | } | 
|---|
| 779 |  | 
|---|
| 780 |  | 
|---|
| 781 | // ------------------------------------------------------------------------ | 
|---|
| 782 | // | 
|---|
| 783 | // Creates a TH1D which contains the projection of the contents of the | 
|---|
| 784 | // MHCamera onto the azimuth angle in the camera. | 
|---|
| 785 | // | 
|---|
| 786 | // If no name is given the newly allocated histogram is removed from | 
|---|
| 787 | // the current directory calling SetDirectory(0) in any other case | 
|---|
| 788 | // the newly created histogram is removed from the current directory | 
|---|
| 789 | // and added to gROOT such the gROOT->FindObject can find the histogram. | 
|---|
| 790 | // | 
|---|
| 791 | // If the standard name "_azi" is given "_azi" is appended to the name | 
|---|
| 792 | // of the MHCamera and the corresponding histogram is searched using | 
|---|
| 793 | // gROOT->FindObject and updated with the present projection. | 
|---|
| 794 | // | 
|---|
| 795 | // It is the responsibility of the user to make sure, that the newly | 
|---|
| 796 | // created histogram is freed correctly. | 
|---|
| 797 | // | 
|---|
| 798 | // Currently the new histogram is restrictred to 60 bins. | 
|---|
| 799 | // Maybe a optimal number can be calulated from the number of | 
|---|
| 800 | // bins on the x-axis of the MHCamera? | 
|---|
| 801 | // | 
|---|
| 802 | // The code was taken mainly from TH2::ProjectX such the interface | 
|---|
| 803 | // is more or less the same than to TH2-projections. | 
|---|
| 804 | // | 
|---|
| 805 | TProfile *MHCamera::AzimuthProfileA(const TArrayI &aidx, const char *name, const Int_t nbins) const | 
|---|
| 806 | { | 
|---|
| 807 | // Create the projection histogram | 
|---|
| 808 | TString pname(name); | 
|---|
| 809 | if (name=="_azi") | 
|---|
| 810 | { | 
|---|
| 811 | pname.Prepend(GetName()); | 
|---|
| 812 | if (aidx.GetSize()>0) | 
|---|
| 813 | { | 
|---|
| 814 | pname += ";"; | 
|---|
| 815 | for (int i=0; i<aidx.GetSize(); i++) | 
|---|
| 816 | pname += aidx[i]; | 
|---|
| 817 | } | 
|---|
| 818 | } | 
|---|
| 819 |  | 
|---|
| 820 | TProfile *h1=0; | 
|---|
| 821 |  | 
|---|
| 822 | //check if histogram with identical name exist | 
|---|
| 823 | TObject *h1obj = gROOT->FindObject(pname); | 
|---|
| 824 | if (h1obj && h1obj->InheritsFrom("TProfile")) { | 
|---|
| 825 | h1 = (TProfile*)h1obj; | 
|---|
| 826 | h1->Reset(); | 
|---|
| 827 | } | 
|---|
| 828 |  | 
|---|
| 829 | if (!h1) | 
|---|
| 830 | { | 
|---|
| 831 |  | 
|---|
| 832 | h1 = new TProfile; | 
|---|
| 833 | h1->UseCurrentStyle(); | 
|---|
| 834 | h1->SetName(pname); | 
|---|
| 835 | h1->SetTitle(GetTitle()); | 
|---|
| 836 | h1->SetDirectory(0); | 
|---|
| 837 | h1->SetXTitle("Azimuth in camera [deg]"); | 
|---|
| 838 | h1->SetYTitle(GetYaxis()->GetTitle()); | 
|---|
| 839 | } | 
|---|
| 840 |  | 
|---|
| 841 | //Double_t min = 0; | 
|---|
| 842 | //Double_t max = 360; | 
|---|
| 843 |  | 
|---|
| 844 | //Int_t newbins=0; | 
|---|
| 845 | //THLimitsFinder::OptimizeLimits(nbins, newbins, min, max, kFALSE); | 
|---|
| 846 |  | 
|---|
| 847 | MBinning bins(nbins, 0, 360); | 
|---|
| 848 | bins.Apply(*h1); | 
|---|
| 849 |  | 
|---|
| 850 | // Fill the projected histogram | 
|---|
| 851 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 852 | { | 
|---|
| 853 | if (IsUsed(idx) && MatchSector(idx, TArrayI(), aidx)) | 
|---|
| 854 | h1->Fill(TMath::ATan2((*fGeomCam)[idx].GetY(),(*fGeomCam)[idx].GetX())*TMath::RadToDeg()+180, | 
|---|
| 855 | GetPixContent(idx)); | 
|---|
| 856 |  | 
|---|
| 857 | } | 
|---|
| 858 |  | 
|---|
| 859 | return h1; | 
|---|
| 860 | } | 
|---|
| 861 |  | 
|---|
| 862 |  | 
|---|
| 863 | // ------------------------------------------------------------------------ | 
|---|
| 864 | // | 
|---|
| 865 | // Resizes the current pad so that the camera is displayed in its | 
|---|
| 866 | // correct aspect ratio | 
|---|
| 867 | // | 
|---|
| 868 | void MHCamera::SetRange() | 
|---|
| 869 | { | 
|---|
| 870 | const Float_t range = fGeomCam->GetMaxRadius()*1.05; | 
|---|
| 871 |  | 
|---|
| 872 | // | 
|---|
| 873 | // Maintain aspect ratio | 
|---|
| 874 | // | 
|---|
| 875 | const float ratio = TestBit(kNoLegend) ? 1 : 1.15; | 
|---|
| 876 |  | 
|---|
| 877 | // | 
|---|
| 878 | // Calculate width and height of the current pad in pixels | 
|---|
| 879 | // | 
|---|
| 880 | Float_t w = gPad->GetWw(); | 
|---|
| 881 | Float_t h = gPad->GetWh()*ratio; | 
|---|
| 882 |  | 
|---|
| 883 | // | 
|---|
| 884 | // This prevents the pad from resizing itself wrongly | 
|---|
| 885 | // | 
|---|
| 886 | if (gPad->GetMother() != gPad) | 
|---|
| 887 | { | 
|---|
| 888 | w *= gPad->GetMother()->GetAbsWNDC(); | 
|---|
| 889 | h *= gPad->GetMother()->GetAbsHNDC(); | 
|---|
| 890 | } | 
|---|
| 891 |  | 
|---|
| 892 | // | 
|---|
| 893 | // Set Range (coordinate system) of pad | 
|---|
| 894 | // | 
|---|
| 895 | gPad->Range(-range, -range, (2*ratio-1)*range, range); | 
|---|
| 896 |  | 
|---|
| 897 | // | 
|---|
| 898 | // Resize Pad to given ratio | 
|---|
| 899 | // | 
|---|
| 900 | if (h<w) | 
|---|
| 901 | gPad->SetPad((1.-h/w)/2, 0, (h/w+1.)/2, 1); | 
|---|
| 902 | else | 
|---|
| 903 | gPad->SetPad(0, (1.-w/h)/2, 1, (w/h+1.)/2); | 
|---|
| 904 | } | 
|---|
| 905 |  | 
|---|
| 906 | // ------------------------------------------------------------------------ | 
|---|
| 907 | // | 
|---|
| 908 | // Updates the pixel colors and paints the pixels | 
|---|
| 909 | // | 
|---|
| 910 | void MHCamera::Update(Bool_t islog, Bool_t isbox, Bool_t iscol, Bool_t issame) | 
|---|
| 911 | { | 
|---|
| 912 | Double_t min = GetMinimum(kFALSE); | 
|---|
| 913 | Double_t max = GetMaximum(kFALSE); | 
|---|
| 914 | if (min==FLT_MAX) | 
|---|
| 915 | { | 
|---|
| 916 | min = 0; | 
|---|
| 917 | max = 1; | 
|---|
| 918 | } | 
|---|
| 919 |  | 
|---|
| 920 | if (min==max) | 
|---|
| 921 | max += 1; | 
|---|
| 922 |  | 
|---|
| 923 | if (!issame) | 
|---|
| 924 | UpdateLegend(min, max, islog); | 
|---|
| 925 |  | 
|---|
| 926 | // Try to estimate the units of the current display. This is only | 
|---|
| 927 | // necessary for 'same' option and allows distorted images of the camera! | 
|---|
| 928 | const Float_t maxr = (1-fGeomCam->GetConvMm2Deg())*fGeomCam->GetMaxRadius()/2; | 
|---|
| 929 | const Float_t conv = !issame || | 
|---|
| 930 | gPad->GetX1()<-maxr || gPad->GetY1()<-maxr || | 
|---|
| 931 | gPad->GetX2()> maxr || gPad->GetY2()>maxr ? 1 : fGeomCam->GetConvMm2Deg(); | 
|---|
| 932 |  | 
|---|
| 933 | MHexagon hex; | 
|---|
| 934 | for (Int_t i=0; i<fNcells-2; i++) | 
|---|
| 935 | { | 
|---|
| 936 | hex.SetFillStyle(issame || (IsTransparent() && !IsUsed(i)) ? 0 : 1001); | 
|---|
| 937 |  | 
|---|
| 938 | if (!issame) | 
|---|
| 939 | { | 
|---|
| 940 | const Bool_t isnan = !TMath::Finite(fArray[i+1]); | 
|---|
| 941 | if (!IsUsed(i) || !iscol || isnan) | 
|---|
| 942 | { | 
|---|
| 943 | hex.SetFillColor(10); | 
|---|
| 944 |  | 
|---|
| 945 | if (isnan) | 
|---|
| 946 | gLog << warn << "MHCamera::Update: " << GetName() << " <" << GetTitle() << "> - Pixel Index #" << i << " contents is not finite..." << endl; | 
|---|
| 947 | } | 
|---|
| 948 | else | 
|---|
| 949 | hex.SetFillColor(GetColor(GetBinContent(i+1), min, max, islog)); | 
|---|
| 950 | } | 
|---|
| 951 |  | 
|---|
| 952 | const MGeomPix &pix = (*fGeomCam)[i]; | 
|---|
| 953 |  | 
|---|
| 954 | Float_t x = pix.GetX()*conv/(fAbberation+1); | 
|---|
| 955 | Float_t y = pix.GetY()*conv/(fAbberation+1); | 
|---|
| 956 | Float_t d = pix.GetD()*conv; | 
|---|
| 957 |  | 
|---|
| 958 | if (!isbox) | 
|---|
| 959 | if (IsUsed(i) || !TestBit(kNoUnused)) | 
|---|
| 960 | hex.PaintHexagon(x, y, d); | 
|---|
| 961 | else | 
|---|
| 962 | if (IsUsed(i) && TMath::Finite(fArray[i+1])) | 
|---|
| 963 | { | 
|---|
| 964 | Float_t size = d*(GetBinContent(i+1)-min)/(max-min); | 
|---|
| 965 | if (size>d) | 
|---|
| 966 | size=d; | 
|---|
| 967 | hex.PaintHexagon(x, y, size); | 
|---|
| 968 | } | 
|---|
| 969 | } | 
|---|
| 970 | } | 
|---|
| 971 |  | 
|---|
| 972 | // ------------------------------------------------------------------------ | 
|---|
| 973 | // | 
|---|
| 974 | // Print minimum and maximum | 
|---|
| 975 | // | 
|---|
| 976 | void MHCamera::Print(Option_t *) const | 
|---|
| 977 | { | 
|---|
| 978 | gLog << all << "Minimum: " << GetMinimum(); | 
|---|
| 979 | if (fMinimum==-1111) | 
|---|
| 980 | gLog << " <autoscaled>"; | 
|---|
| 981 | gLog << endl; | 
|---|
| 982 | gLog << "Maximum: " << GetMaximum(); | 
|---|
| 983 | if (fMaximum==-1111) | 
|---|
| 984 | gLog << " <autoscaled>"; | 
|---|
| 985 | gLog << endl; | 
|---|
| 986 | } | 
|---|
| 987 |  | 
|---|
| 988 | // ------------------------------------------------------------------------ | 
|---|
| 989 | // | 
|---|
| 990 | // Paint the y-axis title | 
|---|
| 991 | // | 
|---|
| 992 | void MHCamera::PaintAxisTitle() | 
|---|
| 993 | { | 
|---|
| 994 | const Float_t range = fGeomCam->GetMaxRadius()*1.05; | 
|---|
| 995 | const Float_t w = (1 + 1.5/sqrt((float)(fNcells-2)))*range; | 
|---|
| 996 |  | 
|---|
| 997 | TLatex *ptitle = new TLatex(w, -.90*range, GetYaxis()->GetTitle()); | 
|---|
| 998 |  | 
|---|
| 999 | ptitle->SetTextSize(0.05); | 
|---|
| 1000 | ptitle->SetTextAlign(21); | 
|---|
| 1001 |  | 
|---|
| 1002 | // box with the histogram title | 
|---|
| 1003 | ptitle->SetTextColor(gStyle->GetTitleTextColor()); | 
|---|
| 1004 | #if ROOT_VERSION_CODE > ROOT_VERSION(3,05,01) | 
|---|
| 1005 | ptitle->SetTextFont(gStyle->GetTitleFont("")); | 
|---|
| 1006 | #endif | 
|---|
| 1007 | ptitle->Paint(); | 
|---|
| 1008 | } | 
|---|
| 1009 |  | 
|---|
| 1010 | // ------------------------------------------------------------------------ | 
|---|
| 1011 | // | 
|---|
| 1012 | // Paints the camera. | 
|---|
| 1013 | // | 
|---|
| 1014 | void MHCamera::Paint(Option_t *o) | 
|---|
| 1015 | { | 
|---|
| 1016 | if (fNcells<=1) | 
|---|
| 1017 | return; | 
|---|
| 1018 |  | 
|---|
| 1019 | TString opt(o); | 
|---|
| 1020 | opt.ToLower(); | 
|---|
| 1021 |  | 
|---|
| 1022 | if (opt.Contains("hist")) | 
|---|
| 1023 | { | 
|---|
| 1024 | opt.ReplaceAll("hist", ""); | 
|---|
| 1025 | opt.ReplaceAll("box", ""); | 
|---|
| 1026 | opt.ReplaceAll("pixelindex", ""); | 
|---|
| 1027 | opt.ReplaceAll("sectorindex", ""); | 
|---|
| 1028 | opt.ReplaceAll("abscontent", ""); | 
|---|
| 1029 | opt.ReplaceAll("content", ""); | 
|---|
| 1030 | opt.ReplaceAll("proj", ""); | 
|---|
| 1031 | opt.ReplaceAll("pal0", ""); | 
|---|
| 1032 | opt.ReplaceAll("pal1", ""); | 
|---|
| 1033 | opt.ReplaceAll("pal2", ""); | 
|---|
| 1034 | opt.ReplaceAll("nopal", ""); | 
|---|
| 1035 | TH1D::Paint(opt); | 
|---|
| 1036 | return; | 
|---|
| 1037 | } | 
|---|
| 1038 |  | 
|---|
| 1039 | if (opt.Contains("proj")) | 
|---|
| 1040 | { | 
|---|
| 1041 | opt.ReplaceAll("proj", ""); | 
|---|
| 1042 | Projection(GetName())->Paint(opt); | 
|---|
| 1043 | return; | 
|---|
| 1044 | } | 
|---|
| 1045 |  | 
|---|
| 1046 | const Bool_t hassame = opt.Contains("same"); | 
|---|
| 1047 | const Bool_t hasbox  = opt.Contains("box"); | 
|---|
| 1048 | const Bool_t hascol  = hasbox ? !opt.Contains("nocol") : kTRUE; | 
|---|
| 1049 |  | 
|---|
| 1050 | if (!hassame) | 
|---|
| 1051 | { | 
|---|
| 1052 | gPad->Clear(); | 
|---|
| 1053 |  | 
|---|
| 1054 | // Maintain aspect ratio | 
|---|
| 1055 | SetRange(); | 
|---|
| 1056 |  | 
|---|
| 1057 | if (GetPainter()) | 
|---|
| 1058 | { | 
|---|
| 1059 | // Paint statistics | 
|---|
| 1060 | if (!TestBit(TH1::kNoStats)) | 
|---|
| 1061 | fPainter->PaintStat(gStyle->GetOptStat(), NULL); | 
|---|
| 1062 |  | 
|---|
| 1063 | // Paint primitives (pixels, color legend, photons, ...) | 
|---|
| 1064 | if (fPainter->InheritsFrom(THistPainter::Class())) | 
|---|
| 1065 | { | 
|---|
| 1066 | static_cast<THistPainter*>(fPainter)->MakeChopt(""); | 
|---|
| 1067 | static_cast<THistPainter*>(fPainter)->PaintTitle(); | 
|---|
| 1068 | } | 
|---|
| 1069 | } | 
|---|
| 1070 | } | 
|---|
| 1071 |  | 
|---|
| 1072 | const Bool_t pal1  = opt.Contains("pal1"); | 
|---|
| 1073 | const Bool_t pal2  = opt.Contains("pal2"); | 
|---|
| 1074 | const Bool_t nopal = opt.Contains("nopal"); | 
|---|
| 1075 |  | 
|---|
| 1076 | if (!pal1 && !pal2 && !nopal) | 
|---|
| 1077 | SetPrettyPalette(); | 
|---|
| 1078 |  | 
|---|
| 1079 | if (pal1) | 
|---|
| 1080 | SetDeepBlueSeaPalette(); | 
|---|
| 1081 |  | 
|---|
| 1082 | if (pal2) | 
|---|
| 1083 | SetInvDeepBlueSeaPalette(); | 
|---|
| 1084 |  | 
|---|
| 1085 | // Update Contents of the pixels and paint legend | 
|---|
| 1086 | Update(gPad->GetLogy(), hasbox, hascol, hassame); | 
|---|
| 1087 |  | 
|---|
| 1088 | if (!hassame) | 
|---|
| 1089 | PaintAxisTitle(); | 
|---|
| 1090 |  | 
|---|
| 1091 | if (opt.Contains("pixelindex")) | 
|---|
| 1092 | { | 
|---|
| 1093 | PaintIndices(0); | 
|---|
| 1094 | return; | 
|---|
| 1095 | } | 
|---|
| 1096 | if (opt.Contains("sectorindex")) | 
|---|
| 1097 | { | 
|---|
| 1098 | PaintIndices(1); | 
|---|
| 1099 | return; | 
|---|
| 1100 | } | 
|---|
| 1101 | if (opt.Contains("abscontent")) | 
|---|
| 1102 | { | 
|---|
| 1103 | PaintIndices(3); | 
|---|
| 1104 | return; | 
|---|
| 1105 | } | 
|---|
| 1106 | if (opt.Contains("content")) | 
|---|
| 1107 | { | 
|---|
| 1108 | PaintIndices(2); | 
|---|
| 1109 | return; | 
|---|
| 1110 | } | 
|---|
| 1111 | if (opt.Contains("pixelentries")) | 
|---|
| 1112 | { | 
|---|
| 1113 | PaintIndices(4); | 
|---|
| 1114 | return; | 
|---|
| 1115 | } | 
|---|
| 1116 | } | 
|---|
| 1117 |  | 
|---|
| 1118 | void MHCamera::SetDrawOption(Option_t *option) | 
|---|
| 1119 | { | 
|---|
| 1120 | // This is a workaround. For some reason MHCamera is | 
|---|
| 1121 | // stored in a TObjLink instead of a TObjOptLink | 
|---|
| 1122 | if (!option || !gPad) | 
|---|
| 1123 | return; | 
|---|
| 1124 |  | 
|---|
| 1125 | TListIter next(gPad->GetListOfPrimitives()); | 
|---|
| 1126 | delete gPad->FindObject("Tframe"); | 
|---|
| 1127 | TObject *obj; | 
|---|
| 1128 | while ((obj = next())) | 
|---|
| 1129 | if (obj == this && (TString)next.GetOption()!=(TString)option) | 
|---|
| 1130 | { | 
|---|
| 1131 | gPad->GetListOfPrimitives()->Remove(this); | 
|---|
| 1132 | gPad->GetListOfPrimitives()->AddFirst(this, option); | 
|---|
| 1133 | return; | 
|---|
| 1134 | } | 
|---|
| 1135 | } | 
|---|
| 1136 |  | 
|---|
| 1137 | // ------------------------------------------------------------------------ | 
|---|
| 1138 | // | 
|---|
| 1139 | // With this function you can change the color palette. For more | 
|---|
| 1140 | // information see TStyle::SetPalette. Only palettes with 50 colors | 
|---|
| 1141 | // are allowed. | 
|---|
| 1142 | // In addition you can use SetPalette(52, 0) to create an inverse | 
|---|
| 1143 | // deep blue sea palette | 
|---|
| 1144 | // | 
|---|
| 1145 | void MHCamera::SetPalette(Int_t ncolors, Int_t *colors) | 
|---|
| 1146 | { | 
|---|
| 1147 | // | 
|---|
| 1148 | // If not enough colors are specified skip this. | 
|---|
| 1149 | // | 
|---|
| 1150 | if (ncolors>1 && ncolors<50) | 
|---|
| 1151 | { | 
|---|
| 1152 | gLog << err << "MHCamera::SetPalette: Only default palettes with 50 colors are allowed... ignored." << endl; | 
|---|
| 1153 | return; | 
|---|
| 1154 | } | 
|---|
| 1155 |  | 
|---|
| 1156 | // | 
|---|
| 1157 | // If ncolors==52 create a reversed deep blue sea palette | 
|---|
| 1158 | // | 
|---|
| 1159 | if (ncolors==52) | 
|---|
| 1160 | { | 
|---|
| 1161 | gStyle->SetPalette(51, NULL); | 
|---|
| 1162 | TArrayI c(kItemsLegend); | 
|---|
| 1163 | for (int i=0; i<kItemsLegend; i++) | 
|---|
| 1164 | c[kItemsLegend-i-1] = gStyle->GetColorPalette(i); | 
|---|
| 1165 | gStyle->SetPalette(kItemsLegend, c.GetArray()); | 
|---|
| 1166 | } | 
|---|
| 1167 | else | 
|---|
| 1168 | gStyle->SetPalette(ncolors, colors); | 
|---|
| 1169 | } | 
|---|
| 1170 |  | 
|---|
| 1171 |  | 
|---|
| 1172 | // ------------------------------------------------------------------------ | 
|---|
| 1173 | // | 
|---|
| 1174 | // Changes the palette of the displayed camera histogram. | 
|---|
| 1175 | // | 
|---|
| 1176 | // Change to the right pad first - otherwise GetDrawOption() might fail. | 
|---|
| 1177 | // | 
|---|
| 1178 | void MHCamera::SetPrettyPalette() | 
|---|
| 1179 | { | 
|---|
| 1180 | TString opt(GetDrawOption()); | 
|---|
| 1181 |  | 
|---|
| 1182 | if (!opt.Contains("hist", TString::kIgnoreCase)) | 
|---|
| 1183 | SetPalette(1, 0); | 
|---|
| 1184 |  | 
|---|
| 1185 | opt.ReplaceAll("pal1", ""); | 
|---|
| 1186 | opt.ReplaceAll("pal2", ""); | 
|---|
| 1187 |  | 
|---|
| 1188 | SetDrawOption(opt); | 
|---|
| 1189 | } | 
|---|
| 1190 |  | 
|---|
| 1191 | // ------------------------------------------------------------------------ | 
|---|
| 1192 | // | 
|---|
| 1193 | // Changes the palette of the displayed camera histogram. | 
|---|
| 1194 | // | 
|---|
| 1195 | // Change to the right pad first - otherwise GetDrawOption() might fail. | 
|---|
| 1196 | // | 
|---|
| 1197 | void MHCamera::SetDeepBlueSeaPalette() | 
|---|
| 1198 | { | 
|---|
| 1199 | TString opt(GetDrawOption()); | 
|---|
| 1200 |  | 
|---|
| 1201 | if (!opt.Contains("hist", TString::kIgnoreCase)) | 
|---|
| 1202 | SetPalette(51, 0); | 
|---|
| 1203 |  | 
|---|
| 1204 | opt.ReplaceAll("pal1", ""); | 
|---|
| 1205 | opt.ReplaceAll("pal2", ""); | 
|---|
| 1206 | opt += "pal1"; | 
|---|
| 1207 |  | 
|---|
| 1208 | SetDrawOption(opt); | 
|---|
| 1209 | } | 
|---|
| 1210 |  | 
|---|
| 1211 | // ------------------------------------------------------------------------ | 
|---|
| 1212 | // | 
|---|
| 1213 | // Changes the palette of the displayed camera histogram. | 
|---|
| 1214 | // | 
|---|
| 1215 | // Change to the right pad first - otherwise GetDrawOption() might fail. | 
|---|
| 1216 | // | 
|---|
| 1217 | void MHCamera::SetInvDeepBlueSeaPalette() | 
|---|
| 1218 | { | 
|---|
| 1219 | TString opt(GetDrawOption()); | 
|---|
| 1220 |  | 
|---|
| 1221 | if (!opt.Contains("hist", TString::kIgnoreCase)) | 
|---|
| 1222 | SetPalette(52, 0); | 
|---|
| 1223 |  | 
|---|
| 1224 | opt.ReplaceAll("pal1", ""); | 
|---|
| 1225 | opt.ReplaceAll("pal2", ""); | 
|---|
| 1226 | opt += "pal2"; | 
|---|
| 1227 |  | 
|---|
| 1228 | SetDrawOption(opt); | 
|---|
| 1229 | } | 
|---|
| 1230 |  | 
|---|
| 1231 | // ------------------------------------------------------------------------ | 
|---|
| 1232 | // | 
|---|
| 1233 | // Paint indices (as text) inside the pixels. Depending of the type- | 
|---|
| 1234 | // argument we paint: | 
|---|
| 1235 | //  0: pixel number | 
|---|
| 1236 | //  1: sector number | 
|---|
| 1237 | //  2: content | 
|---|
| 1238 | // | 
|---|
| 1239 | void MHCamera::PaintIndices(Int_t type) | 
|---|
| 1240 | { | 
|---|
| 1241 | if (fNcells<=1) | 
|---|
| 1242 | return; | 
|---|
| 1243 |  | 
|---|
| 1244 | const Double_t min = GetMinimum(); | 
|---|
| 1245 | const Double_t max = GetMaximum(); | 
|---|
| 1246 |  | 
|---|
| 1247 | if (type==2 && max==min) | 
|---|
| 1248 | return; | 
|---|
| 1249 |  | 
|---|
| 1250 | TText txt; | 
|---|
| 1251 | txt.SetTextFont(122); | 
|---|
| 1252 | txt.SetTextAlign(22);   // centered/centered | 
|---|
| 1253 |  | 
|---|
| 1254 | for (Int_t i=0; i<fNcells-2; i++) | 
|---|
| 1255 | { | 
|---|
| 1256 | const MGeomPix &h = (*fGeomCam)[i]; | 
|---|
| 1257 |  | 
|---|
| 1258 | TString num; | 
|---|
| 1259 | switch (type) | 
|---|
| 1260 | { | 
|---|
| 1261 | case 0: num += i; break; | 
|---|
| 1262 | case 1: num += h.GetSector(); break; | 
|---|
| 1263 | case 2: num += TMath::Nint((fArray[i+1]-min)/(max-min)); break; | 
|---|
| 1264 | case 3: num += TMath::Nint(fArray[i+1]); break; | 
|---|
| 1265 | case 4: num += fBinEntries[i+1]; break; | 
|---|
| 1266 | } | 
|---|
| 1267 |  | 
|---|
| 1268 | // FIXME: Should depend on the color of the pixel... | 
|---|
| 1269 | //(GetColor(GetBinContent(i+1), min, max, 0)); | 
|---|
| 1270 | txt.SetTextColor(kRed); | 
|---|
| 1271 | txt.SetTextSize(0.3*h.GetD()/fGeomCam->GetMaxRadius()/1.05); | 
|---|
| 1272 | txt.PaintText(h.GetX(), h.GetY(), num); | 
|---|
| 1273 | } | 
|---|
| 1274 | } | 
|---|
| 1275 |  | 
|---|
| 1276 | // ------------------------------------------------------------------------ | 
|---|
| 1277 | // | 
|---|
| 1278 | // Call this function to add a MCamEvent on top of the present contents. | 
|---|
| 1279 | // | 
|---|
| 1280 | void MHCamera::AddCamContent(const MCamEvent &event, Int_t type) | 
|---|
| 1281 | { | 
|---|
| 1282 | if (fNcells<=1 || IsFreezed()) | 
|---|
| 1283 | return; | 
|---|
| 1284 |  | 
|---|
| 1285 | // FIXME: Security check missing! | 
|---|
| 1286 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 1287 | { | 
|---|
| 1288 | Double_t val=0; | 
|---|
| 1289 | if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/) | 
|---|
| 1290 | { | 
|---|
| 1291 | SetUsed(idx); | 
|---|
| 1292 | Fill(idx, val); // FIXME: Slow! | 
|---|
| 1293 | } | 
|---|
| 1294 | } | 
|---|
| 1295 | fEntries++; | 
|---|
| 1296 | } | 
|---|
| 1297 |  | 
|---|
| 1298 | // ------------------------------------------------------------------------ | 
|---|
| 1299 | // | 
|---|
| 1300 | // Call this function to add a MCamEvent on top of the present contents. | 
|---|
| 1301 | // | 
|---|
| 1302 | void MHCamera::SetCamError(const MCamEvent &evt, Int_t type) | 
|---|
| 1303 | { | 
|---|
| 1304 |  | 
|---|
| 1305 | if (fNcells<=1 || IsFreezed()) | 
|---|
| 1306 | return; | 
|---|
| 1307 |  | 
|---|
| 1308 | // FIXME: Security check missing! | 
|---|
| 1309 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 1310 | { | 
|---|
| 1311 | Double_t val=0; | 
|---|
| 1312 | if (evt.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/) | 
|---|
| 1313 | SetUsed(idx); | 
|---|
| 1314 |  | 
|---|
| 1315 | SetBinError(idx+1, val); // FIXME: Slow! | 
|---|
| 1316 | } | 
|---|
| 1317 | } | 
|---|
| 1318 |  | 
|---|
| 1319 | Stat_t MHCamera::GetBinContent(Int_t bin) const | 
|---|
| 1320 | { | 
|---|
| 1321 | if (fBuffer) ((TH1D*)this)->BufferEmpty(); | 
|---|
| 1322 | if (bin < 0) bin = 0; | 
|---|
| 1323 | if (bin >= fNcells) bin = fNcells-1; | 
|---|
| 1324 | if (!fArray) return 0; | 
|---|
| 1325 |  | 
|---|
| 1326 | if (!TestBit(kProfile)) | 
|---|
| 1327 | return Stat_t (fArray[bin]); | 
|---|
| 1328 |  | 
|---|
| 1329 | if (fBinEntries.fArray[bin] == 0) return 0; | 
|---|
| 1330 | return fArray[bin]/fBinEntries.fArray[bin]; | 
|---|
| 1331 | } | 
|---|
| 1332 |  | 
|---|
| 1333 | Stat_t MHCamera::GetBinError(Int_t bin) const | 
|---|
| 1334 | { | 
|---|
| 1335 | if (!TestBit(kProfile)) | 
|---|
| 1336 | return TH1D::GetBinError(bin); | 
|---|
| 1337 |  | 
|---|
| 1338 | const UInt_t n = (UInt_t)fBinEntries[bin]; | 
|---|
| 1339 |  | 
|---|
| 1340 | if (n==0) | 
|---|
| 1341 | return 0; | 
|---|
| 1342 |  | 
|---|
| 1343 | const Double_t sqr = fSumw2.fArray[bin] / n; | 
|---|
| 1344 | const Double_t val = fArray[bin]        / n; | 
|---|
| 1345 |  | 
|---|
| 1346 | return sqr>val*val ? TMath::Sqrt(sqr - val*val) / n : 0; | 
|---|
| 1347 |  | 
|---|
| 1348 | /* | 
|---|
| 1349 | Double_t rc = 0; | 
|---|
| 1350 | if (TestBit(kSqrtVariance) && GetEntries()>0) // error on the mean | 
|---|
| 1351 | { | 
|---|
| 1352 | const Double_t error = fSumw2.fArray[bin]/GetEntries(); | 
|---|
| 1353 | const Double_t val   = fArray[bin]/GetEntries(); | 
|---|
| 1354 | rc = val*val>error ? 0 : TMath::Sqrt(error - val*val); | 
|---|
| 1355 | } | 
|---|
| 1356 | else | 
|---|
| 1357 | rc = TH1D::GetBinError(bin); | 
|---|
| 1358 |  | 
|---|
| 1359 | return Profile(rc);*/ | 
|---|
| 1360 | } | 
|---|
| 1361 |  | 
|---|
| 1362 | // ------------------------------------------------------------------------ | 
|---|
| 1363 | // | 
|---|
| 1364 | // Call this function to add a MHCamera on top of the present contents. | 
|---|
| 1365 | // Type: | 
|---|
| 1366 | //  0) bin content | 
|---|
| 1367 | //  1) errors | 
|---|
| 1368 | //  2) rel. errors | 
|---|
| 1369 | // | 
|---|
| 1370 | void MHCamera::AddCamContent(const MHCamera &d, Int_t type) | 
|---|
| 1371 | { | 
|---|
| 1372 | if (fNcells!=d.fNcells || IsFreezed()) | 
|---|
| 1373 | return; | 
|---|
| 1374 |  | 
|---|
| 1375 | // FIXME: Security check missing! | 
|---|
| 1376 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 1377 | if (d.IsUsed(idx)) | 
|---|
| 1378 | SetUsed(idx); | 
|---|
| 1379 |  | 
|---|
| 1380 | switch (type) | 
|---|
| 1381 | { | 
|---|
| 1382 | case 1: | 
|---|
| 1383 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 1384 | Fill(idx, d.GetBinError(idx+1)); | 
|---|
| 1385 | break; | 
|---|
| 1386 | case 2: | 
|---|
| 1387 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 1388 | if (d.GetBinContent(idx+1)!=0) | 
|---|
| 1389 | Fill(idx, TMath::Abs(d.GetBinError(idx+1)/d.GetBinContent(idx+1))); | 
|---|
| 1390 | break; | 
|---|
| 1391 | default: | 
|---|
| 1392 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 1393 | Fill(idx, d.GetBinContent(idx+1)); | 
|---|
| 1394 | break; | 
|---|
| 1395 | } | 
|---|
| 1396 | fEntries++; | 
|---|
| 1397 | } | 
|---|
| 1398 |  | 
|---|
| 1399 | // ------------------------------------------------------------------------ | 
|---|
| 1400 | // | 
|---|
| 1401 | // Call this function to add a TArrayD on top of the present contents. | 
|---|
| 1402 | // | 
|---|
| 1403 | void MHCamera::AddCamContent(const TArrayD &event, const TArrayC *used) | 
|---|
| 1404 | { | 
|---|
| 1405 | if (event.GetSize()!=fNcells-2 || IsFreezed()) | 
|---|
| 1406 | return; | 
|---|
| 1407 |  | 
|---|
| 1408 | if (used && used->GetSize()!=fNcells-2) | 
|---|
| 1409 | return; | 
|---|
| 1410 |  | 
|---|
| 1411 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 1412 | { | 
|---|
| 1413 | Fill(idx, event[idx]); // FIXME: Slow! | 
|---|
| 1414 |  | 
|---|
| 1415 | if (!used || (*used)[idx]) | 
|---|
| 1416 | SetUsed(idx); | 
|---|
| 1417 | } | 
|---|
| 1418 | fEntries++; | 
|---|
| 1419 | } | 
|---|
| 1420 |  | 
|---|
| 1421 | // ------------------------------------------------------------------------ | 
|---|
| 1422 | // | 
|---|
| 1423 | // Call this function to add a MArrayD on top of the present contents. | 
|---|
| 1424 | // | 
|---|
| 1425 | void MHCamera::AddCamContent(const MArrayD &event, const TArrayC *used) | 
|---|
| 1426 | { | 
|---|
| 1427 | if (event.GetSize()!=(UInt_t)(fNcells-2) || IsFreezed()) | 
|---|
| 1428 | return; | 
|---|
| 1429 |  | 
|---|
| 1430 | if (used && used->GetSize()!=fNcells-2) | 
|---|
| 1431 | return; | 
|---|
| 1432 |  | 
|---|
| 1433 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 1434 | { | 
|---|
| 1435 | Fill(idx, event[idx]); // FIXME: Slow! | 
|---|
| 1436 |  | 
|---|
| 1437 | if (!used || (*used)[idx]) | 
|---|
| 1438 | SetUsed(idx); | 
|---|
| 1439 | } | 
|---|
| 1440 | fEntries++; | 
|---|
| 1441 | } | 
|---|
| 1442 |  | 
|---|
| 1443 | // ------------------------------------------------------------------------ | 
|---|
| 1444 | // | 
|---|
| 1445 | // Call this function to add a MCamEvent on top of the present contents. | 
|---|
| 1446 | // 1 is added to each pixel if the contents of MCamEvent>threshold (in case isabove is set to kTRUE == default) | 
|---|
| 1447 | // 1 is added to each pixel if the contents of MCamEvent<threshold (in case isabove is set to kFALSE) | 
|---|
| 1448 | // | 
|---|
| 1449 | // in unused pixel is not counted if it didn't fullfill the condition. | 
|---|
| 1450 | // | 
|---|
| 1451 | void MHCamera::CntCamContent(const MCamEvent &event, Double_t threshold, Int_t type, Bool_t isabove) | 
|---|
| 1452 | { | 
|---|
| 1453 | if (fNcells<=1 || IsFreezed()) | 
|---|
| 1454 | return; | 
|---|
| 1455 |  | 
|---|
| 1456 | // FIXME: Security check missing! | 
|---|
| 1457 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 1458 | { | 
|---|
| 1459 | Double_t val=threshold; | 
|---|
| 1460 | const Bool_t rc = event.GetPixelContent(val, idx, *fGeomCam, type); | 
|---|
| 1461 | if (rc) | 
|---|
| 1462 | SetUsed(idx); | 
|---|
| 1463 |  | 
|---|
| 1464 | const Bool_t cond = | 
|---|
| 1465 | ( isabove && val>threshold) || | 
|---|
| 1466 | (!isabove && val<threshold); | 
|---|
| 1467 |  | 
|---|
| 1468 | Fill(idx, rc && cond ? 1 : 0); | 
|---|
| 1469 | } | 
|---|
| 1470 | fEntries++; | 
|---|
| 1471 | } | 
|---|
| 1472 |  | 
|---|
| 1473 | // ------------------------------------------------------------------------ | 
|---|
| 1474 | // | 
|---|
| 1475 | // Call this function to add a MCamEvent on top of the present contents. | 
|---|
| 1476 | // - the contents of the pixels in event are added to each pixel | 
|---|
| 1477 | //   if the pixel of thresevt<threshold (in case isabove is set | 
|---|
| 1478 | //   to kTRUE == default) | 
|---|
| 1479 | // - the contents of the pixels in event are  added to each pixel | 
|---|
| 1480 | //   if the pixel of thresevt<threshold (in case isabove is set | 
|---|
| 1481 | //   to kFALSE) | 
|---|
| 1482 | // | 
|---|
| 1483 | // in unused pixel is not counted if it didn't fullfill the condition. | 
|---|
| 1484 | // | 
|---|
| 1485 | void MHCamera::CntCamContent(const MCamEvent &event, Int_t type1, const MCamEvent &thresevt, Int_t type2, Double_t threshold, Bool_t isabove) | 
|---|
| 1486 | { | 
|---|
| 1487 | if (fNcells<=1 || IsFreezed()) | 
|---|
| 1488 | return; | 
|---|
| 1489 |  | 
|---|
| 1490 | // FIXME: Security check missing! | 
|---|
| 1491 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 1492 | { | 
|---|
| 1493 | Double_t th=0; | 
|---|
| 1494 | if (!thresevt.GetPixelContent(th, idx, *fGeomCam, type2)) | 
|---|
| 1495 | continue; | 
|---|
| 1496 |  | 
|---|
| 1497 | if ((isabove && th>threshold) || (!isabove && th<threshold)) | 
|---|
| 1498 | continue; | 
|---|
| 1499 |  | 
|---|
| 1500 | Double_t val=th; | 
|---|
| 1501 | if (event.GetPixelContent(val, idx, *fGeomCam, type1)) | 
|---|
| 1502 | { | 
|---|
| 1503 | SetUsed(idx); | 
|---|
| 1504 | Fill(idx, val); | 
|---|
| 1505 | } | 
|---|
| 1506 | } | 
|---|
| 1507 | fEntries++; | 
|---|
| 1508 | } | 
|---|
| 1509 |  | 
|---|
| 1510 | // ------------------------------------------------------------------------ | 
|---|
| 1511 | // | 
|---|
| 1512 | // Call this function to add a MCamEvent on top of the present contents. | 
|---|
| 1513 | // 1 is added to each pixel if the contents of MCamEvent>threshold (in case isabove is set to kTRUE == default) | 
|---|
| 1514 | // 1 is added to each pixel if the contents of MCamEvent<threshold (in case isabove is set to kFALSE) | 
|---|
| 1515 | // | 
|---|
| 1516 | // in unused pixel is not counted if it didn't fullfill the condition. | 
|---|
| 1517 | // | 
|---|
| 1518 | void MHCamera::CntCamContent(const MCamEvent &event, TArrayD threshold, Int_t type, Bool_t isabove) | 
|---|
| 1519 | { | 
|---|
| 1520 | if (fNcells<=1 || IsFreezed()) | 
|---|
| 1521 | return; | 
|---|
| 1522 |  | 
|---|
| 1523 | // FIXME: Security check missing! | 
|---|
| 1524 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 1525 | { | 
|---|
| 1526 | Double_t val=threshold[idx]; | 
|---|
| 1527 | if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/) | 
|---|
| 1528 | { | 
|---|
| 1529 | SetUsed(idx); | 
|---|
| 1530 |  | 
|---|
| 1531 | if (val>threshold[idx] && isabove) | 
|---|
| 1532 | Fill(idx); | 
|---|
| 1533 | if (val<threshold[idx] && !isabove) | 
|---|
| 1534 | Fill(idx); | 
|---|
| 1535 | } | 
|---|
| 1536 | } | 
|---|
| 1537 | fEntries++; | 
|---|
| 1538 | } | 
|---|
| 1539 |  | 
|---|
| 1540 | // ------------------------------------------------------------------------ | 
|---|
| 1541 | // | 
|---|
| 1542 | // Call this function to add a TArrayD on top of the present contents. | 
|---|
| 1543 | // 1 is added to each pixel if the contents of MCamEvent>threshold | 
|---|
| 1544 | // | 
|---|
| 1545 | void MHCamera::CntCamContent(const TArrayD &event, Double_t threshold, Bool_t ispos) | 
|---|
| 1546 | { | 
|---|
| 1547 | if (event.GetSize()!=fNcells-2 || IsFreezed()) | 
|---|
| 1548 | return; | 
|---|
| 1549 |  | 
|---|
| 1550 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 1551 | { | 
|---|
| 1552 | if (event[idx]>threshold) | 
|---|
| 1553 | Fill(idx); | 
|---|
| 1554 |  | 
|---|
| 1555 | if (!ispos || fArray[idx+1]>0) | 
|---|
| 1556 | SetUsed(idx); | 
|---|
| 1557 | } | 
|---|
| 1558 | fEntries++; | 
|---|
| 1559 | } | 
|---|
| 1560 |  | 
|---|
| 1561 | // ------------------------------------------------------------------------ | 
|---|
| 1562 | // | 
|---|
| 1563 | // Fill the pixels with random contents. | 
|---|
| 1564 | // | 
|---|
| 1565 | void MHCamera::FillRandom() | 
|---|
| 1566 | { | 
|---|
| 1567 | if (fNcells<=1 || IsFreezed()) | 
|---|
| 1568 | return; | 
|---|
| 1569 |  | 
|---|
| 1570 | Reset(); | 
|---|
| 1571 |  | 
|---|
| 1572 | // FIXME: Security check missing! | 
|---|
| 1573 | for (Int_t idx=0; idx<fNcells-2; idx++) | 
|---|
| 1574 | { | 
|---|
| 1575 | Fill(idx, gRandom->Uniform()*fGeomCam->GetPixRatio(idx)); | 
|---|
| 1576 | SetUsed(idx); | 
|---|
| 1577 | } | 
|---|
| 1578 | fEntries=1; | 
|---|
| 1579 | } | 
|---|
| 1580 |  | 
|---|
| 1581 |  | 
|---|
| 1582 | // ------------------------------------------------------------------------ | 
|---|
| 1583 | // | 
|---|
| 1584 | // The array must be in increasing order, eg: 2.5, 3.7, 4.9 | 
|---|
| 1585 | // The values in each bin are replaced by the interval in which the value | 
|---|
| 1586 | // fits. In the example we have four intervals | 
|---|
| 1587 | // (<2.5, 2.5-3.7, 3.7-4.9, >4.9). Maximum and minimum are set | 
|---|
| 1588 | // accordingly. | 
|---|
| 1589 | // | 
|---|
| 1590 | void MHCamera::SetLevels(const TArrayF &arr) | 
|---|
| 1591 | { | 
|---|
| 1592 | if (fNcells<=1) | 
|---|
| 1593 | return; | 
|---|
| 1594 |  | 
|---|
| 1595 | for (Int_t i=0; i<fNcells-2; i++) | 
|---|
| 1596 | { | 
|---|
| 1597 | if (!IsUsed(i)) | 
|---|
| 1598 | continue; | 
|---|
| 1599 |  | 
|---|
| 1600 | Int_t j = arr.GetSize(); | 
|---|
| 1601 | while (j && fArray[i+1]<arr[j-1]) | 
|---|
| 1602 | j--; | 
|---|
| 1603 |  | 
|---|
| 1604 | fArray[i+1] = j; | 
|---|
| 1605 | } | 
|---|
| 1606 | SetMaximum(arr.GetSize()); | 
|---|
| 1607 | SetMinimum(0); | 
|---|
| 1608 | } | 
|---|
| 1609 |  | 
|---|
| 1610 | // ------------------------------------------------------------------------ | 
|---|
| 1611 | // | 
|---|
| 1612 | // Reset the all pixel colors to a default value | 
|---|
| 1613 | // | 
|---|
| 1614 | void MHCamera::Reset(Option_t *opt) | 
|---|
| 1615 | { | 
|---|
| 1616 | if (fNcells<=1 || IsFreezed()) | 
|---|
| 1617 | return; | 
|---|
| 1618 |  | 
|---|
| 1619 | TH1::Reset(opt); | 
|---|
| 1620 |  | 
|---|
| 1621 | for (Int_t i=0; i<fNcells-2; i++) | 
|---|
| 1622 | { | 
|---|
| 1623 | fArray[i+1]=0; | 
|---|
| 1624 | fBinEntries[i]=0; | 
|---|
| 1625 | ResetUsed(i); | 
|---|
| 1626 | } | 
|---|
| 1627 |  | 
|---|
| 1628 | fArray[0]         = 0; | 
|---|
| 1629 | fArray[fNcells-1] = 0; | 
|---|
| 1630 | } | 
|---|
| 1631 |  | 
|---|
| 1632 | // ------------------------------------------------------------------------ | 
|---|
| 1633 | // | 
|---|
| 1634 | //  Here we calculate the color index for the current value. | 
|---|
| 1635 | //  The color index is defined with the class TStyle and the | 
|---|
| 1636 | //  Color palette inside. We use the command gStyle->SetPalette(1,0) | 
|---|
| 1637 | //  for the display. So we have to convert the value "wert" into | 
|---|
| 1638 | //  a color index that fits the color palette. | 
|---|
| 1639 | //  The range of the color palette is defined by the values fMinPhe | 
|---|
| 1640 | //  and fMaxRange. Between this values we have 50 color index, starting | 
|---|
| 1641 | //  with 0 up to 49. | 
|---|
| 1642 | // | 
|---|
| 1643 | Int_t MHCamera::GetColor(Float_t val, Float_t min, Float_t max, Bool_t islog) | 
|---|
| 1644 | { | 
|---|
| 1645 | // | 
|---|
| 1646 | //   first treat the over- and under-flows | 
|---|
| 1647 | // | 
|---|
| 1648 | const Int_t maxcolidx = kItemsLegend-1; | 
|---|
| 1649 |  | 
|---|
| 1650 | if (!TMath::Finite(val)) // FIXME: gLog! | 
|---|
| 1651 | return maxcolidx/2; | 
|---|
| 1652 |  | 
|---|
| 1653 | if (val >= max) | 
|---|
| 1654 | return gStyle->GetColorPalette(maxcolidx); | 
|---|
| 1655 |  | 
|---|
| 1656 | if (val <= min) | 
|---|
| 1657 | return gStyle->GetColorPalette(0); | 
|---|
| 1658 |  | 
|---|
| 1659 | // | 
|---|
| 1660 | // calculate the color index | 
|---|
| 1661 | // | 
|---|
| 1662 | Float_t ratio; | 
|---|
| 1663 | if (islog && min>0) | 
|---|
| 1664 | ratio = log10(val/min) / log10(max/min); | 
|---|
| 1665 | else | 
|---|
| 1666 | ratio = (val-min) / (max-min); | 
|---|
| 1667 |  | 
|---|
| 1668 | const Int_t colidx = (Int_t)(ratio*maxcolidx + .5); | 
|---|
| 1669 | return gStyle->GetColorPalette(colidx); | 
|---|
| 1670 | } | 
|---|
| 1671 |  | 
|---|
| 1672 | TPaveStats *MHCamera::GetStatisticBox() | 
|---|
| 1673 | { | 
|---|
| 1674 | TObject *obj = 0; | 
|---|
| 1675 |  | 
|---|
| 1676 | TIter Next(fFunctions); | 
|---|
| 1677 | while ((obj = Next())) | 
|---|
| 1678 | if (obj->InheritsFrom(TPaveStats::Class())) | 
|---|
| 1679 | return static_cast<TPaveStats*>(obj); | 
|---|
| 1680 |  | 
|---|
| 1681 | return NULL; | 
|---|
| 1682 | } | 
|---|
| 1683 |  | 
|---|
| 1684 | // ------------------------------------------------------------------------ | 
|---|
| 1685 | // | 
|---|
| 1686 | //  Change the text on the legend according to the range of the Display | 
|---|
| 1687 | // | 
|---|
| 1688 | void MHCamera::UpdateLegend(Float_t min, Float_t max, Bool_t islog) | 
|---|
| 1689 | { | 
|---|
| 1690 | const Float_t range = fGeomCam->GetMaxRadius()*1.05; | 
|---|
| 1691 |  | 
|---|
| 1692 | if (!TestBit(kNoScale)) | 
|---|
| 1693 | { | 
|---|
| 1694 | TArrow arr; | 
|---|
| 1695 | arr.PaintArrow(-range*.9, -range*.9, -range*.6, -range*.9, 0.025); | 
|---|
| 1696 | arr.PaintArrow(-range*.9, -range*.9, -range*.9, -range*.6, 0.025); | 
|---|
| 1697 |  | 
|---|
| 1698 | TString text; | 
|---|
| 1699 | text += (int)(range*.3); | 
|---|
| 1700 | text += "mm"; | 
|---|
| 1701 |  | 
|---|
| 1702 | TText newtxt2; | 
|---|
| 1703 | newtxt2.SetTextSize(0.04); | 
|---|
| 1704 | newtxt2.PaintText(-range*.85, -range*.85, text); | 
|---|
| 1705 |  | 
|---|
| 1706 | text = ""; | 
|---|
| 1707 | text += Form("%.2f", (float)((int)(range*.3*fGeomCam->GetConvMm2Deg()*10))/10); | 
|---|
| 1708 | text += "\\circ"; | 
|---|
| 1709 | text = text.Strip(TString::kLeading); | 
|---|
| 1710 |  | 
|---|
| 1711 | TLatex latex; | 
|---|
| 1712 | latex.PaintLatex(-range*.85, -range*.75, 0, 0.04, text); | 
|---|
| 1713 | } | 
|---|
| 1714 |  | 
|---|
| 1715 | if (!TestBit(kNoLegend)) | 
|---|
| 1716 | { | 
|---|
| 1717 | TPaveStats *stats = GetStatisticBox(); | 
|---|
| 1718 |  | 
|---|
| 1719 | const Float_t hndc   = 0.92 - (stats ? stats->GetY1NDC() : 1); | 
|---|
| 1720 | const Float_t H      = (0.75-hndc)*range; | 
|---|
| 1721 | const Float_t offset = hndc*range; | 
|---|
| 1722 |  | 
|---|
| 1723 | const Float_t h = 2./kItemsLegend; | 
|---|
| 1724 | const Float_t w = range/sqrt((float)(fNcells-2)); | 
|---|
| 1725 |  | 
|---|
| 1726 | TBox newbox; | 
|---|
| 1727 | TText newtxt; | 
|---|
| 1728 | newtxt.SetTextSize(0.03); | 
|---|
| 1729 | newtxt.SetTextAlign(12); | 
|---|
| 1730 | #if ROOT_VERSION_CODE > ROOT_VERSION(3,01,06) | 
|---|
| 1731 | newtxt.SetBit(/*kNoContextMenu|*/kCannotPick); | 
|---|
| 1732 | newbox.SetBit(/*kNoContextMenu|*/kCannotPick); | 
|---|
| 1733 | #endif | 
|---|
| 1734 |  | 
|---|
| 1735 | const Float_t step   = (islog && min>0 ? log10(max/min) : max-min) / kItemsLegend; | 
|---|
| 1736 | const Int_t   firsts = step*3 < 1e-8 ? 8 : (Int_t)floor(log10(step*3)); | 
|---|
| 1737 | const TString opt    = Form("%%.%if", firsts>0 ? 0 : TMath::Abs(firsts)); | 
|---|
| 1738 |  | 
|---|
| 1739 | for (Int_t i=0; i<kItemsLegend+1; i+=3) | 
|---|
| 1740 | { | 
|---|
| 1741 | Float_t val; | 
|---|
| 1742 | if (islog && min>0) | 
|---|
| 1743 | val = pow(10, step*i) * min; | 
|---|
| 1744 | else | 
|---|
| 1745 | val = min + step*i; | 
|---|
| 1746 |  | 
|---|
| 1747 | //const bool dispexp = max-min>1.5 && fabs(val)>0.1 && fabs(val)<1e6; | 
|---|
| 1748 | newtxt.PaintText(range+1.5*w, H*(i*h-1)-offset, Form(opt, val)); | 
|---|
| 1749 | } | 
|---|
| 1750 |  | 
|---|
| 1751 | for (Int_t i=0; i<kItemsLegend; i++) | 
|---|
| 1752 | { | 
|---|
| 1753 | newbox.SetFillColor(gStyle->GetColorPalette(i)); | 
|---|
| 1754 | newbox.PaintBox(range, H*(i*h-1)-offset, range+w, H*((i+1)*h-1)-offset); | 
|---|
| 1755 | } | 
|---|
| 1756 | } | 
|---|
| 1757 | } | 
|---|
| 1758 |  | 
|---|
| 1759 | // ------------------------------------------------------------------------ | 
|---|
| 1760 | // | 
|---|
| 1761 | // Save primitive as a C++ statement(s) on output stream out | 
|---|
| 1762 | // | 
|---|
| 1763 | void MHCamera::SavePrimitive(ostream &out, Option_t *opt) | 
|---|
| 1764 | { | 
|---|
| 1765 | gLog << err << "MHCamera::SavePrimitive: Must be rewritten!" << endl; | 
|---|
| 1766 | /* | 
|---|
| 1767 | if (!gROOT->ClassSaved(TCanvas::Class())) | 
|---|
| 1768 | fDrawingPad->SavePrimitive(out, opt); | 
|---|
| 1769 |  | 
|---|
| 1770 | out << "   " << fDrawingPad->GetName() << "->SetWindowSize("; | 
|---|
| 1771 | out << fDrawingPad->GetWw() << "," << fDrawingPad->GetWh() << ");" << endl; | 
|---|
| 1772 | */ | 
|---|
| 1773 | } | 
|---|
| 1774 |  | 
|---|
| 1775 | void MHCamera::SavePrimitive(ofstream &out, Option_t *) | 
|---|
| 1776 | { | 
|---|
| 1777 | MHCamera::SavePrimitive(static_cast<ostream&>(out), ""); | 
|---|
| 1778 | } | 
|---|
| 1779 |  | 
|---|
| 1780 | // ------------------------------------------------------------------------ | 
|---|
| 1781 | // | 
|---|
| 1782 | // compute the distance of a point (px,py) to the Camera | 
|---|
| 1783 | // this functions needed for graphical primitives, that | 
|---|
| 1784 | // means without this function you are not able to interact | 
|---|
| 1785 | // with the graphical primitive with the mouse!!! | 
|---|
| 1786 | // | 
|---|
| 1787 | // All calcutations are done in pixel coordinates | 
|---|
| 1788 | // | 
|---|
| 1789 | Int_t MHCamera::DistancetoPrimitive(Int_t px, Int_t py) | 
|---|
| 1790 | { | 
|---|
| 1791 | if (fNcells<=1) | 
|---|
| 1792 | return 999999; | 
|---|
| 1793 |  | 
|---|
| 1794 | TPaveStats *box = (TPaveStats*)gPad->GetPrimitive("stats"); | 
|---|
| 1795 | if (box) | 
|---|
| 1796 | { | 
|---|
| 1797 | const Double_t w = box->GetY2NDC()-box->GetY1NDC(); | 
|---|
| 1798 | box->SetX1NDC(gStyle->GetStatX()-gStyle->GetStatW()); | 
|---|
| 1799 | box->SetY1NDC(gStyle->GetStatY()-w); | 
|---|
| 1800 | box->SetX2NDC(gStyle->GetStatX()); | 
|---|
| 1801 | box->SetY2NDC(gStyle->GetStatY()); | 
|---|
| 1802 | } | 
|---|
| 1803 |  | 
|---|
| 1804 | if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase)) | 
|---|
| 1805 | return TH1D::DistancetoPrimitive(px, py); | 
|---|
| 1806 |  | 
|---|
| 1807 | const Bool_t issame = TString(GetDrawOption()).Contains("same", TString::kIgnoreCase); | 
|---|
| 1808 |  | 
|---|
| 1809 | const Float_t maxr = (1-fGeomCam->GetConvMm2Deg())*fGeomCam->GetMaxRadius()/2; | 
|---|
| 1810 | const Float_t conv = !issame || | 
|---|
| 1811 | gPad->GetX1()<-maxr || gPad->GetY1()<-maxr || | 
|---|
| 1812 | gPad->GetX2()> maxr || gPad->GetY2()>maxr ? 1 : fGeomCam->GetConvMm2Deg(); | 
|---|
| 1813 |  | 
|---|
| 1814 | if (GetPixelIndex(px, py, conv)>=0) | 
|---|
| 1815 | return 0; | 
|---|
| 1816 |  | 
|---|
| 1817 | if (!box) | 
|---|
| 1818 | return 999999; | 
|---|
| 1819 |  | 
|---|
| 1820 | const Int_t dist = box->DistancetoPrimitive(px, py); | 
|---|
| 1821 | if (dist > TPad::GetMaxPickDistance()) | 
|---|
| 1822 | return 999999; | 
|---|
| 1823 |  | 
|---|
| 1824 | gPad->SetSelected(box); | 
|---|
| 1825 | return dist; | 
|---|
| 1826 | } | 
|---|
| 1827 |  | 
|---|
| 1828 | // ------------------------------------------------------------------------ | 
|---|
| 1829 | // | 
|---|
| 1830 | // | 
|---|
| 1831 | Int_t MHCamera::GetPixelIndex(Int_t px, Int_t py, Float_t conv) const | 
|---|
| 1832 | { | 
|---|
| 1833 | if (fNcells<=1) | 
|---|
| 1834 | return -1; | 
|---|
| 1835 |  | 
|---|
| 1836 | Int_t i; | 
|---|
| 1837 | for (i=0; i<fNcells-2; i++) | 
|---|
| 1838 | { | 
|---|
| 1839 | MHexagon hex((*fGeomCam)[i]); | 
|---|
| 1840 | if (hex.DistancetoPrimitive(px, py, conv)>0) | 
|---|
| 1841 | continue; | 
|---|
| 1842 |  | 
|---|
| 1843 | return i; | 
|---|
| 1844 | } | 
|---|
| 1845 | return -1; | 
|---|
| 1846 | } | 
|---|
| 1847 |  | 
|---|
| 1848 | // ------------------------------------------------------------------------ | 
|---|
| 1849 | // | 
|---|
| 1850 | // Returns string containing info about the object at position (px,py). | 
|---|
| 1851 | // Returned string will be re-used (lock in MT environment). | 
|---|
| 1852 | // | 
|---|
| 1853 | char *MHCamera::GetObjectInfo(Int_t px, Int_t py) const | 
|---|
| 1854 | { | 
|---|
| 1855 | if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase)) | 
|---|
| 1856 | return TH1D::GetObjectInfo(px, py); | 
|---|
| 1857 |  | 
|---|
| 1858 | static char info[128]; | 
|---|
| 1859 |  | 
|---|
| 1860 | const Int_t idx=GetPixelIndex(px, py); | 
|---|
| 1861 |  | 
|---|
| 1862 | if (idx<0) | 
|---|
| 1863 | return TObject::GetObjectInfo(px, py); | 
|---|
| 1864 |  | 
|---|
| 1865 | sprintf(info, "Software Pixel Idx: %d (Hardware Id=%d) c=%.1f <%s>", | 
|---|
| 1866 | idx, idx+1, GetBinContent(idx+1), IsUsed(idx)?"on":"off"); | 
|---|
| 1867 | return info; | 
|---|
| 1868 | } | 
|---|
| 1869 |  | 
|---|
| 1870 | // ------------------------------------------------------------------------ | 
|---|
| 1871 | // | 
|---|
| 1872 | // Add a MCamEvent which should be displayed when the user clicks on a | 
|---|
| 1873 | // pixel. | 
|---|
| 1874 | // Warning: The object MUST inherit from TObject AND MCamEvent | 
|---|
| 1875 | // | 
|---|
| 1876 | void MHCamera::AddNotify(TObject *obj) | 
|---|
| 1877 | { | 
|---|
| 1878 | // Make sure, that the object derives from MCamEvent! | 
|---|
| 1879 | MCamEvent *evt = dynamic_cast<MCamEvent*>(obj); | 
|---|
| 1880 | if (!evt) | 
|---|
| 1881 | { | 
|---|
| 1882 | gLog << err << "ERROR: MHCamera::AddNotify - TObject doesn't inherit from MCamEvent... ignored." << endl; | 
|---|
| 1883 | return; | 
|---|
| 1884 | } | 
|---|
| 1885 |  | 
|---|
| 1886 | // Make sure, that it is deleted from the list too, if the obj is deleted | 
|---|
| 1887 | obj->SetBit(kMustCleanup); | 
|---|
| 1888 |  | 
|---|
| 1889 | // Add object to list | 
|---|
| 1890 | fNotify->Add(obj); | 
|---|
| 1891 | } | 
|---|
| 1892 |  | 
|---|
| 1893 | // ------------------------------------------------------------------------ | 
|---|
| 1894 | // | 
|---|
| 1895 | // Execute a mouse event on the camera | 
|---|
| 1896 | // | 
|---|
| 1897 | void MHCamera::ExecuteEvent(Int_t event, Int_t px, Int_t py) | 
|---|
| 1898 | { | 
|---|
| 1899 | if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase)) | 
|---|
| 1900 | { | 
|---|
| 1901 | TH1D::ExecuteEvent(event, px, py); | 
|---|
| 1902 | return; | 
|---|
| 1903 | } | 
|---|
| 1904 | //if (event==kMouseMotion && fStatusBar) | 
|---|
| 1905 | //    fStatusBar->SetText(GetObjectInfo(px, py), 0); | 
|---|
| 1906 | if (event!=kButton1Down) | 
|---|
| 1907 | return; | 
|---|
| 1908 |  | 
|---|
| 1909 | const Int_t idx = GetPixelIndex(px, py); | 
|---|
| 1910 | if (idx<0) | 
|---|
| 1911 | return; | 
|---|
| 1912 |  | 
|---|
| 1913 | gLog << all << GetTitle() << " <" << GetName() << ">" << dec << endl; | 
|---|
| 1914 | gLog << "Software Pixel Idx: " << idx << endl; | 
|---|
| 1915 | gLog << "Hardware Pixel Id:  " << idx+1 << endl; | 
|---|
| 1916 | gLog << "Contents:           " << GetBinContent(idx+1); | 
|---|
| 1917 | if (GetBinError(idx+1)>0) | 
|---|
| 1918 | gLog << " +/- " << GetBinError(idx+1); | 
|---|
| 1919 | gLog << "  <" << (IsUsed(idx)?"on":"off") << ">  n=" << fBinEntries[idx+1] << endl; | 
|---|
| 1920 |  | 
|---|
| 1921 | if (fNotify && fNotify->GetSize()>0) | 
|---|
| 1922 | { | 
|---|
| 1923 | // FIXME: Is there a simpler and more convinient way? | 
|---|
| 1924 |  | 
|---|
| 1925 | // The name which is created here depends on the instance of | 
|---|
| 1926 | // MHCamera and on the pad on which it is drawn --> The name | 
|---|
| 1927 | // is unique. For ExecuteEvent gPad is always correctly set. | 
|---|
| 1928 | const TString name = Form("%p;%p;PixelContent", this, gPad); | 
|---|
| 1929 |  | 
|---|
| 1930 | TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name); | 
|---|
| 1931 | if (old) | 
|---|
| 1932 | old->cd(); | 
|---|
| 1933 | else | 
|---|
| 1934 | new TCanvas(name); | 
|---|
| 1935 |  | 
|---|
| 1936 | /* | 
|---|
| 1937 | TIter Next(gPad->GetListOfPrimitives()); | 
|---|
| 1938 | TObject *o; | 
|---|
| 1939 | while (o=Next()) cout << o << ": " << o->GetName() << " " << o->IsA()->GetName() << endl; | 
|---|
| 1940 | */ | 
|---|
| 1941 |  | 
|---|
| 1942 | // FIXME: Make sure, that the old histograms are really deleted. | 
|---|
| 1943 | //        Are they already deleted? | 
|---|
| 1944 |  | 
|---|
| 1945 | // The dynamic_cast is necessary here: We cannot use ForEach | 
|---|
| 1946 | TIter Next(fNotify); | 
|---|
| 1947 | MCamEvent *evt; | 
|---|
| 1948 | while ((evt=dynamic_cast<MCamEvent*>(Next()))) | 
|---|
| 1949 | evt->DrawPixelContent(idx); | 
|---|
| 1950 |  | 
|---|
| 1951 | gPad->Modified(); | 
|---|
| 1952 | gPad->Update(); | 
|---|
| 1953 | } | 
|---|
| 1954 | } | 
|---|
| 1955 |  | 
|---|
| 1956 | UInt_t MHCamera::GetNumPixels() const | 
|---|
| 1957 | { | 
|---|
| 1958 | return fGeomCam ? fGeomCam->GetNumPixels() : 0; | 
|---|
| 1959 | } | 
|---|
| 1960 |  | 
|---|
| 1961 | TH1 *MHCamera::DrawCopy() const | 
|---|
| 1962 | { | 
|---|
| 1963 | gPad=NULL; | 
|---|
| 1964 | return TH1D::DrawCopy(fName+";cpy"); | 
|---|
| 1965 | } | 
|---|
| 1966 |  | 
|---|
| 1967 | // -------------------------------------------------------------------------- | 
|---|
| 1968 | // | 
|---|
| 1969 | // Draw a projection of MHCamera onto the y-axis values. Depending on the | 
|---|
| 1970 | // variable fit, the following fits are performed: | 
|---|
| 1971 | // | 
|---|
| 1972 | // 0: No fit, simply draw the projection | 
|---|
| 1973 | // 1: Single Gauss (for distributions flat-fielded over the whole camera) | 
|---|
| 1974 | // 2: Double Gauss (for distributions different for inner and outer pixels) | 
|---|
| 1975 | // 3: Triple Gauss (for distributions with inner, outer pixels and outliers) | 
|---|
| 1976 | // 4: flat         (for the probability distributions) | 
|---|
| 1977 | // (1-4:) Moreover, sectors 6,1 and 2 of the camera and sectors 3,4 and 5 are | 
|---|
| 1978 | //        drawn separately, for inner and outer pixels. | 
|---|
| 1979 | // 5: Fit Inner and Outer pixels separately by a single Gaussian | 
|---|
| 1980 | //                 (only for MAGIC cameras) | 
|---|
| 1981 | // 6: Fit Inner and Outer pixels separately by a single Gaussian and display | 
|---|
| 1982 | //                 additionally the two camera halfs separately (for MAGIC camera) | 
|---|
| 1983 | // 7: Single Gauss with TLegend to show the meaning of the colours | 
|---|
| 1984 | // | 
|---|
| 1985 | void MHCamera::DrawProjection(Int_t fit) const | 
|---|
| 1986 | { | 
|---|
| 1987 | TArrayI inner(1); | 
|---|
| 1988 | inner[0] = 0; | 
|---|
| 1989 |  | 
|---|
| 1990 | TArrayI outer(1); | 
|---|
| 1991 | outer[0] = 1; | 
|---|
| 1992 |  | 
|---|
| 1993 | if (fit==5 || fit==6) | 
|---|
| 1994 | { | 
|---|
| 1995 | if (GetGeomCam().InheritsFrom("MGeomCamMagic")) | 
|---|
| 1996 | { | 
|---|
| 1997 | TArrayI s0(6); | 
|---|
| 1998 | s0[0] = 6; | 
|---|
| 1999 | s0[1] = 1; | 
|---|
| 2000 | s0[2] = 2; | 
|---|
| 2001 | s0[3] = 3; | 
|---|
| 2002 | s0[4] = 4; | 
|---|
| 2003 | s0[5] = 5; | 
|---|
| 2004 |  | 
|---|
| 2005 | TArrayI s1(3); | 
|---|
| 2006 | s1[0] = 6; | 
|---|
| 2007 | s1[1] = 1; | 
|---|
| 2008 | s1[2] = 2; | 
|---|
| 2009 |  | 
|---|
| 2010 | TArrayI s2(3); | 
|---|
| 2011 | s2[0] = 3; | 
|---|
| 2012 | s2[1] = 4; | 
|---|
| 2013 | s2[2] = 5; | 
|---|
| 2014 |  | 
|---|
| 2015 | gPad->Clear(); | 
|---|
| 2016 | TVirtualPad *pad = gPad; | 
|---|
| 2017 | pad->Divide(2,1); | 
|---|
| 2018 |  | 
|---|
| 2019 | TH1D *inout[2]; | 
|---|
| 2020 | inout[0] = ProjectionS(s0, inner, "Inner"); | 
|---|
| 2021 | inout[1] = ProjectionS(s0, outer, "Outer"); | 
|---|
| 2022 |  | 
|---|
| 2023 | inout[0]->SetDirectory(NULL); | 
|---|
| 2024 | inout[1]->SetDirectory(NULL); | 
|---|
| 2025 |  | 
|---|
| 2026 | for (int i=0; i<2; i++) | 
|---|
| 2027 | { | 
|---|
| 2028 | pad->cd(i+1); | 
|---|
| 2029 | gPad->SetBorderMode(0); | 
|---|
| 2030 |  | 
|---|
| 2031 | inout[i]->SetLineColor(kRed+i); | 
|---|
| 2032 | inout[i]->SetBit(kCanDelete); | 
|---|
| 2033 | inout[i]->Draw(); | 
|---|
| 2034 | inout[i]->Fit("gaus","Q"); | 
|---|
| 2035 |  | 
|---|
| 2036 | if (fit == 6) | 
|---|
| 2037 | { | 
|---|
| 2038 | TH1D *half[2]; | 
|---|
| 2039 | half[0] = ProjectionS(s1, i==0 ? inner : outer , "Sector 6-1-2"); | 
|---|
| 2040 | half[1] = ProjectionS(s2, i==0 ? inner : outer , "Sector 3-4-5"); | 
|---|
| 2041 |  | 
|---|
| 2042 | for (int j=0; j<2; j++) | 
|---|
| 2043 | { | 
|---|
| 2044 | half[j]->SetLineColor(kRed+i+2*j+1); | 
|---|
| 2045 | half[j]->SetDirectory(NULL); | 
|---|
| 2046 | half[j]->SetBit(kCanDelete); | 
|---|
| 2047 | half[j]->Draw("same"); | 
|---|
| 2048 | } | 
|---|
| 2049 | } | 
|---|
| 2050 |  | 
|---|
| 2051 | } | 
|---|
| 2052 | } | 
|---|
| 2053 | return; | 
|---|
| 2054 | } | 
|---|
| 2055 |  | 
|---|
| 2056 | TH1D *obj2 = (TH1D*)Projection(GetName()); | 
|---|
| 2057 | obj2->SetDirectory(0); | 
|---|
| 2058 | obj2->Draw(); | 
|---|
| 2059 | obj2->SetBit(kCanDelete); | 
|---|
| 2060 |  | 
|---|
| 2061 | if (fit == 0) | 
|---|
| 2062 | return; | 
|---|
| 2063 |  | 
|---|
| 2064 | if (GetGeomCam().InheritsFrom("MGeomCamMagic")) | 
|---|
| 2065 | { | 
|---|
| 2066 | TArrayI s0(3); | 
|---|
| 2067 | s0[0] = 6; | 
|---|
| 2068 | s0[1] = 1; | 
|---|
| 2069 | s0[2] = 2; | 
|---|
| 2070 |  | 
|---|
| 2071 | TArrayI s1(3); | 
|---|
| 2072 | s1[0] = 3; | 
|---|
| 2073 | s1[1] = 4; | 
|---|
| 2074 | s1[2] = 5; | 
|---|
| 2075 |  | 
|---|
| 2076 | TH1D *halfInOut[4]; | 
|---|
| 2077 |  | 
|---|
| 2078 | // Just to get the right (maximum) binning | 
|---|
| 2079 | halfInOut[0] = ProjectionS(s0, inner, "Sector 6-1-2 Inner"); | 
|---|
| 2080 | halfInOut[1] = ProjectionS(s1, inner, "Sector 3-4-5 Inner"); | 
|---|
| 2081 | halfInOut[2] = ProjectionS(s0, outer, "Sector 6-1-2 Outer"); | 
|---|
| 2082 | halfInOut[3] = ProjectionS(s1, outer, "Sector 3-4-5 Outer"); | 
|---|
| 2083 |  | 
|---|
| 2084 | TLegend *leg = new TLegend(0.05,0.65,0.35,0.9); | 
|---|
| 2085 |  | 
|---|
| 2086 | for (int i=0; i<4; i++) | 
|---|
| 2087 | { | 
|---|
| 2088 | halfInOut[i]->SetLineColor(kRed+i); | 
|---|
| 2089 | halfInOut[i]->SetDirectory(0); | 
|---|
| 2090 | halfInOut[i]->SetBit(kCanDelete); | 
|---|
| 2091 | halfInOut[i]->Draw("same"); | 
|---|
| 2092 | leg->AddEntry(halfInOut[i],halfInOut[i]->GetTitle(),"l"); | 
|---|
| 2093 | } | 
|---|
| 2094 |  | 
|---|
| 2095 | if (fit==7) | 
|---|
| 2096 | leg->Draw(); | 
|---|
| 2097 |  | 
|---|
| 2098 | gPad->Modified(); | 
|---|
| 2099 | gPad->Update(); | 
|---|
| 2100 | } | 
|---|
| 2101 |  | 
|---|
| 2102 | const Double_t min   = obj2->GetBinCenter(obj2->GetXaxis()->GetFirst()); | 
|---|
| 2103 | const Double_t max   = obj2->GetBinCenter(obj2->GetXaxis()->GetLast()); | 
|---|
| 2104 | const Double_t integ = obj2->Integral("width")/2.5; | 
|---|
| 2105 | const Double_t mean  = obj2->GetMean(); | 
|---|
| 2106 | const Double_t rms   = obj2->GetRMS(); | 
|---|
| 2107 | const Double_t width = max-min; | 
|---|
| 2108 |  | 
|---|
| 2109 | const TString dgausformula = "([0]-[3])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])" | 
|---|
| 2110 | "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])"; | 
|---|
| 2111 |  | 
|---|
| 2112 | const TString tgausformula = "([0]-[3]-[6])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])" | 
|---|
| 2113 | "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])" | 
|---|
| 2114 | "+[6]/[8]*exp(-0.5*(x-[7])*(x-[7])/[8]/[8])"; | 
|---|
| 2115 |  | 
|---|
| 2116 | TF1 *f=0; | 
|---|
| 2117 | switch (fit) | 
|---|
| 2118 | { | 
|---|
| 2119 | case 1: | 
|---|
| 2120 | f = new TF1("sgaus", "gaus(0)", min, max); | 
|---|
| 2121 | f->SetLineColor(kYellow); | 
|---|
| 2122 | f->SetBit(kCanDelete); | 
|---|
| 2123 | f->SetParNames("Area", "#mu", "#sigma"); | 
|---|
| 2124 | f->SetParameters(integ/rms, mean, rms); | 
|---|
| 2125 | f->SetParLimits(0, 0,   integ); | 
|---|
| 2126 | f->SetParLimits(1, min, max); | 
|---|
| 2127 | f->SetParLimits(2, 0,   width/1.5); | 
|---|
| 2128 |  | 
|---|
| 2129 | obj2->Fit(f, "QLR"); | 
|---|
| 2130 | break; | 
|---|
| 2131 |  | 
|---|
| 2132 | case 2: | 
|---|
| 2133 | f = new TF1("dgaus",dgausformula.Data(),min,max); | 
|---|
| 2134 | f->SetLineColor(kYellow); | 
|---|
| 2135 | f->SetBit(kCanDelete); | 
|---|
| 2136 | f->SetParNames("A_{tot}", "#mu1", "#sigma1", "A2", "#mu2", "#sigma2"); | 
|---|
| 2137 | f->SetParameters(integ,(min+mean)/2.,width/4., | 
|---|
| 2138 | integ/width/2.,(max+mean)/2.,width/4.); | 
|---|
| 2139 | // The left-sided Gauss | 
|---|
| 2140 | f->SetParLimits(0,integ-1.5      , integ+1.5); | 
|---|
| 2141 | f->SetParLimits(1,min+(width/10.), mean); | 
|---|
| 2142 | f->SetParLimits(2,0              , width/2.); | 
|---|
| 2143 | // The right-sided Gauss | 
|---|
| 2144 | f->SetParLimits(3,0   , integ); | 
|---|
| 2145 | f->SetParLimits(4,mean, max-(width/10.)); | 
|---|
| 2146 | f->SetParLimits(5,0   , width/2.); | 
|---|
| 2147 | obj2->Fit(f,"QLRM"); | 
|---|
| 2148 | break; | 
|---|
| 2149 |  | 
|---|
| 2150 | case 3: | 
|---|
| 2151 | f = new TF1("tgaus",tgausformula.Data(),min,max); | 
|---|
| 2152 | f->SetLineColor(kYellow); | 
|---|
| 2153 | f->SetBit(kCanDelete); | 
|---|
| 2154 | f->SetParNames("A_{tot}","#mu_{1}","#sigma_{1}", | 
|---|
| 2155 | "A_{2}","#mu_{2}","#sigma_{2}", | 
|---|
| 2156 | "A_{3}","#mu_{3}","#sigma_{3}"); | 
|---|
| 2157 | f->SetParameters(integ,(min+mean)/2,width/4., | 
|---|
| 2158 | integ/width/3.,(max+mean)/2.,width/4., | 
|---|
| 2159 | integ/width/3.,mean,width/2.); | 
|---|
| 2160 | // The left-sided Gauss | 
|---|
| 2161 | f->SetParLimits(0,integ-1.5,integ+1.5); | 
|---|
| 2162 | f->SetParLimits(1,min+(width/10.),mean); | 
|---|
| 2163 | f->SetParLimits(2,width/15.,width/2.); | 
|---|
| 2164 | // The right-sided Gauss | 
|---|
| 2165 | f->SetParLimits(3,0.,integ); | 
|---|
| 2166 | f->SetParLimits(4,mean,max-(width/10.)); | 
|---|
| 2167 | f->SetParLimits(5,width/15.,width/2.); | 
|---|
| 2168 | // The Gauss describing the outliers | 
|---|
| 2169 | f->SetParLimits(6,0.,integ); | 
|---|
| 2170 | f->SetParLimits(7,min,max); | 
|---|
| 2171 | f->SetParLimits(8,width/4.,width/1.5); | 
|---|
| 2172 | obj2->Fit(f,"QLRM"); | 
|---|
| 2173 | break; | 
|---|
| 2174 |  | 
|---|
| 2175 | case 4: | 
|---|
| 2176 | obj2->Fit("pol0", "Q"); | 
|---|
| 2177 | obj2->GetFunction("pol0")->SetLineColor(kYellow); | 
|---|
| 2178 | break; | 
|---|
| 2179 |  | 
|---|
| 2180 | case 9: | 
|---|
| 2181 | break; | 
|---|
| 2182 |  | 
|---|
| 2183 | default: | 
|---|
| 2184 | obj2->Fit("gaus", "Q"); | 
|---|
| 2185 | obj2->GetFunction("gaus")->SetLineColor(kYellow); | 
|---|
| 2186 | break; | 
|---|
| 2187 | } | 
|---|
| 2188 | } | 
|---|
| 2189 |  | 
|---|
| 2190 | // -------------------------------------------------------------------------- | 
|---|
| 2191 | // | 
|---|
| 2192 | // Draw a projection of MHCamera vs. the radius from the central pixel. | 
|---|
| 2193 | // | 
|---|
| 2194 | // The inner and outer pixels are drawn separately, both fitted by a polynomial | 
|---|
| 2195 | // of grade 1. | 
|---|
| 2196 | // | 
|---|
| 2197 | void MHCamera::DrawRadialProfile() const | 
|---|
| 2198 | { | 
|---|
| 2199 | TProfile *obj2 = (TProfile*)RadialProfile(GetName()); | 
|---|
| 2200 | obj2->SetDirectory(0); | 
|---|
| 2201 | obj2->Draw(); | 
|---|
| 2202 | obj2->SetBit(kCanDelete); | 
|---|
| 2203 |  | 
|---|
| 2204 | if (GetGeomCam().InheritsFrom("MGeomCamMagic")) | 
|---|
| 2205 | { | 
|---|
| 2206 | TArrayI s0(6); | 
|---|
| 2207 | s0[0] = 1; | 
|---|
| 2208 | s0[1] = 2; | 
|---|
| 2209 | s0[2] = 3; | 
|---|
| 2210 | s0[3] = 4; | 
|---|
| 2211 | s0[4] = 5; | 
|---|
| 2212 | s0[5] = 6; | 
|---|
| 2213 |  | 
|---|
| 2214 | TArrayI inner(1); | 
|---|
| 2215 | inner[0] = 0; | 
|---|
| 2216 |  | 
|---|
| 2217 | TArrayI outer(1); | 
|---|
| 2218 | outer[0] = 1; | 
|---|
| 2219 |  | 
|---|
| 2220 | // Just to get the right (maximum) binning | 
|---|
| 2221 | TProfile *half[2]; | 
|---|
| 2222 | half[0] = RadialProfileS(s0, inner,Form("%sInner",GetName())); | 
|---|
| 2223 | half[1] = RadialProfileS(s0, outer,Form("%sOuter",GetName())); | 
|---|
| 2224 |  | 
|---|
| 2225 | for (Int_t i=0; i<2; i++) | 
|---|
| 2226 | { | 
|---|
| 2227 | Double_t min = GetGeomCam().GetMinRadius(i); | 
|---|
| 2228 | Double_t max = GetGeomCam().GetMaxRadius(i); | 
|---|
| 2229 |  | 
|---|
| 2230 | half[i]->SetLineColor(kRed+i); | 
|---|
| 2231 | half[i]->SetDirectory(0); | 
|---|
| 2232 | half[i]->SetBit(kCanDelete); | 
|---|
| 2233 | half[i]->Draw("same"); | 
|---|
| 2234 | half[i]->Fit("pol1","Q","",min,max); | 
|---|
| 2235 | half[i]->GetFunction("pol1")->SetLineColor(kRed+i); | 
|---|
| 2236 | half[i]->GetFunction("pol1")->SetLineWidth(1); | 
|---|
| 2237 | } | 
|---|
| 2238 | } | 
|---|
| 2239 | } | 
|---|
| 2240 |  | 
|---|
| 2241 | // -------------------------------------------------------------------------- | 
|---|
| 2242 | // | 
|---|
| 2243 | // Draw a projection of MHCamera vs. the azimuth angle inside the camera. | 
|---|
| 2244 | // | 
|---|
| 2245 | // The inner and outer pixels are drawn separately. | 
|---|
| 2246 | // The general azimuth profile is fitted by a straight line | 
|---|
| 2247 | // | 
|---|
| 2248 | void MHCamera::DrawAzimuthProfile() const | 
|---|
| 2249 | { | 
|---|
| 2250 | TProfile *obj2 = (TProfile*)AzimuthProfile(GetName()); | 
|---|
| 2251 | obj2->SetDirectory(0); | 
|---|
| 2252 | obj2->Draw(); | 
|---|
| 2253 | obj2->SetBit(kCanDelete); | 
|---|
| 2254 | obj2->Fit("pol0","Q",""); | 
|---|
| 2255 | obj2->GetFunction("pol0")->SetLineWidth(1); | 
|---|
| 2256 |  | 
|---|
| 2257 | if (GetGeomCam().InheritsFrom("MGeomCamMagic")) | 
|---|
| 2258 | { | 
|---|
| 2259 | TArrayI inner(1); | 
|---|
| 2260 | inner[0] = 0; | 
|---|
| 2261 |  | 
|---|
| 2262 | TArrayI outer(1); | 
|---|
| 2263 | outer[0] = 1; | 
|---|
| 2264 |  | 
|---|
| 2265 | // Just to get the right (maximum) binning | 
|---|
| 2266 | TProfile *half[2]; | 
|---|
| 2267 | half[0] = AzimuthProfileA(inner,Form("%sInner",GetName())); | 
|---|
| 2268 | half[1] = AzimuthProfileA(outer,Form("%sOuter",GetName())); | 
|---|
| 2269 |  | 
|---|
| 2270 | for (Int_t i=0; i<2; i++) | 
|---|
| 2271 | { | 
|---|
| 2272 | half[i]->SetLineColor(kRed+i); | 
|---|
| 2273 | half[i]->SetDirectory(0); | 
|---|
| 2274 | half[i]->SetBit(kCanDelete); | 
|---|
| 2275 | half[i]->SetMarkerSize(0.5); | 
|---|
| 2276 | half[i]->Draw("same"); | 
|---|
| 2277 | } | 
|---|
| 2278 | } | 
|---|
| 2279 | } | 
|---|
| 2280 |  | 
|---|
| 2281 | // -------------------------------------------------------------------------- | 
|---|
| 2282 | // | 
|---|
| 2283 | // Draw the MHCamera into the MStatusDisplay: | 
|---|
| 2284 | // | 
|---|
| 2285 | // 1) Draw it as histogram (MHCamera::DrawCopy("hist") | 
|---|
| 2286 | // 2) Draw it as a camera, with MHCamera::SetPrettyPalette() set. | 
|---|
| 2287 | // 3) If "rad" is not zero, draw its values vs. the radius from the camera center. | 
|---|
| 2288 | //    (DrawRadialProfile()) | 
|---|
| 2289 | // 4) Depending on the variable "fit", draw the values projection on the y-axis | 
|---|
| 2290 | //    (DrawProjection()): | 
|---|
| 2291 | //    0: don't draw | 
|---|
| 2292 | //    1: Draw fit to Single Gauss (for distributions flat-fielded over the whole camera) | 
|---|
| 2293 | //    2: Draw and fit to Double Gauss (for distributions different for inner and outer pixels) | 
|---|
| 2294 | //    3: Draw and fit to Triple Gauss (for distributions with inner, outer pixels and outliers) | 
|---|
| 2295 | //    4: Draw and fit to Polynomial grade 0: (for the probability distributions) | 
|---|
| 2296 | //    >4: Draw and don;t fit. | 
|---|
| 2297 | // | 
|---|
| 2298 | void MHCamera::CamDraw(TCanvas &c, const Int_t x, const Int_t y, | 
|---|
| 2299 | const Int_t fit, const Int_t rad, const Int_t azi, | 
|---|
| 2300 | TObject *notify) | 
|---|
| 2301 | { | 
|---|
| 2302 | c.cd(x); | 
|---|
| 2303 | gPad->SetBorderMode(0); | 
|---|
| 2304 | gPad->SetTicks(); | 
|---|
| 2305 | MHCamera *obj1=(MHCamera*)DrawCopy("hist"); | 
|---|
| 2306 | obj1->SetDirectory(NULL); | 
|---|
| 2307 |  | 
|---|
| 2308 | if (notify) | 
|---|
| 2309 | obj1->AddNotify(notify); | 
|---|
| 2310 |  | 
|---|
| 2311 | c.cd(x+y); | 
|---|
| 2312 | gPad->SetBorderMode(0); | 
|---|
| 2313 | obj1->SetPrettyPalette(); | 
|---|
| 2314 | obj1->Draw(); | 
|---|
| 2315 |  | 
|---|
| 2316 | Int_t cnt = 2; | 
|---|
| 2317 |  | 
|---|
| 2318 | if (rad) | 
|---|
| 2319 | { | 
|---|
| 2320 | c.cd(x+2*y); | 
|---|
| 2321 | gPad->SetBorderMode(0); | 
|---|
| 2322 | gPad->SetTicks(); | 
|---|
| 2323 | DrawRadialProfile(); | 
|---|
| 2324 | cnt++; | 
|---|
| 2325 | } | 
|---|
| 2326 |  | 
|---|
| 2327 | if (azi) | 
|---|
| 2328 | { | 
|---|
| 2329 | c.cd(x+cnt*y); | 
|---|
| 2330 | gPad->SetBorderMode(0); | 
|---|
| 2331 | gPad->SetTicks(); | 
|---|
| 2332 | DrawAzimuthProfile(); | 
|---|
| 2333 | cnt++; | 
|---|
| 2334 | } | 
|---|
| 2335 |  | 
|---|
| 2336 | if (!fit) | 
|---|
| 2337 | return; | 
|---|
| 2338 |  | 
|---|
| 2339 | c.cd(x + cnt*y); | 
|---|
| 2340 | gPad->SetBorderMode(0); | 
|---|
| 2341 | gPad->SetTicks(); | 
|---|
| 2342 | DrawProjection(fit); | 
|---|
| 2343 | } | 
|---|