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

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