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

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