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

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