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

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