source: trunk/Mars/mhist/MHCamera.cc@ 19994

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