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

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