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

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