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

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