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

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