source: tags/Mars-V0.9.2/mhist/MHCamera.cc

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