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