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

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