source: trunk/MagicSoft/Mars/mhist/MHCamera.cc@ 9241

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