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

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