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

Last change on this file since 3500 was 3500, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 38.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!
21! Copyright: MAGIC Software Development, 2000-2003
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27//
28// MHCamera
29//
30// Camera Display, based on a TH1D. Pleas be carefull using the
31// underlaying TH1D.
32//
33// To change the scale to a logarithmic scale SetLogy() of the Pad.
34//
35// Be carefull: Entries in this context means Entries/bin or Events
36//
37// FIXME? Maybe MHCamera can take the fLog object from MGeomCam?
38//
39////////////////////////////////////////////////////////////////////////////
40#include "MHCamera.h"
41
42#include <fstream>
43#include <iostream>
44
45#include <TBox.h>
46#include <TArrow.h>
47#include <TLatex.h>
48#include <TStyle.h>
49#include <TCanvas.h>
50#include <TArrayF.h>
51#include <TRandom.h>
52#include <TPaveText.h>
53#include <TPaveStats.h>
54#include <TClonesArray.h>
55#include <THistPainter.h>
56#include <THLimitsFinder.h>
57
58#include "MLog.h"
59#include "MLogManip.h"
60
61#include "MH.h"
62#include "MHexagon.h"
63
64#include "MGeomPix.h"
65#include "MGeomCam.h"
66
67#include "MCamEvent.h"
68
69
70#define kItemsLegend 48 // see SetPalette(1,0)
71
72ClassImp(MHCamera);
73
74using namespace std;
75
76void MHCamera::Init()
77{
78 UseCurrentStyle();
79
80 SetDirectory(NULL);
81
82 SetLineColor(kGreen);
83 SetMarkerStyle(kFullDotMedium);
84 SetXTitle("Pixel Index");
85
86 fNotify = new TList;
87 fNotify->SetBit(kMustCleanup);
88 gROOT->GetListOfCleanups()->Add(fNotify);
89
90#if ROOT_VERSION_CODE < ROOT_VERSION(3,01,06)
91 SetPalette(1, 0);
92#else
93 SetInvDeepBlueSeaPalette();
94#endif
95}
96
97// ------------------------------------------------------------------------
98//
99// Default Constructor. To be used by the root system ONLY.
100//
101MHCamera::MHCamera() : TH1D(), fGeomCam(NULL), fColors(kItemsLegend)
102{
103 Init();
104}
105
106// ------------------------------------------------------------------------
107//
108// Constructor. Makes a clone of MGeomCam. Removed the TH1D from the
109// current directory. Calls Sumw2(). Set the histogram line color
110// (for error bars) to Green and the marker style to kFullDotMedium.
111//
112MHCamera::MHCamera(const MGeomCam &geom, const char *name, const char *title)
113: fGeomCam(NULL), fColors(kItemsLegend)
114{
115 //fGeomCam = (MGeomCam*)geom.Clone();
116 SetGeometry(geom, name, title);
117 Init();
118
119 //
120 // root 3.02
121 // * base/inc/TObject.h:
122 // register BIT(8) as kNoContextMenu. If an object has this bit set it will
123 // not get an automatic context menu when clicked with the right mouse button.
124}
125
126void MHCamera::SetGeometry(const MGeomCam &geom, const char *name, const char *title)
127{
128 SetNameTitle(name, title);
129
130 TAxis &x = *GetXaxis();
131
132 SetBins(geom.GetNumPixels(), 0, 1);
133 x.Set(geom.GetNumPixels(), -0.5, geom.GetNumPixels()-0.5);
134
135 //SetBins(geom.GetNumPixels(), -0.5, geom.GetNumPixels()-0.5);
136 //Rebuild();
137
138 Sumw2(); // necessary?
139
140 if (fGeomCam)
141 delete fGeomCam;
142 fGeomCam = (MGeomCam*)geom.Clone();
143
144 fUsed.Set(geom.GetNumPixels());
145 for (Int_t i=0; i<fNcells-2; i++)
146 ResetUsed(i);
147}
148
149// ------------------------------------------------------------------------
150//
151// Destructor. Deletes the cloned fGeomCam and the notification list.
152//
153MHCamera::~MHCamera()
154{
155 if (fGeomCam)
156 delete fGeomCam;
157 if (fNotify)
158 delete fNotify;
159}
160
161// ------------------------------------------------------------------------
162//
163// Return kTRUE for sector<0. Otherwise return kTRUE only if the specified
164// sector idx matches the sector of the pixel with index idx.
165//
166Bool_t MHCamera::MatchSector(Int_t idx, const TArrayI &sector, const TArrayI &aidx) const
167{
168 const MGeomPix &pix = (*fGeomCam)[idx];
169 return FindVal(sector, pix.GetSector()) && FindVal(aidx, pix.GetAidx());
170}
171
172// ------------------------------------------------------------------------
173//
174// Taken from TH1D::Fill(). Uses the argument directly as bin index.
175// Doesn't increment the number of entries.
176//
177// -*-*-*-*-*-*-*-*Increment bin with abscissa X by 1*-*-*-*-*-*-*-*-*-*-*
178// ==================================
179//
180// if x is less than the low-edge of the first bin, the Underflow bin is incremented
181// if x is greater than the upper edge of last bin, the Overflow bin is incremented
182//
183// If the storage of the sum of squares of weights has been triggered,
184// via the function Sumw2, then the sum of the squares of weights is incremented
185// by 1 in the bin corresponding to x.
186//
187// -*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
188Int_t MHCamera::Fill(Axis_t x)
189{
190#if ROOT_VERSION_CODE > ROOT_VERSION(3,05,00)
191 if (fBuffer) return BufferFill(x,1);
192#endif
193 const Int_t bin = (Int_t)x+1;
194 AddBinContent(bin);
195 if (fSumw2.fN)
196 fSumw2.fArray[bin]++;
197
198 if (bin<=0 || bin>fNcells-2)
199 return -1;
200
201 fTsumw++;
202 fTsumw2++;
203 fTsumwx += x;
204 fTsumwx2 += x*x;
205 return bin;
206}
207
208// ------------------------------------------------------------------------
209//
210// Taken from TH1D::Fill(). Uses the argument directly as bin index.
211// Doesn't increment the number of entries.
212//
213// -*-*-*-*-*-*Increment bin with abscissa X with a weight w*-*-*-*-*-*-*-*
214// =============================================
215//
216// if x is less than the low-edge of the first bin, the Underflow bin is incremented
217// if x is greater than the upper edge of last bin, the Overflow bin is incremented
218//
219// If the storage of the sum of squares of weights has been triggered,
220// via the function Sumw2, then the sum of the squares of weights is incremented
221// by w^2 in the bin corresponding to x.
222//
223// -*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
224Int_t MHCamera::Fill(Axis_t x, Stat_t w)
225{
226#if ROOT_VERSION_CODE > ROOT_VERSION(3,05,00)
227 if (fBuffer) return BufferFill(x,w);
228#endif
229 const Int_t bin = (Int_t)x+1;
230 AddBinContent(bin, w);
231 if (fSumw2.fN)
232 fSumw2.fArray[bin] += w*w;
233
234 if (bin<=0 || bin>fNcells-2)
235 return -1;
236
237 const Stat_t z = (w > 0 ? w : -w);
238 fTsumw += z;
239 fTsumw2 += z*z;
240 fTsumwx += z*x;
241 fTsumwx2 += z*x*x;
242 return bin;
243}
244
245// ------------------------------------------------------------------------
246//
247// Use x and y in millimeters
248//
249Int_t MHCamera::Fill(Axis_t x, Axis_t y, Stat_t w)
250{
251 if (fNcells<=1 || IsFreezed())
252 return -1;
253
254 for (Int_t idx=0; idx<fNcells-2; idx++)
255 {
256 MHexagon hex((*fGeomCam)[idx]);
257 if (hex.DistanceToPrimitive(x, y)>0)
258 continue;
259
260 SetUsed(idx);
261 return Fill(idx, w);
262 }
263 return -1;
264}
265
266// ------------------------------------------------------------------------
267//
268// Return the mean value of all entries which are used if all=kFALSE and
269// of all entries if all=kTRUE if sector<0. If sector>=0 only
270// entries with match the given sector are taken into account.
271//
272Stat_t MHCamera::GetMeanSectors(const TArrayI &sector, const TArrayI &aidx, Bool_t all) const
273{
274 if (fNcells<=1)
275 return 0;
276
277 Int_t n=0;
278
279 Stat_t mean = 0;
280 for (int i=0; i<fNcells-2; i++)
281 {
282 if ((all || IsUsed(i)) && MatchSector(i, sector, aidx))
283 {
284 mean += fArray[i+1];
285 n++;
286 }
287 }
288
289 return Profile(mean/n);
290}
291
292// ------------------------------------------------------------------------
293//
294// Return the variance of all entries which are used if all=kFALSE and
295// of all entries if all=kTRUE if sector<0. If sector>=0 only
296// entries with match the given sector are taken into account.
297//
298Stat_t MHCamera::GetRmsSectors(const TArrayI &sector, const TArrayI &aidx, Bool_t all) const
299{
300 if (fNcells<=1)
301 return -1;
302
303 Int_t n=0;
304
305 Stat_t sum = 0;
306 Stat_t sq = 0;
307 for (int i=0; i<fNcells-2; i++)
308 {
309 if ((all || IsUsed(i)) && MatchSector(i, sector, aidx))
310 {
311 sum += fArray[i+1];
312 sq += fArray[i+1]*fArray[i+1];
313 n++;
314 }
315 }
316
317 sum /= n;
318 sq /= n;
319
320 return Profile(TMath::Sqrt(sq-sum*sum));
321}
322
323// ------------------------------------------------------------------------
324//
325// Return the minimum contents of all pixels (if all is set, otherwise
326// only of all 'used' pixels), fMinimum if fMinimum set. If sector>=0
327// only pixels with matching sector number are taken into account.
328//
329Double_t MHCamera::GetMinimumSectors(const TArrayI &sector, const TArrayI &aidx, Bool_t all) const
330{
331 if (fMinimum != -1111)
332 return fMinimum;
333
334 if (fNcells<=1)
335 return 0;
336
337 Double_t minimum=FLT_MAX;
338
339 for (Int_t idx=0; idx<fNcells-2; idx++)
340 if (MatchSector(idx, sector, aidx) && (all || IsUsed(idx)) && fArray[idx+1]<minimum)
341 minimum = fArray[idx+1];
342
343 return Profile(minimum);
344}
345
346// ------------------------------------------------------------------------
347//
348// Return the maximum contents of all pixels (if all is set, otherwise
349// only of all 'used' pixels), fMaximum if fMaximum set. If sector>=0
350// only pixels with matching sector number are taken into account.
351//
352Double_t MHCamera::GetMaximumSectors(const TArrayI &sector, const TArrayI &aidx, Bool_t all) const
353{
354 if (fMaximum!=-1111)
355 return fMaximum;
356
357 if (fNcells<=1)
358 return 1;
359
360 Double_t maximum=-FLT_MAX;
361 for (Int_t idx=0; idx<fNcells-2; idx++)
362 if (MatchSector(idx, sector, aidx) && (all || IsUsed(idx)) && fArray[idx+1]>maximum)
363 maximum = fArray[idx+1];
364
365 return Profile(maximum);
366}
367
368// ------------------------------------------------------------------------
369//
370// Call this function to draw the camera layout into your canvas.
371// Setup a drawing canvas. Add this object and all child objects
372// (hexagons, etc) to the current pad. If no pad exists a new one is
373// created. (To access the 'real' pad containing the camera you have
374// to do a cd(1) in the current layer.
375//
376// To draw a camera into its own pad do something like:
377//
378// MGeomCamMagic m;
379// MHCamera *d=new MHCamera(m);
380//
381// TCanvas *c = new TCanvas;
382// c->Divide(2,1);
383// c->cd(1);
384//
385// d->FillRandom();
386// d->Draw();
387// d->SetBit(kCanDelete);
388//
389// There are several drawing options:
390// 'hist' Draw as a standard TH1 histogram (value vs. pixel index)
391// 'box' Draw hexagons which size is in respect to its contents
392// 'nocol' Leave the 'boxed' hexagons empty
393// 'pixelindex' Display the pixel index in each pixel
394// 'sectorindex' Display the sector index in each pixel
395// 'content' Display the relative content aligned to GetMaximum() and
396// GeMinimum() ((val-min)/(max-min))
397// 'proj' Display the y-projection of the histogram
398// 'same' Draw trandparent pixels on top of an existing pad. This
399// makes it possible to draw the camera image on top of an
400// existing TH2, but also allows for distorted camera images
401//
402void MHCamera::Draw(Option_t *option)
403{
404 const Bool_t hassame = TString(option).Contains("same", TString::kIgnoreCase) && gPad;
405
406 // root 3.02:
407 // gPad->SetFixedAspectRatio()
408 const Color_t col = gPad ? gPad->GetFillColor() : 16;
409 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas("CamDisplay", "Mars Camera Display", 656, 600);
410
411 if (!hassame)
412 {
413 pad->SetBorderMode(0);
414 pad->SetFillColor(col);
415
416 //
417 // Create an own pad for the MHCamera-Object which can be
418 // resized in paint to keep the correct aspect ratio
419 //
420 pad->Divide(1, 1, 0, 0, col);
421 pad->cd(1);
422 gPad->SetBorderMode(0);
423 }
424
425 AppendPad(option);
426 //fGeomCam->AppendPad();
427
428 //
429 // Do not change gPad. The user should not see, that Draw
430 // changes gPad...
431 //
432 if (!hassame)
433 pad->cd();
434}
435
436// ------------------------------------------------------------------------
437//
438// This is TObject::DrawClone but completely ignores
439// gROOT->GetSelectedPad(). tbretz had trouble with this in the past.
440// If this makes trouble please write a bug report.
441//
442TObject *MHCamera::DrawClone(Option_t *option) const
443{
444 // Draw a clone of this object in the current pad
445
446 //TVirtualPad *pad = gROOT->GetSelectedPad();
447 TVirtualPad *padsav = gPad;
448 //if (pad) pad->cd();
449
450 TObject *newobj = Clone();
451
452 if (!newobj)
453 return 0;
454
455 /*
456 if (pad) {
457 if (strlen(option)) pad->GetListOfPrimitives()->Add(newobj,option);
458 else pad->GetListOfPrimitives()->Add(newobj,GetDrawOption());
459 pad->Modified(kTRUE);
460 pad->Update();
461 if (padsav) padsav->cd();
462 return newobj;
463 }
464 */
465
466 TString opt(option);
467 opt.ToLower();
468
469 newobj->Draw(opt.IsNull() ? GetDrawOption() : option);
470
471 if (padsav)
472 padsav->cd();
473
474 return newobj;
475}
476
477// ------------------------------------------------------------------------
478//
479// Creates a TH1D which contains the projection of the contents of the
480// MHCamera onto the y-axis. The maximum and minimum are calculated
481// such that a slighly wider range than (GetMinimum(), GetMaximum()) is
482// displayed using THLimitsFinder::OptimizeLimits.
483//
484// If no name is given the newly allocated histogram is removed from
485// the current directory calling SetDirectory(0) in any other case
486// the newly created histogram is removed from the current directory
487// and added to gROOT such the gROOT->FindObject can find the histogram.
488//
489// If the standard name "_py" is given "_py" is appended to the name
490// of the MHCamera and the corresponding histogram is searched using
491// gROOT->FindObject and updated with the present projection.
492//
493// It is the responsibility of the user to make sure, that the newly
494// created histogram is freed correctly.
495//
496// Currently the new histogram is restrictred to 50 bins.
497// Maybe a optimal number can be calulated from the number of
498// bins on the x-axis of the MHCamera?
499//
500// The code was taken mainly from TH2::ProjectX such the interface
501// is more or less the same than to TH2-projections.
502//
503// If sector>=0 only entries with matching sector index are taken
504// into account.
505//
506TH1D *MHCamera::ProjectionS(const TArrayI &sector, const TArrayI &aidx, const char *name) const
507{
508 Int_t nbins = 50;
509
510 // Create the projection histogram
511 TString pname(name);
512 if (name=="_py")
513 {
514 pname.Prepend(GetName());
515 if (sector.GetSize()>0)
516 {
517 pname += ";";
518 for (int i=0; i<sector.GetSize(); i++)
519 pname += sector[i];
520 }
521 if (aidx.GetSize()>0)
522 {
523 pname += ";";
524 for (int i=0; i<aidx.GetSize(); i++)
525 pname += aidx[i];
526 }
527 }
528
529 TH1D *h1=0;
530
531 //check if histogram with identical name exist
532 TObject *h1obj = gROOT->FindObject(pname);
533 if (h1obj && h1obj->InheritsFrom("TH1D")) {
534 h1 = (TH1D*)h1obj;
535 h1->Reset();
536 }
537
538 if (!h1)
539 {
540 Double_t min = GetMinimumSectors(sector, aidx);
541 Double_t max = GetMaximumSectors(sector, aidx);
542
543 Int_t newbins=0;
544
545 THLimitsFinder::OptimizeLimits(nbins, newbins, min, max, kFALSE);
546
547 h1 = new TH1D(pname, GetTitle(), nbins, min, max);
548 h1->SetDirectory(pname.IsNull() ? NULL : gROOT);
549 h1->SetXTitle(GetYaxis()->GetTitle());
550 h1->SetYTitle("Counts");
551 //h1->Sumw2();
552 }
553
554 // Fill the projected histogram
555 for (Int_t idx=0; idx<fNcells-2; idx++)
556 if (IsUsed(idx) && MatchSector(idx, sector, aidx))
557 h1->Fill(GetBinContent(idx+1));
558
559 return h1;
560}
561
562// ------------------------------------------------------------------------
563//
564// Resizes the current pad so that the camera is displayed in its
565// correct aspect ratio
566//
567void MHCamera::SetRange()
568{
569 const Float_t range = fGeomCam->GetMaxRadius()*1.05;
570
571 //
572 // Maintain aspect ratio
573 //
574 const float ratio = TestBit(kNoLegend) ? 1 : 1.15;
575
576 //
577 // Calculate width and height of the current pad in pixels
578 //
579 Float_t w = gPad->GetWw();
580 Float_t h = gPad->GetWh()*ratio;
581
582 //
583 // This prevents the pad from resizing itself wrongly
584 //
585 if (gPad->GetMother() != gPad)
586 {
587 w *= gPad->GetMother()->GetAbsWNDC();
588 h *= gPad->GetMother()->GetAbsHNDC();
589 }
590
591 //
592 // Set Range (coordinate system) of pad
593 //
594 gPad->Range(-range, -range, (2*ratio-1)*range, range);
595
596 //
597 // Resize Pad to given ratio
598 //
599 if (h<w)
600 gPad->SetPad((1.-h/w)/2, 0, (h/w+1.)/2, 1);
601 else
602 gPad->SetPad(0, (1.-w/h)/2, 1, (w/h+1.)/2);
603}
604
605// ------------------------------------------------------------------------
606//
607// Updates the pixel colors and paints the pixels
608//
609void MHCamera::Update(Bool_t islog, Bool_t isbox, Bool_t iscol, Bool_t issame)
610{
611 Double_t min = GetMinimum(kFALSE);
612 Double_t max = GetMaximum(kFALSE);
613 if (min==FLT_MAX)
614 {
615 min = 0;
616 max = 1;
617 }
618
619 if (min==max)
620 max += 1;
621
622 if (!issame)
623 UpdateLegend(min, max, islog);
624
625 // Try to estimate the units of the current display. This is only
626 // necessary for 'same' option and allows distorted images of the camera!
627 const Float_t maxr = (1-fGeomCam->GetConvMm2Deg())*fGeomCam->GetMaxRadius()/2;
628 const Float_t conv = !issame ||
629 gPad->GetX1()<-maxr || gPad->GetY1()<-maxr ||
630 gPad->GetX2()> maxr || gPad->GetY2()>maxr ? 1 : fGeomCam->GetConvMm2Deg();
631
632 MHexagon hex;
633 for (Int_t i=0; i<fNcells-2; i++)
634 {
635 hex.SetFillStyle(issame ? 4000 : 1001);
636
637 if (!issame)
638 {
639 if (IsUsed(i) && iscol)
640 {
641 if (TMath::IsNaN(fArray[i+1]))
642 gLog << warn << "MHCamera::Update: " << GetName() << " <" << GetTitle() << "> - Pixel Index #" << i << " contents is NaN (Not a Number)..." << endl;
643
644 hex.SetFillColor(GetColor(GetBinContent(i+1), min, max, islog));
645 }
646 else
647 hex.SetFillColor(10);
648 }
649
650 const MGeomPix &pix = (*fGeomCam)[i];
651
652 const Float_t x = pix.GetX()*conv;
653 const Float_t y = pix.GetY()*conv;
654 const Float_t d = pix.GetD()*conv;
655
656 if (!isbox)
657 hex.PaintHexagon(x, y, d);
658 else
659 if (IsUsed(i) && !TMath::IsNaN(fArray[i+1]))
660 {
661 Float_t size = d*(GetBinContent(i+1)-min)/(max-min);
662 if (size>d)
663 size=d;
664 hex.PaintHexagon(x, y, size);
665 }
666 }
667}
668
669// ------------------------------------------------------------------------
670//
671// Print minimum and maximum
672//
673void MHCamera::Print(Option_t *) const
674{
675 gLog << all << "Minimum: " << GetMinimum();
676 if (fMinimum==-1111)
677 gLog << " <autoscaled>";
678 gLog << endl;
679 gLog << "Maximum: " << GetMaximum();
680 if (fMaximum==-1111)
681 gLog << " <autoscaled>";
682 gLog << endl;
683}
684
685// ------------------------------------------------------------------------
686//
687// Paint the y-axis title
688//
689void MHCamera::PaintAxisTitle()
690{
691 const Float_t range = fGeomCam->GetMaxRadius()*1.05;
692 const Float_t w = (1 + 1.5/sqrt((float)(fNcells-2)))*range;
693
694 TLatex *ptitle = new TLatex(w, -.90*range, GetYaxis()->GetTitle());
695
696 ptitle->SetTextSize(0.05);
697 ptitle->SetTextAlign(21);
698
699 // box with the histogram title
700 ptitle->SetTextColor(gStyle->GetTitleTextColor());
701#if ROOT_VERSION_CODE > ROOT_VERSION(3,05,01)
702 ptitle->SetTextFont(gStyle->GetTitleFont(""));
703#endif
704 ptitle->Paint();
705}
706
707// ------------------------------------------------------------------------
708//
709// Paints the camera.
710//
711void MHCamera::Paint(Option_t *o)
712{
713 if (fNcells<=1)
714 return;
715
716 TString opt(o);
717 opt.ToLower();
718
719 if (opt.Contains("hist"))
720 {
721 opt.ReplaceAll("hist", "");
722 TH1D::Paint(opt);
723 return;
724 }
725
726 if (opt.Contains("proj"))
727 {
728 opt.ReplaceAll("proj", "");
729 Projection(GetName())->Paint(opt);
730 return;
731 }
732
733 const Bool_t hassame = opt.Contains("same");
734 const Bool_t hasbox = opt.Contains("box");
735 const Bool_t hascol = hasbox ? !opt.Contains("nocol") : kTRUE;
736
737 if (!hassame)
738 {
739 gPad->Clear();
740
741 // Maintain aspect ratio
742 SetRange();
743
744 if (GetPainter())
745 {
746 // Paint statistics
747 if (!TestBit(TH1::kNoStats))
748 fPainter->PaintStat(gStyle->GetOptStat(), NULL);
749
750 // Paint primitives (pixels, color legend, photons, ...)
751 if (fPainter->InheritsFrom(THistPainter::Class()))
752 {
753 static_cast<THistPainter*>(fPainter)->MakeChopt("");
754 static_cast<THistPainter*>(fPainter)->PaintTitle();
755 }
756 }
757 }
758
759 // Update Contents of the pixels and paint legend
760 Update(gPad->GetLogy(), hasbox, hascol, hassame);
761
762 if (!hassame)
763 PaintAxisTitle();
764
765 if (opt.Contains("pixelindex"))
766 PaintIndices(0);
767 if (opt.Contains("sectorindex"))
768 PaintIndices(1);
769 if (opt.Contains("content"))
770 PaintIndices(2);
771}
772
773// ------------------------------------------------------------------------
774//
775// With this function you can change the color palette. For more
776// information see TStyle::SetPalette. Only palettes with 50 colors
777// are allowed.
778// In addition you can use SetPalette(52, 0) to create an inverse
779// deep blue sea palette
780//
781void MHCamera::SetPalette(Int_t ncolors, Int_t *colors)
782{
783 //
784 // If not enough colors are specified skip this.
785 //
786 if (ncolors>1 && ncolors<50)
787 {
788 gLog << err << "MHCamera::SetPalette: Only default palettes with 50 colors are allowed... ignored." << endl;
789 return;
790 }
791
792 //
793 // If ncolors==52 create a reversed deep blue sea palette
794 //
795 if (ncolors==52)
796 {
797 gStyle->SetPalette(51, NULL);
798 TArrayI c(kItemsLegend);
799 for (int i=0; i<kItemsLegend; i++)
800 c[kItemsLegend-i-1] = gStyle->GetColorPalette(i);
801 gStyle->SetPalette(kItemsLegend, c.GetArray());
802 }
803 else
804 gStyle->SetPalette(ncolors, colors);
805
806 fColors.Set(kItemsLegend);
807 for (int i=0; i<kItemsLegend; i++)
808 fColors[i] = gStyle->GetColorPalette(i);
809}
810
811
812void MHCamera::SetPrettyPalette()
813{
814 if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
815 SetPalette(1, 0);
816}
817
818void MHCamera::SetDeepBlueSeaPalette()
819{
820 if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
821 SetPalette(51, 0);
822}
823
824void MHCamera::SetInvDeepBlueSeaPalette()
825{
826 if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
827 SetPalette(52, 0);
828}
829
830void MHCamera::PaintIndices(Int_t type)
831{
832 if (fNcells<=1)
833 return;
834
835 const Double_t min = GetMinimum();
836 const Double_t max = GetMaximum();
837
838 if (type==2 && max==min)
839 return;
840
841 TText txt;
842 txt.SetTextFont(122);
843 txt.SetTextAlign(22); // centered/centered
844
845 for (Int_t i=0; i<fNcells-2; i++)
846 {
847 const MGeomPix &h = (*fGeomCam)[i];
848
849 TString num;
850 switch (type)
851 {
852 case 0: num += i; break;
853 case 1: num += h.GetSector(); break;
854 case 2: num += (Int_t)((fArray[i+1]-min)/(max-min)); break;
855 }
856
857 // FIXME: Should depend on the color of the pixel...
858 //(GetColor(GetBinContent(i+1), min, max, 0));
859 txt.SetTextColor(kRed);
860 txt.SetTextSize(0.3*h.GetD()/fGeomCam->GetMaxRadius()/1.05);
861 txt.PaintText(h.GetX(), h.GetY(), num);
862 }
863}
864
865// ------------------------------------------------------------------------
866//
867// Call this function to add a MCamEvent on top of the present contents.
868//
869void MHCamera::AddCamContent(const MCamEvent &event, Int_t type)
870{
871 if (fNcells<=1 || IsFreezed())
872 return;
873
874 // FIXME: Security check missing!
875 for (Int_t idx=0; idx<fNcells-2; idx++)
876 {
877 Double_t val=0;
878 if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
879 SetUsed(idx);
880
881 Fill(idx, val); // FIXME: Slow!
882 }
883 fEntries++;
884}
885
886// ------------------------------------------------------------------------
887//
888// Call this function to add a MCamEvent on top of the present contents.
889//
890void MHCamera::SetCamError(const MCamEvent &evt, Int_t type)
891{
892
893 if (fNcells<=1 || IsFreezed())
894 return;
895
896 // FIXME: Security check missing!
897 for (Int_t idx=0; idx<fNcells-2; idx++)
898 {
899 Double_t val=0;
900 if (evt.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
901 SetUsed(idx);
902
903 SetBinError(idx+1, val); // FIXME: Slow!
904 }
905}
906
907Stat_t MHCamera::GetBinError(Int_t bin) const
908{
909 Double_t rc = 0;
910 if (TestBit(kVariance) && GetEntries()>0) // error on the mean
911 {
912 const Double_t error = fSumw2.fArray[bin]/GetEntries();
913 const Double_t val = fArray[bin]/GetEntries();
914 rc = TMath::Sqrt(error - val*val);
915 }
916 else
917 {
918 rc = TH1D::GetBinError(bin);
919 rc = Profile(rc);
920 }
921
922 return rc;
923}
924
925// ------------------------------------------------------------------------
926//
927// Call this function to add a MHCamera on top of the present contents.
928// Type:
929// 0) bin content
930// 1) errors
931// 2) rel. errors
932//
933void MHCamera::AddCamContent(const MHCamera &d, Int_t type)
934{
935 if (fNcells!=d.fNcells || IsFreezed())
936 return;
937
938 // FIXME: Security check missing!
939 for (Int_t idx=0; idx<fNcells-2; idx++)
940 if (d.IsUsed(idx))
941 SetUsed(idx);
942
943 switch (type)
944 {
945 case 1:
946 for (Int_t idx=0; idx<fNcells-2; idx++)
947 Fill(idx, d.GetBinError(idx+1));
948 break;
949 case 2:
950 for (Int_t idx=0; idx<fNcells-2; idx++)
951 if (d.GetBinContent(idx+1)!=0)
952 Fill(idx, TMath::Abs(d.GetBinError(idx+1)/d.GetBinContent(idx+1)));
953 break;
954 default:
955 for (Int_t idx=0; idx<fNcells-2; idx++)
956 Fill(idx, d.GetBinContent(idx+1));
957 break;
958 }
959 fEntries++;
960}
961
962// ------------------------------------------------------------------------
963//
964// Call this function to add a TArrayD on top of the present contents.
965//
966void MHCamera::AddCamContent(const TArrayD &event, const TArrayC *used)
967{
968 if (event.GetSize()!=fNcells-2 || IsFreezed())
969 return;
970
971 if (used && used->GetSize()!=fNcells-2)
972 return;
973
974 for (Int_t idx=0; idx<fNcells-2; idx++)
975 {
976 Fill(idx, const_cast<TArrayD&>(event)[idx]); // FIXME: Slow!
977
978 if (!used || (*const_cast<TArrayC*>(used))[idx])
979 SetUsed(idx);
980 }
981 fEntries++;
982}
983
984// ------------------------------------------------------------------------
985//
986// Call this function to add a MCamEvent on top of the present contents.
987// 1 is added to each pixel if the contents of MCamEvent>threshold
988//
989void MHCamera::CntCamContent(const MCamEvent &event, Double_t threshold, Int_t type)
990{
991 if (fNcells<=1 || IsFreezed())
992 return;
993
994 // FIXME: Security check missing!
995 for (Int_t idx=0; idx<fNcells-2; idx++)
996 {
997 Double_t val=threshold;
998 if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
999 SetUsed(idx);
1000
1001 if (val>threshold)
1002 Fill(idx);
1003 }
1004 fEntries++;
1005}
1006
1007// ------------------------------------------------------------------------
1008//
1009// Call this function to add a TArrayD on top of the present contents.
1010// 1 is added to each pixel if the contents of MCamEvent>threshold
1011//
1012void MHCamera::CntCamContent(const TArrayD &event, Double_t threshold, Bool_t ispos)
1013{
1014 if (event.GetSize()!=fNcells-2 || IsFreezed())
1015 return;
1016
1017 for (Int_t idx=0; idx<fNcells-2; idx++)
1018 {
1019 if (const_cast<TArrayD&>(event)[idx]>threshold)
1020 Fill(idx);
1021
1022 if (!ispos || fArray[idx+1]>0)
1023 SetUsed(idx);
1024 }
1025 fEntries++;
1026}
1027
1028// ------------------------------------------------------------------------
1029//
1030// Fill the pixels with random contents.
1031//
1032void MHCamera::FillRandom()
1033{
1034 if (fNcells<=1 || IsFreezed())
1035 return;
1036
1037 Reset();
1038
1039 // FIXME: Security check missing!
1040 for (Int_t idx=0; idx<fNcells-2; idx++)
1041 {
1042 Fill(idx, gRandom->Uniform()*fGeomCam->GetPixRatio(idx));
1043 SetUsed(idx);
1044 }
1045 fEntries=1;
1046}
1047
1048
1049// ------------------------------------------------------------------------
1050//
1051// The array must be in increasing order, eg: 2.5, 3.7, 4.9
1052// The values in each bin are replaced by the interval in which the value
1053// fits. In the example we have four intervals
1054// (<2.5, 2.5-3.7, 3.7-4.9, >4.9). Maximum and minimum are set
1055// accordingly.
1056//
1057void MHCamera::SetLevels(const TArrayF &arr)
1058{
1059 if (fNcells<=1)
1060 return;
1061
1062 for (Int_t i=0; i<fNcells-2; i++)
1063 {
1064 if (!IsUsed(i))
1065 continue;
1066
1067 Int_t j = arr.GetSize();
1068 while (j && fArray[i+1]<arr[j-1])
1069 j--;
1070
1071 fArray[i+1] = j;
1072 }
1073 SetMaximum(arr.GetSize());
1074 SetMinimum(0);
1075}
1076
1077// ------------------------------------------------------------------------
1078//
1079// Reset the all pixel colors to a default value
1080//
1081void MHCamera::Reset(Option_t *opt)
1082{
1083 if (fNcells<=1 || IsFreezed())
1084 return;
1085
1086 TH1::Reset(opt);
1087
1088 for (Int_t i=0; i<fNcells-2; i++)
1089 {
1090 fArray[i+1]=0;
1091 ResetUsed(i);
1092 }
1093 fArray[0] = 0;
1094 fArray[fNcells-1] = 0;
1095}
1096
1097// ------------------------------------------------------------------------
1098//
1099// Here we calculate the color index for the current value.
1100// The color index is defined with the class TStyle and the
1101// Color palette inside. We use the command gStyle->SetPalette(1,0)
1102// for the display. So we have to convert the value "wert" into
1103// a color index that fits the color palette.
1104// The range of the color palette is defined by the values fMinPhe
1105// and fMaxRange. Between this values we have 50 color index, starting
1106// with 0 up to 49.
1107//
1108Int_t MHCamera::GetColor(Float_t val, Float_t min, Float_t max, Bool_t islog)
1109{
1110 if (TMath::IsNaN(val)) // FIXME: gLog!
1111 return 10;
1112
1113 //
1114 // first treat the over- and under-flows
1115 //
1116 const Int_t maxcolidx = kItemsLegend-1;
1117
1118 if (val >= max)
1119 return fColors[maxcolidx];
1120
1121 if (val <= min)
1122 return fColors[0];
1123
1124 //
1125 // calculate the color index
1126 //
1127 Float_t ratio;
1128 if (islog && min>0)
1129 ratio = log10(val/min) / log10(max/min);
1130 else
1131 ratio = (val-min) / (max-min);
1132
1133 const Int_t colidx = (Int_t)(ratio*maxcolidx + .5);
1134 return fColors[colidx];
1135}
1136
1137TPaveStats *MHCamera::GetStatisticBox()
1138{
1139 TObject *obj = 0;
1140
1141 TIter Next(fFunctions);
1142 while ((obj = Next()))
1143 if (obj->InheritsFrom(TPaveStats::Class()))
1144 return static_cast<TPaveStats*>(obj);
1145
1146 return NULL;
1147}
1148
1149// ------------------------------------------------------------------------
1150//
1151// Change the text on the legend according to the range of the Display
1152//
1153void MHCamera::UpdateLegend(Float_t min, Float_t max, Bool_t islog)
1154{
1155 const Float_t range = fGeomCam->GetMaxRadius()*1.05;
1156
1157 TArrow arr;
1158 arr.PaintArrow(-range*.9, -range*.9, -range*.6, -range*.9, 0.025);
1159 arr.PaintArrow(-range*.9, -range*.9, -range*.9, -range*.6, 0.025);
1160
1161 TString text;
1162 text += (int)(range*.3);
1163 text += "mm";
1164
1165 TText newtxt2;
1166 newtxt2.SetTextSize(0.04);
1167 newtxt2.PaintText(-range*.85, -range*.85, text);
1168
1169 text = "";
1170 text += (float)((int)(range*.3*fGeomCam->GetConvMm2Deg()*10))/10;
1171 text += "\\circ";
1172 text = text.Strip(TString::kLeading);
1173
1174 TLatex latex;
1175 latex.PaintLatex(-range*.85, -range*.75, 0, 0.04, text);
1176
1177 if (TestBit(kNoLegend))
1178 return;
1179
1180 TPaveStats *stats = GetStatisticBox();
1181
1182 const Float_t hndc = 0.92 - (stats ? stats->GetY1NDC() : 1);
1183 const Float_t H = (0.75-hndc)*range;
1184 const Float_t offset = hndc*range;
1185
1186 const Float_t h = 2./kItemsLegend;
1187 const Float_t w = range/sqrt((float)(fNcells-2));
1188
1189 TBox newbox;
1190 TText newtxt;
1191 newtxt.SetTextSize(0.03);
1192 newtxt.SetTextAlign(12);
1193#if ROOT_VERSION_CODE > ROOT_VERSION(3,01,06)
1194 newtxt.SetBit(/*kNoContextMenu|*/kCannotPick);
1195 newbox.SetBit(/*kNoContextMenu|*/kCannotPick);
1196#endif
1197
1198 const Float_t step = (islog && min>0 ? log10(max/min) : max-min) / kItemsLegend;
1199 const Int_t firsts = step*3 < 1e-8 ? 8 : (Int_t)floor(log10(step*3));
1200 const TString opt = Form("%%.%if", firsts>0 ? 0 : TMath::Abs(firsts));
1201
1202 for (Int_t i=0; i<kItemsLegend+1; i+=3)
1203 {
1204 Float_t val;
1205 if (islog && min>0)
1206 val = pow(10, step*i) * min;
1207 else
1208 val = min + step*i;
1209
1210 //const bool dispexp = max-min>1.5 && fabs(val)>0.1 && fabs(val)<1e6;
1211 newtxt.PaintText(range+1.5*w, H*(i*h-1)-offset, Form(opt, val));
1212 }
1213
1214 for (Int_t i=0; i<kItemsLegend; i++)
1215 {
1216 newbox.SetFillColor(fColors[i]);
1217 newbox.PaintBox(range, H*(i*h-1)-offset, range+w, H*((i+1)*h-1)-offset);
1218 }
1219}
1220
1221// ------------------------------------------------------------------------
1222//
1223// Save primitive as a C++ statement(s) on output stream out
1224//
1225void MHCamera::SavePrimitive(ofstream &out, Option_t *opt)
1226{
1227 gLog << err << "MHCamera::SavePrimitive: Must be rewritten!" << endl;
1228 /*
1229 if (!gROOT->ClassSaved(TCanvas::Class()))
1230 fDrawingPad->SavePrimitive(out, opt);
1231
1232 out << " " << fDrawingPad->GetName() << "->SetWindowSize(";
1233 out << fDrawingPad->GetWw() << "," << fDrawingPad->GetWh() << ");" << endl;
1234 */
1235}
1236
1237// ------------------------------------------------------------------------
1238//
1239// compute the distance of a point (px,py) to the Camera
1240// this functions needed for graphical primitives, that
1241// means without this function you are not able to interact
1242// with the graphical primitive with the mouse!!!
1243//
1244// All calcutations are done in pixel coordinates
1245//
1246Int_t MHCamera::DistancetoPrimitive(Int_t px, Int_t py)
1247{
1248 if (fNcells<=1)
1249 return 999999;
1250
1251 const Int_t kMaxDiff = 7;
1252
1253 TPaveStats *box = (TPaveStats*)gPad->GetPrimitive("stats");
1254 if (box)
1255 {
1256 const Double_t w = box->GetY2NDC()-box->GetY1NDC();
1257 box->SetX1NDC(gStyle->GetStatX()-gStyle->GetStatW());
1258 box->SetY1NDC(gStyle->GetStatY()-w);
1259 box->SetX2NDC(gStyle->GetStatX());
1260 box->SetY2NDC(gStyle->GetStatY());
1261 }
1262
1263 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1264 return TH1D::DistancetoPrimitive(px, py);
1265
1266 for (Int_t i=0; i<fNcells-2; i++)
1267 {
1268 MHexagon hex((*fGeomCam)[i]);
1269 if (hex.DistancetoPrimitive(px, py)==0)
1270 return 0;
1271 }
1272
1273 if (!box)
1274 return 999999;
1275
1276 const Int_t dist = box->DistancetoPrimitive(px, py);
1277 if (dist > kMaxDiff)
1278 return 999999;
1279
1280 gPad->SetSelected(box);
1281 return dist;
1282}
1283
1284// ------------------------------------------------------------------------
1285//
1286//
1287Int_t MHCamera::GetPixelIndex(Int_t px, Int_t py) const
1288{
1289 if (fNcells<=1)
1290 return -1;
1291
1292 Int_t i;
1293 for (i=0; i<fNcells-2; i++)
1294 {
1295 MHexagon hex((*fGeomCam)[i]);
1296 if (hex.DistancetoPrimitive(px, py)>0)
1297 continue;
1298
1299 return i;
1300 }
1301 return -1;
1302}
1303
1304// ------------------------------------------------------------------------
1305//
1306// Returns string containing info about the object at position (px,py).
1307// Returned string will be re-used (lock in MT environment).
1308//
1309char *MHCamera::GetObjectInfo(Int_t px, Int_t py) const
1310{
1311 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1312 return TH1D::GetObjectInfo(px, py);
1313
1314 static char info[128];
1315
1316 const Int_t idx=GetPixelIndex(px, py);
1317
1318 if (idx<0)
1319 return TObject::GetObjectInfo(px, py);
1320
1321 sprintf(info, "Software Pixel Index: %d (Hardware Id=%d)", idx, idx+1);
1322 return info;
1323}
1324
1325// ------------------------------------------------------------------------
1326//
1327// Add a MCamEvent which should be displayed when the user clicks on a
1328// pixel.
1329// Warning: The object MUST inherit from TObject AND MCamEvent
1330//
1331void MHCamera::AddNotify(TObject *obj)
1332{
1333 // Make sure, that the object derives from MCamEvent!
1334 MCamEvent *evt = dynamic_cast<MCamEvent*>(obj);
1335 if (!evt)
1336 {
1337 gLog << err << "ERROR: MHCamera::AddNotify - TObject doesn't inherit from MCamEvent... ignored." << endl;
1338 return;
1339 }
1340
1341 // Make sure, that it is deleted from the list too, if the obj is deleted
1342 obj->SetBit(kMustCleanup);
1343
1344 // Add object to list
1345 fNotify->Add(obj);
1346}
1347
1348// ------------------------------------------------------------------------
1349//
1350// Execute a mouse event on the camera
1351//
1352void MHCamera::ExecuteEvent(Int_t event, Int_t px, Int_t py)
1353{
1354 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1355 {
1356 TH1D::ExecuteEvent(event, px, py);
1357 return;
1358 }
1359 //if (event==kMouseMotion && fStatusBar)
1360 // fStatusBar->SetText(GetObjectInfo(px, py), 0);
1361 if (event!=kButton1Down)
1362 return;
1363
1364 const Int_t idx = GetPixelIndex(px, py);
1365 if (idx<0)
1366 return;
1367
1368 gLog << all << GetTitle() << " <" << GetName() << ">" << endl;
1369 gLog << "Software Pixel Idx: " << idx << endl;
1370 gLog << "Hardware Pixel Id: " << idx+1 << endl;
1371 gLog << "Contents: " << GetBinContent(idx+1);
1372 if (GetBinError(idx+1)>0)
1373 gLog << " +/- " << GetBinError(idx+1);
1374 gLog << " <" << (IsUsed(idx)?"on":"off") << ">" << endl;
1375
1376 if (fNotify && fNotify->GetSize()>0)
1377 {
1378 // FIXME: Is there a simpler and more convinient way?
1379
1380 // The name which is created here depends on the instance of
1381 // MHCamera and on the pad on which it is drawn --> The name
1382 // is unique. For ExecuteEvent gPad is always correctly set.
1383 const TString name = Form("%p;%p;PixelContent", this, gPad);
1384
1385 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
1386 if (old)
1387 old->cd();
1388 else
1389 new TCanvas(name);
1390
1391 /*
1392 TIter Next(gPad->GetListOfPrimitives());
1393 TObject *o;
1394 while (o=Next()) cout << o << ": " << o->GetName() << " " << o->IsA()->GetName() << endl;
1395 */
1396
1397 // FIXME: Make sure, that the old histograms are really deleted.
1398 // Are they already deleted?
1399
1400 // The dynamic_cast is necessary here: We cannot use ForEach
1401 TIter Next(fNotify);
1402 MCamEvent *evt;
1403 while ((evt=dynamic_cast<MCamEvent*>(Next())))
1404 evt->DrawPixelContent(idx);
1405
1406 gPad->Modified();
1407 gPad->Update();
1408 }
1409}
1410
1411UInt_t MHCamera::GetNumPixels() const
1412{
1413 return fGeomCam->GetNumPixels();
1414}
1415
1416TH1 *MHCamera::DrawCopy() const
1417{
1418 gPad=NULL;
1419 return TH1D::DrawCopy(fName+";cpy");
1420}
Note: See TracBrowser for help on using the repository browser.