source: trunk/MagicSoft/Mars/mcalib/MCalibrationIntensityCam.cc@ 7109

Last change on this file since 7109 was 7005, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 13.2 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): Markus Gaug 11/2003 <mailto:markus@ifae.es>
19!
20! Copyright: MAGIC Software Development, 2000-2004
21!
22!
23\* ======================================================================== */
24/////////////////////////////////////////////////////////////////////////////
25//
26// MCalibrationIntensityCam
27//
28// Base class for intensity calibration results
29//
30// Contains TOrdCollections for the following objects:
31// - fCams: Array of classes derived from MCalibrationCam, one entry
32// per calibration camera result.
33// - fHists: Array of classes derived from MHCalibrationPix, one entry
34// per calibration camera result and area index
35//
36// See also: MCalibrationIntensityChargeCam, MCalibrationIntensityQECam,
37// MCalibrationIntensityRelTimeCam,
38// MCalibrationCam, MCalibrationPix,
39// MCalibrationQECam, MCalibrationQEPix,
40// MHCalibrationChargePix, MHCalibrationChargeCam
41// MCalibrationChargeBlindPix, MCalibrationChargePINDiode
42//
43// ClassVersion 2:
44// + fHists
45//
46/////////////////////////////////////////////////////////////////////////////
47#include "MCalibrationIntensityCam.h"
48
49#include <TOrdCollection.h>
50
51#include "MLog.h"
52#include "MLogManip.h"
53
54#include "MGeomCam.h"
55#include "MHCalibrationCam.h"
56
57ClassImp(MCalibrationIntensityCam);
58
59using namespace std;
60
61// --------------------------------------------------------------------------
62//
63// Default constructor.
64//
65// Set the following pointer to NULL:
66// - fCams
67// - fHists
68//
69MCalibrationIntensityCam::MCalibrationIntensityCam(const char *name, const char *title)
70{
71
72 fName = name ? name : "MCalibrationIntensityCam";
73 fTitle = title ? title : "Base container for the Intensity Calibration";
74
75 fCams = new TOrdCollection;
76 fCams->SetOwner();
77
78 fHists = new TOrdCollection;
79 fHists->SetOwner();
80}
81
82
83// --------------------------------------------------------------------------
84//
85// Deletes the histograms if they exist
86//
87MCalibrationIntensityCam::~MCalibrationIntensityCam()
88{
89 delete fCams;
90 delete fHists;
91}
92
93// --------------------------------------------------------------------------
94//
95// Add a new MHCalibrationCam to fHists
96//
97void MCalibrationIntensityCam::AddHist( const MHCalibrationCam *cam)
98{
99 const Int_t size = fHists->GetSize();
100 fHists->AddAt((TObject*)cam,size);
101
102 if (size != GetSize()-1)
103 *fLog << warn << "Histogram Cams and Calibration Cams size mismatch! " << endl;
104}
105
106// --------------------------------------------------------------------------
107//
108// Add a new MCalibrationCam to fCams, give it the name "name" and initialize
109// it with geom.
110//
111void MCalibrationIntensityCam::AddToList( const char* name, const MGeomCam &geom)
112{
113 InitSize(GetSize()+1);
114 GetCam()->SetName(name);
115 GetCam()->Init(geom);
116}
117
118
119
120// --------------------------------------------------------------------------
121//
122// Copy 'constructor'
123//
124void MCalibrationIntensityCam::Copy(TObject& object) const
125{
126 MCalibrationIntensityCam &calib = (MCalibrationIntensityCam&)object;
127
128 MParContainer::Copy(calib);
129
130 calib.fOffsets = fOffsets;
131 calib.fSlopes = fSlopes;
132
133 const UInt_t n = GetSize();
134 if (n != 0)
135 {
136 calib.InitSize(n);
137 for (UInt_t i=0; i<n; i++)
138 {
139 GetCam(i)->Copy(*(calib.GetCam(i)));
140 GetHist(i)->Copy(*(calib.GetHist(i)));
141 }
142 }
143}
144
145// -----------------------------------------------------
146//
147// Calls Clear() for all entries fCams
148//
149void MCalibrationIntensityCam::Clear(Option_t *o)
150{
151 fCams->ForEach(MCalibrationCam, Clear)();
152 fHists->ForEach(MHCalibrationCam, Clear)();
153}
154
155// -----------------------------------------------------
156//
157// Calls Print(o) for all entries fCams
158//
159void MCalibrationIntensityCam::Print(Option_t *o) const
160{
161 fCams->ForEach(MCalibrationCam, Print)(o);
162 fHists->ForEach(MHCalibrationCam, Print)(o);
163}
164
165// -----------------------------------------------------
166//
167// Not yet installed...
168//
169void MCalibrationIntensityCam::DrawHiLoFits()
170{
171
172 /*
173 if (!fOffsets)
174 fOffsets = new TH1D("pp","Offsets of the HiGain LoGain Fit",100,-600.,400.);
175 if (!fSlopes)
176 fSlopes = new TH1D("mm","Slopes of the HiGain LoGain Fit",100,-2.,2.);
177 if (!fOffvsSlope)
178 fOffvsSlope = new TH2D("aa","Slopes vs Offsets of the HiGain LoGain Fit",100,-600.,400.,100,-2.,2.);
179
180 TIter Next(fPixels);
181 MCalibrationPix *pix;
182 MHCalibrationPixel *hist;
183 while ((pix=(MCalibrationPix*)Next()))
184 {
185 hist = pix->GetHist();
186 hist->FitHiGainvsLoGain();
187 fOffsets->Fill(hist->GetOffset(),1.);
188 fSlopes->Fill(hist->GetSlope(),1.);
189 fOffvsSlope->Fill(hist->GetOffset(),hist->GetSlope(),1.);
190 }
191
192 TCanvas *c1 = new TCanvas();
193
194 c1->Divide(1,3);
195 c1->cd(1);
196 fOffsets->Draw();
197 gPad->Modified();
198 gPad->Update();
199
200 c1->cd(2);
201 fSlopes->Draw();
202 gPad->Modified();
203 gPad->Update();
204
205 c1->cd(3);
206 fOffvsSlope->Draw("col1");
207 gPad->Modified();
208 gPad->Update();
209 */
210}
211
212// -------------------------------------------------------------------
213//
214// Initialize the objects inside the TOrdCollection using the
215// virtual function Add().
216//
217// InitSize can only increase the size, but not shrink.
218//
219// It can be called more than one time. New Containers are
220// added only from the current size to the argument i.
221//
222void MCalibrationIntensityCam::InitSize(const UInt_t i)
223{
224
225 const UInt_t save = GetSize();
226
227 if (i==save)
228 return;
229
230 if (i>save)
231 Add(save,i);
232}
233
234// -------------------------------------------------------------------
235//
236// Add MCalibrationCams in the ranges from - to. In order to initialize
237// from MCalibrationCam derived containers, overwrite this function
238//
239void MCalibrationIntensityCam::Add(const UInt_t from, const UInt_t to)
240{
241 for (UInt_t i=from; i<to; i++)
242 fCams->AddAt(new MCalibrationCam,i);
243}
244
245// -------------------------------------------------------------------
246//
247// If size is still 0, Intialize a first Cam.
248// Calls Init(geom) for all fCams
249//
250void MCalibrationIntensityCam::Init(const MGeomCam &geom)
251{
252 if (GetSize() == 0)
253 InitSize(1);
254
255 fCams->ForEach(MCalibrationCam,Init)(geom);
256}
257
258
259// --------------------------------------------------------------------------
260//
261// Returns the current size of the TOrdCollection fCams
262// independently if the MCalibrationCam is filled with values or not.
263//
264const Int_t MCalibrationIntensityCam::GetSize() const
265{
266 return fCams->GetSize();
267}
268
269// --------------------------------------------------------------------------
270//
271// Get i-th pixel from current camera
272//
273MCalibrationPix &MCalibrationIntensityCam::operator[](UInt_t i)
274{
275 return (*GetCam())[i];
276}
277
278// --------------------------------------------------------------------------
279//
280// Get i-th pixel from current camera
281//
282const MCalibrationPix &MCalibrationIntensityCam::operator[](UInt_t i) const
283{
284 return (*GetCam())[i];
285}
286
287
288// --------------------------------------------------------------------------
289//
290// Returns the current size of the TOrdCollection fAverageAreas of the current camera.
291//
292const Int_t MCalibrationIntensityCam::GetAverageAreas() const
293{
294 return GetCam()->GetAverageAreas();
295}
296
297// --------------------------------------------------------------------------
298//
299// Get i-th High Gain pixel Area from the current camera
300//
301MCalibrationPix &MCalibrationIntensityCam::GetAverageArea(UInt_t i)
302{
303 return GetCam()->GetAverageArea(i);
304}
305
306// --------------------------------------------------------------------------
307//
308// Get i-th High Gain pixel Area from the current camera
309//
310const MCalibrationPix &MCalibrationIntensityCam::GetAverageArea(UInt_t i) const
311{
312 return GetCam()->GetAverageArea(i);
313}
314
315// --------------------------------------------------------------------------
316//
317// Get i-th High Gain pixel Area from the current camera
318//
319MBadPixelsPix &MCalibrationIntensityCam::GetAverageBadArea(UInt_t i)
320{
321 return GetCam()->GetAverageBadArea(i);
322}
323
324// --------------------------------------------------------------------------
325//
326// Get i-th High Gain pixel Area from the current camera
327//
328const MBadPixelsPix &MCalibrationIntensityCam::GetAverageBadArea(UInt_t i) const
329{
330 return GetCam()->GetAverageBadArea(i);
331}
332
333// --------------------------------------------------------------------------
334//
335// Returns the current size of the TOrdCollection fAverageSectors or the current camera
336//
337const Int_t MCalibrationIntensityCam::GetAverageSectors() const
338{
339 return GetCam()->GetAverageSectors();
340}
341
342// --------------------------------------------------------------------------
343//
344// Get i-th High Gain Sector from the current camera
345//
346MCalibrationPix &MCalibrationIntensityCam::GetAverageSector(UInt_t i)
347{
348 return GetCam()->GetAverageSector(i);
349}
350
351// --------------------------------------------------------------------------
352//
353// Get i-th High Gain Sector from the current camera
354//
355const MCalibrationPix &MCalibrationIntensityCam::GetAverageSector(UInt_t i) const
356{
357 return GetCam()->GetAverageSector(i);
358}
359
360// --------------------------------------------------------------------------
361//
362// Get i-th High Gain Sector from the current camera
363//
364MBadPixelsPix &MCalibrationIntensityCam::GetAverageBadSector(UInt_t i)
365{
366 return GetCam()->GetAverageBadSector(i);
367}
368
369// --------------------------------------------------------------------------
370//
371// Get i-th High Gain Sector from the current camera
372//
373const MBadPixelsPix &MCalibrationIntensityCam::GetAverageBadSector(UInt_t i) const
374{
375 return GetCam()->GetAverageBadSector(i);
376}
377
378
379// --------------------------------------------------------------------------
380//
381// Get i-th camera
382//
383MCalibrationCam *MCalibrationIntensityCam::GetCam(Int_t i)
384{
385 return static_cast<MCalibrationCam*>(i==-1 ? fCams->Last() : fCams->At(i));
386}
387
388// --------------------------------------------------------------------------
389//
390// Get i-th camera
391//
392const MCalibrationCam *MCalibrationIntensityCam::GetCam(Int_t i) const
393{
394 return static_cast<MCalibrationCam*>(i==-1 ? fCams->Last() : fCams->At(i));
395}
396
397// --------------------------------------------------------------------------
398//
399// Get camera with name 'name'
400//
401MCalibrationCam *MCalibrationIntensityCam::GetCam(const char *name )
402{
403 return static_cast<MCalibrationCam*>(fCams->FindObject(name));
404}
405
406// --------------------------------------------------------------------------
407//
408// Get camera with name 'name'
409//
410const MCalibrationCam *MCalibrationIntensityCam::GetCam(const char *name ) const
411{
412 return static_cast<MCalibrationCam*>(fCams->FindObject(name));
413}
414
415// --------------------------------------------------------------------------
416//
417// Calls GetPixelContent for the current entry in fCams
418//
419Bool_t MCalibrationIntensityCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
420{
421 return GetCam()->GetPixelContent(val,idx,cam,type);
422}
423
424// --------------------------------------------------------------------------
425//
426// Calls DrawPixelContent for the current entry in fCams
427//
428void MCalibrationIntensityCam::DrawPixelContent( Int_t num ) const
429{
430 return GetCam()->DrawPixelContent(num);
431}
432
433Int_t MCalibrationIntensityCam::CountNumEntries(const MCalibrationCam::PulserColor_t col) const
434{
435
436 Int_t size = 0;
437
438 if (col == MCalibrationCam::kNONE)
439 return GetSize();
440 else
441 for (Int_t i=0;i<GetSize();i++)
442 {
443 const MCalibrationCam *cam = GetCam(i);
444 if (cam->GetPulserColor() == col)
445 size++;
446 }
447
448 return size;
449}
450
451// --------------------------------------------------------------------------
452//
453// Get i-th histogram class
454//
455MHCalibrationCam *MCalibrationIntensityCam::GetHist(Int_t i)
456{
457 return static_cast<MHCalibrationCam*>(i==-1 ? fHists->Last() : fHists->At(i));
458}
459
460// --------------------------------------------------------------------------
461//
462// Get i-th histogram class
463//
464const MHCalibrationCam *MCalibrationIntensityCam::GetHist(Int_t i) const
465{
466 return static_cast<MHCalibrationCam*>(i==-1 ? fHists->Last() : fHists->At(i));
467}
468
469// --------------------------------------------------------------------------
470//
471// Get histogram class with name 'name'
472//
473MHCalibrationCam *MCalibrationIntensityCam::GetHist(const char *name )
474{
475 return static_cast<MHCalibrationCam*>(fHists->FindObject(name));
476}
477
478// --------------------------------------------------------------------------
479//
480// Get histogram class with name 'name'
481//
482const MHCalibrationCam *MCalibrationIntensityCam::GetHist(const char *name ) const
483{
484 return static_cast<MHCalibrationCam*>(fHists->FindObject(name));
485}
Note: See TracBrowser for help on using the repository browser.