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

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