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

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