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

Last change on this file since 4301 was 4301, checked in by gaug, 20 years ago
*** empty log message ***
File size: 57.8 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
544{
545 Int_t nbins = 50;
546
547 // Create the projection histogram
548 TString pname(name);
549 if (name=="_py")
550 {
551 pname.Prepend(GetName());
552 if (sector.GetSize()>0)
553 {
554 pname += ";";
555 for (int i=0; i<sector.GetSize(); i++)
556 pname += sector[i];
557 }
558 if (aidx.GetSize()>0)
559 {
560 pname += ";";
561 for (int i=0; i<aidx.GetSize(); i++)
562 pname += aidx[i];
563 }
564 }
565
566 TH1D *h1=0;
567
568 //check if histogram with identical name exist
569 TObject *h1obj = gROOT->FindObject(pname);
570 if (h1obj && h1obj->InheritsFrom("TH1D")) {
571 h1 = (TH1D*)h1obj;
572 h1->Reset();
573 }
574
575 if (!h1)
576 {
577 Double_t min = GetMinimumSectors(sector, aidx);
578 Double_t max = GetMaximumSectors(sector, aidx);
579
580 Int_t newbins=0;
581
582 THLimitsFinder::OptimizeLimits(nbins, newbins, min, max, kFALSE);
583
584 h1 = new TH1D(pname, GetTitle(), nbins, min, max);
585 h1->SetDirectory(pname.IsNull() ? NULL : gROOT);
586 h1->SetXTitle(GetYaxis()->GetTitle());
587 h1->SetYTitle("Counts");
588 //h1->Sumw2();
589 }
590
591 // Fill the projected histogram
592 for (Int_t idx=0; idx<fNcells-2; idx++)
593 if (IsUsed(idx) && MatchSector(idx, sector, aidx))
594 h1->Fill(GetBinContent(idx+1));
595
596 return h1;
597}
598
599// ------------------------------------------------------------------------
600//
601// Creates a TH1D which contains the projection of the contents of the
602// MHCamera onto the radius from the camera center.
603// The maximum and minimum are calculated
604// such that a slighly wider range than (GetMinimum(), GetMaximum()) is
605// displayed using THLimitsFinder::OptimizeLimits.
606//
607// If no name is given the newly allocated histogram is removed from
608// the current directory calling SetDirectory(0) in any other case
609// the newly created histogram is removed from the current directory
610// and added to gROOT such the gROOT->FindObject can find the histogram.
611//
612// If the standard name "_rad" is given "_rad" is appended to the name
613// of the MHCamera and the corresponding histogram is searched using
614// gROOT->FindObject and updated with the present projection.
615//
616// It is the responsibility of the user to make sure, that the newly
617// created histogram is freed correctly.
618//
619// Currently the new histogram is restrictred to 50 bins.
620// Maybe a optimal number can be calulated from the number of
621// bins on the x-axis of the MHCamera?
622//
623// The code was taken mainly from TH2::ProjectX such the interface
624// is more or less the same than to TH2-projections.
625//
626// If sector>=0 only entries with matching sector index are taken
627// into account.
628//
629TProfile *MHCamera::RadialProfileS(const TArrayI &sector, const TArrayI &aidx, const char *name, const Int_t nbins) const
630{
631
632 // Create the projection histogram
633 TString pname(name);
634 if (name=="_rad")
635 {
636 pname.Prepend(GetName());
637 if (sector.GetSize()>0)
638 {
639 pname += ";";
640 for (int i=0; i<sector.GetSize(); i++)
641 pname += sector[i];
642 }
643 if (aidx.GetSize()>0)
644 {
645 pname += ";";
646 for (int i=0; i<aidx.GetSize(); i++)
647 pname += aidx[i];
648 }
649 }
650
651 TProfile *h1=0;
652
653 //check if histogram with identical name exist
654 TObject *h1obj = gROOT->FindObject(pname);
655 if (h1obj && h1obj->InheritsFrom("TProfile")) {
656 h1 = (TProfile*)h1obj;
657 h1->Reset();
658 }
659
660 if (!h1)
661 {
662
663 Double_t min = 0.;
664 Double_t max = fGeomCam->GetMaxRadius();
665
666 Int_t newbins=0;
667
668 THLimitsFinder::OptimizeLimits(nbins, newbins, min, max, kFALSE);
669
670 h1 = new TProfile(pname, GetTitle(), nbins, min, max);
671 h1->SetDirectory(pname.IsNull() ? NULL : gROOT);
672 h1->SetXTitle("Radius from camera center [mm]");
673 h1->SetYTitle(GetYaxis()->GetTitle());
674 }
675
676 // Fill the projected histogram
677 for (Int_t idx=0; idx<fNcells-2; idx++)
678 if (IsUsed(idx) && MatchSector(idx, sector, aidx))
679 h1->Fill(TMath::Hypot((*fGeomCam)[idx].GetX(),(*fGeomCam)[idx].GetY()),
680 GetBinContent(idx+1));
681 return h1;
682}
683
684
685// ------------------------------------------------------------------------
686//
687// Creates a TH1D which contains the projection of the contents of the
688// MHCamera onto the azimuth angle in the camera.
689//
690// If no name is given the newly allocated histogram is removed from
691// the current directory calling SetDirectory(0) in any other case
692// the newly created histogram is removed from the current directory
693// and added to gROOT such the gROOT->FindObject can find the histogram.
694//
695// If the standard name "_azi" is given "_azi" is appended to the name
696// of the MHCamera and the corresponding histogram is searched using
697// gROOT->FindObject and updated with the present projection.
698//
699// It is the responsibility of the user to make sure, that the newly
700// created histogram is freed correctly.
701//
702// Currently the new histogram is restrictred to 60 bins.
703// Maybe a optimal number can be calulated from the number of
704// bins on the x-axis of the MHCamera?
705//
706// The code was taken mainly from TH2::ProjectX such the interface
707// is more or less the same than to TH2-projections.
708//
709TProfile *MHCamera::AzimuthProfileA(const TArrayI &aidx, const char *name, const Int_t nbins) const
710{
711
712 // Create the projection histogram
713 TString pname(name);
714 if (name=="_azi")
715 {
716 pname.Prepend(GetName());
717 if (aidx.GetSize()>0)
718 {
719 pname += ";";
720 for (int i=0; i<aidx.GetSize(); i++)
721 pname += aidx[i];
722 }
723 }
724
725 TProfile *h1=0;
726
727 //check if histogram with identical name exist
728 TObject *h1obj = gROOT->FindObject(pname);
729 if (h1obj && h1obj->InheritsFrom("TProfile")) {
730 h1 = (TProfile*)h1obj;
731 h1->Reset();
732 }
733
734 if (!h1)
735 {
736
737 Double_t min = -0.5;
738 Double_t max = 359.5;
739
740 Int_t newbins=0;
741
742 THLimitsFinder::OptimizeLimits(nbins, newbins, min, max, kFALSE);
743
744 h1 = new TProfile(pname, GetTitle(), nbins, min, max);
745 h1->SetDirectory(pname.IsNull() ? NULL : gROOT);
746 h1->SetXTitle("Azimuth in camera [deg]");
747 h1->SetYTitle(GetYaxis()->GetTitle());
748 }
749
750 // Fill the projected histogram
751 for (Int_t idx=0; idx<fNcells-2; idx++)
752 {
753 if (IsUsed(idx) && MatchSector(idx, TArrayI(), aidx))
754 h1->Fill(TMath::ATan2((*fGeomCam)[idx].GetY(),(*fGeomCam)[idx].GetX())*180./TMath::Pi()+180.,
755 GetPixContent(idx));
756
757 }
758
759 return h1;
760}
761
762
763// ------------------------------------------------------------------------
764//
765// Resizes the current pad so that the camera is displayed in its
766// correct aspect ratio
767//
768void MHCamera::SetRange()
769{
770 const Float_t range = fGeomCam->GetMaxRadius()*1.05;
771
772 //
773 // Maintain aspect ratio
774 //
775 const float ratio = TestBit(kNoLegend) ? 1 : 1.15;
776
777 //
778 // Calculate width and height of the current pad in pixels
779 //
780 Float_t w = gPad->GetWw();
781 Float_t h = gPad->GetWh()*ratio;
782
783 //
784 // This prevents the pad from resizing itself wrongly
785 //
786 if (gPad->GetMother() != gPad)
787 {
788 w *= gPad->GetMother()->GetAbsWNDC();
789 h *= gPad->GetMother()->GetAbsHNDC();
790 }
791
792 //
793 // Set Range (coordinate system) of pad
794 //
795 gPad->Range(-range, -range, (2*ratio-1)*range, range);
796
797 //
798 // Resize Pad to given ratio
799 //
800 if (h<w)
801 gPad->SetPad((1.-h/w)/2, 0, (h/w+1.)/2, 1);
802 else
803 gPad->SetPad(0, (1.-w/h)/2, 1, (w/h+1.)/2);
804}
805
806// ------------------------------------------------------------------------
807//
808// Updates the pixel colors and paints the pixels
809//
810void MHCamera::Update(Bool_t islog, Bool_t isbox, Bool_t iscol, Bool_t issame)
811{
812 Double_t min = GetMinimum(kFALSE);
813 Double_t max = GetMaximum(kFALSE);
814 if (min==FLT_MAX)
815 {
816 min = 0;
817 max = 1;
818 }
819
820 if (min==max)
821 max += 1;
822
823 if (!issame)
824 UpdateLegend(min, max, islog);
825
826 // Try to estimate the units of the current display. This is only
827 // necessary for 'same' option and allows distorted images of the camera!
828 const Float_t maxr = (1-fGeomCam->GetConvMm2Deg())*fGeomCam->GetMaxRadius()/2;
829 const Float_t conv = !issame ||
830 gPad->GetX1()<-maxr || gPad->GetY1()<-maxr ||
831 gPad->GetX2()> maxr || gPad->GetY2()>maxr ? 1 : fGeomCam->GetConvMm2Deg();
832
833 MHexagon hex;
834 for (Int_t i=0; i<fNcells-2; i++)
835 {
836 hex.SetFillStyle(issame ? 0 : 1001);
837
838 if (!issame)
839 {
840 if (IsUsed(i) && iscol)
841 {
842 if (TMath::IsNaN(fArray[i+1]))
843 gLog << warn << "MHCamera::Update: " << GetName() << " <" << GetTitle() << "> - Pixel Index #" << i << " contents is NaN (Not a Number)..." << endl;
844
845 hex.SetFillColor(GetColor(GetBinContent(i+1), min, max, islog));
846 }
847 else
848 hex.SetFillColor(10);
849 }
850
851 const MGeomPix &pix = (*fGeomCam)[i];
852
853 Float_t x = pix.GetX()*conv/(fAbberation+1);
854 Float_t y = pix.GetY()*conv/(fAbberation+1);
855 Float_t d = pix.GetD()*conv;
856
857 if (!isbox)
858 hex.PaintHexagon(x, y, d);
859 else
860 if (IsUsed(i) && !TMath::IsNaN(fArray[i+1]))
861 {
862 Float_t size = d*(GetBinContent(i+1)-min)/(max-min);
863 if (size>d)
864 size=d;
865 hex.PaintHexagon(x, y, size);
866 }
867 }
868}
869
870// ------------------------------------------------------------------------
871//
872// Print minimum and maximum
873//
874void MHCamera::Print(Option_t *) const
875{
876 gLog << all << "Minimum: " << GetMinimum();
877 if (fMinimum==-1111)
878 gLog << " <autoscaled>";
879 gLog << endl;
880 gLog << "Maximum: " << GetMaximum();
881 if (fMaximum==-1111)
882 gLog << " <autoscaled>";
883 gLog << endl;
884}
885
886// ------------------------------------------------------------------------
887//
888// Paint the y-axis title
889//
890void MHCamera::PaintAxisTitle()
891{
892 const Float_t range = fGeomCam->GetMaxRadius()*1.05;
893 const Float_t w = (1 + 1.5/sqrt((float)(fNcells-2)))*range;
894
895 TLatex *ptitle = new TLatex(w, -.90*range, GetYaxis()->GetTitle());
896
897 ptitle->SetTextSize(0.05);
898 ptitle->SetTextAlign(21);
899
900 // box with the histogram title
901 ptitle->SetTextColor(gStyle->GetTitleTextColor());
902#if ROOT_VERSION_CODE > ROOT_VERSION(3,05,01)
903 ptitle->SetTextFont(gStyle->GetTitleFont(""));
904#endif
905 ptitle->Paint();
906}
907
908// ------------------------------------------------------------------------
909//
910// Paints the camera.
911//
912void MHCamera::Paint(Option_t *o)
913{
914 if (fNcells<=1)
915 return;
916
917 TString opt(o);
918 opt.ToLower();
919
920 if (opt.Contains("hist"))
921 {
922 opt.ReplaceAll("hist", "");
923 TH1D::Paint(opt);
924 return;
925 }
926
927 if (opt.Contains("proj"))
928 {
929 opt.ReplaceAll("proj", "");
930 Projection(GetName())->Paint(opt);
931 return;
932 }
933
934 const Bool_t hassame = opt.Contains("same");
935 const Bool_t hasbox = opt.Contains("box");
936 const Bool_t hascol = hasbox ? !opt.Contains("nocol") : kTRUE;
937
938 if (!hassame)
939 {
940 gPad->Clear();
941
942 // Maintain aspect ratio
943 SetRange();
944
945 if (GetPainter())
946 {
947 // Paint statistics
948 if (!TestBit(TH1::kNoStats))
949 fPainter->PaintStat(gStyle->GetOptStat(), NULL);
950
951 // Paint primitives (pixels, color legend, photons, ...)
952 if (fPainter->InheritsFrom(THistPainter::Class()))
953 {
954 static_cast<THistPainter*>(fPainter)->MakeChopt("");
955 static_cast<THistPainter*>(fPainter)->PaintTitle();
956 }
957 }
958 }
959
960 // Update Contents of the pixels and paint legend
961 Update(gPad->GetLogy(), hasbox, hascol, hassame);
962
963 if (!hassame)
964 PaintAxisTitle();
965
966 if (opt.Contains("pixelindex"))
967 PaintIndices(0);
968 if (opt.Contains("sectorindex"))
969 PaintIndices(1);
970 if (opt.Contains("content"))
971 PaintIndices(2);
972}
973
974// ------------------------------------------------------------------------
975//
976// With this function you can change the color palette. For more
977// information see TStyle::SetPalette. Only palettes with 50 colors
978// are allowed.
979// In addition you can use SetPalette(52, 0) to create an inverse
980// deep blue sea palette
981//
982void MHCamera::SetPalette(Int_t ncolors, Int_t *colors)
983{
984 //
985 // If not enough colors are specified skip this.
986 //
987 if (ncolors>1 && ncolors<50)
988 {
989 gLog << err << "MHCamera::SetPalette: Only default palettes with 50 colors are allowed... ignored." << endl;
990 return;
991 }
992
993 //
994 // If ncolors==52 create a reversed deep blue sea palette
995 //
996 if (ncolors==52)
997 {
998 gStyle->SetPalette(51, NULL);
999 TArrayI c(kItemsLegend);
1000 for (int i=0; i<kItemsLegend; i++)
1001 c[kItemsLegend-i-1] = gStyle->GetColorPalette(i);
1002 gStyle->SetPalette(kItemsLegend, c.GetArray());
1003 }
1004 else
1005 gStyle->SetPalette(ncolors, colors);
1006
1007 fColors.Set(kItemsLegend);
1008 for (int i=0; i<kItemsLegend; i++)
1009 fColors[i] = gStyle->GetColorPalette(i);
1010}
1011
1012
1013// ------------------------------------------------------------------------
1014//
1015// Changes the palette of the displayed camera histogram.
1016//
1017// Change to the right pad first - otherwise GetDrawOption() might fail.
1018//
1019void MHCamera::SetPrettyPalette()
1020{
1021 if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1022 SetPalette(1, 0);
1023}
1024
1025// ------------------------------------------------------------------------
1026//
1027// Changes the palette of the displayed camera histogram.
1028//
1029// Change to the right pad first - otherwise GetDrawOption() might fail.
1030//
1031void MHCamera::SetDeepBlueSeaPalette()
1032{
1033 if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1034 SetPalette(51, 0);
1035}
1036
1037// ------------------------------------------------------------------------
1038//
1039// Changes the palette of the displayed camera histogram.
1040//
1041// Change to the right pad first - otherwise GetDrawOption() might fail.
1042//
1043void MHCamera::SetInvDeepBlueSeaPalette()
1044{
1045 if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1046 SetPalette(52, 0);
1047}
1048
1049// ------------------------------------------------------------------------
1050//
1051// Paint indices (as text) inside the pixels. Depending of the type-
1052// argument we paint:
1053// 0: pixel number
1054// 1: sector number
1055// 2: content
1056//
1057void MHCamera::PaintIndices(Int_t type)
1058{
1059 if (fNcells<=1)
1060 return;
1061
1062 const Double_t min = GetMinimum();
1063 const Double_t max = GetMaximum();
1064
1065 if (type==2 && max==min)
1066 return;
1067
1068 TText txt;
1069 txt.SetTextFont(122);
1070 txt.SetTextAlign(22); // centered/centered
1071
1072 for (Int_t i=0; i<fNcells-2; i++)
1073 {
1074 const MGeomPix &h = (*fGeomCam)[i];
1075
1076 TString num;
1077 switch (type)
1078 {
1079 case 0: num += i; break;
1080 case 1: num += h.GetSector(); break;
1081 case 2: num += (Int_t)((fArray[i+1]-min)/(max-min)); break;
1082 }
1083
1084 // FIXME: Should depend on the color of the pixel...
1085 //(GetColor(GetBinContent(i+1), min, max, 0));
1086 txt.SetTextColor(kRed);
1087 txt.SetTextSize(0.3*h.GetD()/fGeomCam->GetMaxRadius()/1.05);
1088 txt.PaintText(h.GetX(), h.GetY(), num);
1089 }
1090}
1091
1092// ------------------------------------------------------------------------
1093//
1094// Call this function to add a MCamEvent on top of the present contents.
1095//
1096void MHCamera::AddCamContent(const MCamEvent &event, Int_t type)
1097{
1098 if (fNcells<=1 || IsFreezed())
1099 return;
1100
1101 // FIXME: Security check missing!
1102 for (Int_t idx=0; idx<fNcells-2; idx++)
1103 {
1104 Double_t val=0;
1105 if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
1106 SetUsed(idx);
1107
1108 Fill(idx, val); // FIXME: Slow!
1109 }
1110 fEntries++;
1111}
1112
1113// ------------------------------------------------------------------------
1114//
1115// Call this function to add a MCamEvent on top of the present contents.
1116//
1117void MHCamera::SetCamError(const MCamEvent &evt, Int_t type)
1118{
1119
1120 if (fNcells<=1 || IsFreezed())
1121 return;
1122
1123 // FIXME: Security check missing!
1124 for (Int_t idx=0; idx<fNcells-2; idx++)
1125 {
1126 Double_t val=0;
1127 if (evt.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
1128 SetUsed(idx);
1129
1130 SetBinError(idx+1, val); // FIXME: Slow!
1131 }
1132}
1133
1134Stat_t MHCamera::GetBinError(Int_t bin) const
1135{
1136 Double_t rc = 0;
1137 if (TestBit(kVariance) && GetEntries()>0) // error on the mean
1138 {
1139 const Double_t error = fSumw2.fArray[bin]/GetEntries();
1140 const Double_t val = fArray[bin]/GetEntries();
1141 rc = TMath::Sqrt(error - val*val);
1142 }
1143 else
1144 {
1145 rc = TH1D::GetBinError(bin);
1146 rc = Profile(rc);
1147 }
1148
1149 return rc;
1150}
1151
1152// ------------------------------------------------------------------------
1153//
1154// Call this function to add a MHCamera on top of the present contents.
1155// Type:
1156// 0) bin content
1157// 1) errors
1158// 2) rel. errors
1159//
1160void MHCamera::AddCamContent(const MHCamera &d, Int_t type)
1161{
1162 if (fNcells!=d.fNcells || IsFreezed())
1163 return;
1164
1165 // FIXME: Security check missing!
1166 for (Int_t idx=0; idx<fNcells-2; idx++)
1167 if (d.IsUsed(idx))
1168 SetUsed(idx);
1169
1170 switch (type)
1171 {
1172 case 1:
1173 for (Int_t idx=0; idx<fNcells-2; idx++)
1174 Fill(idx, d.GetBinError(idx+1));
1175 break;
1176 case 2:
1177 for (Int_t idx=0; idx<fNcells-2; idx++)
1178 if (d.GetBinContent(idx+1)!=0)
1179 Fill(idx, TMath::Abs(d.GetBinError(idx+1)/d.GetBinContent(idx+1)));
1180 break;
1181 default:
1182 for (Int_t idx=0; idx<fNcells-2; idx++)
1183 Fill(idx, d.GetBinContent(idx+1));
1184 break;
1185 }
1186 fEntries++;
1187}
1188
1189// ------------------------------------------------------------------------
1190//
1191// Call this function to add a TArrayD on top of the present contents.
1192//
1193void MHCamera::AddCamContent(const TArrayD &event, const TArrayC *used)
1194{
1195 if (event.GetSize()!=fNcells-2 || IsFreezed())
1196 return;
1197
1198 if (used && used->GetSize()!=fNcells-2)
1199 return;
1200
1201 for (Int_t idx=0; idx<fNcells-2; idx++)
1202 {
1203 Fill(idx, const_cast<TArrayD&>(event)[idx]); // FIXME: Slow!
1204
1205 if (!used || (*const_cast<TArrayC*>(used))[idx])
1206 SetUsed(idx);
1207 }
1208 fEntries++;
1209}
1210
1211// ------------------------------------------------------------------------
1212//
1213// Call this function to add a MCamEvent on top of the present contents.
1214// 1 is added to each pixel if the contents of MCamEvent>threshold
1215//
1216void MHCamera::CntCamContent(const MCamEvent &event, Double_t threshold, Int_t type)
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)
1229 Fill(idx);
1230 }
1231 fEntries++;
1232}
1233
1234// ------------------------------------------------------------------------
1235//
1236// Call this function to add a TArrayD on top of the present contents.
1237// 1 is added to each pixel if the contents of MCamEvent>threshold
1238//
1239void MHCamera::CntCamContent(const TArrayD &event, Double_t threshold, Bool_t ispos)
1240{
1241 if (event.GetSize()!=fNcells-2 || IsFreezed())
1242 return;
1243
1244 for (Int_t idx=0; idx<fNcells-2; idx++)
1245 {
1246 if (const_cast<TArrayD&>(event)[idx]>threshold)
1247 Fill(idx);
1248
1249 if (!ispos || fArray[idx+1]>0)
1250 SetUsed(idx);
1251 }
1252 fEntries++;
1253}
1254
1255// ------------------------------------------------------------------------
1256//
1257// Fill the pixels with random contents.
1258//
1259void MHCamera::FillRandom()
1260{
1261 if (fNcells<=1 || IsFreezed())
1262 return;
1263
1264 Reset();
1265
1266 // FIXME: Security check missing!
1267 for (Int_t idx=0; idx<fNcells-2; idx++)
1268 {
1269 Fill(idx, gRandom->Uniform()*fGeomCam->GetPixRatio(idx));
1270 SetUsed(idx);
1271 }
1272 fEntries=1;
1273}
1274
1275
1276// ------------------------------------------------------------------------
1277//
1278// The array must be in increasing order, eg: 2.5, 3.7, 4.9
1279// The values in each bin are replaced by the interval in which the value
1280// fits. In the example we have four intervals
1281// (<2.5, 2.5-3.7, 3.7-4.9, >4.9). Maximum and minimum are set
1282// accordingly.
1283//
1284void MHCamera::SetLevels(const TArrayF &arr)
1285{
1286 if (fNcells<=1)
1287 return;
1288
1289 for (Int_t i=0; i<fNcells-2; i++)
1290 {
1291 if (!IsUsed(i))
1292 continue;
1293
1294 Int_t j = arr.GetSize();
1295 while (j && fArray[i+1]<arr[j-1])
1296 j--;
1297
1298 fArray[i+1] = j;
1299 }
1300 SetMaximum(arr.GetSize());
1301 SetMinimum(0);
1302}
1303
1304// ------------------------------------------------------------------------
1305//
1306// Reset the all pixel colors to a default value
1307//
1308void MHCamera::Reset(Option_t *opt)
1309{
1310 if (fNcells<=1 || IsFreezed())
1311 return;
1312
1313 TH1::Reset(opt);
1314
1315 for (Int_t i=0; i<fNcells-2; i++)
1316 {
1317 fArray[i+1]=0;
1318 ResetUsed(i);
1319 }
1320 fArray[0] = 0;
1321 fArray[fNcells-1] = 0;
1322}
1323
1324// ------------------------------------------------------------------------
1325//
1326// Here we calculate the color index for the current value.
1327// The color index is defined with the class TStyle and the
1328// Color palette inside. We use the command gStyle->SetPalette(1,0)
1329// for the display. So we have to convert the value "wert" into
1330// a color index that fits the color palette.
1331// The range of the color palette is defined by the values fMinPhe
1332// and fMaxRange. Between this values we have 50 color index, starting
1333// with 0 up to 49.
1334//
1335Int_t MHCamera::GetColor(Float_t val, Float_t min, Float_t max, Bool_t islog)
1336{
1337 if (TMath::IsNaN(val)) // FIXME: gLog!
1338 return 10;
1339
1340 //
1341 // first treat the over- and under-flows
1342 //
1343 const Int_t maxcolidx = kItemsLegend-1;
1344
1345 if (val >= max)
1346 return fColors[maxcolidx];
1347
1348 if (val <= min)
1349 return fColors[0];
1350
1351 //
1352 // calculate the color index
1353 //
1354 Float_t ratio;
1355 if (islog && min>0)
1356 ratio = log10(val/min) / log10(max/min);
1357 else
1358 ratio = (val-min) / (max-min);
1359
1360 const Int_t colidx = (Int_t)(ratio*maxcolidx + .5);
1361 return fColors[colidx];
1362}
1363
1364TPaveStats *MHCamera::GetStatisticBox()
1365{
1366 TObject *obj = 0;
1367
1368 TIter Next(fFunctions);
1369 while ((obj = Next()))
1370 if (obj->InheritsFrom(TPaveStats::Class()))
1371 return static_cast<TPaveStats*>(obj);
1372
1373 return NULL;
1374}
1375
1376// ------------------------------------------------------------------------
1377//
1378// Change the text on the legend according to the range of the Display
1379//
1380void MHCamera::UpdateLegend(Float_t min, Float_t max, Bool_t islog)
1381{
1382 const Float_t range = fGeomCam->GetMaxRadius()*1.05;
1383
1384 TArrow arr;
1385 arr.PaintArrow(-range*.9, -range*.9, -range*.6, -range*.9, 0.025);
1386 arr.PaintArrow(-range*.9, -range*.9, -range*.9, -range*.6, 0.025);
1387
1388 TString text;
1389 text += (int)(range*.3);
1390 text += "mm";
1391
1392 TText newtxt2;
1393 newtxt2.SetTextSize(0.04);
1394 newtxt2.PaintText(-range*.85, -range*.85, text);
1395
1396 text = "";
1397 text += (float)((int)(range*.3*fGeomCam->GetConvMm2Deg()*10))/10;
1398 text += "\\circ";
1399 text = text.Strip(TString::kLeading);
1400
1401 TLatex latex;
1402 latex.PaintLatex(-range*.85, -range*.75, 0, 0.04, text);
1403
1404 if (TestBit(kNoLegend))
1405 return;
1406
1407 TPaveStats *stats = GetStatisticBox();
1408
1409 const Float_t hndc = 0.92 - (stats ? stats->GetY1NDC() : 1);
1410 const Float_t H = (0.75-hndc)*range;
1411 const Float_t offset = hndc*range;
1412
1413 const Float_t h = 2./kItemsLegend;
1414 const Float_t w = range/sqrt((float)(fNcells-2));
1415
1416 TBox newbox;
1417 TText newtxt;
1418 newtxt.SetTextSize(0.03);
1419 newtxt.SetTextAlign(12);
1420#if ROOT_VERSION_CODE > ROOT_VERSION(3,01,06)
1421 newtxt.SetBit(/*kNoContextMenu|*/kCannotPick);
1422 newbox.SetBit(/*kNoContextMenu|*/kCannotPick);
1423#endif
1424
1425 const Float_t step = (islog && min>0 ? log10(max/min) : max-min) / kItemsLegend;
1426 const Int_t firsts = step*3 < 1e-8 ? 8 : (Int_t)floor(log10(step*3));
1427 const TString opt = Form("%%.%if", firsts>0 ? 0 : TMath::Abs(firsts));
1428
1429 for (Int_t i=0; i<kItemsLegend+1; i+=3)
1430 {
1431 Float_t val;
1432 if (islog && min>0)
1433 val = pow(10, step*i) * min;
1434 else
1435 val = min + step*i;
1436
1437 //const bool dispexp = max-min>1.5 && fabs(val)>0.1 && fabs(val)<1e6;
1438 newtxt.PaintText(range+1.5*w, H*(i*h-1)-offset, Form(opt, val));
1439 }
1440
1441 for (Int_t i=0; i<kItemsLegend; i++)
1442 {
1443 newbox.SetFillColor(fColors[i]);
1444 newbox.PaintBox(range, H*(i*h-1)-offset, range+w, H*((i+1)*h-1)-offset);
1445 }
1446}
1447
1448// ------------------------------------------------------------------------
1449//
1450// Save primitive as a C++ statement(s) on output stream out
1451//
1452void MHCamera::SavePrimitive(ofstream &out, Option_t *opt)
1453{
1454 gLog << err << "MHCamera::SavePrimitive: Must be rewritten!" << endl;
1455 /*
1456 if (!gROOT->ClassSaved(TCanvas::Class()))
1457 fDrawingPad->SavePrimitive(out, opt);
1458
1459 out << " " << fDrawingPad->GetName() << "->SetWindowSize(";
1460 out << fDrawingPad->GetWw() << "," << fDrawingPad->GetWh() << ");" << endl;
1461 */
1462}
1463
1464// ------------------------------------------------------------------------
1465//
1466// compute the distance of a point (px,py) to the Camera
1467// this functions needed for graphical primitives, that
1468// means without this function you are not able to interact
1469// with the graphical primitive with the mouse!!!
1470//
1471// All calcutations are done in pixel coordinates
1472//
1473Int_t MHCamera::DistancetoPrimitive(Int_t px, Int_t py)
1474{
1475 if (fNcells<=1)
1476 return 999999;
1477
1478 TPaveStats *box = (TPaveStats*)gPad->GetPrimitive("stats");
1479 if (box)
1480 {
1481 const Double_t w = box->GetY2NDC()-box->GetY1NDC();
1482 box->SetX1NDC(gStyle->GetStatX()-gStyle->GetStatW());
1483 box->SetY1NDC(gStyle->GetStatY()-w);
1484 box->SetX2NDC(gStyle->GetStatX());
1485 box->SetY2NDC(gStyle->GetStatY());
1486 }
1487
1488 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1489 return TH1D::DistancetoPrimitive(px, py);
1490
1491 const Bool_t issame = TString(GetDrawOption()).Contains("same", TString::kIgnoreCase);
1492
1493 const Float_t maxr = (1-fGeomCam->GetConvMm2Deg())*fGeomCam->GetMaxRadius()/2;
1494 const Float_t conv = !issame ||
1495 gPad->GetX1()<-maxr || gPad->GetY1()<-maxr ||
1496 gPad->GetX2()> maxr || gPad->GetY2()>maxr ? 1 : fGeomCam->GetConvMm2Deg();
1497
1498 if (GetPixelIndex(px, py, conv)>=0)
1499 return 0;
1500
1501 if (!box)
1502 return 999999;
1503
1504 const Int_t dist = box->DistancetoPrimitive(px, py);
1505 if (dist > TPad::GetMaxPickDistance())
1506 return 999999;
1507
1508 gPad->SetSelected(box);
1509 return dist;
1510}
1511
1512// ------------------------------------------------------------------------
1513//
1514//
1515Int_t MHCamera::GetPixelIndex(Int_t px, Int_t py, Float_t conv) const
1516{
1517 if (fNcells<=1)
1518 return -1;
1519
1520 Int_t i;
1521 for (i=0; i<fNcells-2; i++)
1522 {
1523 MHexagon hex((*fGeomCam)[i]);
1524 if (hex.DistancetoPrimitive(px, py, conv)>0)
1525 continue;
1526
1527 return i;
1528 }
1529 return -1;
1530}
1531
1532// ------------------------------------------------------------------------
1533//
1534// Returns string containing info about the object at position (px,py).
1535// Returned string will be re-used (lock in MT environment).
1536//
1537char *MHCamera::GetObjectInfo(Int_t px, Int_t py) const
1538{
1539 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1540 return TH1D::GetObjectInfo(px, py);
1541
1542 static char info[128];
1543
1544 const Int_t idx=GetPixelIndex(px, py);
1545
1546 if (idx<0)
1547 return TObject::GetObjectInfo(px, py);
1548
1549 sprintf(info, "Software Pixel Index: %d (Hardware Id=%d)", idx, idx+1);
1550 return info;
1551}
1552
1553// ------------------------------------------------------------------------
1554//
1555// Add a MCamEvent which should be displayed when the user clicks on a
1556// pixel.
1557// Warning: The object MUST inherit from TObject AND MCamEvent
1558//
1559void MHCamera::AddNotify(TObject *obj)
1560{
1561 // Make sure, that the object derives from MCamEvent!
1562 MCamEvent *evt = dynamic_cast<MCamEvent*>(obj);
1563 if (!evt)
1564 {
1565 gLog << err << "ERROR: MHCamera::AddNotify - TObject doesn't inherit from MCamEvent... ignored." << endl;
1566 return;
1567 }
1568
1569 // Make sure, that it is deleted from the list too, if the obj is deleted
1570 obj->SetBit(kMustCleanup);
1571
1572 // Add object to list
1573 fNotify->Add(obj);
1574}
1575
1576// ------------------------------------------------------------------------
1577//
1578// Execute a mouse event on the camera
1579//
1580void MHCamera::ExecuteEvent(Int_t event, Int_t px, Int_t py)
1581{
1582 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1583 {
1584 TH1D::ExecuteEvent(event, px, py);
1585 return;
1586 }
1587 //if (event==kMouseMotion && fStatusBar)
1588 // fStatusBar->SetText(GetObjectInfo(px, py), 0);
1589 if (event!=kButton1Down)
1590 return;
1591
1592 const Int_t idx = GetPixelIndex(px, py);
1593 if (idx<0)
1594 return;
1595
1596 gLog << all << GetTitle() << " <" << GetName() << ">" << endl;
1597 gLog << "Software Pixel Idx: " << idx << endl;
1598 gLog << "Hardware Pixel Id: " << idx+1 << endl;
1599 gLog << "Contents: " << GetBinContent(idx+1);
1600 if (GetBinError(idx+1)>0)
1601 gLog << " +/- " << GetBinError(idx+1);
1602 gLog << " <" << (IsUsed(idx)?"on":"off") << ">" << endl;
1603
1604 if (fNotify && fNotify->GetSize()>0)
1605 {
1606 // FIXME: Is there a simpler and more convinient way?
1607
1608 // The name which is created here depends on the instance of
1609 // MHCamera and on the pad on which it is drawn --> The name
1610 // is unique. For ExecuteEvent gPad is always correctly set.
1611 const TString name = Form("%p;%p;PixelContent", this, gPad);
1612
1613 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
1614 if (old)
1615 old->cd();
1616 else
1617 new TCanvas(name);
1618
1619 /*
1620 TIter Next(gPad->GetListOfPrimitives());
1621 TObject *o;
1622 while (o=Next()) cout << o << ": " << o->GetName() << " " << o->IsA()->GetName() << endl;
1623 */
1624
1625 // FIXME: Make sure, that the old histograms are really deleted.
1626 // Are they already deleted?
1627
1628 // The dynamic_cast is necessary here: We cannot use ForEach
1629 TIter Next(fNotify);
1630 MCamEvent *evt;
1631 while ((evt=dynamic_cast<MCamEvent*>(Next())))
1632 evt->DrawPixelContent(idx);
1633
1634 gPad->Modified();
1635 gPad->Update();
1636 }
1637}
1638
1639UInt_t MHCamera::GetNumPixels() const
1640{
1641 return fGeomCam->GetNumPixels();
1642}
1643
1644TH1 *MHCamera::DrawCopy() const
1645{
1646 gPad=NULL;
1647 return TH1D::DrawCopy(fName+";cpy");
1648}
1649
1650
1651// --------------------------------------------------------------------------
1652//
1653// Draw a projection of MHCamera onto the y-axis values. Depending on the
1654// variable fit, the following fits are performed:
1655//
1656// 0: No fit, simply draw the projection
1657// 1: Single Gauss (for distributions flat-fielded over the whole camera)
1658// 2: Double Gauss (for distributions different for inner and outer pixels)
1659// 3: Triple Gauss (for distributions with inner, outer pixels and outliers)
1660// 4: flat (for the probability distributions)
1661// (1-4:) Moreover, sectors 6,1 and 2 of the camera and sectors 3,4 and 5 are
1662// drawn separately, for inner and outer pixels.
1663// 5: Fit Inner and Outer pixels separately by a single Gaussian
1664// (only for MAGIC cameras)
1665// 6: Fit Inner and Outer pixels separately by a single Gaussian and display
1666// additionally the two camera halfs separately (for MAGIC camera)
1667// 7: Single Gauss with TLegend to show the meaning of the colours
1668//
1669//
1670void MHCamera::DrawProjection(Int_t fit) const
1671{
1672
1673 TArrayI inner(1);
1674 inner[0] = 0;
1675
1676 TArrayI outer(1);
1677 outer[0] = 1;
1678
1679 if (fit==5 || fit==6)
1680 {
1681
1682 if (GetGeomCam().InheritsFrom("MGeomCamMagic"))
1683 {
1684 TArrayI s0(6);
1685 s0[0] = 6;
1686 s0[1] = 1;
1687 s0[2] = 2;
1688 s0[3] = 3;
1689 s0[4] = 4;
1690 s0[5] = 5;
1691
1692 TArrayI s1(3);
1693 s1[0] = 6;
1694 s1[1] = 1;
1695 s1[2] = 2;
1696
1697 TArrayI s2(3);
1698 s2[0] = 3;
1699 s2[1] = 4;
1700 s2[2] = 5;
1701
1702 gPad->Clear();
1703 TVirtualPad *pad = gPad;
1704 pad->Divide(2,1);
1705
1706 TH1D *inout[2];
1707 inout[0] = ProjectionS(s0, inner, "Inner");
1708 inout[1] = ProjectionS(s0, outer, "Outer");
1709
1710 inout[0]->SetDirectory(NULL);
1711 inout[1]->SetDirectory(NULL);
1712
1713 for (int i=0; i<2; i++)
1714 {
1715 pad->cd(i+1);
1716 inout[i]->SetLineColor(kRed+i);
1717 inout[i]->SetBit(kCanDelete);
1718 inout[i]->Draw();
1719 inout[i]->Fit("gaus","Q");
1720
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 TLegend *leg = new TLegend(0.05,0.65,0.35,0.9);
1781
1782 for (int i=0; i<4; i++)
1783 {
1784 halfInOut[i]->SetLineColor(kRed+i);
1785 halfInOut[i]->SetDirectory(0);
1786 halfInOut[i]->SetBit(kCanDelete);
1787 halfInOut[i]->Draw("same");
1788 leg->AddEntry(halfInOut[i],halfInOut[i]->GetTitle(),"l");
1789 }
1790
1791 if (fit==7)
1792 leg->Draw();
1793
1794 gPad->Modified();
1795 gPad->Update();
1796 }
1797
1798 const Double_t min = obj2->GetBinCenter(obj2->GetXaxis()->GetFirst());
1799 const Double_t max = obj2->GetBinCenter(obj2->GetXaxis()->GetLast());
1800 const Double_t integ = obj2->Integral("width")/2.5;
1801 const Double_t mean = obj2->GetMean();
1802 const Double_t rms = obj2->GetRMS();
1803 const Double_t width = max-min;
1804
1805 const TString dgausformula = "([0]-[3])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
1806 "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])";
1807
1808 const TString tgausformula = "([0]-[3]-[6])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
1809 "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])"
1810 "+[6]/[8]*exp(-0.5*(x-[7])*(x-[7])/[8]/[8])";
1811 TF1 *f=0;
1812 switch (fit)
1813 {
1814 case 1:
1815 f = new TF1("sgaus", "gaus(0)", min, max);
1816 f->SetLineColor(kYellow);
1817 f->SetBit(kCanDelete);
1818 f->SetParNames("Area", "#mu", "#sigma");
1819 f->SetParameters(integ/rms, mean, rms);
1820 f->SetParLimits(0, 0, integ);
1821 f->SetParLimits(1, min, max);
1822 f->SetParLimits(2, 0, width/1.5);
1823
1824 obj2->Fit(f, "QLR");
1825 break;
1826
1827 case 2:
1828 f = new TF1("dgaus",dgausformula.Data(),min,max);
1829 f->SetLineColor(kYellow);
1830 f->SetBit(kCanDelete);
1831 f->SetParNames("A_{tot}", "#mu1", "#sigma1", "A2", "#mu2", "#sigma2");
1832 f->SetParameters(integ,(min+mean)/2.,width/4.,
1833 integ/width/2.,(max+mean)/2.,width/4.);
1834 // The left-sided Gauss
1835 f->SetParLimits(0,integ-1.5 , integ+1.5);
1836 f->SetParLimits(1,min+(width/10.), mean);
1837 f->SetParLimits(2,0 , width/2.);
1838 // The right-sided Gauss
1839 f->SetParLimits(3,0 , integ);
1840 f->SetParLimits(4,mean, max-(width/10.));
1841 f->SetParLimits(5,0 , width/2.);
1842 obj2->Fit(f,"QLRM");
1843 break;
1844
1845 case 3:
1846 f = new TF1("tgaus",tgausformula.Data(),min,max);
1847 f->SetLineColor(kYellow);
1848 f->SetBit(kCanDelete);
1849 f->SetParNames("A_{tot}","#mu_{1}","#sigma_{1}",
1850 "A_{2}","#mu_{2}","#sigma_{2}",
1851 "A_{3}","#mu_{3}","#sigma_{3}");
1852 f->SetParameters(integ,(min+mean)/2,width/4.,
1853 integ/width/3.,(max+mean)/2.,width/4.,
1854 integ/width/3.,mean,width/2.);
1855 // The left-sided Gauss
1856 f->SetParLimits(0,integ-1.5,integ+1.5);
1857 f->SetParLimits(1,min+(width/10.),mean);
1858 f->SetParLimits(2,width/15.,width/2.);
1859 // The right-sided Gauss
1860 f->SetParLimits(3,0.,integ);
1861 f->SetParLimits(4,mean,max-(width/10.));
1862 f->SetParLimits(5,width/15.,width/2.);
1863 // The Gauss describing the outliers
1864 f->SetParLimits(6,0.,integ);
1865 f->SetParLimits(7,min,max);
1866 f->SetParLimits(8,width/4.,width/1.5);
1867 obj2->Fit(f,"QLRM");
1868 break;
1869
1870 case 4:
1871 obj2->Fit("pol0", "Q");
1872 obj2->GetFunction("pol0")->SetLineColor(kYellow);
1873 break;
1874
1875 case 9:
1876 break;
1877
1878 default:
1879 obj2->Fit("gaus", "Q");
1880 obj2->GetFunction("gaus")->SetLineColor(kYellow);
1881 break;
1882 }
1883}
1884
1885
1886// --------------------------------------------------------------------------
1887//
1888// Draw a projection of MHCamera vs. the radius from the central pixel.
1889//
1890// The inner and outer pixels are drawn separately, both fitted by a polynomial
1891// of grade 1.
1892//
1893void MHCamera::DrawRadialProfile() const
1894{
1895
1896 TProfile *obj2 = (TProfile*)RadialProfile(GetName());
1897 obj2->SetDirectory(0);
1898 obj2->Draw();
1899 obj2->SetBit(kCanDelete);
1900
1901 if (GetGeomCam().InheritsFrom("MGeomCamMagic"))
1902 {
1903
1904 TArrayI s0(6);
1905 s0[0] = 1;
1906 s0[1] = 2;
1907 s0[2] = 3;
1908 s0[3] = 4;
1909 s0[4] = 5;
1910 s0[5] = 6;
1911
1912 TArrayI inner(1);
1913 inner[0] = 0;
1914
1915 TArrayI outer(1);
1916 outer[0] = 1;
1917
1918 // Just to get the right (maximum) binning
1919 TProfile *half[2];
1920 half[0] = RadialProfileS(s0, inner,Form("%s%s",GetName(),"Inner"));
1921 half[1] = RadialProfileS(s0, outer,Form("%s%s",GetName(),"Outer"));
1922
1923 for (Int_t i=0; i<2; i++)
1924 {
1925 Double_t min = GetGeomCam().GetMinRadius(i);
1926 Double_t max = GetGeomCam().GetMaxRadius(i);
1927
1928 half[i]->SetLineColor(kRed+i);
1929 half[i]->SetDirectory(0);
1930 half[i]->SetBit(kCanDelete);
1931 half[i]->Draw("same");
1932 half[i]->Fit("pol1","Q","",min,max);
1933 half[i]->GetFunction("pol1")->SetLineColor(kRed+i);
1934 half[i]->GetFunction("pol1")->SetLineWidth(1);
1935 }
1936 }
1937}
1938
1939// --------------------------------------------------------------------------
1940//
1941// Draw a projection of MHCamera vs. the azimuth angle inside the camera.
1942//
1943// The inner and outer pixels are drawn separately.
1944// The general azimuth profile is fitted by a straight line
1945//
1946void MHCamera::DrawAzimuthProfile() const
1947{
1948
1949 TProfile *obj2 = (TProfile*)AzimuthProfile(GetName());
1950 obj2->SetDirectory(0);
1951 obj2->Draw();
1952 obj2->SetBit(kCanDelete);
1953 obj2->Fit("pol0","Q","");
1954 obj2->GetFunction("pol0")->SetLineWidth(1);
1955
1956 if (GetGeomCam().InheritsFrom("MGeomCamMagic"))
1957 {
1958
1959 TArrayI inner(1);
1960 inner[0] = 0;
1961
1962 TArrayI outer(1);
1963 outer[0] = 1;
1964
1965 // Just to get the right (maximum) binning
1966 TProfile *half[2];
1967 half[0] = AzimuthProfileA(inner,Form("%s%s",GetName(),"Inner"));
1968 half[1] = AzimuthProfileA(outer,Form("%s%s",GetName(),"Outer"));
1969
1970 for (Int_t i=0; i<2; i++)
1971 {
1972 half[i]->SetLineColor(kRed+i);
1973 half[i]->SetDirectory(0);
1974 half[i]->SetBit(kCanDelete);
1975 half[i]->SetMarkerSize(0.5);
1976 half[i]->Draw("same");
1977 }
1978 }
1979}
1980
1981// --------------------------------------------------------------------------
1982//
1983// Draw the MHCamera into the MStatusDisplay:
1984//
1985// 1) Draw it as histogram (MHCamera::DrawCopy("hist")
1986// 2) Draw it as a camera, with MHCamera::SetPrettyPalette() set.
1987// 3) If "rad" is not zero, draw its values vs. the radius from the camera center.
1988// (DrawRadialProfile())
1989// 4) Depending on the variable "fit", draw the values projection on the y-axis
1990// (DrawProjection()):
1991// 0: don't draw
1992// 1: Draw fit to Single Gauss (for distributions flat-fielded over the whole camera)
1993// 2: Draw and fit to Double Gauss (for distributions different for inner and outer pixels)
1994// 3: Draw and fit to Triple Gauss (for distributions with inner, outer pixels and outliers)
1995// 4: Draw and fit to Polynomial grade 0: (for the probability distributions)
1996// >4: Draw and don;t fit.
1997//
1998void MHCamera::CamDraw(TCanvas &c, const Int_t x, const Int_t y,
1999 const Int_t fit, const Int_t rad, const Int_t azi,
2000 TObject *notify)
2001{
2002
2003 c.cd(x);
2004 gPad->SetBorderMode(0);
2005 gPad->SetTicks();
2006 MHCamera *obj1=(MHCamera*)DrawCopy("hist");
2007 obj1->SetDirectory(NULL);
2008
2009 if (notify)
2010 obj1->AddNotify(notify);
2011
2012 c.cd(x+y);
2013 gPad->SetBorderMode(0);
2014 obj1->SetPrettyPalette();
2015 obj1->Draw();
2016
2017 Int_t cnt = 2;
2018
2019 if (rad)
2020 {
2021 c.cd(x+2*y);
2022 gPad->SetBorderMode(0);
2023 gPad->SetTicks();
2024 DrawRadialProfile();
2025 cnt++;
2026 }
2027
2028 if (azi)
2029 {
2030 c.cd(x+cnt*y);
2031 gPad->SetBorderMode(0);
2032 gPad->SetTicks();
2033 DrawAzimuthProfile();
2034 cnt++;
2035 }
2036
2037 if (!fit)
2038 return;
2039
2040 c.cd(x + cnt*y);
2041 gPad->SetBorderMode(0);
2042 gPad->SetTicks();
2043 DrawProjection(fit);
2044}
2045
2046
Note: See TracBrowser for help on using the repository browser.