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

Last change on this file since 4280 was 4280, checked in by gaug, 20 years ago
*** empty log message ***
File size: 57.7 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz, 05/2002 <mailto:tbretz@astro.uni-wuerzburg.de>
19! Author(s): Harald Kornmayer, 1/2001
20! Author(s): Markus Gaug, 03/2004 <mailto:markus@ifae.es>
21!
22! Copyright: MAGIC Software Development, 2000-2004
23!
24!
25\* ======================================================================== */
26
27/////////////////////////////////////////////////////////////////////////////
28//
29// MHCamera
30//
31// Camera Display, based on a TH1D. Pleas be carefull using the
32// underlaying TH1D.
33//
34// To change the scale to a logarithmic scale SetLogy() of the Pad.
35//
36// You can correct for the abberation. Assuming that the distance
37// between the mean position of the light distribution and the position
38// of a perfect reflection on a perfect mirror in the distance r on
39// the camera plane is dr it is d = a*dr while a is the abberation
40// constant (for the MAGIC mirror it is roughly 0.0713). You can
41// set this constant by calling SetAbberation(a) which results in a
42// 'corrected' display (all outer pixels are shifted towards the center
43// of the camera to correct for this abberation)
44//
45// Be carefull: Entries in this context means Entries/bin or Events
46//
47// FIXME? Maybe MHCamera can take the fLog object from MGeomCam?
48//
49////////////////////////////////////////////////////////////////////////////
50#include "MHCamera.h"
51
52#include <fstream>
53#include <iostream>
54
55#include <TBox.h>
56#include <TArrow.h>
57#include <TLatex.h>
58#include <TStyle.h>
59#include <TCanvas.h>
60#include <TArrayF.h>
61#include <TRandom.h>
62#include <TPaveText.h>
63#include <TPaveStats.h>
64#include <TClonesArray.h>
65#include <THistPainter.h>
66#include <THLimitsFinder.h>
67#include <TProfile.h>
68#include <TH1.h>
69#include <TF1.h>
70#include <TCanvas.h>
71
72#include "MLog.h"
73#include "MLogManip.h"
74
75#include "MH.h"
76#include "MHexagon.h"
77
78#include "MGeomPix.h"
79#include "MGeomCam.h"
80
81#include "MCamEvent.h"
82
83#define kItemsLegend 48 // see SetPalette(1,0)
84
85ClassImp(MHCamera);
86
87using namespace std;
88
89void MHCamera::Init()
90{
91 UseCurrentStyle();
92
93 SetDirectory(NULL);
94
95 SetLineColor(kGreen);
96 SetMarkerStyle(kFullDotMedium);
97 SetXTitle("Pixel Index");
98
99 fNotify = new TList;
100 fNotify->SetBit(kMustCleanup);
101 gROOT->GetListOfCleanups()->Add(fNotify);
102
103 TVirtualPad *save = gPad;
104 gPad = 0;
105#if ROOT_VERSION_CODE < ROOT_VERSION(3,01,06)
106 SetPalette(1, 0);
107#else
108 SetInvDeepBlueSeaPalette();
109#endif
110 gPad = save;
111}
112
113// ------------------------------------------------------------------------
114//
115// Default Constructor. To be used by the root system ONLY.
116//
117MHCamera::MHCamera() : TH1D(), fGeomCam(NULL), fColors(kItemsLegend), fAbberation(0)
118{
119 Init();
120}
121
122// ------------------------------------------------------------------------
123//
124// Constructor. Makes a clone of MGeomCam. Removed the TH1D from the
125// current directory. Calls Sumw2(). Set the histogram line color
126// (for error bars) to Green and the marker style to kFullDotMedium.
127//
128MHCamera::MHCamera(const MGeomCam &geom, const char *name, const char *title)
129: fGeomCam(NULL), fColors(kItemsLegend), fAbberation(0)
130{
131 //fGeomCam = (MGeomCam*)geom.Clone();
132 SetGeometry(geom, name, title);
133 Init();
134
135 //
136 // root 3.02
137 // * base/inc/TObject.h:
138 // register BIT(8) as kNoContextMenu. If an object has this bit set it will
139 // not get an automatic context menu when clicked with the right mouse button.
140}
141
142void MHCamera::SetGeometry(const MGeomCam &geom, const char *name, const char *title)
143{
144 SetNameTitle(name, title);
145
146 TAxis &x = *GetXaxis();
147
148 SetBins(geom.GetNumPixels(), 0, 1);
149 x.Set(geom.GetNumPixels(), -0.5, geom.GetNumPixels()-0.5);
150
151 //SetBins(geom.GetNumPixels(), -0.5, geom.GetNumPixels()-0.5);
152 //Rebuild();
153
154 Sumw2(); // necessary?
155
156 if (fGeomCam)
157 delete fGeomCam;
158 fGeomCam = (MGeomCam*)geom.Clone();
159
160 fUsed.Set(geom.GetNumPixels());
161 for (Int_t i=0; i<fNcells-2; i++)
162 ResetUsed(i);
163}
164
165// ------------------------------------------------------------------------
166//
167// Destructor. Deletes the cloned fGeomCam and the notification list.
168//
169MHCamera::~MHCamera()
170{
171 if (fGeomCam)
172 delete fGeomCam;
173 if (fNotify)
174 delete fNotify;
175}
176
177// ------------------------------------------------------------------------
178//
179// Return kTRUE for sector<0. Otherwise return kTRUE only if the specified
180// sector idx matches the sector of the pixel with index idx.
181//
182Bool_t MHCamera::MatchSector(Int_t idx, const TArrayI &sector, const TArrayI &aidx) const
183{
184 const MGeomPix &pix = (*fGeomCam)[idx];
185 return FindVal(sector, pix.GetSector()) && FindVal(aidx, pix.GetAidx());
186}
187
188// ------------------------------------------------------------------------
189//
190// Taken from TH1D::Fill(). Uses the argument directly as bin index.
191// Doesn't increment the number of entries.
192//
193// -*-*-*-*-*-*-*-*Increment bin with abscissa X by 1*-*-*-*-*-*-*-*-*-*-*
194// ==================================
195//
196// if x is less than the low-edge of the first bin, the Underflow bin is incremented
197// if x is greater than the upper edge of last bin, the Overflow bin is incremented
198//
199// If the storage of the sum of squares of weights has been triggered,
200// via the function Sumw2, then the sum of the squares of weights is incremented
201// by 1 in the bin corresponding to x.
202//
203// -*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
204Int_t MHCamera::Fill(Axis_t x)
205{
206#if ROOT_VERSION_CODE > ROOT_VERSION(3,05,00)
207 if (fBuffer) return BufferFill(x,1);
208#endif
209 const Int_t bin = (Int_t)x+1;
210 AddBinContent(bin);
211 if (fSumw2.fN)
212 fSumw2.fArray[bin]++;
213
214 if (bin<=0 || bin>fNcells-2)
215 return -1;
216
217 fTsumw++;
218 fTsumw2++;
219 fTsumwx += x;
220 fTsumwx2 += x*x;
221 return bin;
222}
223
224// ------------------------------------------------------------------------
225//
226// Taken from TH1D::Fill(). Uses the argument directly as bin index.
227// Doesn't increment the number of entries.
228//
229// -*-*-*-*-*-*Increment bin with abscissa X with a weight w*-*-*-*-*-*-*-*
230// =============================================
231//
232// if x is less than the low-edge of the first bin, the Underflow bin is incremented
233// if x is greater than the upper edge of last bin, the Overflow bin is incremented
234//
235// If the storage of the sum of squares of weights has been triggered,
236// via the function Sumw2, then the sum of the squares of weights is incremented
237// by w^2 in the bin corresponding to x.
238//
239// -*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
240Int_t MHCamera::Fill(Axis_t x, Stat_t w)
241{
242#if ROOT_VERSION_CODE > ROOT_VERSION(3,05,00)
243 if (fBuffer) return BufferFill(x,w);
244#endif
245 const Int_t bin = (Int_t)x+1;
246 AddBinContent(bin, w);
247 if (fSumw2.fN)
248 fSumw2.fArray[bin] += w*w;
249
250 if (bin<=0 || bin>fNcells-2)
251 return -1;
252
253 const Stat_t z = (w > 0 ? w : -w);
254 fTsumw += z;
255 fTsumw2 += z*z;
256 fTsumwx += z*x;
257 fTsumwx2 += z*x*x;
258 return bin;
259}
260
261// ------------------------------------------------------------------------
262//
263// Use x and y in millimeters
264//
265Int_t MHCamera::Fill(Axis_t x, Axis_t y, Stat_t w)
266{
267 if (fNcells<=1 || IsFreezed())
268 return -1;
269
270 for (Int_t idx=0; idx<fNcells-2; idx++)
271 {
272 MHexagon hex((*fGeomCam)[idx]);
273 if (hex.DistanceToPrimitive(x, y)>0)
274 continue;
275
276 SetUsed(idx);
277 return Fill(idx, w);
278 }
279 return -1;
280}
281
282// ------------------------------------------------------------------------
283//
284// Call this if you want to change the display status (displayed or not)
285// for all pixels. val==0 means that the pixel is not displayed.
286//
287void MHCamera::SetUsed(const TArrayC &arr)
288{
289 if (fNcells-2 != arr.GetSize())
290 {
291 gLog << warn << "WARNING - MHCamera::SetUsed: array size mismatch... ignored." << endl;
292 return;
293 }
294
295 for (Int_t idx=0; idx<fNcells-2; idx++)
296 arr[idx] ? SetUsed(idx) : ResetUsed(idx);
297}
298
299// ------------------------------------------------------------------------
300//
301// Return the mean value of all entries which are used if all=kFALSE and
302// of all entries if all=kTRUE if sector<0. If sector>=0 only
303// entries with match the given sector are taken into account.
304//
305Stat_t MHCamera::GetMeanSectors(const TArrayI &sector, const TArrayI &aidx, Bool_t all) const
306{
307 if (fNcells<=1)
308 return 0;
309
310 Int_t n=0;
311
312 Stat_t mean = 0;
313 for (int i=0; i<fNcells-2; i++)
314 {
315 if ((all || IsUsed(i)) && MatchSector(i, sector, aidx))
316 {
317 mean += fArray[i+1];
318 n++;
319 }
320 }
321
322 return n==0 ? 0 : Profile(mean/n);
323}
324
325// ------------------------------------------------------------------------
326//
327// Return the variance of all entries which are used if all=kFALSE and
328// of all entries if all=kTRUE if sector<0. If sector>=0 only
329// entries with match the given sector are taken into account.
330//
331Stat_t MHCamera::GetRmsSectors(const TArrayI &sector, const TArrayI &aidx, Bool_t all) const
332{
333 if (fNcells<=1)
334 return -1;
335
336 Int_t n=0;
337
338 Stat_t sum = 0;
339 Stat_t sq = 0;
340 for (int i=0; i<fNcells-2; i++)
341 {
342 if ((all || IsUsed(i)) && MatchSector(i, sector, aidx))
343 {
344 sum += fArray[i+1];
345 sq += fArray[i+1]*fArray[i+1];
346 n++;
347 }
348 }
349
350 if (n==0)
351 return 0;
352
353 sum /= n;
354 sq /= n;
355
356 return Profile(TMath::Sqrt(sq-sum*sum));
357}
358
359// ------------------------------------------------------------------------
360//
361// Return the minimum contents of all pixels (if all is set, otherwise
362// only of all 'used' pixels), fMinimum if fMinimum set. If sector>=0
363// only pixels with matching sector number are taken into account.
364//
365Double_t MHCamera::GetMinimumSectors(const TArrayI &sector, const TArrayI &aidx, Bool_t all) const
366{
367 if (fMinimum != -1111)
368 return fMinimum;
369
370 if (fNcells<=1)
371 return 0;
372
373 Double_t minimum=FLT_MAX;
374
375 for (Int_t idx=0; idx<fNcells-2; idx++)
376 if (MatchSector(idx, sector, aidx) && (all || IsUsed(idx)) && fArray[idx+1]<minimum)
377 minimum = fArray[idx+1];
378
379 return Profile(minimum);
380}
381
382// ------------------------------------------------------------------------
383//
384// Return the maximum contents of all pixels (if all is set, otherwise
385// only of all 'used' pixels), fMaximum if fMaximum set. If sector>=0
386// only pixels with matching sector number are taken into account.
387//
388Double_t MHCamera::GetMaximumSectors(const TArrayI &sector, const TArrayI &aidx, Bool_t all) const
389{
390 if (fMaximum!=-1111)
391 return fMaximum;
392
393 if (fNcells<=1)
394 return 1;
395
396 Double_t maximum=-FLT_MAX;
397 for (Int_t idx=0; idx<fNcells-2; idx++)
398 if (MatchSector(idx, sector, aidx) && (all || IsUsed(idx)) && fArray[idx+1]>maximum)
399 maximum = fArray[idx+1];
400
401 return Profile(maximum);
402}
403
404// ------------------------------------------------------------------------
405//
406// Call this function to draw the camera layout into your canvas.
407// Setup a drawing canvas. Add this object and all child objects
408// (hexagons, etc) to the current pad. If no pad exists a new one is
409// created. (To access the 'real' pad containing the camera you have
410// to do a cd(1) in the current layer.
411//
412// To draw a camera into its own pad do something like:
413//
414// MGeomCamMagic m;
415// MHCamera *d=new MHCamera(m);
416//
417// TCanvas *c = new TCanvas;
418// c->Divide(2,1);
419// c->cd(1);
420//
421// d->FillRandom();
422// d->Draw();
423// d->SetBit(kCanDelete);
424//
425// There are several drawing options:
426// 'hist' Draw as a standard TH1 histogram (value vs. pixel index)
427// 'box' Draw hexagons which size is in respect to its contents
428// 'nocol' Leave the 'boxed' hexagons empty
429// 'pixelindex' Display the pixel index in each pixel
430// 'sectorindex' Display the sector index in each pixel
431// 'content' Display the relative content aligned to GetMaximum() and
432// GeMinimum() ((val-min)/(max-min))
433// 'proj' Display the y-projection of the histogram
434// 'same' Draw trandparent pixels on top of an existing pad. This
435// makes it possible to draw the camera image on top of an
436// existing TH2, but also allows for distorted camera images
437//
438void MHCamera::Draw(Option_t *option)
439{
440 const Bool_t hassame = TString(option).Contains("same", TString::kIgnoreCase) && gPad;
441
442 // root 3.02:
443 // gPad->SetFixedAspectRatio()
444 const Color_t col = gPad ? gPad->GetFillColor() : 16;
445 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas("CamDisplay", "Mars Camera Display", 656, 600);
446
447 if (!hassame)
448 {
449 pad->SetBorderMode(0);
450 pad->SetFillColor(col);
451
452 //
453 // Create an own pad for the MHCamera-Object which can be
454 // resized in paint to keep the correct aspect ratio
455 //
456 pad->Divide(1, 1, 0, 0, col);
457 pad->cd(1);
458 gPad->SetBorderMode(0);
459 }
460
461 AppendPad(option);
462 //fGeomCam->AppendPad();
463
464 //
465 // Do not change gPad. The user should not see, that Draw
466 // changes gPad...
467 //
468 if (!hassame)
469 pad->cd();
470}
471
472// ------------------------------------------------------------------------
473//
474// This is TObject::DrawClone but completely ignores
475// gROOT->GetSelectedPad(). tbretz had trouble with this in the past.
476// If this makes trouble please write a bug report.
477//
478TObject *MHCamera::DrawClone(Option_t *option) const
479{
480 // Draw a clone of this object in the current pad
481
482 //TVirtualPad *pad = gROOT->GetSelectedPad();
483 TVirtualPad *padsav = gPad;
484 //if (pad) pad->cd();
485
486 TObject *newobj = Clone();
487
488 if (!newobj)
489 return 0;
490
491 /*
492 if (pad) {
493 if (strlen(option)) pad->GetListOfPrimitives()->Add(newobj,option);
494 else pad->GetListOfPrimitives()->Add(newobj,GetDrawOption());
495 pad->Modified(kTRUE);
496 pad->Update();
497 if (padsav) padsav->cd();
498 return newobj;
499 }
500 */
501
502 TString opt(option);
503 opt.ToLower();
504
505 newobj->Draw(opt.IsNull() ? GetDrawOption() : option);
506
507 if (padsav)
508 padsav->cd();
509
510 return newobj;
511}
512
513// ------------------------------------------------------------------------
514//
515// Creates a TH1D which contains the projection of the contents of the
516// MHCamera onto the y-axis. The maximum and minimum are calculated
517// such that a slighly wider range than (GetMinimum(), GetMaximum()) is
518// displayed using THLimitsFinder::OptimizeLimits.
519//
520// If no name is given the newly allocated histogram is removed from
521// the current directory calling SetDirectory(0) in any other case
522// the newly created histogram is removed from the current directory
523// and added to gROOT such the gROOT->FindObject can find the histogram.
524//
525// If the standard name "_py" is given "_py" is appended to the name
526// of the MHCamera and the corresponding histogram is searched using
527// gROOT->FindObject and updated with the present projection.
528//
529// It is the responsibility of the user to make sure, that the newly
530// created histogram is freed correctly.
531//
532// Currently the new histogram is restrictred to 50 bins.
533// Maybe a optimal number can be calulated from the number of
534// bins on the x-axis of the MHCamera?
535//
536// The code was taken mainly from TH2::ProjectX such the interface
537// is more or less the same than to TH2-projections.
538//
539// If sector>=0 only entries with matching sector index are taken
540// into account.
541//
542TH1D *MHCamera::ProjectionS(const TArrayI &sector, const TArrayI &aidx, const char *name) const
543{
544 Int_t nbins = 50;
545
546 // Create the projection histogram
547 TString pname(name);
548 if (name=="_py")
549 {
550 pname.Prepend(GetName());
551 if (sector.GetSize()>0)
552 {
553 pname += ";";
554 for (int i=0; i<sector.GetSize(); i++)
555 pname += sector[i];
556 }
557 if (aidx.GetSize()>0)
558 {
559 pname += ";";
560 for (int i=0; i<aidx.GetSize(); i++)
561 pname += aidx[i];
562 }
563 }
564
565 TH1D *h1=0;
566
567 //check if histogram with identical name exist
568 TObject *h1obj = gROOT->FindObject(pname);
569 if (h1obj && h1obj->InheritsFrom("TH1D")) {
570 h1 = (TH1D*)h1obj;
571 h1->Reset();
572 }
573
574 if (!h1)
575 {
576 Double_t min = GetMinimumSectors(sector, aidx);
577 Double_t max = GetMaximumSectors(sector, aidx);
578
579 Int_t newbins=0;
580
581 THLimitsFinder::OptimizeLimits(nbins, newbins, min, max, kFALSE);
582
583 h1 = new TH1D(pname, GetTitle(), nbins, min, max);
584 h1->SetDirectory(pname.IsNull() ? NULL : gROOT);
585 h1->SetXTitle(GetYaxis()->GetTitle());
586 h1->SetYTitle("Counts");
587 //h1->Sumw2();
588 }
589
590 // Fill the projected histogram
591 for (Int_t idx=0; idx<fNcells-2; idx++)
592 if (IsUsed(idx) && MatchSector(idx, sector, aidx))
593 h1->Fill(GetBinContent(idx+1));
594
595 return h1;
596}
597
598// ------------------------------------------------------------------------
599//
600// Creates a TH1D which contains the projection of the contents of the
601// MHCamera onto the radius from the camera center.
602// The maximum and minimum are calculated
603// such that a slighly wider range than (GetMinimum(), GetMaximum()) is
604// displayed using THLimitsFinder::OptimizeLimits.
605//
606// If no name is given the newly allocated histogram is removed from
607// the current directory calling SetDirectory(0) in any other case
608// the newly created histogram is removed from the current directory
609// and added to gROOT such the gROOT->FindObject can find the histogram.
610//
611// If the standard name "_rad" is given "_rad" is appended to the name
612// of the MHCamera and the corresponding histogram is searched using
613// gROOT->FindObject and updated with the present projection.
614//
615// It is the responsibility of the user to make sure, that the newly
616// created histogram is freed correctly.
617//
618// Currently the new histogram is restrictred to 50 bins.
619// Maybe a optimal number can be calulated from the number of
620// bins on the x-axis of the MHCamera?
621//
622// The code was taken mainly from TH2::ProjectX such the interface
623// is more or less the same than to TH2-projections.
624//
625// If sector>=0 only entries with matching sector index are taken
626// into account.
627//
628TProfile *MHCamera::RadialProfileS(const TArrayI &sector, const TArrayI &aidx, const char *name, const Int_t nbins) const
629{
630
631 // Create the projection histogram
632 TString pname(name);
633 if (name=="_rad")
634 {
635 pname.Prepend(GetName());
636 if (sector.GetSize()>0)
637 {
638 pname += ";";
639 for (int i=0; i<sector.GetSize(); i++)
640 pname += sector[i];
641 }
642 if (aidx.GetSize()>0)
643 {
644 pname += ";";
645 for (int i=0; i<aidx.GetSize(); i++)
646 pname += aidx[i];
647 }
648 }
649
650 TProfile *h1=0;
651
652 //check if histogram with identical name exist
653 TObject *h1obj = gROOT->FindObject(pname);
654 if (h1obj && h1obj->InheritsFrom("TProfile")) {
655 h1 = (TProfile*)h1obj;
656 h1->Reset();
657 }
658
659 if (!h1)
660 {
661
662 Double_t min = 0.;
663 Double_t max = fGeomCam->GetMaxRadius();
664
665 Int_t newbins=0;
666
667 THLimitsFinder::OptimizeLimits(nbins, newbins, min, max, kFALSE);
668
669 h1 = new TProfile(pname, GetTitle(), nbins, min, max);
670 h1->SetDirectory(pname.IsNull() ? NULL : gROOT);
671 h1->SetXTitle("Radius from camera center [mm]");
672 h1->SetYTitle(GetYaxis()->GetTitle());
673 }
674
675 // Fill the projected histogram
676 for (Int_t idx=0; idx<fNcells-2; idx++)
677 if (IsUsed(idx) && MatchSector(idx, sector, aidx))
678 h1->Fill(TMath::Hypot((*fGeomCam)[idx].GetX(),(*fGeomCam)[idx].GetY()),
679 GetBinContent(idx+1));
680 return h1;
681}
682
683
684// ------------------------------------------------------------------------
685//
686// Creates a TH1D which contains the projection of the contents of the
687// MHCamera onto the azimuth angle in the camera.
688//
689// If no name is given the newly allocated histogram is removed from
690// the current directory calling SetDirectory(0) in any other case
691// the newly created histogram is removed from the current directory
692// and added to gROOT such the gROOT->FindObject can find the histogram.
693//
694// If the standard name "_azi" is given "_azi" is appended to the name
695// of the MHCamera and the corresponding histogram is searched using
696// gROOT->FindObject and updated with the present projection.
697//
698// It is the responsibility of the user to make sure, that the newly
699// created histogram is freed correctly.
700//
701// Currently the new histogram is restrictred to 60 bins.
702// Maybe a optimal number can be calulated from the number of
703// bins on the x-axis of the MHCamera?
704//
705// The code was taken mainly from TH2::ProjectX such the interface
706// is more or less the same than to TH2-projections.
707//
708TProfile *MHCamera::AzimuthProfileA(const TArrayI &aidx, const char *name, const Int_t nbins) const
709{
710
711 // Create the projection histogram
712 TString pname(name);
713 if (name=="_azi")
714 {
715 pname.Prepend(GetName());
716 if (aidx.GetSize()>0)
717 {
718 pname += ";";
719 for (int i=0; i<aidx.GetSize(); i++)
720 pname += aidx[i];
721 }
722 }
723
724 TProfile *h1=0;
725
726 //check if histogram with identical name exist
727 TObject *h1obj = gROOT->FindObject(pname);
728 if (h1obj && h1obj->InheritsFrom("TProfile")) {
729 h1 = (TProfile*)h1obj;
730 h1->Reset();
731 }
732
733 if (!h1)
734 {
735
736 Double_t min = -0.5;
737 Double_t max = 359.5;
738
739 Int_t newbins=0;
740
741 THLimitsFinder::OptimizeLimits(nbins, newbins, min, max, kFALSE);
742
743 h1 = new TProfile(pname, GetTitle(), nbins, min, max);
744 h1->SetDirectory(pname.IsNull() ? NULL : gROOT);
745 h1->SetXTitle("Azimuth in camera [deg]");
746 h1->SetYTitle(GetYaxis()->GetTitle());
747 }
748
749 // Fill the projected histogram
750 for (Int_t idx=0; idx<fNcells-2; idx++)
751 {
752 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
1214//
1215void MHCamera::CntCamContent(const MCamEvent &event, Double_t threshold, Int_t type)
1216{
1217 if (fNcells<=1 || IsFreezed())
1218 return;
1219
1220 // FIXME: Security check missing!
1221 for (Int_t idx=0; idx<fNcells-2; idx++)
1222 {
1223 Double_t val=threshold;
1224 if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
1225 SetUsed(idx);
1226
1227 if (val>threshold)
1228 Fill(idx);
1229 }
1230 fEntries++;
1231}
1232
1233// ------------------------------------------------------------------------
1234//
1235// Call this function to add a TArrayD on top of the present contents.
1236// 1 is added to each pixel if the contents of MCamEvent>threshold
1237//
1238void MHCamera::CntCamContent(const TArrayD &event, Double_t threshold, Bool_t ispos)
1239{
1240 if (event.GetSize()!=fNcells-2 || IsFreezed())
1241 return;
1242
1243 for (Int_t idx=0; idx<fNcells-2; idx++)
1244 {
1245 if (const_cast<TArrayD&>(event)[idx]>threshold)
1246 Fill(idx);
1247
1248 if (!ispos || fArray[idx+1]>0)
1249 SetUsed(idx);
1250 }
1251 fEntries++;
1252}
1253
1254// ------------------------------------------------------------------------
1255//
1256// Fill the pixels with random contents.
1257//
1258void MHCamera::FillRandom()
1259{
1260 if (fNcells<=1 || IsFreezed())
1261 return;
1262
1263 Reset();
1264
1265 // FIXME: Security check missing!
1266 for (Int_t idx=0; idx<fNcells-2; idx++)
1267 {
1268 Fill(idx, gRandom->Uniform()*fGeomCam->GetPixRatio(idx));
1269 SetUsed(idx);
1270 }
1271 fEntries=1;
1272}
1273
1274
1275// ------------------------------------------------------------------------
1276//
1277// The array must be in increasing order, eg: 2.5, 3.7, 4.9
1278// The values in each bin are replaced by the interval in which the value
1279// fits. In the example we have four intervals
1280// (<2.5, 2.5-3.7, 3.7-4.9, >4.9). Maximum and minimum are set
1281// accordingly.
1282//
1283void MHCamera::SetLevels(const TArrayF &arr)
1284{
1285 if (fNcells<=1)
1286 return;
1287
1288 for (Int_t i=0; i<fNcells-2; i++)
1289 {
1290 if (!IsUsed(i))
1291 continue;
1292
1293 Int_t j = arr.GetSize();
1294 while (j && fArray[i+1]<arr[j-1])
1295 j--;
1296
1297 fArray[i+1] = j;
1298 }
1299 SetMaximum(arr.GetSize());
1300 SetMinimum(0);
1301}
1302
1303// ------------------------------------------------------------------------
1304//
1305// Reset the all pixel colors to a default value
1306//
1307void MHCamera::Reset(Option_t *opt)
1308{
1309 if (fNcells<=1 || IsFreezed())
1310 return;
1311
1312 TH1::Reset(opt);
1313
1314 for (Int_t i=0; i<fNcells-2; i++)
1315 {
1316 fArray[i+1]=0;
1317 ResetUsed(i);
1318 }
1319 fArray[0] = 0;
1320 fArray[fNcells-1] = 0;
1321}
1322
1323// ------------------------------------------------------------------------
1324//
1325// Here we calculate the color index for the current value.
1326// The color index is defined with the class TStyle and the
1327// Color palette inside. We use the command gStyle->SetPalette(1,0)
1328// for the display. So we have to convert the value "wert" into
1329// a color index that fits the color palette.
1330// The range of the color palette is defined by the values fMinPhe
1331// and fMaxRange. Between this values we have 50 color index, starting
1332// with 0 up to 49.
1333//
1334Int_t MHCamera::GetColor(Float_t val, Float_t min, Float_t max, Bool_t islog)
1335{
1336 if (TMath::IsNaN(val)) // FIXME: gLog!
1337 return 10;
1338
1339 //
1340 // first treat the over- and under-flows
1341 //
1342 const Int_t maxcolidx = kItemsLegend-1;
1343
1344 if (val >= max)
1345 return fColors[maxcolidx];
1346
1347 if (val <= min)
1348 return fColors[0];
1349
1350 //
1351 // calculate the color index
1352 //
1353 Float_t ratio;
1354 if (islog && min>0)
1355 ratio = log10(val/min) / log10(max/min);
1356 else
1357 ratio = (val-min) / (max-min);
1358
1359 const Int_t colidx = (Int_t)(ratio*maxcolidx + .5);
1360 return fColors[colidx];
1361}
1362
1363TPaveStats *MHCamera::GetStatisticBox()
1364{
1365 TObject *obj = 0;
1366
1367 TIter Next(fFunctions);
1368 while ((obj = Next()))
1369 if (obj->InheritsFrom(TPaveStats::Class()))
1370 return static_cast<TPaveStats*>(obj);
1371
1372 return NULL;
1373}
1374
1375// ------------------------------------------------------------------------
1376//
1377// Change the text on the legend according to the range of the Display
1378//
1379void MHCamera::UpdateLegend(Float_t min, Float_t max, Bool_t islog)
1380{
1381 const Float_t range = fGeomCam->GetMaxRadius()*1.05;
1382
1383 TArrow arr;
1384 arr.PaintArrow(-range*.9, -range*.9, -range*.6, -range*.9, 0.025);
1385 arr.PaintArrow(-range*.9, -range*.9, -range*.9, -range*.6, 0.025);
1386
1387 TString text;
1388 text += (int)(range*.3);
1389 text += "mm";
1390
1391 TText newtxt2;
1392 newtxt2.SetTextSize(0.04);
1393 newtxt2.PaintText(-range*.85, -range*.85, text);
1394
1395 text = "";
1396 text += (float)((int)(range*.3*fGeomCam->GetConvMm2Deg()*10))/10;
1397 text += "\\circ";
1398 text = text.Strip(TString::kLeading);
1399
1400 TLatex latex;
1401 latex.PaintLatex(-range*.85, -range*.75, 0, 0.04, text);
1402
1403 if (TestBit(kNoLegend))
1404 return;
1405
1406 TPaveStats *stats = GetStatisticBox();
1407
1408 const Float_t hndc = 0.92 - (stats ? stats->GetY1NDC() : 1);
1409 const Float_t H = (0.75-hndc)*range;
1410 const Float_t offset = hndc*range;
1411
1412 const Float_t h = 2./kItemsLegend;
1413 const Float_t w = range/sqrt((float)(fNcells-2));
1414
1415 TBox newbox;
1416 TText newtxt;
1417 newtxt.SetTextSize(0.03);
1418 newtxt.SetTextAlign(12);
1419#if ROOT_VERSION_CODE > ROOT_VERSION(3,01,06)
1420 newtxt.SetBit(/*kNoContextMenu|*/kCannotPick);
1421 newbox.SetBit(/*kNoContextMenu|*/kCannotPick);
1422#endif
1423
1424 const Float_t step = (islog && min>0 ? log10(max/min) : max-min) / kItemsLegend;
1425 const Int_t firsts = step*3 < 1e-8 ? 8 : (Int_t)floor(log10(step*3));
1426 const TString opt = Form("%%.%if", firsts>0 ? 0 : TMath::Abs(firsts));
1427
1428 for (Int_t i=0; i<kItemsLegend+1; i+=3)
1429 {
1430 Float_t val;
1431 if (islog && min>0)
1432 val = pow(10, step*i) * min;
1433 else
1434 val = min + step*i;
1435
1436 //const bool dispexp = max-min>1.5 && fabs(val)>0.1 && fabs(val)<1e6;
1437 newtxt.PaintText(range+1.5*w, H*(i*h-1)-offset, Form(opt, val));
1438 }
1439
1440 for (Int_t i=0; i<kItemsLegend; i++)
1441 {
1442 newbox.SetFillColor(fColors[i]);
1443 newbox.PaintBox(range, H*(i*h-1)-offset, range+w, H*((i+1)*h-1)-offset);
1444 }
1445}
1446
1447// ------------------------------------------------------------------------
1448//
1449// Save primitive as a C++ statement(s) on output stream out
1450//
1451void MHCamera::SavePrimitive(ofstream &out, Option_t *opt)
1452{
1453 gLog << err << "MHCamera::SavePrimitive: Must be rewritten!" << endl;
1454 /*
1455 if (!gROOT->ClassSaved(TCanvas::Class()))
1456 fDrawingPad->SavePrimitive(out, opt);
1457
1458 out << " " << fDrawingPad->GetName() << "->SetWindowSize(";
1459 out << fDrawingPad->GetWw() << "," << fDrawingPad->GetWh() << ");" << endl;
1460 */
1461}
1462
1463// ------------------------------------------------------------------------
1464//
1465// compute the distance of a point (px,py) to the Camera
1466// this functions needed for graphical primitives, that
1467// means without this function you are not able to interact
1468// with the graphical primitive with the mouse!!!
1469//
1470// All calcutations are done in pixel coordinates
1471//
1472Int_t MHCamera::DistancetoPrimitive(Int_t px, Int_t py)
1473{
1474 if (fNcells<=1)
1475 return 999999;
1476
1477 TPaveStats *box = (TPaveStats*)gPad->GetPrimitive("stats");
1478 if (box)
1479 {
1480 const Double_t w = box->GetY2NDC()-box->GetY1NDC();
1481 box->SetX1NDC(gStyle->GetStatX()-gStyle->GetStatW());
1482 box->SetY1NDC(gStyle->GetStatY()-w);
1483 box->SetX2NDC(gStyle->GetStatX());
1484 box->SetY2NDC(gStyle->GetStatY());
1485 }
1486
1487 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1488 return TH1D::DistancetoPrimitive(px, py);
1489
1490 const Bool_t issame = TString(GetDrawOption()).Contains("same", TString::kIgnoreCase);
1491
1492 const Float_t maxr = (1-fGeomCam->GetConvMm2Deg())*fGeomCam->GetMaxRadius()/2;
1493 const Float_t conv = !issame ||
1494 gPad->GetX1()<-maxr || gPad->GetY1()<-maxr ||
1495 gPad->GetX2()> maxr || gPad->GetY2()>maxr ? 1 : fGeomCam->GetConvMm2Deg();
1496
1497 if (GetPixelIndex(px, py, conv)>=0)
1498 return 0;
1499
1500 if (!box)
1501 return 999999;
1502
1503 const Int_t dist = box->DistancetoPrimitive(px, py);
1504 if (dist > TPad::GetMaxPickDistance())
1505 return 999999;
1506
1507 gPad->SetSelected(box);
1508 return dist;
1509}
1510
1511// ------------------------------------------------------------------------
1512//
1513//
1514Int_t MHCamera::GetPixelIndex(Int_t px, Int_t py, Float_t conv) const
1515{
1516 if (fNcells<=1)
1517 return -1;
1518
1519 Int_t i;
1520 for (i=0; i<fNcells-2; i++)
1521 {
1522 MHexagon hex((*fGeomCam)[i]);
1523 if (hex.DistancetoPrimitive(px, py, conv)>0)
1524 continue;
1525
1526 return i;
1527 }
1528 return -1;
1529}
1530
1531// ------------------------------------------------------------------------
1532//
1533// Returns string containing info about the object at position (px,py).
1534// Returned string will be re-used (lock in MT environment).
1535//
1536char *MHCamera::GetObjectInfo(Int_t px, Int_t py) const
1537{
1538 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1539 return TH1D::GetObjectInfo(px, py);
1540
1541 static char info[128];
1542
1543 const Int_t idx=GetPixelIndex(px, py);
1544
1545 if (idx<0)
1546 return TObject::GetObjectInfo(px, py);
1547
1548 sprintf(info, "Software Pixel Index: %d (Hardware Id=%d)", idx, idx+1);
1549 return info;
1550}
1551
1552// ------------------------------------------------------------------------
1553//
1554// Add a MCamEvent which should be displayed when the user clicks on a
1555// pixel.
1556// Warning: The object MUST inherit from TObject AND MCamEvent
1557//
1558void MHCamera::AddNotify(TObject *obj)
1559{
1560 // Make sure, that the object derives from MCamEvent!
1561 MCamEvent *evt = dynamic_cast<MCamEvent*>(obj);
1562 if (!evt)
1563 {
1564 gLog << err << "ERROR: MHCamera::AddNotify - TObject doesn't inherit from MCamEvent... ignored." << endl;
1565 return;
1566 }
1567
1568 // Make sure, that it is deleted from the list too, if the obj is deleted
1569 obj->SetBit(kMustCleanup);
1570
1571 // Add object to list
1572 fNotify->Add(obj);
1573}
1574
1575// ------------------------------------------------------------------------
1576//
1577// Execute a mouse event on the camera
1578//
1579void MHCamera::ExecuteEvent(Int_t event, Int_t px, Int_t py)
1580{
1581 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1582 {
1583 TH1D::ExecuteEvent(event, px, py);
1584 return;
1585 }
1586 //if (event==kMouseMotion && fStatusBar)
1587 // fStatusBar->SetText(GetObjectInfo(px, py), 0);
1588 if (event!=kButton1Down)
1589 return;
1590
1591 const Int_t idx = GetPixelIndex(px, py);
1592 if (idx<0)
1593 return;
1594
1595 gLog << all << GetTitle() << " <" << GetName() << ">" << endl;
1596 gLog << "Software Pixel Idx: " << idx << endl;
1597 gLog << "Hardware Pixel Id: " << idx+1 << endl;
1598 gLog << "Contents: " << GetBinContent(idx+1);
1599 if (GetBinError(idx+1)>0)
1600 gLog << " +/- " << GetBinError(idx+1);
1601 gLog << " <" << (IsUsed(idx)?"on":"off") << ">" << endl;
1602
1603 if (fNotify && fNotify->GetSize()>0)
1604 {
1605 // FIXME: Is there a simpler and more convinient way?
1606
1607 // The name which is created here depends on the instance of
1608 // MHCamera and on the pad on which it is drawn --> The name
1609 // is unique. For ExecuteEvent gPad is always correctly set.
1610 const TString name = Form("%p;%p;PixelContent", this, gPad);
1611
1612 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
1613 if (old)
1614 old->cd();
1615 else
1616 new TCanvas(name);
1617
1618 /*
1619 TIter Next(gPad->GetListOfPrimitives());
1620 TObject *o;
1621 while (o=Next()) cout << o << ": " << o->GetName() << " " << o->IsA()->GetName() << endl;
1622 */
1623
1624 // FIXME: Make sure, that the old histograms are really deleted.
1625 // Are they already deleted?
1626
1627 // The dynamic_cast is necessary here: We cannot use ForEach
1628 TIter Next(fNotify);
1629 MCamEvent *evt;
1630 while ((evt=dynamic_cast<MCamEvent*>(Next())))
1631 evt->DrawPixelContent(idx);
1632
1633 gPad->Modified();
1634 gPad->Update();
1635 }
1636}
1637
1638UInt_t MHCamera::GetNumPixels() const
1639{
1640 return fGeomCam->GetNumPixels();
1641}
1642
1643TH1 *MHCamera::DrawCopy() const
1644{
1645 gPad=NULL;
1646 return TH1D::DrawCopy(fName+";cpy");
1647}
1648
1649
1650// --------------------------------------------------------------------------
1651//
1652// Draw a projection of MHCamera onto the y-axis values. Depending on the
1653// variable fit, the following fits are performed:
1654//
1655// 0: No fit, simply draw the projection
1656// 1: Single Gauss (for distributions flat-fielded over the whole camera)
1657// 2: Double Gauss (for distributions different for inner and outer pixels)
1658// 3: Triple Gauss (for distributions with inner, outer pixels and outliers)
1659// 4: flat (for the probability distributions)
1660// (1-4:) Moreover, sectors 6,1 and 2 of the camera and sectors 3,4 and 5 are
1661// drawn separately, for inner and outer pixels.
1662// 5: Fit Inner and Outer pixels separately by a single Gaussian
1663// (only for MAGIC cameras)
1664// 6: Fit Inner and Outer pixels separately by a single Gaussian and display
1665// additionally the two camera halfs separately (for MAGIC camera)
1666//
1667//
1668void MHCamera::DrawProjection(Int_t fit) const
1669{
1670
1671 TArrayI inner(1);
1672 inner[0] = 0;
1673
1674 TArrayI outer(1);
1675 outer[0] = 1;
1676
1677 if (fit==5 || fit==6)
1678 {
1679
1680 if (GetGeomCam().InheritsFrom("MGeomCamMagic"))
1681 {
1682 TArrayI s0(6);
1683 s0[0] = 6;
1684 s0[1] = 1;
1685 s0[2] = 2;
1686 s0[3] = 3;
1687 s0[4] = 4;
1688 s0[5] = 5;
1689
1690 TArrayI s1(3);
1691 s1[0] = 6;
1692 s1[1] = 1;
1693 s1[2] = 2;
1694
1695 TArrayI s2(3);
1696 s2[0] = 3;
1697 s2[1] = 4;
1698 s2[2] = 5;
1699
1700 gPad->Clear();
1701 TVirtualPad *pad = gPad;
1702 pad->Divide(2,1);
1703
1704 TH1D *inout[2];
1705 inout[0] = ProjectionS(s0, inner, "Inner");
1706 inout[1] = ProjectionS(s0, outer, "Outer");
1707
1708 inout[0]->SetDirectory(NULL);
1709 inout[1]->SetDirectory(NULL);
1710
1711 for (int i=0; i<2; i++)
1712 {
1713 pad->cd(i+1);
1714 inout[i]->SetLineColor(kRed+i);
1715 inout[i]->SetBit(kCanDelete);
1716 inout[i]->Draw();
1717 inout[i]->Fit("gaus","Q");
1718 /*
1719 gPad->Modified();
1720 gPad->Update();
1721 gPad->SaveAs(Form("%s%s%i%s",GetName(),"AreaIdx",i,".eps"));
1722 */
1723 if (fit == 6)
1724 {
1725 TH1D *half[2];
1726 half[0] = ProjectionS(s1, i==0 ? inner : outer , "Sector 6-1-2");
1727 half[1] = ProjectionS(s2, i==0 ? inner : outer , "Sector 3-4-5");
1728
1729 for (int j=0; j<2; j++)
1730 {
1731 half[j]->SetLineColor(kRed+i+2*j+1);
1732 half[j]->SetDirectory(NULL);
1733 half[j]->SetBit(kCanDelete);
1734 half[j]->Draw("same");
1735 }
1736 }
1737
1738 }
1739
1740 gLog << all << GetName()
1741 << Form("%s%5.3f%s%3.2f"," Inner Pixels: ",
1742 inout[0]->GetFunction("gaus")->GetParameter(1),"+-",
1743 inout[0]->GetFunction("gaus")->GetParameter(2));
1744 gLog << Form("%s%5.3f%s%3.2f"," Outer Pixels: ",
1745 inout[1]->GetFunction("gaus")->GetParameter(1),"+-",
1746 inout[1]->GetFunction("gaus")->GetParameter(2)) << endl;
1747
1748 }
1749 return;
1750 }
1751
1752 TH1D *obj2 = (TH1D*)Projection(GetName());
1753 obj2->SetDirectory(0);
1754 obj2->Draw();
1755 obj2->SetBit(kCanDelete);
1756
1757 if (fit == 0)
1758 return;
1759
1760 if (GetGeomCam().InheritsFrom("MGeomCamMagic"))
1761 {
1762 TArrayI s0(3);
1763 s0[0] = 6;
1764 s0[1] = 1;
1765 s0[2] = 2;
1766
1767 TArrayI s1(3);
1768 s1[0] = 3;
1769 s1[1] = 4;
1770 s1[2] = 5;
1771
1772
1773 TH1D *halfInOut[4];
1774
1775 // Just to get the right (maximum) binning
1776 halfInOut[0] = ProjectionS(s0, inner, "Sector 6-1-2 Inner");
1777 halfInOut[1] = ProjectionS(s1, inner, "Sector 3-4-5 Inner");
1778 halfInOut[2] = ProjectionS(s0, outer, "Sector 6-1-2 Outer");
1779 halfInOut[3] = ProjectionS(s1, outer, "Sector 3-4-5 Outer");
1780
1781 for (int i=0; i<4; i++)
1782 {
1783 halfInOut[i]->SetLineColor(kRed+i);
1784 halfInOut[i]->SetDirectory(0);
1785 halfInOut[i]->SetBit(kCanDelete);
1786 halfInOut[i]->Draw("same");
1787 }
1788 }
1789
1790 const Double_t min = obj2->GetBinCenter(obj2->GetXaxis()->GetFirst());
1791 const Double_t max = obj2->GetBinCenter(obj2->GetXaxis()->GetLast());
1792 const Double_t integ = obj2->Integral("width")/2.5;
1793 const Double_t mean = obj2->GetMean();
1794 const Double_t rms = obj2->GetRMS();
1795 const Double_t width = max-min;
1796
1797 const TString dgausformula = "([0]-[3])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
1798 "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])";
1799
1800 const TString tgausformula = "([0]-[3]-[6])/[2]*exp(-0.5*(x-[1])*(x-[1])/[2]/[2])"
1801 "+[3]/[5]*exp(-0.5*(x-[4])*(x-[4])/[5]/[5])"
1802 "+[6]/[8]*exp(-0.5*(x-[7])*(x-[7])/[8]/[8])";
1803 TF1 *f=0;
1804 switch (fit)
1805 {
1806 case 1:
1807 f = new TF1("sgaus", "gaus(0)", min, max);
1808 f->SetLineColor(kYellow);
1809 f->SetBit(kCanDelete);
1810 f->SetParNames("Area", "#mu", "#sigma");
1811 f->SetParameters(integ/rms, mean, rms);
1812 f->SetParLimits(0, 0, integ);
1813 f->SetParLimits(1, min, max);
1814 f->SetParLimits(2, 0, width/1.5);
1815
1816 obj2->Fit(f, "QLR");
1817 break;
1818
1819 case 2:
1820 f = new TF1("dgaus",dgausformula.Data(),min,max);
1821 f->SetLineColor(kYellow);
1822 f->SetBit(kCanDelete);
1823 f->SetParNames("A_{tot}", "#mu1", "#sigma1", "A2", "#mu2", "#sigma2");
1824 f->SetParameters(integ,(min+mean)/2.,width/4.,
1825 integ/width/2.,(max+mean)/2.,width/4.);
1826 // The left-sided Gauss
1827 f->SetParLimits(0,integ-1.5 , integ+1.5);
1828 f->SetParLimits(1,min+(width/10.), mean);
1829 f->SetParLimits(2,0 , width/2.);
1830 // The right-sided Gauss
1831 f->SetParLimits(3,0 , integ);
1832 f->SetParLimits(4,mean, max-(width/10.));
1833 f->SetParLimits(5,0 , width/2.);
1834 obj2->Fit(f,"QLRM");
1835 break;
1836
1837 case 3:
1838 f = new TF1("tgaus",tgausformula.Data(),min,max);
1839 f->SetLineColor(kYellow);
1840 f->SetBit(kCanDelete);
1841 f->SetParNames("A_{tot}","#mu_{1}","#sigma_{1}",
1842 "A_{2}","#mu_{2}","#sigma_{2}",
1843 "A_{3}","#mu_{3}","#sigma_{3}");
1844 f->SetParameters(integ,(min+mean)/2,width/4.,
1845 integ/width/3.,(max+mean)/2.,width/4.,
1846 integ/width/3.,mean,width/2.);
1847 // The left-sided Gauss
1848 f->SetParLimits(0,integ-1.5,integ+1.5);
1849 f->SetParLimits(1,min+(width/10.),mean);
1850 f->SetParLimits(2,width/15.,width/2.);
1851 // The right-sided Gauss
1852 f->SetParLimits(3,0.,integ);
1853 f->SetParLimits(4,mean,max-(width/10.));
1854 f->SetParLimits(5,width/15.,width/2.);
1855 // The Gauss describing the outliers
1856 f->SetParLimits(6,0.,integ);
1857 f->SetParLimits(7,min,max);
1858 f->SetParLimits(8,width/4.,width/1.5);
1859 obj2->Fit(f,"QLRM");
1860 break;
1861
1862 case 4:
1863 obj2->Fit("pol0", "Q");
1864 obj2->GetFunction("pol0")->SetLineColor(kYellow);
1865 break;
1866
1867 case 9:
1868 break;
1869
1870 default:
1871 obj2->Fit("gaus", "Q");
1872 obj2->GetFunction("gaus")->SetLineColor(kYellow);
1873 break;
1874 }
1875}
1876
1877
1878// --------------------------------------------------------------------------
1879//
1880// Draw a projection of MHCamera vs. the radius from the central pixel.
1881//
1882// The inner and outer pixels are drawn separately, both fitted by a polynomial
1883// of grade 1.
1884//
1885void MHCamera::DrawRadialProfile() const
1886{
1887
1888 TProfile *obj2 = (TProfile*)RadialProfile(GetName());
1889 obj2->SetDirectory(0);
1890 obj2->Draw();
1891 obj2->SetBit(kCanDelete);
1892
1893 if (GetGeomCam().InheritsFrom("MGeomCamMagic"))
1894 {
1895
1896 TArrayI s0(6);
1897 s0[0] = 1;
1898 s0[1] = 2;
1899 s0[2] = 3;
1900 s0[3] = 4;
1901 s0[4] = 5;
1902 s0[5] = 6;
1903
1904 TArrayI inner(1);
1905 inner[0] = 0;
1906
1907 TArrayI outer(1);
1908 outer[0] = 1;
1909
1910 // Just to get the right (maximum) binning
1911 TProfile *half[2];
1912 half[0] = RadialProfileS(s0, inner,Form("%s%s",GetName(),"Inner"));
1913 half[1] = RadialProfileS(s0, outer,Form("%s%s",GetName(),"Outer"));
1914
1915 for (Int_t i=0; i<2; i++)
1916 {
1917 Double_t min = GetGeomCam().GetMinRadius(i);
1918 Double_t max = GetGeomCam().GetMaxRadius(i);
1919
1920 half[i]->SetLineColor(kRed+i);
1921 half[i]->SetDirectory(0);
1922 half[i]->SetBit(kCanDelete);
1923 half[i]->Draw("same");
1924 half[i]->Fit("pol1","Q","",min,max);
1925 half[i]->GetFunction("pol1")->SetLineColor(kRed+i);
1926 half[i]->GetFunction("pol1")->SetLineWidth(1);
1927 }
1928 }
1929}
1930
1931// --------------------------------------------------------------------------
1932//
1933// Draw a projection of MHCamera vs. the azimuth angle inside the camera.
1934//
1935// The inner and outer pixels are drawn separately.
1936// The general azimuth profile is fitted by a straight line
1937//
1938void MHCamera::DrawAzimuthProfile() const
1939{
1940
1941 TProfile *obj2 = (TProfile*)AzimuthProfile(GetName());
1942 obj2->SetDirectory(0);
1943 obj2->SetLineColor(kRed);
1944 obj2->Draw();
1945 obj2->SetBit(kCanDelete);
1946 obj2->Fit("pol0","Q","");
1947 obj2->GetFunction("pol0")->SetLineColor(kRed);
1948 obj2->GetFunction("pol0")->SetLineWidth(1);
1949
1950 if (GetGeomCam().InheritsFrom("MGeomCamMagic"))
1951 {
1952
1953 TArrayI inner(1);
1954 inner[0] = 0;
1955
1956 TArrayI outer(1);
1957 outer[0] = 1;
1958
1959 // Just to get the right (maximum) binning
1960 TProfile *half[2];
1961 half[0] = AzimuthProfileA(inner,Form("%s%s",GetName(),"Inner"));
1962 half[1] = AzimuthProfileA(outer,Form("%s%s",GetName(),"Outer"));
1963
1964 for (Int_t i=0; i<2; i++)
1965 {
1966 half[i]->SetLineColor(kRed+i+1);
1967 half[i]->SetDirectory(0);
1968 half[i]->SetBit(kCanDelete);
1969 half[i]->SetMarkerSize(0.5);
1970 half[i]->Draw("same");
1971 }
1972 }
1973}
1974
1975// --------------------------------------------------------------------------
1976//
1977// Draw the MHCamera into the MStatusDisplay:
1978//
1979// 1) Draw it as histogram (MHCamera::DrawCopy("hist")
1980// 2) Draw it as a camera, with MHCamera::SetPrettyPalette() set.
1981// 3) If "rad" is not zero, draw its values vs. the radius from the camera center.
1982// (DrawRadialProfile())
1983// 4) Depending on the variable "fit", draw the values projection on the y-axis
1984// (DrawProjection()):
1985// 0: don't draw
1986// 1: Draw fit to Single Gauss (for distributions flat-fielded over the whole camera)
1987// 2: Draw and fit to Double Gauss (for distributions different for inner and outer pixels)
1988// 3: Draw and fit to Triple Gauss (for distributions with inner, outer pixels and outliers)
1989// 4: Draw and fit to Polynomial grade 0: (for the probability distributions)
1990// >4: Draw and don;t fit.
1991//
1992void MHCamera::CamDraw(TCanvas &c, const Int_t x, const Int_t y,
1993 const Int_t fit, const Int_t rad, const Int_t azi,
1994 TObject *notify)
1995{
1996
1997 c.cd(x);
1998 gPad->SetBorderMode(0);
1999 gPad->SetTicks();
2000 MHCamera *obj1=(MHCamera*)DrawCopy("hist");
2001 obj1->SetDirectory(NULL);
2002
2003 if (notify)
2004 obj1->AddNotify(notify);
2005
2006 c.cd(x+y);
2007 gPad->SetBorderMode(0);
2008 obj1->SetPrettyPalette();
2009 obj1->Draw();
2010
2011 Int_t cnt = 2;
2012
2013 if (rad)
2014 {
2015 c.cd(x+2*y);
2016 gPad->SetBorderMode(0);
2017 gPad->SetTicks();
2018 DrawRadialProfile();
2019 cnt++;
2020 }
2021
2022 if (azi)
2023 {
2024 c.cd(x+cnt*y);
2025 gPad->SetBorderMode(0);
2026 gPad->SetTicks();
2027 DrawAzimuthProfile();
2028 cnt++;
2029 }
2030
2031 if (!fit)
2032 return;
2033
2034 c.cd(x + cnt*y);
2035 gPad->SetBorderMode(0);
2036 gPad->SetTicks();
2037 DrawProjection(fit);
2038}
2039
2040
Note: See TracBrowser for help on using the repository browser.