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

Last change on this file since 13142 was 13001, checked in by tbretz, 13 years ago
Added AddCamDifference and weights to some Add-functions.
File size: 73.8 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::AddMeanShift(const MCamEvent &event, Int_t type, Stat_t w)
1362{
1363 if (fNcells<=1 || IsFreezed())
1364 return;
1365
1366 const Double_t mean = event.GetCameraMean(*fGeomCam, type);
1367
1368 // FIXME: Security check missing!
1369 for (Int_t idx=0; idx<fNcells-2; idx++)
1370 {
1371 Double_t val=0;
1372 if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
1373 {
1374 SetUsed(idx);
1375 Fill(idx, (val-mean)*w); // FIXME: Slow!
1376 }
1377 }
1378 fEntries++;
1379}
1380
1381// ------------------------------------------------------------------------
1382//
1383// Call this function to add a MCamEvent on top of the present contents.
1384//
1385void MHCamera::AddMedianShift(const MCamEvent &event, Int_t type, Stat_t w)
1386{
1387 if (fNcells<=1 || IsFreezed())
1388 return;
1389
1390 const Double_t median = event.GetCameraMedian(*fGeomCam, type);
1391
1392 // FIXME: Security check missing!
1393 for (Int_t idx=0; idx<fNcells-2; idx++)
1394 {
1395 Double_t val=0;
1396 if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
1397 {
1398 SetUsed(idx);
1399 Fill(idx, (val-median)*w); // FIXME: Slow!
1400 }
1401 }
1402 fEntries++;
1403}
1404
1405// ------------------------------------------------------------------------
1406//
1407// Call this function to add a MCamEvent on top of the present contents.
1408//
1409void MHCamera::AddCamContent(const MCamEvent &event, Int_t type, Stat_t w)
1410{
1411 if (fNcells<=1 || IsFreezed())
1412 return;
1413
1414 // FIXME: Security check missing!
1415 for (Int_t idx=0; idx<fNcells-2; idx++)
1416 {
1417 Double_t val=0;
1418 if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
1419 {
1420 SetUsed(idx);
1421 Fill(idx, val*w); // FIXME: Slow!
1422 }
1423 }
1424 fEntries++;
1425}
1426
1427// ------------------------------------------------------------------------
1428//
1429void MHCamera::AddCamDifference(const MCamEvent &event1, const MCamEvent &event2, Int_t type)
1430{
1431 if (fNcells<=1 || IsFreezed())
1432 return;
1433
1434 // FIXME: Security check missing!
1435 for (Int_t idx=0; idx<fNcells-2; idx++)
1436 {
1437 Double_t val1=0;
1438 if (!event1.GetPixelContent(val1, idx, *fGeomCam, type))
1439 continue;
1440
1441 Double_t val2=0;
1442 if (!event2.GetPixelContent(val2, idx, *fGeomCam, type))
1443 continue;
1444
1445 SetUsed(idx);
1446 Fill(idx, val1-val2); // FIXME: Slow!
1447 }
1448 fEntries++;
1449}
1450
1451// ------------------------------------------------------------------------
1452//
1453// Call this function to add a MCamEvent on top of the present contents.
1454//
1455void MHCamera::SetCamError(const MCamEvent &evt, Int_t type)
1456{
1457
1458 if (fNcells<=1 || IsFreezed())
1459 return;
1460
1461 // FIXME: Security check missing!
1462 for (Int_t idx=0; idx<fNcells-2; idx++)
1463 {
1464 Double_t val=0;
1465 if (evt.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
1466 SetUsed(idx);
1467
1468 SetBinError(idx+1, val); // FIXME: Slow!
1469 }
1470}
1471
1472Stat_t MHCamera::GetBinContent(Int_t bin) const
1473{
1474 if (fBuffer) ((TH1D*)this)->BufferEmpty();
1475 if (bin < 0) bin = 0;
1476 if (bin >= fNcells) bin = fNcells-1;
1477 if (!fArray) return 0;
1478
1479 if (!TestBit(kProfile))
1480 return Stat_t (fArray[bin]);
1481
1482 if (fBinEntries.fArray[bin] == 0) return 0;
1483 return fArray[bin]/fBinEntries.fArray[bin];
1484}
1485
1486// ------------------------------------------------------------------------
1487//
1488// In the case the kProfile flag is set the spread of the bin is returned.
1489// If you want to have the mean error instead set the kErrorMean bit via
1490// SetBit(kErrorMean) first.
1491//
1492Stat_t MHCamera::GetBinError(Int_t bin) const
1493{
1494 if (!TestBit(kProfile))
1495 return TH1D::GetBinError(bin);
1496
1497 const UInt_t n = (UInt_t)fBinEntries[bin];
1498
1499 if (n==0)
1500 return 0;
1501
1502 const Double_t sqr = fSumw2.fArray[bin] / n;
1503 const Double_t val = fArray[bin] / n;
1504
1505 const Double_t spread = sqr>val*val ? TMath::Sqrt(sqr - val*val) : 0;
1506
1507 return TestBit(kErrorMean) ? spread/TMath::Sqrt(n) : spread;
1508
1509 /*
1510 Double_t rc = 0;
1511 if (TestBit(kSqrtVariance) && GetEntries()>0) // error on the mean
1512 {
1513 const Double_t error = fSumw2.fArray[bin]/GetEntries();
1514 const Double_t val = fArray[bin]/GetEntries();
1515 rc = val*val>error ? 0 : TMath::Sqrt(error - val*val);
1516 }
1517 else
1518 rc = TH1D::GetBinError(bin);
1519
1520 return Profile(rc);*/
1521}
1522
1523// ------------------------------------------------------------------------
1524//
1525// Call this function to add a MHCamera on top of the present contents.
1526// Type:
1527// 0) bin content
1528// 1) errors
1529// 2) rel. errors
1530//
1531void MHCamera::AddCamContent(const MHCamera &d, Int_t type)
1532{
1533 if (fNcells!=d.fNcells || IsFreezed())
1534 return;
1535
1536 // FIXME: Security check missing!
1537 for (Int_t idx=0; idx<fNcells-2; idx++)
1538 if (d.IsUsed(idx))
1539 SetUsed(idx);
1540
1541 switch (type)
1542 {
1543 case 1:
1544 // Under-/Overflow bins not handled!
1545 for (Int_t idx=0; idx<fNcells-2; idx++)
1546 if (d.IsUsed(idx))
1547 Fill(idx, d.GetBinError(idx+1));
1548 fEntries++;
1549 break;
1550 case 2:
1551 // Under-/Overflow bins not handled!
1552 for (Int_t idx=0; idx<fNcells-2; idx++)
1553 if (d.GetBinContent(idx+1)!=0 && d.IsUsed(idx))
1554 Fill(idx, TMath::Abs(d.GetBinError(idx+1)/d.GetBinContent(idx+1)));
1555 fEntries++;
1556 break;
1557 default:
1558 if (TestBit(kProfile)!=d.TestBit(kProfile))
1559 gLog << warn << "WARNING - You have tried to call AddCamContent for two different kind of histograms (kProfile set or not)." << endl;
1560
1561 // environment
1562 fEntries += d.fEntries;
1563 fTsumw += d.fTsumw;
1564 fTsumw2 += d.fTsumw2;
1565 fTsumwx += d.fTsumwx;
1566 fTsumwx2 += d.fTsumwx2;
1567 // Bin contents
1568 for (Int_t idx=1; idx<fNcells-1; idx++)
1569 {
1570 if (!d.IsUsed(idx-1))
1571 continue;
1572
1573 fArray[idx] += d.fArray[idx];
1574 fBinEntries[idx] += d.fBinEntries[idx];
1575 fSumw2.fArray[idx] += d.fSumw2.fArray[idx];
1576 }
1577 // Underflow bin
1578 fArray[0] += d.fArray[0];
1579 fBinEntries[0] += d.fBinEntries[0];
1580 fSumw2.fArray[0] += d.fSumw2.fArray[0];
1581 // Overflow bin
1582 fArray[fNcells-1] += d.fArray[fNcells-1];
1583 fBinEntries[fNcells-1] += d.fBinEntries[fNcells-1];
1584 fSumw2.fArray[fNcells-1] += d.fSumw2.fArray[fNcells-1];
1585 break;
1586/* default:
1587 if (TestBit(kProfile)!=d.TestBit(kProfile))
1588 gLog << warn << "WARNING - You have tried to call AddCamContent for two different kind of histograms (kProfile set or not)." << endl;
1589
1590 for (Int_t idx=0; idx<fNcells-2; idx++)
1591 Fill(idx, d.GetBinContent(idx+1));
1592 break;*/
1593 }
1594 fEntries++;
1595}
1596
1597// ------------------------------------------------------------------------
1598//
1599// Call this function to add a TArrayD on top of the present contents.
1600//
1601void MHCamera::AddCamContent(const TArrayD &event, const TArrayC *used)
1602{
1603 if (event.GetSize()!=fNcells-2 || IsFreezed())
1604 return;
1605
1606 if (used && used->GetSize()!=fNcells-2)
1607 return;
1608
1609 for (Int_t idx=0; idx<fNcells-2; idx++)
1610 {
1611 Fill(idx, event[idx]); // FIXME: Slow!
1612
1613 if (!used || (*used)[idx])
1614 SetUsed(idx);
1615 }
1616 fEntries++;
1617}
1618
1619// ------------------------------------------------------------------------
1620//
1621// Call this function to add a MArrayD on top of the present contents.
1622//
1623void MHCamera::AddCamContent(const MArrayD &event, const TArrayC *used)
1624{
1625 if (event.GetSize()!=(UInt_t)(fNcells-2) || IsFreezed())
1626 return;
1627
1628 if (used && used->GetSize()!=fNcells-2)
1629 return;
1630
1631 for (Int_t idx=0; idx<fNcells-2; idx++)
1632 {
1633 Fill(idx, event[idx]); // FIXME: Slow!
1634
1635 if (!used || (*used)[idx])
1636 SetUsed(idx);
1637 }
1638 fEntries++;
1639}
1640
1641// ------------------------------------------------------------------------
1642//
1643// Call this function to add a MCamEvent on top of the present contents.
1644// 1 is added to each pixel if the contents of MCamEvent>threshold
1645// (in case isabove is set to kTRUE == default)
1646// 1 is added to each pixel if the contents of MCamEvent<threshold
1647// (in case isabove is set to kFALSE)
1648//
1649// in unused pixel is not counted if it didn't fullfill the condition.
1650//
1651void MHCamera::CntCamContent(const MCamEvent &event, Double_t threshold, Int_t type, Bool_t isabove)
1652{
1653 if (fNcells<=1 || IsFreezed())
1654 return;
1655
1656 // FIXME: Security check missing!
1657 for (Int_t idx=0; idx<fNcells-2; idx++)
1658 {
1659 Double_t val=threshold;
1660 const Bool_t rc = event.GetPixelContent(val, idx, *fGeomCam, type);
1661 if (rc)
1662 SetUsed(idx);
1663
1664 const Bool_t cond =
1665 ( isabove && val>threshold) ||
1666 (!isabove && val<threshold);
1667
1668 Fill(idx, rc && cond ? 1 : 0);
1669 }
1670 fEntries++;
1671}
1672
1673// ------------------------------------------------------------------------
1674//
1675// Call this function to add a MCamEvent on top of the present contents.
1676// - the contents of the pixels in event are added to each pixel
1677// if the pixel of thresevt<threshold (in case isabove is set
1678// to kTRUE == default)
1679// - the contents of the pixels in event are added to each pixel
1680// if the pixel of thresevt<threshold (in case isabove is set
1681// to kFALSE)
1682//
1683// in unused pixel is not counted if it didn't fullfill the condition.
1684//
1685void MHCamera::CntCamContent(const MCamEvent &event, Int_t type1, const MCamEvent &thresevt, Int_t type2, Double_t threshold, Bool_t isabove)
1686{
1687 if (fNcells<=1 || IsFreezed())
1688 return;
1689
1690 // FIXME: Security check missing!
1691 for (Int_t idx=0; idx<fNcells-2; idx++)
1692 {
1693 Double_t th=0;
1694 if (!thresevt.GetPixelContent(th, idx, *fGeomCam, type2))
1695 continue;
1696
1697 if ((isabove && th>threshold) || (!isabove && th<threshold))
1698 continue;
1699
1700 Double_t val=th;
1701 if (event.GetPixelContent(val, idx, *fGeomCam, type1))
1702 {
1703 SetUsed(idx);
1704 Fill(idx, val);
1705 }
1706 }
1707 fEntries++;
1708}
1709
1710// ------------------------------------------------------------------------
1711//
1712// Call this function to add a MCamEvent on top of the present contents.
1713// 1 is added to each pixel if the contents of MCamEvent>threshold
1714// (in case isabove is set to kTRUE == default)
1715// 1 is added to each pixel if the contents of MCamEvent<threshold
1716// (in case isabove is set to kFALSE)
1717//
1718// in unused pixel is not counted if it didn't fullfill the condition.
1719//
1720void MHCamera::CntCamContent(const MCamEvent &event, TArrayD threshold, Int_t type, Bool_t isabove)
1721{
1722 if (fNcells<=1 || IsFreezed())
1723 return;
1724
1725 // FIXME: Security check missing!
1726 for (Int_t idx=0; idx<fNcells-2; idx++)
1727 {
1728 Double_t val=threshold[idx];
1729 if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
1730 {
1731 SetUsed(idx);
1732
1733 if (val>threshold[idx] && isabove)
1734 Fill(idx);
1735 if (val<threshold[idx] && !isabove)
1736 Fill(idx);
1737 }
1738 }
1739 fEntries++;
1740}
1741
1742// ------------------------------------------------------------------------
1743//
1744// Call this function to add a TArrayD on top of the present contents.
1745// 1 is added to each pixel if the contents of MCamEvent>threshold
1746//
1747void MHCamera::CntCamContent(const TArrayD &event, Double_t threshold, Bool_t ispos)
1748{
1749 if (event.GetSize()!=fNcells-2 || IsFreezed())
1750 return;
1751
1752 for (Int_t idx=0; idx<fNcells-2; idx++)
1753 {
1754 if (event[idx]>threshold)
1755 Fill(idx);
1756
1757 if (!ispos || fArray[idx+1]>0)
1758 SetUsed(idx);
1759 }
1760 fEntries++;
1761}
1762
1763// ------------------------------------------------------------------------
1764//
1765// Call this function to add a MCamEvent on top of the present contents.
1766// 1 is added to each pixel if the contents of MCamEvent>threshold
1767// (in case isabove is set to kTRUE == default)
1768// 1 is added to each pixel if the contents of MCamEvent<threshold
1769// (in case isabove is set to kFALSE)
1770//
1771// in unused pixel is not counted if it didn't fullfill the condition.
1772//
1773void MHCamera::SetMaxCamContent(const MCamEvent &event, Int_t type)
1774{
1775 if (fNcells<=1 || IsFreezed())
1776 return;
1777
1778 // FIXME: Security check missing!
1779 for (Int_t idx=0; idx<fNcells-2; idx++)
1780 {
1781 Double_t val=0;
1782 const Bool_t rc = event.GetPixelContent(val, idx, *fGeomCam, type);
1783 if (!rc)
1784 continue;
1785
1786 if (!IsUsed(idx))
1787 {
1788 fArray[idx+1] = val;
1789 SetUsed(idx);
1790 fBinEntries.fArray[idx+1]=1;
1791 }
1792 else
1793 if (val>fArray[idx+1])
1794 fArray[idx+1] = val;
1795 }
1796 fEntries++;
1797}
1798
1799// ------------------------------------------------------------------------
1800//
1801// Call this function to add a MCamEvent on top of the present contents.
1802// 1 is added to each pixel if the contents of MCamEvent>threshold
1803// (in case isabove is set to kTRUE == default)
1804// 1 is added to each pixel if the contents of MCamEvent<threshold
1805// (in case isabove is set to kFALSE)
1806//
1807// in unused pixel is not counted if it didn't fullfill the condition.
1808//
1809void MHCamera::SetMinCamContent(const MCamEvent &event, Int_t type)
1810{
1811 if (fNcells<=1 || IsFreezed())
1812 return;
1813
1814 // FIXME: Security check missing!
1815 for (Int_t idx=0; idx<fNcells-2; idx++)
1816 {
1817 Double_t val=0;
1818 const Bool_t rc = event.GetPixelContent(val, idx, *fGeomCam, type);
1819 if (!rc)
1820 continue;
1821
1822 if (!IsUsed(idx))
1823 {
1824 fArray[idx+1] = val;
1825 SetUsed(idx);
1826 fBinEntries.fArray[idx+1]=1;
1827 }
1828 else
1829 if (val<fArray[idx+1])
1830 fArray[idx+1] = val;
1831 }
1832 fEntries++;
1833}
1834
1835// ------------------------------------------------------------------------
1836//
1837// Fill the pixels with random contents.
1838//
1839void MHCamera::FillRandom()
1840{
1841 if (fNcells<=1 || IsFreezed())
1842 return;
1843
1844 Reset();
1845
1846 // FIXME: Security check missing!
1847 for (Int_t idx=0; idx<fNcells-2; idx++)
1848 {
1849 Fill(idx, gRandom->Uniform()*fGeomCam->GetPixRatio(idx));
1850 SetUsed(idx);
1851 }
1852 fEntries=1;
1853}
1854
1855
1856// ------------------------------------------------------------------------
1857//
1858// The array must be in increasing order, eg: 2.5, 3.7, 4.9
1859// The values in each bin are replaced by the interval in which the value
1860// fits. In the example we have four intervals
1861// (<2.5, 2.5-3.7, 3.7-4.9, >4.9). Maximum and minimum are set
1862// accordingly.
1863//
1864void MHCamera::SetLevels(const TArrayF &arr)
1865{
1866 if (fNcells<=1)
1867 return;
1868
1869 for (Int_t i=0; i<fNcells-2; i++)
1870 {
1871 if (!IsUsed(i))
1872 continue;
1873
1874 Int_t j = arr.GetSize();
1875 while (j && fArray[i+1]<arr[j-1])
1876 j--;
1877
1878 fArray[i+1] = j;
1879 }
1880 SetMaximum(arr.GetSize());
1881 SetMinimum(0);
1882}
1883
1884// ------------------------------------------------------------------------
1885//
1886// Reset the all pixel colors to a default value
1887//
1888void MHCamera::Reset(Option_t *opt)
1889{
1890 if (fNcells<=1 || IsFreezed())
1891 return;
1892
1893 TH1::Reset(opt);
1894
1895 fUsed.Reset();
1896 fBinEntries.Reset();
1897
1898 for (Int_t i=0; i<fNcells; i++)
1899 fArray[i] = 0;
1900}
1901
1902// ------------------------------------------------------------------------
1903//
1904// Here we calculate the color index for the current value.
1905// The color index is defined with the class TStyle and the
1906// Color palette inside. We use the command gStyle->SetPalette(1,0)
1907// for the display. So we have to convert the value "wert" into
1908// a color index that fits the color palette.
1909// The range of the color palette is defined by the values fMinPhe
1910// and fMaxRange. Between this values we have 50 color index, starting
1911// with 0 up to 49.
1912//
1913Int_t MHCamera::GetColor(Float_t val, Float_t min, Float_t max, Bool_t islog)
1914{
1915 //
1916 // first treat the over- and under-flows
1917 //
1918 const Int_t ncol = GetContour()==0?50:GetContour();
1919
1920 const Int_t maxcolidx = ncol-1;
1921
1922 if (!TMath::Finite(val)) // FIXME: gLog!
1923 return maxcolidx/2;
1924
1925 if (val >= max)
1926 return gStyle->GetColorPalette(maxcolidx);
1927
1928 if (val <= min)
1929 return gStyle->GetColorPalette(0);
1930
1931 //
1932 // calculate the color index
1933 //
1934 Float_t ratio;
1935 if (islog && min>0)
1936 ratio = log10(val/min) / log10(max/min);
1937 else
1938 ratio = (val-min) / (max-min);
1939
1940 const Int_t colidx = TMath::FloorNint(ratio*ncol);
1941 return gStyle->GetColorPalette(colidx);
1942}
1943
1944TPaveStats *MHCamera::GetStatisticBox()
1945{
1946 TObject *obj = 0;
1947
1948 TIter Next(fFunctions);
1949 while ((obj = Next()))
1950 if (obj->InheritsFrom(TPaveStats::Class()))
1951 return static_cast<TPaveStats*>(obj);
1952
1953 return NULL;
1954}
1955
1956// ------------------------------------------------------------------------
1957//
1958// Change the text on the legend according to the range of the Display
1959//
1960void MHCamera::UpdateLegend(Float_t min, Float_t max, Bool_t islog)
1961{
1962 const Float_t range = fGeomCam->GetMaxRadius();
1963
1964 if (!TestBit(kNoScale))
1965 {
1966 TArrow arr;
1967 arr.PaintArrow(-range*.99, -range*.99, -range*.60, -range*.99, 0.010);
1968 arr.PaintArrow(-range*.99, -range*.99, -range*.99, -range*.60, 0.010);
1969
1970 TString text;
1971 text += (int)(range*.3);
1972 text += "mm";
1973
1974 TText newtxt2;
1975 newtxt2.SetTextSize(0.04);
1976 newtxt2.PaintText(-range*.95, -range*.95, text);
1977
1978 text = MString::Format("%.2f\\circ", range*(0.99-0.60)*fGeomCam->GetConvMm2Deg());
1979 text = text.Strip(TString::kLeading);
1980
1981 TLatex latex;
1982 latex.PaintLatex(-range*.95, -range*.85, 0, 0.04, text);
1983 }
1984
1985 if (!TestBit(kNoLegend))
1986 {
1987 TPaveStats *stats = GetStatisticBox();
1988
1989 const Float_t hndc = 0.88 - (stats ? stats->GetY1NDC() : 1);
1990 const Float_t H = (0.80-hndc)*range;
1991 const Float_t offset = hndc*range;
1992
1993 const Int_t ncol = GetContour()==0 ? 50 : GetContour();
1994
1995 const Float_t h = 2./ncol;
1996 const Float_t w = range/sqrt((float)(fNcells-2));
1997
1998 TBox newbox;
1999 TText newtxt;
2000 newtxt.SetTextSize(0.03);
2001 newtxt.SetTextAlign(12);
2002#if ROOT_VERSION_CODE > ROOT_VERSION(3,01,06)
2003 newtxt.SetBit(/*kNoContextMenu|*/kCannotPick);
2004 newbox.SetBit(/*kNoContextMenu|*/kCannotPick);
2005#endif
2006
2007 const Float_t step = (islog && min>0 ? log10(max/min) : max-min);
2008 const Int_t firsts = step/48*3 < 1e-8 ? 8 : (Int_t)floor(log10(step/48*3));
2009 const TString opt = MString::Format("%%.%if", firsts>0 ? 0 : TMath::Abs(firsts));
2010/*
2011 for (Int_t i=0; i<ncol+1; i+=3)
2012 {
2013 Float_t val;
2014 if (islog && min>0)
2015 val = pow(10, step*i) * min;
2016 else
2017 val = min + step*i;
2018
2019 //const bool dispexp = max-min>1.5 && fabs(val)>0.1 && fabs(val)<1e6;
2020 newtxt.PaintText(range+1.5*w, H*(i*h-1)-offset, Form(opt, val));
2021 }
2022 */
2023 const MBinning bins(25, min, max, islog&&min>0?"log":"lin");
2024
2025 for (Int_t i=0; i<=25; i++)
2026 newtxt.PaintText(range+1.5*w, H*(i*ncol/25*h-1)-offset, MString::Format(opt, bins[i]));
2027
2028 for (Int_t i=0; i<ncol; i++)
2029 {
2030 newbox.SetFillColor(gStyle->GetColorPalette(i));
2031 newbox.PaintBox(range, H*(i*h-1)-offset, range+w, H*((i+1)*h-1)-offset);
2032 }
2033 }
2034}
2035
2036// ------------------------------------------------------------------------
2037//
2038// Save primitive as a C++ statement(s) on output stream out
2039//
2040void MHCamera::SavePrimitive(ostream &out, Option_t *opt)
2041{
2042 gLog << err << "MHCamera::SavePrimitive: Must be rewritten!" << endl;
2043 /*
2044 if (!gROOT->ClassSaved(TCanvas::Class()))
2045 fDrawingPad->SavePrimitive(out, opt);
2046
2047 out << " " << fDrawingPad->GetName() << "->SetWindowSize(";
2048 out << fDrawingPad->GetWw() << "," << fDrawingPad->GetWh() << ");" << endl;
2049 */
2050}
2051
2052void MHCamera::SavePrimitive(ofstream &out, Option_t *)
2053{
2054 MHCamera::SavePrimitive(static_cast<ostream&>(out), "");
2055}
2056
2057// ------------------------------------------------------------------------
2058//
2059// compute the distance of a point (px,py) to the Camera
2060// this functions needed for graphical primitives, that
2061// means without this function you are not able to interact
2062// with the graphical primitive with the mouse!!!
2063//
2064// All calcutations are done in pixel coordinates
2065//
2066Int_t MHCamera::DistancetoPrimitive(Int_t px, Int_t py)
2067{
2068 if (fNcells<=1)
2069 return 999999;
2070
2071 TPaveStats *box = (TPaveStats*)gPad->GetPrimitive("stats");
2072 if (box)
2073 {
2074 const Double_t w = box->GetY2NDC()-box->GetY1NDC();
2075 box->SetX1NDC(gStyle->GetStatX()-gStyle->GetStatW());
2076 box->SetY1NDC(gStyle->GetStatY()-w);
2077 box->SetX2NDC(gStyle->GetStatX());
2078 box->SetY2NDC(gStyle->GetStatY());
2079 }
2080
2081 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
2082 return TH1D::DistancetoPrimitive(px, py);
2083
2084 const Bool_t issame = TString(GetDrawOption()).Contains("same", TString::kIgnoreCase);
2085
2086 const Float_t maxr = (1-fGeomCam->GetConvMm2Deg())*fGeomCam->GetMaxRadius()/2;
2087 const Float_t conv = !issame ||
2088 gPad->GetX1()<-maxr || gPad->GetY1()<-maxr ||
2089 gPad->GetX2()> maxr || gPad->GetY2()>maxr ? 1 : fGeomCam->GetConvMm2Deg();
2090
2091 if (GetPixelIndex(px, py, conv)>=0)
2092 return 0;
2093
2094 if (!box)
2095 return 999999;
2096
2097 const Int_t dist = box->DistancetoPrimitive(px, py);
2098 if (dist > TPad::GetMaxPickDistance())
2099 return 999999;
2100
2101 gPad->SetSelected(box);
2102 return dist;
2103}
2104
2105// ------------------------------------------------------------------------
2106//
2107//
2108Int_t MHCamera::GetPixelIndex(Int_t px, Int_t py, Float_t conv) const
2109{
2110 if (fNcells<=1)
2111 return -1;
2112
2113 for (Int_t i=0; i<fNcells-2; i++)
2114 if ((*fGeomCam)[i].DistancetoPrimitive(TMath::Nint(px*conv), TMath::Nint(py*conv))<=0)
2115 return i;
2116
2117 return -1;
2118}
2119
2120// ------------------------------------------------------------------------
2121//
2122// Returns string containing info about the object at position (px,py).
2123// Returned string will be re-used (lock in MT environment).
2124//
2125char *MHCamera::GetObjectInfo(Int_t px, Int_t py) const
2126{
2127 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
2128#if ROOT_VERSION_CODE < ROOT_VERSION(5,22,00)
2129 return MH::GetObjectInfoH(px, py, *this);
2130#else
2131 return TH1D::GetObjectInfo(px, py);
2132#endif
2133
2134 static char info[128];
2135
2136 const Int_t idx=GetPixelIndex(px, py);
2137
2138 if (idx<0)
2139 return TObject::GetObjectInfo(px, py);
2140
2141 sprintf(info, "Software Pixel Idx: %d (Hardware Id=%d) c=%.1f <%s>",
2142 idx, idx+1, GetBinContent(idx+1), IsUsed(idx)?"on":"off");
2143 return info;
2144}
2145
2146// ------------------------------------------------------------------------
2147//
2148// Add a MCamEvent which should be displayed when the user clicks on a
2149// pixel.
2150// Warning: The object MUST inherit from TObject AND MCamEvent
2151//
2152void MHCamera::AddNotify(TObject *obj)
2153{
2154 // Make sure, that the object derives from MCamEvent!
2155 MCamEvent *evt = dynamic_cast<MCamEvent*>(obj);
2156 if (!evt)
2157 {
2158 gLog << err << "ERROR: MHCamera::AddNotify - TObject doesn't inherit from MCamEvent... ignored." << endl;
2159 return;
2160 }
2161
2162 // Make sure, that it is deleted from the list too, if the obj is deleted
2163 obj->SetBit(kMustCleanup);
2164
2165 // Add object to list
2166 fNotify->Add(obj);
2167}
2168
2169// ------------------------------------------------------------------------
2170//
2171// Execute a mouse event on the camera
2172//
2173void MHCamera::ExecuteEvent(Int_t event, Int_t px, Int_t py)
2174{
2175 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
2176 {
2177 TH1D::ExecuteEvent(event, px, py);
2178 return;
2179 }
2180 //if (event==kMouseMotion && fStatusBar)
2181 // fStatusBar->SetText(GetObjectInfo(px, py), 0);
2182 if (event!=kButton1Down)
2183 return;
2184
2185 const Int_t idx = GetPixelIndex(px, py);
2186 if (idx<0)
2187 return;
2188
2189 gLog << all << GetTitle() << " <" << GetName() << ">" << dec << endl;
2190 gLog << "Software Pixel Idx: " << idx << endl;
2191 gLog << "Hardware Pixel Id: " << idx+1 << endl;
2192 gLog << "Contents: " << GetBinContent(idx+1);
2193 if (GetBinError(idx+1)>0)
2194 gLog << " +/- " << GetBinError(idx+1);
2195 gLog << " <" << (IsUsed(idx)?"on":"off") << "> n=" << fBinEntries[idx+1] << endl;
2196
2197 if (fNotify && fNotify->GetSize()>0)
2198 {
2199 // FIXME: Is there a simpler and more convinient way?
2200
2201 // The name which is created here depends on the instance of
2202 // MHCamera and on the pad on which it is drawn --> The name
2203 // is unique. For ExecuteEvent gPad is always correctly set.
2204 const TString name = MString::Format("%p;%p;PixelContent", this, gPad);
2205
2206 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
2207 if (old)
2208 old->cd();
2209 else
2210 new TCanvas(name);
2211
2212 /*
2213 TIter Next(gPad->GetListOfPrimitives());
2214 TObject *o;
2215 while (o=Next()) cout << o << ": " << o->GetName() << " " << o->IsA()->GetName() << endl;
2216 */
2217
2218 // FIXME: Make sure, that the old histograms are really deleted.
2219 // Are they already deleted?
2220
2221 // The dynamic_cast is necessary here: We cannot use ForEach
2222 TIter Next(fNotify);
2223 MCamEvent *evt;
2224 while ((evt=dynamic_cast<MCamEvent*>(Next())))
2225 evt->DrawPixelContent(idx);
2226
2227 gPad->Modified();
2228 gPad->Update();
2229 }
2230}
2231
2232UInt_t MHCamera::GetNumPixels() const
2233{
2234 return fGeomCam ? fGeomCam->GetNumPixels() : 0;
2235}
2236
2237TH1 *MHCamera::DrawCopy() const
2238{
2239 gPad=NULL;
2240 return TH1D::DrawCopy(fName+";cpy");
2241}
2242
2243// --------------------------------------------------------------------------
2244//
2245// Draw a projection of MHCamera onto the y-axis values. Depending on the
2246// variable fit, the following fits are performed:
2247//
2248// 0: No fit, simply draw the projection
2249// 1: Single Gauss (for distributions flat-fielded over the whole camera)
2250// 2: Double Gauss (for distributions different for inner and outer pixels)
2251// 3: Triple Gauss (for distributions with inner, outer pixels and outliers)
2252// 4: flat (for the probability distributions)
2253// (1-4:) Moreover, sectors 6,1 and 2 of the camera and sectors 3,4 and 5 are
2254// drawn separately, for inner and outer pixels.
2255// 5: Fit Inner and Outer pixels separately by a single Gaussian
2256// (only for MAGIC cameras)
2257// 6: Fit Inner and Outer pixels separately by a single Gaussian and display
2258// additionally the two camera halfs separately (for MAGIC camera)
2259// 7: Single Gauss with TLegend to show the meaning of the colours
2260//
2261void MHCamera::DrawProjection(Int_t fit) const
2262{
2263 if (fit==5 || fit==6)
2264 {
2265 const UInt_t n = fGeomCam->GetNumAreas();
2266
2267 TVirtualPad *pad = gPad;
2268 pad->Divide(n, 1, 1e-5, 1e-5);;
2269
2270 for (UInt_t i=0; i<n; i++)
2271 {
2272 pad->cd(i+1);
2273 gPad->SetBorderMode(0);
2274 gPad->SetRightMargin(0.025);
2275 gPad->SetTicks();
2276
2277 TH1D &h = *ProjectionS(TArrayI(), TArrayI(1, (Int_t*)&i),
2278 MString::Format("%s_%d", GetName(), i));
2279 h.SetTitle(MString::Format("%s %d", GetTitle(), i));
2280 h.SetDirectory(NULL);
2281 h.SetBit(kCanDelete);
2282 h.Draw();
2283
2284 TAxis *xaxe = h.GetXaxis();
2285 TAxis *yaxe = h.GetYaxis();
2286
2287 xaxe->CenterTitle();
2288 yaxe->CenterTitle();
2289 xaxe->SetTitleSize(0.06);
2290 yaxe->SetTitleSize(0.06);
2291 xaxe->SetTitleOffset(0.8);
2292 yaxe->SetTitleOffset(0.85);
2293 xaxe->SetLabelSize(0.05);
2294 yaxe->SetLabelSize(0.05);
2295 if (i>0)
2296 yaxe->SetTitle("");
2297
2298 if (fit==5)
2299 continue;
2300
2301 h.Fit("gaus", "Q");
2302
2303 TF1 *f = h.GetFunction("gaus");
2304 if (f)
2305 {
2306 f->SetLineWidth(2);
2307 f->SetLineColor(kBlue);
2308 }
2309 }
2310 return;
2311 }
2312
2313 TH1D *obj2 = (TH1D*)Projection(GetName());
2314 obj2->SetDirectory(0);
2315 obj2->Draw();
2316 obj2->SetBit(kCanDelete);
2317
2318 if (fit==0)
2319 return;
2320
2321 const Double_t xmin = obj2->GetBinCenter(obj2->GetXaxis()->GetFirst());
2322 const Double_t xmax = obj2->GetBinCenter(obj2->GetXaxis()->GetLast());
2323 const Double_t integ = obj2->Integral("width")/2.5;
2324 const Double_t max = obj2->GetMaximum();
2325 const Double_t mean = obj2->GetMean();
2326 const Double_t rms = obj2->GetRMS();
2327 const Double_t width = xmax-xmin;
2328
2329 const char *dgausformula = "([0]-[3])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
2330 "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])";
2331
2332 const char *tgausformula = "([0]-[3]-[6])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
2333 "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])"
2334 "+[6]/[8]*exp(-0.5*(x-[7])*(x-[7])/[8]/[8])";
2335
2336 TF1 *f=0;
2337 switch (fit)
2338 {
2339 case 1: // Single gauss
2340 f = new TF1("sgaus", "gaus", xmin, xmax);
2341 f->SetLineColor(kBlue);
2342 f->SetBit(kCanDelete);
2343 f->SetParNames("Max", "#mu", "#sigma");
2344 f->SetParameters(max, mean, rms);
2345 f->SetParLimits(0, 0, max*2);
2346 f->SetParLimits(1, xmin, xmax);
2347 f->SetParLimits(2, 0, width/1.5);
2348
2349 obj2->Fit(f, "QLR");
2350 break;
2351
2352 case 2: // Double gauss
2353 f = new TF1("dgaus", dgausformula, xmin, xmax);
2354 f->SetLineColor(kBlue);
2355 f->SetBit(kCanDelete);
2356 f->SetParNames("A_{tot}", "#mu1", "#sigma1", "A2", "#mu2", "#sigma2");
2357 f->SetParameters(integ, (xmin+mean)/2., width/4.,
2358 integ/width/2., (xmax+mean)/2., width/4.);
2359 // The left-sided Gauss
2360 f->SetParLimits(0, integ-1.5, integ+1.5);
2361 f->SetParLimits(1, xmin+(width/10.), mean);
2362 f->SetParLimits(2, 0, width/2.);
2363 // The right-sided Gauss
2364 f->SetParLimits(3, 0, integ);
2365 f->SetParLimits(4, mean, xmax-(width/10.));
2366 f->SetParLimits(5, 0, width/2.);
2367 obj2->Fit(f,"QLRM");
2368 break;
2369
2370 case 3: // Triple gauss
2371 f = new TF1("tgaus", tgausformula, xmin,xmax);
2372 f->SetLineColor(kBlue);
2373 f->SetBit(kCanDelete);
2374 f->SetParNames("A_{tot}","#mu_{1}","#sigma_{1}",
2375 "A_{2}","#mu_{2}","#sigma_{2}",
2376 "A_{3}","#mu_{3}","#sigma_{3}");
2377 f->SetParameters(integ, (xmin+mean)/2, width/4.,
2378 integ/width/3., (xmax+mean)/2., width/4.,
2379 integ/width/3., mean,width/2.);
2380 // The left-sided Gauss
2381 f->SetParLimits(0, integ-1.5, integ+1.5);
2382 f->SetParLimits(1, xmin+(width/10.), mean);
2383 f->SetParLimits(2, width/15., width/2.);
2384 // The right-sided Gauss
2385 f->SetParLimits(3, 0., integ);
2386 f->SetParLimits(4, mean, xmax-(width/10.));
2387 f->SetParLimits(5, width/15., width/2.);
2388 // The Gauss describing the outliers
2389 f->SetParLimits(6, 0., integ);
2390 f->SetParLimits(7, xmin, xmax);
2391 f->SetParLimits(8, width/4., width/1.5);
2392 obj2->Fit(f,"QLRM");
2393 break;
2394
2395 case 4:
2396 obj2->Fit("pol0", "Q");
2397 if (obj2->GetFunction("pol0"))
2398 obj2->GetFunction("pol0")->SetLineColor(kBlue);
2399 break;
2400
2401 case 9:
2402 break;
2403
2404 default:
2405 obj2->Fit("gaus", "Q");
2406 if (obj2->GetFunction("gaus"))
2407 obj2->GetFunction("gaus")->SetLineColor(kBlue);
2408 break;
2409 }
2410}
2411
2412// --------------------------------------------------------------------------
2413//
2414// Draw a projection of MHCamera vs. the radius from the central pixel.
2415//
2416// The inner and outer pixels are drawn separately, both fitted by a polynomial
2417// of grade 1.
2418//
2419void MHCamera::DrawRadialProfile() const
2420{
2421 TProfile *obj2 = (TProfile*)RadialProfile(GetName());
2422 obj2->SetDirectory(0);
2423 obj2->Draw();
2424 obj2->SetBit(kCanDelete);
2425/*
2426 if (GetGeomCam().InheritsFrom("MGeomCamMagic"))
2427 {
2428 TArrayI s0(6);
2429 s0[0] = 1;
2430 s0[1] = 2;
2431 s0[2] = 3;
2432 s0[3] = 4;
2433 s0[4] = 5;
2434 s0[5] = 6;
2435
2436 TArrayI inner(1);
2437 inner[0] = 0;
2438
2439 TArrayI outer(1);
2440 outer[0] = 1;
2441
2442 // Just to get the right (maximum) binning
2443 TProfile *half[2];
2444 half[0] = RadialProfileS(s0, inner, MString::Format("%sInner",GetName()));
2445 half[1] = RadialProfileS(s0, outer, MString::Format("%sOuter",GetName()));
2446
2447 for (Int_t i=0; i<2; i++)
2448 {
2449 Double_t min = GetGeomCam().GetMinRadius(i);
2450 Double_t max = GetGeomCam().GetMaxRadius(i);
2451
2452 half[i]->SetLineColor(kRed+i);
2453 half[i]->SetDirectory(0);
2454 half[i]->SetBit(kCanDelete);
2455 half[i]->Draw("same");
2456 half[i]->Fit("pol1","Q","",min,max);
2457 half[i]->GetFunction("pol1")->SetLineColor(kRed+i);
2458 half[i]->GetFunction("pol1")->SetLineWidth(1);
2459 }
2460 }
2461*/
2462}
2463
2464// --------------------------------------------------------------------------
2465//
2466// Draw a projection of MHCamera vs. the azimuth angle inside the camera.
2467//
2468// The inner and outer pixels are drawn separately.
2469// The general azimuth profile is fitted by a straight line
2470//
2471void MHCamera::DrawAzimuthProfile() const
2472{
2473 TProfile *obj2 = (TProfile*)AzimuthProfile(GetName());
2474 obj2->SetDirectory(0);
2475 obj2->Draw();
2476 obj2->SetBit(kCanDelete);
2477 obj2->Fit("pol0","Q","");
2478 if (obj2->GetFunction("pol0"))
2479 obj2->GetFunction("pol0")->SetLineWidth(1);
2480/*
2481 if (GetGeomCam().InheritsFrom("MGeomCamMagic"))
2482 {
2483 TArrayI inner(1);
2484 inner[0] = 0;
2485
2486 TArrayI outer(1);
2487 outer[0] = 1;
2488
2489 // Just to get the right (maximum) binning
2490 TProfile *half[2];
2491 half[0] = AzimuthProfileA(inner, MString::Format("%sInner",GetName()));
2492 half[1] = AzimuthProfileA(outer, MString::Format("%sOuter",GetName()));
2493
2494 for (Int_t i=0; i<2; i++)
2495 {
2496 half[i]->SetLineColor(kRed+i);
2497 half[i]->SetDirectory(0);
2498 half[i]->SetBit(kCanDelete);
2499 half[i]->SetMarkerSize(0.5);
2500 half[i]->Draw("same");
2501 }
2502 }
2503*/
2504}
2505
2506// --------------------------------------------------------------------------
2507//
2508// Draw the MHCamera into the MStatusDisplay:
2509//
2510// 1) Draw it as histogram (MHCamera::DrawCopy("hist")
2511// 2) Draw it as a camera, with MHCamera::SetPrettyPalette() set.
2512// 3) If "rad" is not zero, draw its values vs. the radius from the camera center.
2513// (DrawRadialProfile())
2514// 4) Depending on the variable "fit", draw the values projection on the y-axis
2515// (DrawProjection()):
2516// 0: don't draw
2517// 1: Draw fit to Single Gauss (for distributions flat-fielded over the whole camera)
2518// 2: Draw and fit to Double Gauss (for distributions different for inner and outer pixels)
2519// 3: Draw and fit to Triple Gauss (for distributions with inner, outer pixels and outliers)
2520// 4: Draw and fit to Polynomial grade 0: (for the probability distributions)
2521// >4: Draw and don;t fit.
2522//
2523void MHCamera::CamDraw(TCanvas &c, const Int_t x, const Int_t y,
2524 const Int_t fit, const Int_t rad, const Int_t azi,
2525 TObject *notify)
2526{
2527 c.cd(x);
2528 gPad->SetBorderMode(0);
2529 gPad->SetRightMargin(0.02);
2530 gPad->SetTicks();
2531 MHCamera *obj1=(MHCamera*)DrawCopy("hist");
2532 obj1->SetDirectory(NULL);
2533 obj1->SetStats(kFALSE);
2534
2535 if (notify)
2536 obj1->AddNotify(notify);
2537
2538 c.cd(x+y);
2539 gPad->SetBorderMode(0);
2540 obj1->SetPrettyPalette();
2541 obj1->Draw();
2542
2543 Int_t cnt = 2;
2544
2545 if (rad)
2546 {
2547 c.cd(x+2*y);
2548 gPad->SetBorderMode(0);
2549 gPad->SetTicks();
2550 DrawRadialProfile();
2551 cnt++;
2552 }
2553
2554 if (azi)
2555 {
2556 c.cd(x+cnt*y);
2557 gPad->SetBorderMode(0);
2558 gPad->SetTicks();
2559 DrawAzimuthProfile();
2560 cnt++;
2561 }
2562
2563 if (fit<0)
2564 return;
2565
2566 c.cd(x + cnt*y);
2567 gPad->SetBorderMode(0);
2568 gPad->SetRightMargin(0.025);
2569 gPad->SetTicks();
2570 DrawProjection(fit);
2571}
Note: See TracBrowser for help on using the repository browser.