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

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