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

Last change on this file since 2478 was 2476, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 29.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////////////////////////////////////////////////////////////////////////////
36#include "MHCamera.h"
37
38#include <fstream>
39#include <iostream>
40
41#include <TBox.h>
42#include <TArrow.h>
43#include <TLatex.h>
44#include <TStyle.h>
45#include <TCanvas.h>
46#include <TArrayF.h>
47#include <TRandom.h>
48#include <TPaveText.h>
49#include <TPaveStats.h>
50#include <TClonesArray.h>
51#include <THistPainter.h>
52
53#include "MLog.h"
54#include "MLogManip.h"
55
56#include "MH.h"
57#include "MHexagon.h"
58
59#include "MGeomPix.h"
60#include "MGeomCam.h"
61
62#include "MCerPhotPix.h"
63#include "MCerPhotEvt.h"
64
65#include "MPedestalPix.h"
66#include "MPedestalCam.h"
67
68#include "MCurrents.h"
69#include "MCamEvent.h"
70
71#include "MImgCleanStd.h"
72
73
74#define kItemsLegend 48 // see SetPalette(1,0)
75
76ClassImp(MHCamera);
77
78using namespace std;
79
80void MHCamera::Init()
81{
82 UseCurrentStyle();
83
84 SetDirectory(NULL);
85
86 SetLineColor(kGreen);
87 SetMarkerStyle(kFullDotMedium);
88 SetXTitle("Pixel Index");
89
90 fNotify = new TList;
91
92#if ROOT_VERSION_CODE < ROOT_VERSION(3,01,06)
93 SetPalette(1, 0);
94#else
95 SetPalette(51, 0);
96#endif
97}
98
99// ------------------------------------------------------------------------
100//
101// Default Constructor. To be used by the root system ONLY.
102//
103MHCamera::MHCamera() : TH1D(), fGeomCam(NULL), fColors(kItemsLegend)
104{
105 Init();
106}
107
108// ------------------------------------------------------------------------
109//
110// Constructor. Makes a clone of MGeomCam. Removed the TH1D from the
111// current directory. Calls Sumw2(). Set the histogram line color
112// (for error bars) to Green and the marker style to kFullDotMedium.
113//
114MHCamera::MHCamera(const MGeomCam &geom, const char *name, const char *title)
115: fGeomCam(NULL), /*TH1D(name, title, geom.GetNumPixels(), -0.5, geom.GetNumPixels()-0.5),
116fUsed(geom.GetNumPixels()),*/ fColors(kItemsLegend)
117{
118 //fGeomCam = (MGeomCam*)geom.Clone();
119 SetGeometry(geom, name, title);
120 Init();
121
122 //
123 // root 3.02
124 // * base/inc/TObject.h:
125 // register BIT(8) as kNoContextMenu. If an object has this bit set it will
126 // not get an automatic context menu when clicked with the right mouse button.
127}
128
129void MHCamera::SetGeometry(const MGeomCam &geom, const char *name, const char *title)
130{
131 SetNameTitle(name, title);
132
133 TAxis &x = *GetXaxis();
134
135 SetBins(geom.GetNumPixels(), 0, 1);
136 x.Set(geom.GetNumPixels(), -0.5, geom.GetNumPixels()-0.5);
137
138 //SetBins(geom.GetNumPixels(), -0.5, geom.GetNumPixels()-0.5);
139 //Rebuild();
140
141 Sumw2(); // necessary?
142
143 if (fGeomCam)
144 delete fGeomCam;
145 fGeomCam = (MGeomCam*)geom.Clone();
146
147 fUsed.Set(geom.GetNumPixels());
148 for (Int_t i=0; i<fNcells-2; i++)
149 ResetUsed(i);
150
151}
152
153// ------------------------------------------------------------------------
154//
155// Destructor. Deletes the cloned fGeomCam and the notification list.
156//
157MHCamera::~MHCamera()
158{
159 if (fGeomCam)
160 delete fGeomCam;
161 if (fNotify)
162 delete fNotify;
163}
164
165// ------------------------------------------------------------------------
166//
167// Taken from TH1D::Fill(). Uses the argument directly as bin index.
168// Doesn't increment the number of entries.
169//
170// -*-*-*-*-*-*-*-*Increment bin with abscissa X by 1*-*-*-*-*-*-*-*-*-*-*
171// ==================================
172//
173// if x is less than the low-edge of the first bin, the Underflow bin is incremented
174// if x is greater than the upper edge of last bin, the Overflow bin is incremented
175//
176// If the storage of the sum of squares of weights has been triggered,
177// via the function Sumw2, then the sum of the squares of weights is incremented
178// by 1 in the bin corresponding to x.
179//
180// -*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
181Int_t MHCamera::Fill(Axis_t x)
182{
183
184#if ROOT_VERSION_CODE > ROOT_VERSION(3,05,00)
185 if (fBuffer) return BufferFill(x,1);
186#endif
187 const Int_t bin = (Int_t)x+1;
188 AddBinContent(bin);
189 if (fSumw2.fN)
190 fSumw2.fArray[bin]++;
191
192 if (bin<=0 || bin>fNcells-2)
193 return -1;
194
195 fTsumw++;
196 fTsumw2++;
197 fTsumwx += x;
198 fTsumwx2 += x*x;
199 return bin;
200}
201
202// ------------------------------------------------------------------------
203//
204// Taken from TH1D::Fill(). Uses the argument directly as bin index.
205// Doesn't increment the number of entries.
206//
207// -*-*-*-*-*-*Increment bin with abscissa X with a weight w*-*-*-*-*-*-*-*
208// =============================================
209//
210// if x is less than the low-edge of the first bin, the Underflow bin is incremented
211// if x is greater than the upper edge of last bin, the Overflow bin is incremented
212//
213// If the storage of the sum of squares of weights has been triggered,
214// via the function Sumw2, then the sum of the squares of weights is incremented
215// by w^2 in the bin corresponding to x.
216//
217// -*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
218Int_t MHCamera::Fill(Axis_t x, Stat_t w)
219{
220#if ROOT_VERSION_CODE > ROOT_VERSION(3,05,00)
221 if (fBuffer) return BufferFill(x,w);
222#endif
223 const Int_t bin = (Int_t)x+1;
224 AddBinContent(bin, w);
225 if (fSumw2.fN)
226 fSumw2.fArray[bin] += w*w;
227
228 if (bin<=0 || bin>fNcells-2)
229 return -1;
230
231 const Stat_t z = (w > 0 ? w : -w);
232 fTsumw += z;
233 fTsumw2 += z*z;
234 fTsumwx += z*x;
235 fTsumwx2 += z*x*x;
236 return bin;
237}
238
239// ------------------------------------------------------------------------
240//
241// Use x and y in millimeters
242//
243Int_t MHCamera::Fill(Axis_t x, Axis_t y, Stat_t w)
244{
245 if (fNcells<=1 || IsFreezed())
246 return -1;
247
248 for (Int_t idx=0; idx<fNcells-2; idx++)
249 {
250 MHexagon hex((*fGeomCam)[idx]);
251 if (hex.DistanceToPrimitive(x, y)>0)
252 continue;
253
254 SetUsed(idx);
255 return Fill(idx, w);
256 }
257 return -1;
258}
259
260Stat_t MHCamera::GetMean(Int_t axis) const
261{
262 if (fNcells<=1)
263 return 0;
264
265 Stat_t mean = 0;
266 for (int i=1; i<fNcells-1; i++)
267 mean += fArray[i];
268
269 return Profile(mean/(fNcells-2));
270}
271
272Stat_t MHCamera::GetRMS(Int_t axis) const
273{
274 if (fNcells<=1)
275 return -1;
276
277 const Int_t n = fNcells-2;
278
279 Stat_t sum = 0;
280 Stat_t sq = 0;
281 for (int i=1; i<n+1; i++)
282 {
283 sum += fArray[i];
284 sq += fArray[i]*fArray[i];
285 }
286
287 sum /= n;
288 sq /= n;
289
290 return Profile(sqrt(sq-sum*sum));
291}
292
293// ------------------------------------------------------------------------
294//
295// Return the minimum contents of all pixels (if all is set, otherwise
296// only of all 'used' pixels), fMinimum if fMinimum set
297//
298Double_t MHCamera::GetMinimum(Bool_t all) const
299{
300 if (fMinimum != -1111)
301 return Profile(fMinimum);
302
303 if (fNcells<=1)
304 return 0;
305
306 Double_t minimum=FLT_MAX;
307
308 if (all)
309 {
310 for (Int_t idx=0; idx<fNcells-2; idx++)
311 if (fArray[idx+1] < minimum)
312 minimum = fArray[idx+1];
313 }
314 else
315 {
316 for (Int_t idx=0; idx<fNcells-2; idx++)
317 if (IsUsed(idx) && fArray[idx+1] < minimum)
318 minimum = fArray[idx+1];
319 }
320 return Profile(minimum);
321}
322
323// ------------------------------------------------------------------------
324//
325// Return the maximum contents of all pixels (if all is set, otherwise
326// only of all 'used' pixels), fMaximum if fMaximum set
327//
328Double_t MHCamera::GetMaximum(Bool_t all) const
329{
330 if (fMaximum!=-1111)
331 return Profile(fMaximum);
332
333 if (fNcells<=1)
334 return 1;
335
336 Double_t maximum=-FLT_MAX;
337 if (all)
338 {
339 for (Int_t idx=0; idx<fNcells-2; idx++)
340 if (fArray[idx+1] > maximum)
341 maximum = fArray[idx+1];
342 }
343 else
344 {
345 for (Int_t idx=0; idx<fNcells-2; idx++)
346 if (IsUsed(idx) && fArray[idx+1] > maximum)
347 maximum = fArray[idx+1];
348 }
349 return Profile(maximum);
350}
351
352// ------------------------------------------------------------------------
353//
354// Call this function to draw the camera layout into your canvas.
355// Setup a drawing canvas. Add this object and all child objects
356// (hexagons, etc) to the current pad. If no pad exists a new one is
357// created.
358//
359// To draw a camera into its own pad do something like:
360//
361// TCanvas *c = new TCanvas;
362// c->Divide(2,1);
363// MGeomCamMagic m;
364// MHCamera *d=new MHCamera(&m);
365// d->FillRandom();
366// c->cd(1);
367// gPad->SetBorderMode(0);
368// gPad->Divide(1,1);
369// gPad->cd(1);
370// d->Draw();
371// d->SetBit(kCanDelete);
372//
373void MHCamera::Draw(Option_t *option)
374{
375 // root 3.02:
376 // gPad->SetFixedAspectRatio()
377 Int_t col = 16;
378
379 if (gPad)
380 col = gPad->GetFillColor();
381
382 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas("CamDisplay", "Mars Camera Display", 656, 600);
383 pad->SetBorderMode(0);
384 pad->SetFillColor(col);
385
386 AppendPad(option);
387}
388
389// ------------------------------------------------------------------------
390//
391// Resizes the current pad so that the camera is displayed in its
392// correct aspect ratio
393//
394void MHCamera::SetRange()
395{
396 const Float_t range = fGeomCam->GetMaxRadius()*1.05;
397
398 //
399 // Maintain aspect ratio
400 //
401 const float ratio = 1.15;
402
403 //
404 // Calculate width and height of the current pad in pixels
405 //
406 Float_t w = gPad->GetWw();
407 Float_t h = gPad->GetWh()*ratio;
408
409 //
410 // This prevents the pad from resizing itself wrongly
411 //
412 if (gPad->GetMother() != gPad)
413 {
414 w *= gPad->GetMother()->GetAbsWNDC();
415 h *= gPad->GetMother()->GetAbsHNDC();
416 }
417
418 //
419 // Set Range (coordinate system) of pad
420 //
421 gPad->Range(-range, -range, (2*ratio-1)*range, range);
422
423 //
424 // Resize Pad to given ratio
425 //
426 if (h<w)
427 gPad->SetPad((1.-h/w)/2, 0, (h/w+1.)/2, 1);
428 else
429 gPad->SetPad(0, (1.-w/h)/2, 1, (w/h+1.)/2);
430}
431
432// ------------------------------------------------------------------------
433//
434// Updates the pixel colors and paints the pixels
435//
436void MHCamera::Update(Bool_t islog, Bool_t isbox, Bool_t iscol)
437{
438 Double_t min = GetMinimum(kFALSE);
439 Double_t max = GetMaximum(kFALSE);
440 if (min==FLT_MAX)
441 {
442 min = 0;
443 max = 1;
444 }
445
446 if (min==max)
447 max += 1;
448
449 UpdateLegend(min, max, islog);
450
451 MHexagon hex;
452 for (Int_t i=0; i<fNcells-2; i++)
453 {
454 if (IsUsed(i) && iscol)
455 {
456 if (TMath::IsNaN(fArray[i+1]))
457 gLog << warn << "MHCamera::Update: " << GetName() << " <" << GetTitle() << "> - Pixel Index #" << i << " contents is NaN (Not a Number)..." << endl;
458
459 hex.SetFillColor(GetColor(GetBinContent(i+1), min, max, islog));
460 }
461 else
462 hex.SetFillColor(10);
463
464 MGeomPix &pix = (*fGeomCam)[i];
465 if (!isbox)
466 hex.PaintHexagon(pix.GetX(), pix.GetY(), pix.GetD());
467 else
468 if (IsUsed(i) && !TMath::IsNaN(fArray[i+1]))
469 {
470 Float_t size = pix.GetD()*(GetBinContent(i+1)-min)/(max-min);
471 if (size>pix.GetD())
472 size=pix.GetD();
473 hex.PaintHexagon(pix.GetX(), pix.GetY(), size);
474 }
475 }
476}
477
478// ------------------------------------------------------------------------
479//
480// Print minimum and maximum
481//
482void MHCamera::Print(Option_t *) const
483{
484 cout << "Minimum: " << GetMinimum();
485 if (fMinimum==-1111)
486 cout << " <autoscaled>";
487 cout << endl;
488 cout << "Maximum: " << GetMaximum();
489 if (fMaximum==-1111)
490 cout << " <autoscaled>";
491 cout << endl;
492}
493
494// ------------------------------------------------------------------------
495//
496// Paint the y-axis title
497//
498void MHCamera::PaintAxisTitle()
499{
500 const Float_t range = fGeomCam->GetMaxRadius()*1.05;
501 const Float_t w = (1 + 1.5/sqrt((float)(fNcells-2)))*range;
502
503 TLatex *ptitle = new TLatex(w, -.90*range, GetYaxis()->GetTitle());
504
505 ptitle->SetTextSize(0.05);
506 ptitle->SetTextAlign(21);
507
508 // box with the histogram title
509 ptitle->SetTextColor(gStyle->GetTitleTextColor());
510#if ROOT_VERSION_CODE > ROOT_VERSION(3,05,01)
511 ptitle->SetTextFont(gStyle->GetTitleFont(""));
512#endif
513 ptitle->Paint();
514}
515
516// ------------------------------------------------------------------------
517//
518// Paints the camera.
519//
520void MHCamera::Paint(Option_t *o)
521{
522 if (fNcells<=1)
523 return;
524
525 TString opt(o);
526 opt.ToLower();
527
528 if (opt.Contains("hist"))
529 {
530 opt.ReplaceAll("hist", "");
531 TH1D::Paint(opt);
532 return;
533 }
534
535 // Maintain aspect ratio
536 SetRange();
537
538 Bool_t isbox = opt.Contains("box");
539 Bool_t iscol = isbox ? !opt.Contains("nocol") : 1;
540
541 GetPainter();
542 if (fPainter)
543 {
544 if (!TestBit(TH1::kNoStats))
545 {
546 fPainter->SetHistogram(this);
547 fPainter->PaintStat(gStyle->GetOptStat(), NULL);
548 }
549
550 // Paint primitives (pixels, color legend, photons, ...)
551 if (fPainter->InheritsFrom(THistPainter::Class()))
552 ((THistPainter*)fPainter)->PaintTitle();
553 }
554
555 // Update Contents of the pixels and paint legend
556 Update(gPad->GetLogy(), isbox, iscol);
557 PaintAxisTitle();
558}
559
560// ------------------------------------------------------------------------
561//
562// With this function you can change the color palette. For more
563// information see TStyle::SetPalette. Only palettes with 50 colors
564// are allowed.
565// In addition you can use SetPalette(52, 0) to create an inverse
566// deep blue sea palette
567//
568void MHCamera::SetPalette(Int_t ncolors, Int_t *colors)
569{
570 //
571 // If not enough colors are specified skip this.
572 //
573 if (ncolors>1 && ncolors<50)
574 {
575 cout << "MHCamera::SetPalette: Only default palettes with 50 colors are allowed... ignored." << endl;
576 return;
577 }
578
579 //
580 // If ncolors==52 create a reversed deep blue sea palette
581 //
582 if (ncolors==52)
583 {
584 gStyle->SetPalette(51, NULL);
585 TArrayI c(kItemsLegend);
586 for (int i=0; i<kItemsLegend; i++)
587 c[kItemsLegend-i-1] = gStyle->GetColorPalette(i);
588 gStyle->SetPalette(kItemsLegend, c.GetArray());
589 }
590 else
591 gStyle->SetPalette(ncolors, colors);
592
593 fColors.Set(kItemsLegend);
594 for (int i=0; i<kItemsLegend; i++)
595 fColors[i] = gStyle->GetColorPalette(i);
596}
597
598
599void MHCamera::SetPrettyPalette()
600{
601 if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
602 SetPalette(1, 0);
603}
604
605void MHCamera::SetDeepBlueSeaPalette()
606{
607 if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
608 SetPalette(51, 0);
609}
610
611void MHCamera::SetInvDeepBlueSeaPalette()
612{
613 if (!TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
614 SetPalette(52, 0);
615}
616
617void MHCamera::DrawPixelIndices()
618{
619 if (fNcells<=1)
620 return;
621
622 // FIXME: Is this correct?
623 for (int i=0; i<kItemsLegend; i++)
624 fColors[i] = 16;
625
626 if (!gPad)
627 Draw();
628
629 TText txt;
630 txt.SetTextFont(122);
631 txt.SetTextAlign(22); // centered/centered
632
633 for (Int_t i=0; i<fNcells-2; i++)
634 {
635 TString num;
636 num += i;
637
638 const MGeomPix &h = (*fGeomCam)[i];
639 TText *nt = txt.DrawText(h.GetX(), h.GetY(), num);
640 nt->SetTextSize(0.3*h.GetD()/fGeomCam->GetMaxRadius()/1.05);
641 }
642}
643
644void MHCamera::DrawSectorIndices()
645{
646 if (fNcells<=1)
647 return;
648
649 for (int i=0; i<kItemsLegend; i++)
650 fColors[i] = 16;
651
652 if (!gPad)
653 Draw();
654
655 TText txt;
656 txt.SetTextFont(122);
657 txt.SetTextAlign(22); // centered/centered
658
659 for (Int_t i=0; i<fNcells-2; i++)
660 {
661 TString num;
662 num += (*fGeomCam)[i].GetSector();
663
664 const MGeomPix &h = (*fGeomCam)[i];
665 TText *nt = txt.DrawText(h.GetX(), h.GetY(), num);
666 nt->SetTextSize(0.3*h.GetD()/fGeomCam->GetMaxRadius()/1.05);
667 }
668}
669
670// ------------------------------------------------------------------------
671//
672// Call this function to add a MCamEvent on top of the present contents.
673//
674void MHCamera::AddCamContent(const MCamEvent &event, Int_t type)
675{
676 if (fNcells<=1 || IsFreezed())
677 return;
678
679 // FIXME: Security check missing!
680 for (Int_t idx=0; idx<fNcells-2; idx++)
681 {
682 Double_t val=0;
683 if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
684 SetUsed(idx);
685
686 Fill(idx, val); // FIXME: Slow!
687 }
688 fEntries++;
689}
690
691// ------------------------------------------------------------------------
692//
693// Call this function to add a MHCamera on top of the present contents.
694// Type:
695// 0) bin content
696// 1) errors
697// 2) rel. errors
698//
699void MHCamera::AddCamContent(const MHCamera &d, Int_t type)
700{
701 if (fNcells!=d.fNcells || IsFreezed())
702 return;
703
704 // FIXME: Security check missing!
705 for (Int_t idx=0; idx<fNcells-2; idx++)
706 if (d.IsUsed(idx))
707 SetUsed(idx);
708
709 switch (type)
710 {
711 case 1:
712 for (Int_t idx=0; idx<fNcells-2; idx++)
713 Fill(idx, d.GetBinError(idx+1));
714 break;
715 case 2:
716 for (Int_t idx=0; idx<fNcells-2; idx++)
717 if (d.GetBinContent(idx+1)!=0)
718 Fill(idx, TMath::Abs(d.GetBinError(idx+1)/d.GetBinContent(idx+1)));
719 break;
720 default:
721 for (Int_t idx=0; idx<fNcells-2; idx++)
722 Fill(idx, d.GetBinContent(idx+1));
723 break;
724 }
725 fEntries++;
726}
727
728// ------------------------------------------------------------------------
729//
730// Call this function to add a TArrayD on top of the present contents.
731//
732void MHCamera::AddCamContent(const TArrayD &event, const TArrayC *used)
733{
734 if (event.GetSize()!=fNcells-2 || IsFreezed())
735 return;
736
737 if (used && used->GetSize()!=fNcells-2)
738 return;
739
740 for (Int_t idx=0; idx<fNcells-2; idx++)
741 {
742 Fill(idx, const_cast<TArrayD&>(event)[idx]); // FIXME: Slow!
743
744 if (!used || (*const_cast<TArrayC*>(used))[idx])
745 SetUsed(idx);
746 }
747 fEntries++;
748}
749
750// ------------------------------------------------------------------------
751//
752// Call this function to add a MCamEvent on top of the present contents.
753// 1 is added to each pixel if the contents of MCamEvent>threshold
754//
755void MHCamera::CntCamContent(const MCamEvent &event, Double_t threshold, Int_t type)
756{
757 if (fNcells<=1 || IsFreezed())
758 return;
759
760 // FIXME: Security check missing!
761 for (Int_t idx=0; idx<fNcells-2; idx++)
762 {
763 Double_t val=threshold;
764 if (event.GetPixelContent(val, idx, *fGeomCam, type)/* && !IsUsed(idx)*/)
765 SetUsed(idx);
766
767 if (val>threshold)
768 Fill(idx);
769 }
770 fEntries++;
771}
772
773// ------------------------------------------------------------------------
774//
775// Call this function to add a TArrayD on top of the present contents.
776// 1 is added to each pixel if the contents of MCamEvent>threshold
777//
778void MHCamera::CntCamContent(const TArrayD &event, Double_t threshold, Bool_t ispos)
779{
780 if (event.GetSize()!=fNcells-2 || IsFreezed())
781 return;
782
783 for (Int_t idx=0; idx<fNcells-2; idx++)
784 {
785 if (const_cast<TArrayD&>(event)[idx]>threshold)
786 Fill(idx);
787
788 if (!ispos || fArray[idx+1]>0)
789 SetUsed(idx);
790 }
791 fEntries++;
792}
793
794// ------------------------------------------------------------------------
795//
796// Fill the pixels with random contents.
797//
798void MHCamera::FillRandom()
799{
800 if (fNcells<=1 || IsFreezed())
801 return;
802
803 Reset();
804
805 // FIXME: Security check missing!
806 for (Int_t idx=0; idx<fNcells-2; idx++)
807 {
808 Fill(idx, gRandom->Uniform()*fGeomCam->GetPixRatio(idx));
809 SetUsed(idx);
810 }
811 fEntries=1;
812}
813
814
815// ------------------------------------------------------------------------
816//
817// The array must be in increasing order, eg: 2.5, 3.7, 4.9
818// The values in each bin are replaced by the interval in which the value
819// fits. In the example we have four intervals
820// (<2.5, 2.5-3.7, 3.7-4.9, >4.9). Maximum and minimum are set
821// accordingly.
822//
823void MHCamera::SetLevels(const TArrayF &arr)
824{
825 if (fNcells<=1)
826 return;
827
828 for (Int_t i=0; i<fNcells-2; i++)
829 {
830 if (!IsUsed(i))
831 continue;
832
833 Int_t j = arr.GetSize();
834 while (j && fArray[i+1]<arr[j-1])
835 j--;
836
837 fArray[i+1] = j;
838 }
839 SetMaximum(arr.GetSize());
840 SetMinimum(0);
841}
842
843// ------------------------------------------------------------------------
844//
845// Reset the all pixel colors to a default value
846//
847void MHCamera::Reset(Option_t *opt)
848{
849 if (fNcells<=1 || IsFreezed())
850 return;
851
852 TH1::Reset(opt);
853
854 for (Int_t i=0; i<fNcells-2; i++)
855 {
856 fArray[i+1]=0;
857 ResetUsed(i);
858 }
859 fArray[0] = 0;
860 fArray[fNcells-1] = 0;
861}
862
863// ------------------------------------------------------------------------
864//
865// Here we calculate the color index for the current value.
866// The color index is defined with the class TStyle and the
867// Color palette inside. We use the command gStyle->SetPalette(1,0)
868// for the display. So we have to convert the value "wert" into
869// a color index that fits the color palette.
870// The range of the color palette is defined by the values fMinPhe
871// and fMaxRange. Between this values we have 50 color index, starting
872// with 0 up to 49.
873//
874Int_t MHCamera::GetColor(Float_t val, Float_t min, Float_t max, Bool_t islog)
875{
876 if (TMath::IsNaN(val)) // FIXME: gLog!
877 return 10;
878
879 //
880 // first treat the over- and under-flows
881 //
882 const Int_t maxcolidx = kItemsLegend-1;
883
884 if (val >= max)
885 return fColors[maxcolidx];
886
887 if (val <= min)
888 return fColors[0];
889
890 //
891 // calculate the color index
892 //
893 Float_t ratio;
894 if (islog && min>0)
895 ratio = log10(val/min) / log10(max/min);
896 else
897 ratio = (val-min) / (max-min);
898
899 const Int_t colidx = (Int_t)(ratio*maxcolidx + .5);
900 return fColors[colidx];
901}
902
903TPaveStats *MHCamera::GetStatisticBox()
904{
905 TObject *obj = 0;
906
907 TIter Next(fFunctions);
908 while ((obj = Next()))
909 if (obj->InheritsFrom(TPaveStats::Class()))
910 return static_cast<TPaveStats*>(obj);
911
912 return NULL;
913}
914
915// ------------------------------------------------------------------------
916//
917// Change the text on the legend according to the range of the Display
918//
919void MHCamera::UpdateLegend(Float_t min, Float_t max, Bool_t islog)
920{
921 TPaveStats *stats = GetStatisticBox();
922
923 const Float_t hndc = 0.92 - (stats ? stats->GetY1NDC() : 1);
924 const Float_t range = fGeomCam->GetMaxRadius()*1.05;
925 const Float_t H = (0.75-hndc)*range;
926 const Float_t offset = hndc*range;
927
928 const Float_t h = 2./kItemsLegend;
929 const Float_t w = range/sqrt((float)(fNcells-2));
930
931 TBox newbox;
932 TText newtxt;
933 newtxt.SetTextSize(0.03);
934 newtxt.SetTextAlign(12);
935#if ROOT_VERSION_CODE > ROOT_VERSION(3,01,06)
936 newtxt.SetBit(/*kNoContextMenu|*/kCannotPick);
937 newbox.SetBit(/*kNoContextMenu|*/kCannotPick);
938#endif
939
940 const Float_t step = (islog && min>0 ? log10(max/min) : max-min) / kItemsLegend;
941 const Int_t firsts = step*3 < 1e-8 ? 8 : (Int_t)floor(log10(step*3));
942 const TString opt = Form("%%.%if", firsts>0 ? 0 : TMath::Abs(firsts));
943
944 for (Int_t i=0; i<kItemsLegend+1; i+=3)
945 {
946 Float_t val;
947 if (islog && min>0)
948 val = pow(10, step*i) * min;
949 else
950 val = min + step*i;
951
952 //const bool dispexp = max-min>1.5 && fabs(val)>0.1 && fabs(val)<1e6;
953 newtxt.PaintText(range+1.5*w, H*(i*h-1)-offset, Form(opt, val));
954 }
955
956 for (Int_t i=0; i<kItemsLegend; i++)
957 {
958 newbox.SetFillColor(fColors[i]);
959 newbox.PaintBox(range, H*(i*h-1)-offset, range+w, H*((i+1)*h-1)-offset);
960 }
961
962 TArrow arr;
963 arr.PaintArrow(-range*.9, -range*.9, -range*.6, -range*.9, 0.025);
964 arr.PaintArrow(-range*.9, -range*.9, -range*.9, -range*.6, 0.025);
965
966 TString text;
967 text += (int)(range*.3);
968 text += "mm";
969
970 TText newtxt2;
971 newtxt2.SetTextSize(0.04);
972 newtxt2.PaintText(-range*.85, -range*.85, text);
973
974 text = "";
975 text += (float)((int)(range*.3*fGeomCam->GetConvMm2Deg()*10))/10;
976 text += "\\circ";
977 text = text.Strip(TString::kLeading);
978
979 TLatex latex;
980 latex.PaintLatex(-range*.85, -range*.75, 0, 0.04, text);
981}
982
983// ------------------------------------------------------------------------
984//
985// Save primitive as a C++ statement(s) on output stream out
986//
987void MHCamera::SavePrimitive(ofstream &out, Option_t *opt)
988{
989 cout << "MHCamera::SavePrimitive: Must be rewritten!" << endl;
990 /*
991 if (!gROOT->ClassSaved(TCanvas::Class()))
992 fDrawingPad->SavePrimitive(out, opt);
993
994 out << " " << fDrawingPad->GetName() << "->SetWindowSize(";
995 out << fDrawingPad->GetWw() << "," << fDrawingPad->GetWh() << ");" << endl;
996 */
997}
998
999// ------------------------------------------------------------------------
1000//
1001// compute the distance of a point (px,py) to the Camera
1002// this functions needed for graphical primitives, that
1003// means without this function you are not able to interact
1004// with the graphical primitive with the mouse!!!
1005//
1006// All calcutations are done in pixel coordinates
1007//
1008Int_t MHCamera::DistancetoPrimitive(Int_t px, Int_t py)
1009{
1010 if (fNcells<=1)
1011 return 999999;
1012
1013 const Int_t kMaxDiff = 7;
1014
1015 TPaveStats *box = (TPaveStats*)gPad->GetPrimitive("stats");
1016 if (box)
1017 {
1018 const Double_t w = box->GetY2NDC()-box->GetY1NDC();
1019 box->SetX1NDC(gStyle->GetStatX()-gStyle->GetStatW());
1020 box->SetY1NDC(gStyle->GetStatY()-w);
1021 box->SetX2NDC(gStyle->GetStatX());
1022 box->SetY2NDC(gStyle->GetStatY());
1023 }
1024
1025 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1026 return TH1D::DistancetoPrimitive(px, py);
1027
1028 for (Int_t i=0; i<fNcells-2; i++)
1029 {
1030 MHexagon hex((*fGeomCam)[i]);
1031 if (hex.DistancetoPrimitive(px, py)==0)
1032 return 0;
1033 }
1034
1035 if (!box)
1036 return 999999;
1037
1038 const Int_t dist = box->DistancetoPrimitive(px, py);
1039 if (dist > kMaxDiff)
1040 return 999999;
1041
1042 gPad->SetSelected(box);
1043 return dist;
1044}
1045
1046// ------------------------------------------------------------------------
1047//
1048// Function introduced (31-01-03) WILL BE REMOVED IN THE FUTURE! DON'T
1049// USE IT!
1050//
1051void MHCamera::SetPix(const Int_t idx, const Int_t color, Float_t min, Float_t max)
1052{
1053 fArray[idx+1] = color;
1054 SetUsed(idx);
1055}
1056
1057// ------------------------------------------------------------------------
1058//
1059//
1060Int_t MHCamera::GetPixelIndex(Int_t px, Int_t py) const
1061{
1062 if (fNcells<=1)
1063 return -1;
1064
1065 Int_t i;
1066 for (i=0; i<fNcells-2; i++)
1067 {
1068 MHexagon hex((*fGeomCam)[i]);
1069 if (hex.DistancetoPrimitive(px, py)>0)
1070 continue;
1071
1072 return i;
1073 }
1074 return -1;
1075}
1076
1077// ------------------------------------------------------------------------
1078//
1079// Returns string containing info about the object at position (px,py).
1080// Returned string will be re-used (lock in MT environment).
1081//
1082char *MHCamera::GetObjectInfo(Int_t px, Int_t py) const
1083{
1084 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1085 return TH1D::GetObjectInfo(px, py);
1086
1087 static char info[128];
1088
1089 const Int_t idx=GetPixelIndex(px, py);
1090
1091 if (idx<0)
1092 return TObject::GetObjectInfo(px, py);
1093
1094 sprintf(info, "Software Pixel Index: %d (Hardware Id=%d)", idx, idx+1);
1095 return info;
1096}
1097
1098// ------------------------------------------------------------------------
1099//
1100// Execute a mouse event on the camera
1101//
1102void MHCamera::ExecuteEvent(Int_t event, Int_t px, Int_t py)
1103{
1104 if (TString(GetDrawOption()).Contains("hist", TString::kIgnoreCase))
1105 {
1106 TH1D::ExecuteEvent(event, px, py);
1107 return;
1108 }
1109 //if (event==kMouseMotion && fStatusBar)
1110 // fStatusBar->SetText(GetObjectInfo(px, py), 0);
1111 if (event!=kButton1Down)
1112 return;
1113
1114 const Int_t idx = GetPixelIndex(px, py);
1115 if (idx<0)
1116 return;
1117
1118 cout << GetTitle() << " <" << GetName() << ">" << endl;
1119 cout << "Software Pixel Index: " << idx << endl;
1120 cout << "Hardware Pixel Id: " << idx+1 << endl;
1121 cout << "Contents: " << GetBinContent(idx+1) << " <";
1122 cout << (IsUsed(idx)?"on":"off");
1123 cout << ">" << endl;
1124
1125 if (fNotify && fNotify->GetSize()>0)
1126 {
1127 // FIXME: Is there a simpler and more convinient way?
1128
1129 // The name which is created here depends on the instance of
1130 // MHCamera and on the pad on which it is drawn --> The name
1131 // is unique. For ExecuteEvent gPad is always correctly set.
1132 const TString name = Form("%p;%p;PixelContent", this, gPad);
1133
1134 TCanvas *old = (TCanvas*)gROOT->GetListOfCanvases()->FindObject(name);
1135 if (old)
1136 old->cd();
1137 else
1138 new TCanvas(name);
1139
1140 /*
1141 TIter Next(gPad->GetListOfPrimitives());
1142 TObject *o;
1143 while (o=Next()) cout << o->GetName() << " " << o->IsA()->GetName() << endl;
1144 */
1145
1146 // FIXME: Make sure, that the old histograms are really deleted.
1147 // Are they already deleted?
1148 fNotify->ForEach(MCamEvent, DrawPixelContent)(idx);
1149 gPad->Modified();
1150 gPad->Update();
1151 }
1152}
1153
1154UInt_t MHCamera::GetNumPixels() const
1155{
1156 return fGeomCam->GetNumPixels();
1157}
1158
1159TH1 *MHCamera::DrawCopy() const
1160{
1161 gPad=NULL;
1162 return TH1D::DrawCopy(fName+";cpy");
1163}
Note: See TracBrowser for help on using the repository browser.