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 02/2004 <mailto:markus@ifae.es>
|
---|
19 | !
|
---|
20 | ! Copyright: MAGIC Software Development, 2000-2004
|
---|
21 | !
|
---|
22 | !
|
---|
23 | \* ======================================================================== */
|
---|
24 | /////////////////////////////////////////////////////////////////////////////
|
---|
25 | //
|
---|
26 | // MHCalibrationCam
|
---|
27 | //
|
---|
28 | // Base class for camera calibration classes. Incorporates the TOrdCollection's:
|
---|
29 | // - fHiGainArray (for calibrated High Gains per pixel)
|
---|
30 | // - fLoGainArray (for calibrated Low Gains per pixel)
|
---|
31 | // - fAverageHiGainAreas (for averaged High Gains events per camera area index)
|
---|
32 | // - fAverageLoGainAreas (for averaged High Gains events per camera area index)
|
---|
33 | // - fAverageHiGainSectors (for averaged High Gains events per camera sector )
|
---|
34 | // - fAverageLoGainSectors (for averaged High Gains events per camera sector )
|
---|
35 | // These TOrdCollection's are called by their default constructors, thus no objects
|
---|
36 | // are created, until the derived class does so.
|
---|
37 | //
|
---|
38 | // The corresponding operators: [],() and the operators GetAverageHiGainArea(),
|
---|
39 | // GetAverageLoGainArea(), GetAverageHiGainSector() and GetAverageLoGainSector()
|
---|
40 | // have to be cast to the corresponding class. It is assumed that all classes
|
---|
41 | // dealing with calibration pixels derive from MHCalibrationPix.
|
---|
42 | //
|
---|
43 | // The following flag can be set:
|
---|
44 | // - SetAverageing() for calculating the event-by-event averages.
|
---|
45 | // - SetDebug() for debug output.
|
---|
46 | // - SetLoGain() for the case that low-gain slices are available, but
|
---|
47 | // MRawRunHeader::GetNumLoGainSlices() gives still 0.
|
---|
48 | // - SetCheckSize() for testing the sizes of the f*Arrays in ReInit() if they
|
---|
49 | // are coinciding with the equiv. sizes in MGeomCam
|
---|
50 | //
|
---|
51 | // The flag kLoGain steers if the low-gain signal is treated at all or not.
|
---|
52 | // The flag kAverageing steers if the event-by-event averages are treated at all.
|
---|
53 | //
|
---|
54 | // Class Version 5:
|
---|
55 | // + Double_t fLowerFitLimitHiGain; // Lower limit for the fit range for the hi-gain hist
|
---|
56 | // + Double_t fUpperFitLimitHiGain; // Upper limit for the fit range for the hi-gain hist
|
---|
57 | // + Double_t fLowerFitLimitLoGain; // Lower limit for the fit range for the lo-gain hist
|
---|
58 | // + Double_t fUpperFitLimitLoGain; // Upper limit for the fit range for the lo-gain hist
|
---|
59 | // + Bool_t fIsHiGainFitRanges; // Are high-gain fit ranges defined?
|
---|
60 | // + Bool_t fIsLoGainFitRanges; // Are low-gain fit ranges defined?
|
---|
61 | //
|
---|
62 | //
|
---|
63 | /////////////////////////////////////////////////////////////////////////////
|
---|
64 | #include "MHCalibrationCam.h"
|
---|
65 | #include "MHCalibrationPix.h"
|
---|
66 |
|
---|
67 | #include <TVirtualPad.h>
|
---|
68 | #include <TCanvas.h>
|
---|
69 | #include <TPad.h>
|
---|
70 | #include <TText.h>
|
---|
71 | #include <TPaveText.h>
|
---|
72 | #include <TOrdCollection.h>
|
---|
73 | #include <TROOT.h>
|
---|
74 |
|
---|
75 | #include "MLog.h"
|
---|
76 | #include "MLogManip.h"
|
---|
77 |
|
---|
78 | #include "MCalibrationIntensityCam.h"
|
---|
79 | #include "MCalibrationCam.h"
|
---|
80 | #include "MCalibrationPix.h"
|
---|
81 |
|
---|
82 | #include "MBadPixelsIntensityCam.h"
|
---|
83 | #include "MBadPixelsCam.h"
|
---|
84 | #include "MBadPixelsPix.h"
|
---|
85 |
|
---|
86 | #include "MGeomCam.h"
|
---|
87 | #include "MGeomPix.h"
|
---|
88 |
|
---|
89 | #include "MParList.h"
|
---|
90 |
|
---|
91 | #include "MRawRunHeader.h"
|
---|
92 |
|
---|
93 | ClassImp(MHCalibrationCam);
|
---|
94 |
|
---|
95 | using namespace std;
|
---|
96 |
|
---|
97 | const Double_t MHCalibrationCam::fgLowerFitLimitHiGain = 0;
|
---|
98 | const Double_t MHCalibrationCam::fgUpperFitLimitHiGain = 0;
|
---|
99 | const Double_t MHCalibrationCam::fgLowerFitLimitLoGain = 0;
|
---|
100 | const Double_t MHCalibrationCam::fgUpperFitLimitLoGain = 0;
|
---|
101 |
|
---|
102 | const Int_t MHCalibrationCam::fgPulserFrequency = 500;
|
---|
103 | const Float_t MHCalibrationCam::fgProbLimit = 0.0001;
|
---|
104 | const Float_t MHCalibrationCam::fgOverflowLimit = 0.005;
|
---|
105 |
|
---|
106 | const TString MHCalibrationCam::gsHistName = "Hist";
|
---|
107 | const TString MHCalibrationCam::gsHistTitle = "";
|
---|
108 | const TString MHCalibrationCam::gsHistXTitle = "";
|
---|
109 | const TString MHCalibrationCam::gsHistYTitle = "Nr. events";
|
---|
110 |
|
---|
111 | // --------------------------------------------------------------------------
|
---|
112 | //
|
---|
113 | // Default Constructor.
|
---|
114 | //
|
---|
115 | // Sets:
|
---|
116 | // - all pointers to NULL
|
---|
117 | //
|
---|
118 | // Initializes and sets owner of:
|
---|
119 | // - fHiGainArray, fLoGainArray
|
---|
120 | // - fAverageHiGainAreas, fAverageLoGainAreas
|
---|
121 | // - fAverageHiGainSectors, fAverageLoGainSectors
|
---|
122 | //
|
---|
123 | // Initializes:
|
---|
124 | // - fPulserFrequency to fgPulserFrequency
|
---|
125 | // - fProbLimit to fgProbLimit
|
---|
126 | // - fOverflowLimit to fgOverflowLimit
|
---|
127 | //
|
---|
128 | // - SetAveregeing (kTRUE);
|
---|
129 | // - SetDebug (kFALSE);
|
---|
130 | // - SetLoGain (kTRUE);
|
---|
131 | //- SetOscillations(kTRUE);
|
---|
132 | //- SetSizeCheck (kTRUE);
|
---|
133 | //- SetInterlaced (kFALSE);
|
---|
134 | //- SetLowerFitLimitHiGain();
|
---|
135 | //- SetUpperFitLimitHiGain();
|
---|
136 | //- SetLowerFitLimitLoGain();
|
---|
137 | //- SetUpperFitLimitLoGain();
|
---|
138 | //
|
---|
139 | MHCalibrationCam::MHCalibrationCam(const char *name, const char *title)
|
---|
140 | : fIsHiGainFitRanges(kFALSE), fIsLoGainFitRanges(kFALSE),
|
---|
141 | fHistName(gsHistName),fHistTitle(gsHistTitle),
|
---|
142 | fHistXTitle(gsHistXTitle),fHistYTitle(gsHistYTitle),
|
---|
143 | fColor(MCalibrationCam::kNONE), fIntensBad(NULL),
|
---|
144 | fBadPixels(NULL), fIntensCam(NULL), fCam(NULL), fGeom(NULL),
|
---|
145 | fRunHeader(NULL)
|
---|
146 | {
|
---|
147 |
|
---|
148 | fHiGainArray = new TOrdCollection;
|
---|
149 | fHiGainArray->SetOwner();
|
---|
150 |
|
---|
151 | fLoGainArray = new TOrdCollection;
|
---|
152 | fLoGainArray->SetOwner();
|
---|
153 |
|
---|
154 | fAverageHiGainAreas = new TOrdCollection;
|
---|
155 | fAverageHiGainAreas->SetOwner();
|
---|
156 |
|
---|
157 | fAverageLoGainAreas = new TOrdCollection;
|
---|
158 | fAverageLoGainAreas->SetOwner();
|
---|
159 |
|
---|
160 | fAverageHiGainSectors = new TOrdCollection;
|
---|
161 | fAverageHiGainSectors->SetOwner();
|
---|
162 |
|
---|
163 | fAverageLoGainSectors = new TOrdCollection;
|
---|
164 | fAverageLoGainSectors->SetOwner();
|
---|
165 |
|
---|
166 | SetPulserFrequency();
|
---|
167 | SetProbLimit();
|
---|
168 | SetOverflowLimit();
|
---|
169 |
|
---|
170 | SetAverageing (kTRUE);
|
---|
171 | SetDebug (kFALSE);
|
---|
172 | SetLoGain (kTRUE);
|
---|
173 | SetOscillations(kTRUE);
|
---|
174 | SetSizeCheck (kTRUE);
|
---|
175 | SetIsReset (kTRUE);
|
---|
176 |
|
---|
177 | SetLowerFitLimitHiGain();
|
---|
178 | SetUpperFitLimitHiGain();
|
---|
179 | SetLowerFitLimitLoGain();
|
---|
180 | SetUpperFitLimitLoGain();
|
---|
181 | }
|
---|
182 |
|
---|
183 | // --------------------------------------------------------------------------
|
---|
184 | //
|
---|
185 | // Deletes the TOrdCollection of:
|
---|
186 | // - fHiGainArray, fLoGainArray
|
---|
187 | // - fAverageHiGainAreas, fAverageLoGainAreas
|
---|
188 | // - fAverageHiGainSectors, fAverageLoGainSectors
|
---|
189 | //
|
---|
190 | MHCalibrationCam::~MHCalibrationCam()
|
---|
191 | {
|
---|
192 |
|
---|
193 | delete fHiGainArray;
|
---|
194 | delete fLoGainArray;
|
---|
195 |
|
---|
196 | delete fAverageHiGainAreas;
|
---|
197 | delete fAverageLoGainAreas;
|
---|
198 |
|
---|
199 | delete fAverageHiGainSectors;
|
---|
200 | delete fAverageLoGainSectors;
|
---|
201 | /*
|
---|
202 |
|
---|
203 | Remove ( fHiGainArray );
|
---|
204 | Remove ( fLoGainArray );
|
---|
205 |
|
---|
206 | Remove ( fAverageHiGainAreas );
|
---|
207 | Remove ( fAverageLoGainAreas );
|
---|
208 |
|
---|
209 | Remove ( fAverageHiGainSectors );
|
---|
210 | Remove ( fAverageLoGainSectors );
|
---|
211 | */
|
---|
212 |
|
---|
213 | }
|
---|
214 |
|
---|
215 | void MHCalibrationCam::Remove(TOrdCollection *col)
|
---|
216 | {
|
---|
217 |
|
---|
218 | if (!col)
|
---|
219 | return;
|
---|
220 |
|
---|
221 | TOrdCollectionIter Next(col);
|
---|
222 |
|
---|
223 | Int_t count = 0;
|
---|
224 |
|
---|
225 | while (MHCalibrationPix *obj = (MHCalibrationPix*)Next())
|
---|
226 | {
|
---|
227 | *fLog << ++count << " " << obj << flush;
|
---|
228 | if (obj && obj->IsOnHeap())
|
---|
229 | {
|
---|
230 | obj->Draw();
|
---|
231 | delete obj;
|
---|
232 | }
|
---|
233 | }
|
---|
234 |
|
---|
235 | delete col;
|
---|
236 | }
|
---|
237 |
|
---|
238 | // --------------------------------------------------------------------------
|
---|
239 | //
|
---|
240 | // Returns size of fHiGainArray
|
---|
241 | //
|
---|
242 | const Int_t MHCalibrationCam::GetSize() const
|
---|
243 | {
|
---|
244 | return fHiGainArray->GetSize();
|
---|
245 | }
|
---|
246 |
|
---|
247 | // --------------------------------------------------------------------------
|
---|
248 | //
|
---|
249 | // Get i-th High Gain pixel (pixel number)
|
---|
250 | //
|
---|
251 | MHCalibrationPix &MHCalibrationCam::operator[](UInt_t i)
|
---|
252 | {
|
---|
253 | return *static_cast<MHCalibrationPix*>(fHiGainArray->At(i));
|
---|
254 | }
|
---|
255 |
|
---|
256 | // --------------------------------------------------------------------------
|
---|
257 | //
|
---|
258 | // Get i-th High Gain pixel (pixel number)
|
---|
259 | //
|
---|
260 | const MHCalibrationPix &MHCalibrationCam::operator[](UInt_t i) const
|
---|
261 | {
|
---|
262 | return *static_cast<MHCalibrationPix*>(fHiGainArray->At(i));
|
---|
263 | }
|
---|
264 |
|
---|
265 | // --------------------------------------------------------------------------
|
---|
266 | //
|
---|
267 | // Get i-th Low Gain pixel (pixel number)
|
---|
268 | //
|
---|
269 | MHCalibrationPix &MHCalibrationCam::operator()(UInt_t i)
|
---|
270 | {
|
---|
271 | return *static_cast<MHCalibrationPix*>(fLoGainArray->At(i));
|
---|
272 | }
|
---|
273 |
|
---|
274 | // --------------------------------------------------------------------------
|
---|
275 | //
|
---|
276 | // Get i-th Low Gain pixel (pixel number)
|
---|
277 | //
|
---|
278 | const MHCalibrationPix &MHCalibrationCam::operator()(UInt_t i) const
|
---|
279 | {
|
---|
280 | return *static_cast<MHCalibrationPix*>(fLoGainArray->At(i));
|
---|
281 | }
|
---|
282 |
|
---|
283 | // --------------------------------------------------------------------------
|
---|
284 | //
|
---|
285 | // Returns the current size of the TOrdCollection fAverageHiGainAreas
|
---|
286 | // independently if the MHCalibrationPix is filled with values or not.
|
---|
287 | //
|
---|
288 | const Int_t MHCalibrationCam::GetAverageAreas() const
|
---|
289 | {
|
---|
290 | return fAverageHiGainAreas->GetSize();
|
---|
291 | }
|
---|
292 |
|
---|
293 | // --------------------------------------------------------------------------
|
---|
294 | //
|
---|
295 | // Get i-th High Gain pixel Area (area number)
|
---|
296 | //
|
---|
297 | MHCalibrationPix &MHCalibrationCam::GetAverageHiGainArea(UInt_t i)
|
---|
298 | {
|
---|
299 | return *static_cast<MHCalibrationPix*>(fAverageHiGainAreas->At(i));
|
---|
300 | }
|
---|
301 |
|
---|
302 | // --------------------------------------------------------------------------
|
---|
303 | //
|
---|
304 | // Get i-th High Gain pixel Area (area number)
|
---|
305 | //
|
---|
306 | const MHCalibrationPix &MHCalibrationCam::GetAverageHiGainArea(UInt_t i) const
|
---|
307 | {
|
---|
308 | return *static_cast<MHCalibrationPix *>(fAverageHiGainAreas->At(i));
|
---|
309 | }
|
---|
310 |
|
---|
311 | // --------------------------------------------------------------------------
|
---|
312 | //
|
---|
313 | // Get i-th Low Gain pixel Area (area number)
|
---|
314 | //
|
---|
315 | MHCalibrationPix &MHCalibrationCam::GetAverageLoGainArea(UInt_t i)
|
---|
316 | {
|
---|
317 | return *static_cast<MHCalibrationPix*>(fAverageLoGainAreas->At(i));
|
---|
318 | }
|
---|
319 |
|
---|
320 | // --------------------------------------------------------------------------
|
---|
321 | //
|
---|
322 | // Get i-th Low Gain pixel Area (area number)
|
---|
323 | //
|
---|
324 | const MHCalibrationPix &MHCalibrationCam::GetAverageLoGainArea(UInt_t i) const
|
---|
325 | {
|
---|
326 | return *static_cast<MHCalibrationPix*>(fAverageLoGainAreas->At(i));
|
---|
327 | }
|
---|
328 |
|
---|
329 | // --------------------------------------------------------------------------
|
---|
330 | //
|
---|
331 | // Returns the current size of the TOrdCollection fAverageHiGainSectors
|
---|
332 | // independently if the MHCalibrationPix is filled with values or not.
|
---|
333 | //
|
---|
334 | const Int_t MHCalibrationCam::GetAverageSectors() const
|
---|
335 | {
|
---|
336 | return fAverageHiGainSectors->GetSize();
|
---|
337 | }
|
---|
338 |
|
---|
339 | // --------------------------------------------------------------------------
|
---|
340 | //
|
---|
341 | // Get i-th High Gain Sector (sector number)
|
---|
342 | //
|
---|
343 | MHCalibrationPix &MHCalibrationCam::GetAverageHiGainSector(UInt_t i)
|
---|
344 | {
|
---|
345 | return *static_cast<MHCalibrationPix*>(fAverageHiGainSectors->At(i));
|
---|
346 | }
|
---|
347 |
|
---|
348 | // --------------------------------------------------------------------------
|
---|
349 | //
|
---|
350 | // Get i-th High Gain Sector (sector number)
|
---|
351 | //
|
---|
352 | const MHCalibrationPix &MHCalibrationCam::GetAverageHiGainSector(UInt_t i) const
|
---|
353 | {
|
---|
354 | return *static_cast<MHCalibrationPix*>(fAverageHiGainSectors->At(i));
|
---|
355 | }
|
---|
356 |
|
---|
357 | // --------------------------------------------------------------------------
|
---|
358 | //
|
---|
359 | // Get i-th Low Gain Sector (sector number)
|
---|
360 | //
|
---|
361 | MHCalibrationPix &MHCalibrationCam::GetAverageLoGainSector(UInt_t i)
|
---|
362 | {
|
---|
363 | return *static_cast<MHCalibrationPix*>(fAverageLoGainSectors->At(i));
|
---|
364 | }
|
---|
365 |
|
---|
366 | // --------------------------------------------------------------------------
|
---|
367 | //
|
---|
368 | // Get i-th Low Gain Sector (sector number)
|
---|
369 | //
|
---|
370 | const MHCalibrationPix &MHCalibrationCam::GetAverageLoGainSector(UInt_t i) const
|
---|
371 | {
|
---|
372 | return *static_cast<MHCalibrationPix*>(fAverageLoGainSectors->At(i));
|
---|
373 | }
|
---|
374 |
|
---|
375 | // --------------------------------------------------------------------------
|
---|
376 | //
|
---|
377 | // Calls ResetHistTitles()
|
---|
378 | //
|
---|
379 | // Calls Reset() for each entry in:
|
---|
380 | // - fHiGainArray, fLoGainArray
|
---|
381 | // - fAverageHiGainAreas, fAverageLoGainAreas
|
---|
382 | // - fAverageHiGainSectors, fAverageLoGainSectors
|
---|
383 | //
|
---|
384 | void MHCalibrationCam::ResetHists()
|
---|
385 | {
|
---|
386 | SetIsReset();
|
---|
387 |
|
---|
388 | ResetHistTitles();
|
---|
389 |
|
---|
390 | if (fHiGainArray)
|
---|
391 | { fHiGainArray->ForEach(MHCalibrationPix,Reset)(); }
|
---|
392 |
|
---|
393 | if (IsAverageing())
|
---|
394 | {
|
---|
395 | if (fAverageHiGainAreas)
|
---|
396 | { fAverageHiGainAreas->ForEach(MHCalibrationPix,Reset)(); }
|
---|
397 | if (fAverageHiGainSectors)
|
---|
398 | { fAverageHiGainSectors->ForEach(MHCalibrationPix,Reset)(); }
|
---|
399 | }
|
---|
400 |
|
---|
401 | if (!IsLoGain())
|
---|
402 | return;
|
---|
403 |
|
---|
404 | if (fLoGainArray)
|
---|
405 | { fLoGainArray->ForEach(MHCalibrationPix,Reset)(); }
|
---|
406 | if (IsAverageing())
|
---|
407 | {
|
---|
408 | if (fAverageLoGainAreas)
|
---|
409 | { fAverageLoGainAreas->ForEach(MHCalibrationPix,Reset)(); }
|
---|
410 | if (fAverageLoGainSectors)
|
---|
411 | { fAverageLoGainSectors->ForEach(MHCalibrationPix,Reset)(); }
|
---|
412 | }
|
---|
413 | }
|
---|
414 |
|
---|
415 | // --------------------------------------------------------------------------
|
---|
416 | //
|
---|
417 | // Resets the histogram titles for each entry in:
|
---|
418 | // - fHiGainArray, fLoGainArray
|
---|
419 | // - fAverageHiGainAreas, fAverageLoGainAreas
|
---|
420 | // - fAverageHiGainSectors, fAverageLoGainSectors
|
---|
421 | //
|
---|
422 | void MHCalibrationCam::ResetHistTitles()
|
---|
423 | {
|
---|
424 |
|
---|
425 | TH1F *h;
|
---|
426 |
|
---|
427 | if (fHiGainArray)
|
---|
428 | for (Int_t i=0;i<fHiGainArray->GetSize(); i++)
|
---|
429 | {
|
---|
430 | MHCalibrationPix &pix = (*this)[i];
|
---|
431 | h = pix.GetHGausHist();
|
---|
432 | h->SetName (Form("H%sHiGainPix%04d",fHistName.Data(),i));
|
---|
433 | h->SetTitle(Form("%s High Gain Pixel %04d Runs: ",fHistTitle.Data(),i));
|
---|
434 | h->SetXTitle(fHistXTitle.Data());
|
---|
435 | h->SetYTitle(fHistYTitle.Data());
|
---|
436 | }
|
---|
437 |
|
---|
438 | if (IsAverageing())
|
---|
439 | {
|
---|
440 | if (fAverageHiGainAreas)
|
---|
441 | for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
|
---|
442 | {
|
---|
443 | MHCalibrationPix &pix = GetAverageHiGainArea(j);
|
---|
444 | h = pix.GetHGausHist();
|
---|
445 | h->SetName (Form("H%sHiGainArea%d",fHistName.Data(),j));
|
---|
446 | h->SetXTitle(fHistXTitle.Data());
|
---|
447 | h->SetYTitle(fHistYTitle.Data());
|
---|
448 | if (fGeom->InheritsFrom("MGeomCamMagic"))
|
---|
449 | h->SetTitle(Form("%s averaged on event-by-event basis %s High Gain Runs: ",
|
---|
450 | fHistTitle.Data(), j==0 ? "Inner Pixels" : "Outer Pixels"));
|
---|
451 | else
|
---|
452 | h->SetTitle(Form("%s averaged on event-by-event basis High Gain Area Idx %d Runs: ",
|
---|
453 | fHistTitle.Data(), j));
|
---|
454 | }
|
---|
455 |
|
---|
456 | if (fAverageHiGainSectors)
|
---|
457 | for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
|
---|
458 | {
|
---|
459 | MHCalibrationPix &pix = GetAverageHiGainSector(j);
|
---|
460 | h = pix.GetHGausHist();
|
---|
461 | h->SetName (Form("H%sHiGainSector%02d",fHistName.Data(),j));
|
---|
462 | h->SetTitle(Form("%s averaged on event-by-event basis High Gain Sector %02d Runs: ",
|
---|
463 | fHistTitle.Data(), j));
|
---|
464 | h->SetXTitle(fHistXTitle.Data());
|
---|
465 | h->SetYTitle(fHistYTitle.Data());
|
---|
466 | }
|
---|
467 | }
|
---|
468 |
|
---|
469 | if (!IsLoGain())
|
---|
470 | return;
|
---|
471 |
|
---|
472 | if (fLoGainArray)
|
---|
473 | for (Int_t i=0;i<fLoGainArray->GetSize(); i++)
|
---|
474 | {
|
---|
475 | MHCalibrationPix &pix = (*this)(i);
|
---|
476 | h = pix.GetHGausHist();
|
---|
477 | h->SetName (Form("H%sLoGainPix%04d",fHistName.Data(),i));
|
---|
478 | h->SetTitle(Form("%s Low Gain Pixel %04d Runs: ", fHistTitle.Data(),i));
|
---|
479 | h->SetXTitle(fHistXTitle.Data());
|
---|
480 | h->SetYTitle(fHistYTitle.Data());
|
---|
481 | }
|
---|
482 |
|
---|
483 | if (IsAverageing())
|
---|
484 | {
|
---|
485 | if (fAverageLoGainAreas)
|
---|
486 | for (Int_t j=0; j<fAverageLoGainAreas->GetSize(); j++)
|
---|
487 | {
|
---|
488 | MHCalibrationPix &pix = GetAverageLoGainArea(j);
|
---|
489 | h = pix.GetHGausHist();
|
---|
490 | h->SetName (Form("H%sLoGainArea%d", fHistName.Data(), j));
|
---|
491 | h->SetXTitle(fHistXTitle.Data());
|
---|
492 | h->SetYTitle(fHistYTitle.Data());
|
---|
493 | if (fGeom->InheritsFrom("MGeomCamMagic"))
|
---|
494 | h->SetTitle(Form("%s averaged on event-by-event basis %s Low Gain Runs: ",
|
---|
495 | fHistTitle.Data(), j==0 ? "Inner Pixels" : "Outer Pixels"));
|
---|
496 | else
|
---|
497 | h->SetTitle(Form("%s averaged on event-by-event basis Low Gain Area Idx %d Runs: ",
|
---|
498 | fHistTitle.Data(), j));
|
---|
499 | }
|
---|
500 |
|
---|
501 | if (fAverageLoGainSectors)
|
---|
502 | for (Int_t j=0; j<fAverageLoGainSectors->GetSize(); j++)
|
---|
503 | {
|
---|
504 | MHCalibrationPix &pix = GetAverageLoGainSector(j);
|
---|
505 | h = pix.GetHGausHist();
|
---|
506 | h->SetName (Form("H%sLoGainSector%02d",fHistName.Data(),j));
|
---|
507 | h->SetTitle(Form("%s averaged on event-by-event basis Low Gain Sector %02d Runs: ",
|
---|
508 | fHistTitle.Data(), j));
|
---|
509 | h->SetXTitle(fHistXTitle.Data());
|
---|
510 | h->SetYTitle(fHistYTitle.Data());
|
---|
511 | }
|
---|
512 | }
|
---|
513 | }
|
---|
514 |
|
---|
515 | // --------------------------------------------------------------------------
|
---|
516 | //
|
---|
517 | // Gets the pointers to:
|
---|
518 | // - MGeomCam
|
---|
519 | //
|
---|
520 | // Calls SetupHists(const MParList *pList)
|
---|
521 | //
|
---|
522 | // Calls Delete-Function of:
|
---|
523 | // - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
|
---|
524 | // - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
|
---|
525 | // - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
|
---|
526 | //
|
---|
527 | Bool_t MHCalibrationCam::SetupFill(const MParList *pList)
|
---|
528 | {
|
---|
529 |
|
---|
530 | fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
|
---|
531 | if (!fGeom)
|
---|
532 | {
|
---|
533 | *fLog << err << GetDescriptor()
|
---|
534 | << ": MGeomCam not found... aborting." << endl;
|
---|
535 | return kFALSE;
|
---|
536 | }
|
---|
537 |
|
---|
538 | fRunHeader = (MRawRunHeader*)pList->FindObject("MRawRunHeader");
|
---|
539 | if (!fRunHeader)
|
---|
540 | {
|
---|
541 | *fLog << warn << GetDescriptor()
|
---|
542 | << ": MRawRunHeader not found... will not store run numbers." << endl;
|
---|
543 | }
|
---|
544 |
|
---|
545 | return SetupHists(pList);
|
---|
546 | }
|
---|
547 |
|
---|
548 |
|
---|
549 | // --------------------------------------------------------------------------
|
---|
550 | //
|
---|
551 | // Searches MRawEvtHeader to find the correct pulser colour
|
---|
552 | //
|
---|
553 | // Gets or creates the pointers to:
|
---|
554 | // - MBadPixelsIntensityCam
|
---|
555 | // - MBadPixelsCam
|
---|
556 | //
|
---|
557 | // Searches pointer to:
|
---|
558 | // - MArrivalTimeCam
|
---|
559 | //
|
---|
560 | // Initializes, if empty to MArrivalTimeCam::GetSize() for:
|
---|
561 | // - MHCalibrationCam::fHiGainArray, MHCalibrationCam::fLoGainArray
|
---|
562 | //
|
---|
563 | // Initializes, if empty to MGeomCam::GetNumAreas() for:
|
---|
564 | // - MHCalibrationCam::fAverageHiGainAreas, MHCalibrationCam::fAverageLoGainAreas
|
---|
565 | //
|
---|
566 | // Initializes, if empty to MGeomCam::GetNumSectors() for:
|
---|
567 | // - MHCalibrationCam::fAverageHiGainSectors, MHCalibrationCam::fAverageLoGainSectors
|
---|
568 | //
|
---|
569 | // Initializes TArray's to MGeomCam::GetNumAreas and MGeomCam::GetNumSectors, respectively
|
---|
570 | // Fills with number of valid pixels (if !MBadPixelsPix::IsBad()):
|
---|
571 | // - MHCalibrationCam::fAverageAreaNum[area index]
|
---|
572 | // - MHCalibrationCam::fAverageSectorNum[area index]
|
---|
573 | //
|
---|
574 | // Calls InitializeHists() for every entry in:
|
---|
575 | // - MHCalibrationCam::fHiGainArray
|
---|
576 | // - MHCalibrationCam::fAverageHiGainAreas
|
---|
577 | // - MHCalibrationCam::fAverageHiGainSectors
|
---|
578 | //
|
---|
579 | // Sets Titles and Names for the Histograms
|
---|
580 | // - MHCalibrationCam::fAverageHiGainAreas
|
---|
581 | // - MHCalibrationCam::fAverageHiGainSectors
|
---|
582 | //
|
---|
583 | // Retrieves the run numbers from MRawRunHeader and stores them in fRunNumbers
|
---|
584 | //
|
---|
585 | Bool_t MHCalibrationCam::ReInit(MParList *pList)
|
---|
586 | {
|
---|
587 |
|
---|
588 | const Int_t npixels = fGeom->GetNumPixels();
|
---|
589 | const Int_t nsectors = fGeom->GetNumSectors();
|
---|
590 | const Int_t nareas = fGeom->GetNumAreas();
|
---|
591 |
|
---|
592 | fIntensBad = (MBadPixelsIntensityCam*)pList->FindObject(AddSerialNumber("MBadPixelsIntensityCam"));
|
---|
593 | if (fIntensBad)
|
---|
594 | *fLog << inf << "Found MBadPixelsIntensityCam ... " << endl;
|
---|
595 | else
|
---|
596 | {
|
---|
597 | fBadPixels = (MBadPixelsCam*)pList->FindObject("MBadPixelsCam");
|
---|
598 | if (!fBadPixels)
|
---|
599 | {
|
---|
600 |
|
---|
601 | fBadPixels = (MBadPixelsCam*)pList->FindCreateObj(AddSerialNumber("MBadPixelsCam"));
|
---|
602 | if (!fBadPixels)
|
---|
603 | return kFALSE;
|
---|
604 |
|
---|
605 | fBadPixels->InitSize(npixels);
|
---|
606 | }
|
---|
607 | }
|
---|
608 |
|
---|
609 | if (IsAverageing())
|
---|
610 | {
|
---|
611 | //
|
---|
612 | // The function TArrayF::Set() already sets all entries to 0.
|
---|
613 | //
|
---|
614 | fAverageAreaNum. Set(nareas);
|
---|
615 | fAverageAreaSat. Set(nareas);
|
---|
616 | fAverageAreaSigma. Set(nareas);
|
---|
617 | fAverageAreaSigmaVar. Set(nareas);
|
---|
618 | fAverageAreaRelSigma. Set(nareas);
|
---|
619 | fAverageAreaRelSigmaVar.Set(nareas);
|
---|
620 | fAverageSectorNum. Set(nsectors);
|
---|
621 |
|
---|
622 | for (Int_t aidx=0; aidx<nareas; aidx++)
|
---|
623 | fAverageAreaNum[aidx] = 0;
|
---|
624 |
|
---|
625 | for (Int_t sector=0; sector<nsectors; sector++)
|
---|
626 | fAverageSectorNum[sector] = 0;
|
---|
627 |
|
---|
628 | for (Int_t i=0; i<npixels; i++)
|
---|
629 | {
|
---|
630 |
|
---|
631 | MBadPixelsPix &bad = fIntensBad ? (*fIntensBad)[i] : (*fBadPixels)[i];
|
---|
632 | if (bad.IsBad())
|
---|
633 | continue;
|
---|
634 |
|
---|
635 | fAverageAreaNum [(*fGeom)[i].GetAidx() ]++;
|
---|
636 | fAverageSectorNum[(*fGeom)[i].GetSector()]++;
|
---|
637 | }
|
---|
638 | }
|
---|
639 |
|
---|
640 | //
|
---|
641 | // Because ReInit has been called, a new run number is added
|
---|
642 | //
|
---|
643 | fRunNumbers.Set(fRunNumbers.GetSize()+1);
|
---|
644 |
|
---|
645 | if (fRunHeader)
|
---|
646 | {
|
---|
647 | fRunNumbers[fRunNumbers.GetSize()-1] = fRunHeader->GetRunNumber();
|
---|
648 | if (IsLoGain())
|
---|
649 | SetLoGain(fRunHeader->GetNumSamplesLoGain());
|
---|
650 | }
|
---|
651 |
|
---|
652 | if (!ReInitHists(pList))
|
---|
653 | return kFALSE;
|
---|
654 |
|
---|
655 | ResetHistTitles();
|
---|
656 |
|
---|
657 | if (!fRunHeader)
|
---|
658 | return kTRUE;
|
---|
659 |
|
---|
660 | for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
|
---|
661 | {
|
---|
662 | TH1F *h = (*this)[i].GetHGausHist();
|
---|
663 | h->SetTitle( Form("%s%08d ", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]));
|
---|
664 | }
|
---|
665 |
|
---|
666 | if (IsLoGain())
|
---|
667 | for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
|
---|
668 | {
|
---|
669 | TH1F *h = (*this)(i).GetHGausHist();
|
---|
670 | h->SetTitle( Form("%s%08d ", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]));
|
---|
671 | }
|
---|
672 |
|
---|
673 | if (!IsAverageing())
|
---|
674 | return kTRUE;
|
---|
675 |
|
---|
676 | for (Int_t j=0; j<nareas; j++)
|
---|
677 | {
|
---|
678 | TH1F *h = GetAverageHiGainArea(j).GetHGausHist();
|
---|
679 | h->SetTitle( Form("%s%08d ", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]));
|
---|
680 | }
|
---|
681 |
|
---|
682 | if (IsLoGain())
|
---|
683 | for (Int_t j=0; j<nareas; j++)
|
---|
684 | {
|
---|
685 | TH1F *h = GetAverageLoGainArea(j).GetHGausHist();
|
---|
686 | h->SetTitle( Form("%s%08d ", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]));
|
---|
687 | }
|
---|
688 |
|
---|
689 | for (Int_t j=0; j<nsectors; j++)
|
---|
690 | {
|
---|
691 | TH1F *h = GetAverageHiGainSector(j).GetHGausHist();
|
---|
692 | h->SetTitle( Form("%s%08d ", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]));
|
---|
693 | }
|
---|
694 |
|
---|
695 | if (IsLoGain())
|
---|
696 | for (Int_t j=0; j<nsectors; j++)
|
---|
697 | {
|
---|
698 | TH1F *h = GetAverageLoGainSector(j).GetHGausHist();
|
---|
699 | h->SetTitle( Form("%s%08d ", h->GetTitle(),fRunNumbers[fRunNumbers.GetSize()-1]));
|
---|
700 | }
|
---|
701 |
|
---|
702 | return kTRUE;
|
---|
703 | }
|
---|
704 |
|
---|
705 | //--------------------------------------------------------------------------------------
|
---|
706 | //
|
---|
707 | // Initializes the High Gain Arrays:
|
---|
708 | //
|
---|
709 | // - For every entry in the expanded arrays:
|
---|
710 | // * Initialize an MHCalibrationPix
|
---|
711 | // * Set Binning from fNbins, fFirst and fLast
|
---|
712 | // * Set Histgram names and titles from fHistName and fHistTitle
|
---|
713 | // * Set X-axis and Y-axis titles with fHistXTitle and fHistYTitle
|
---|
714 | // * Call InitHists
|
---|
715 | //
|
---|
716 | void MHCalibrationCam::InitHiGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
|
---|
717 | {
|
---|
718 |
|
---|
719 | if (fHiGainArray->GetSize()==0)
|
---|
720 | {
|
---|
721 | for (Int_t i=0; i<npixels; i++)
|
---|
722 | {
|
---|
723 | fHiGainArray->AddAt(new MHCalibrationPix(Form("%sHiGainPix%04d",fHistName.Data(),i),
|
---|
724 | Form("%s High Gain Pixel %4d",fHistTitle.Data(),i)),i);
|
---|
725 |
|
---|
726 | MHCalibrationPix &pix = (*this)[i];
|
---|
727 | pix.SetNbins(fNbins);
|
---|
728 | pix.SetFirst(fFirst);
|
---|
729 | pix.SetLast (fLast);
|
---|
730 | pix.SetProbLimit(fProbLimit);
|
---|
731 |
|
---|
732 | MBadPixelsPix &bad = fIntensBad ? (*fIntensBad)[i] : (*fBadPixels)[i];
|
---|
733 | InitHists(pix,bad,i);
|
---|
734 |
|
---|
735 | if (fCam)
|
---|
736 | (*fCam)[i].SetPixId(i);
|
---|
737 | }
|
---|
738 | }
|
---|
739 |
|
---|
740 | if (!IsAverageing())
|
---|
741 | return;
|
---|
742 |
|
---|
743 | if (fAverageHiGainAreas->GetSize()==0)
|
---|
744 | {
|
---|
745 | for (Int_t j=0; j<nareas; j++)
|
---|
746 | {
|
---|
747 | fAverageHiGainAreas->AddAt(new MHCalibrationPix(Form("%sHiGainArea%d",fHistName.Data(),j),
|
---|
748 | Form("%s High Gain Area Idx %d",fHistTitle.Data(),j)),j);
|
---|
749 |
|
---|
750 | MHCalibrationPix &pix = GetAverageHiGainArea(j);
|
---|
751 |
|
---|
752 | pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
|
---|
753 | pix.SetFirst(fFirst);
|
---|
754 | pix.SetLast (fLast);
|
---|
755 |
|
---|
756 | if (fGeom && fGeom->InheritsFrom("MGeomCamMagic"))
|
---|
757 | {
|
---|
758 | pix.InitBins();
|
---|
759 | pix.SetEventFrequency(fPulserFrequency);
|
---|
760 | }
|
---|
761 | else
|
---|
762 | InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadArea(j) : fCam->GetAverageBadArea(j),j);
|
---|
763 | }
|
---|
764 | }
|
---|
765 |
|
---|
766 | if (fAverageHiGainSectors->GetSize()==0)
|
---|
767 | {
|
---|
768 | for (Int_t j=0; j<nsectors; j++)
|
---|
769 | {
|
---|
770 | fAverageHiGainSectors->AddAt(new MHCalibrationPix(Form("%sHiGainSector%02d",fHistName.Data(),j),
|
---|
771 | Form("%s High Gain Sector %02d",fHistTitle.Data(),j)),j);
|
---|
772 | MHCalibrationPix &pix = GetAverageHiGainSector(j);
|
---|
773 |
|
---|
774 | pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nsectors));
|
---|
775 | pix.SetFirst(fFirst);
|
---|
776 | pix.SetLast (fLast);
|
---|
777 |
|
---|
778 | InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadSector(j) : fCam->GetAverageBadSector(j),j);
|
---|
779 | }
|
---|
780 | }
|
---|
781 | }
|
---|
782 |
|
---|
783 | //--------------------------------------------------------------------------------------
|
---|
784 | //
|
---|
785 | // Return, if IsLoGain() is kFALSE
|
---|
786 | //
|
---|
787 | // Initializes the Low Gain Arrays:
|
---|
788 | //
|
---|
789 | // - For every entry in the expanded arrays:
|
---|
790 | // * Initialize an MHCalibrationPix
|
---|
791 | // * Set Binning from fNbins, fFirst and fLast
|
---|
792 | // * Set Histgram names and titles from fHistName and fHistTitle
|
---|
793 | // * Set X-axis and Y-axis titles with fHistXTitle and fHistYTitle
|
---|
794 | // * Call InitHists
|
---|
795 | //
|
---|
796 | void MHCalibrationCam::InitLoGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
|
---|
797 | {
|
---|
798 |
|
---|
799 | if (!IsLoGain())
|
---|
800 | return;
|
---|
801 |
|
---|
802 | if (fLoGainArray->GetSize()==0)
|
---|
803 | {
|
---|
804 | for (Int_t i=0; i<npixels; i++)
|
---|
805 | {
|
---|
806 | fLoGainArray->AddAt(new MHCalibrationPix(Form("%sLoGainPix%04d",fHistName.Data(),i),
|
---|
807 | Form("%s Low Gain Pixel%04d",fHistTitle.Data(),i)),i);
|
---|
808 |
|
---|
809 | MHCalibrationPix &pix = (*this)(i);
|
---|
810 | pix.SetNbins(fNbins);
|
---|
811 | pix.SetFirst(fFirst);
|
---|
812 | pix.SetLast (fLast);
|
---|
813 | pix.SetProbLimit(fProbLimit);
|
---|
814 |
|
---|
815 | MBadPixelsPix &bad = fIntensBad ? (*fIntensBad)[i] : (*fBadPixels)[i];
|
---|
816 | InitHists(pix,bad,i);
|
---|
817 | }
|
---|
818 | }
|
---|
819 |
|
---|
820 | if (!IsAverageing())
|
---|
821 | return;
|
---|
822 |
|
---|
823 | if (fAverageLoGainAreas->GetSize()==0)
|
---|
824 | {
|
---|
825 | for (Int_t j=0; j<nareas; j++)
|
---|
826 | {
|
---|
827 | fAverageLoGainAreas->AddAt(new MHCalibrationPix(Form("%sLoGainArea%d",fHistName.Data(),j),
|
---|
828 | Form("%s Low Gain Area Idx %d",fHistTitle.Data(),j)),j);
|
---|
829 |
|
---|
830 | MHCalibrationPix &pix = GetAverageLoGainArea(j);
|
---|
831 |
|
---|
832 | pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
|
---|
833 | pix.SetFirst(fFirst);
|
---|
834 | pix.SetLast (fLast);
|
---|
835 |
|
---|
836 | if (fGeom && fGeom->InheritsFrom("MGeomCamMagic"))
|
---|
837 | {
|
---|
838 | pix.InitBins();
|
---|
839 | pix.SetEventFrequency(fPulserFrequency);
|
---|
840 | }
|
---|
841 | else
|
---|
842 | InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadArea(j) : fCam->GetAverageBadArea(j),j);
|
---|
843 | }
|
---|
844 | }
|
---|
845 |
|
---|
846 | if (fAverageLoGainSectors->GetSize()==0)
|
---|
847 | {
|
---|
848 | for (Int_t j=0; j<nsectors; j++)
|
---|
849 | {
|
---|
850 | fAverageLoGainSectors->AddAt(new MHCalibrationPix(Form("%sLoGainSector%02d",fHistName.Data(),j),
|
---|
851 | Form("%s Low Gain Sector %02d",fHistTitle.Data(),j)),j);
|
---|
852 | MHCalibrationPix &pix = GetAverageLoGainSector(j);
|
---|
853 |
|
---|
854 | pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nsectors));
|
---|
855 | pix.SetFirst(fFirst);
|
---|
856 | pix.SetLast (fLast);
|
---|
857 |
|
---|
858 | InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadSector(j) : fCam->GetAverageBadSector(j),j);
|
---|
859 | }
|
---|
860 | }
|
---|
861 | }
|
---|
862 |
|
---|
863 | //--------------------------------------------------------------------------------
|
---|
864 | //
|
---|
865 | // Retrieves from MGeomCam:
|
---|
866 | // - number of pixels
|
---|
867 | // - number of pixel areas
|
---|
868 | // - number of sectors
|
---|
869 | //
|
---|
870 | // Return kFALSE, if sizes of the TOrdCollections do not match npixels, nareas or nsectors
|
---|
871 | //
|
---|
872 | // Call FillHists()
|
---|
873 | //
|
---|
874 | Bool_t MHCalibrationCam::Fill(const MParContainer *par, const Stat_t w)
|
---|
875 | {
|
---|
876 | SetIsReset(kFALSE);
|
---|
877 |
|
---|
878 | if (!IsSizeCheck())
|
---|
879 | return FillHists(par,w);
|
---|
880 |
|
---|
881 | const Int_t npixels = fGeom->GetNumPixels();
|
---|
882 | const Int_t nareas = fGeom->GetNumAreas();
|
---|
883 | const Int_t nsectors = fGeom->GetNumSectors();
|
---|
884 |
|
---|
885 | //
|
---|
886 | // Hi-Gain OrdCollections
|
---|
887 | //
|
---|
888 | if (fHiGainArray->GetSize() != npixels)
|
---|
889 | {
|
---|
890 | *fLog << err << "ERROR - Size mismatch in number of pixels... abort." << endl;
|
---|
891 | return kFALSE;
|
---|
892 | }
|
---|
893 |
|
---|
894 | if (IsLoGain())
|
---|
895 | {
|
---|
896 | if (fLoGainArray->GetSize() != npixels)
|
---|
897 | {
|
---|
898 | *fLog << err << "ERROR - Size mismatch in number of pixels... abort." << endl;
|
---|
899 | return kFALSE;
|
---|
900 | }
|
---|
901 | }
|
---|
902 |
|
---|
903 | if (!IsAverageing())
|
---|
904 | return FillHists(par,w);
|
---|
905 |
|
---|
906 | if (fAverageHiGainAreas->GetSize() != nareas)
|
---|
907 | {
|
---|
908 | *fLog << err << "ERROR - Size mismatch in number of areas ... abort." << endl;
|
---|
909 | return kFALSE;
|
---|
910 | }
|
---|
911 |
|
---|
912 | if (fAverageHiGainSectors->GetSize() != nsectors)
|
---|
913 | {
|
---|
914 | *fLog << err << "ERROR - Size mismatch in number of sectors ... abort." << endl;
|
---|
915 | return kFALSE;
|
---|
916 | }
|
---|
917 |
|
---|
918 | if (IsLoGain())
|
---|
919 | {
|
---|
920 |
|
---|
921 | if (fAverageLoGainAreas->GetSize() != nareas)
|
---|
922 | {
|
---|
923 | *fLog << err << "ERROR - Size mismatch in number of areas ... abort." << endl;
|
---|
924 | return kFALSE;
|
---|
925 | }
|
---|
926 |
|
---|
927 | if (fAverageLoGainSectors->GetSize() != nsectors)
|
---|
928 | {
|
---|
929 | *fLog << err << "ERROR - Size mismatch in number of sectors ... abort." << endl;
|
---|
930 | return kFALSE;
|
---|
931 | }
|
---|
932 | }
|
---|
933 |
|
---|
934 | return FillHists(par, w);
|
---|
935 | }
|
---|
936 |
|
---|
937 | // --------------------------------------------------------------------------
|
---|
938 | //
|
---|
939 | // 0) Ask if fHiGainArray and fLoGainArray have been initialized,
|
---|
940 | // otherwise return kFALSE.
|
---|
941 | // 1) FinalizeHists()
|
---|
942 | // 2) FinalizeBadPixels()
|
---|
943 | // 3) CalcAverageSigma()
|
---|
944 | //
|
---|
945 | Bool_t MHCalibrationCam::Finalize()
|
---|
946 | {
|
---|
947 | if (IsReset())
|
---|
948 | return kTRUE;
|
---|
949 |
|
---|
950 | if (GetNumExecutions() < 2)
|
---|
951 | return kTRUE;
|
---|
952 |
|
---|
953 | if (fHiGainArray->GetSize() == 0 && fLoGainArray->GetSize() == 0)
|
---|
954 | {
|
---|
955 | *fLog << err << GetDescriptor()
|
---|
956 | << ": ERROR - Both (HiGain and LoGain) histogram arrays have not been initialized... abort." << endl;
|
---|
957 | return kFALSE;
|
---|
958 | }
|
---|
959 |
|
---|
960 | for (Int_t i=0; i<fAverageHiGainAreas->GetSize(); i++)
|
---|
961 | {
|
---|
962 | TH1F *h = GetAverageHiGainArea(i).GetHGausHist();
|
---|
963 | switch (fColor)
|
---|
964 | {
|
---|
965 | case MCalibrationCam::kNONE:
|
---|
966 | break;
|
---|
967 | case MCalibrationCam::kBLUE:
|
---|
968 | h->SetTitle( Form("%s%s", h->GetTitle(),"BLUE "));
|
---|
969 | break;
|
---|
970 | case MCalibrationCam::kGREEN:
|
---|
971 | h->SetTitle( Form("%s%s", h->GetTitle(),"GREEN "));
|
---|
972 | break;
|
---|
973 | case MCalibrationCam::kUV:
|
---|
974 | h->SetTitle( Form("%s%s", h->GetTitle(),"UV "));
|
---|
975 | break;
|
---|
976 | case MCalibrationCam::kCT1:
|
---|
977 | h->SetTitle( Form("%s%s", h->GetTitle(),"CT1-Pulser "));
|
---|
978 | break;
|
---|
979 | }
|
---|
980 | }
|
---|
981 |
|
---|
982 | for (Int_t i=0; i<fAverageLoGainAreas->GetSize(); i++)
|
---|
983 | {
|
---|
984 | TH1F *h = GetAverageLoGainArea(i).GetHGausHist();
|
---|
985 | switch (fColor)
|
---|
986 | {
|
---|
987 | case MCalibrationCam::kNONE:
|
---|
988 | break;
|
---|
989 | case MCalibrationCam::kBLUE:
|
---|
990 | h->SetTitle( Form("%s%s", h->GetTitle(),"BLUE "));
|
---|
991 | break;
|
---|
992 | case MCalibrationCam::kGREEN:
|
---|
993 | h->SetTitle( Form("%s%s", h->GetTitle(),"GREEN "));
|
---|
994 | break;
|
---|
995 | case MCalibrationCam::kUV:
|
---|
996 | h->SetTitle( Form("%s%s", h->GetTitle(),"UV "));
|
---|
997 | break;
|
---|
998 | case MCalibrationCam::kCT1:
|
---|
999 | h->SetTitle( Form("%s%s", h->GetTitle(),"CT1-Pulser "));
|
---|
1000 | break;
|
---|
1001 | }
|
---|
1002 | }
|
---|
1003 |
|
---|
1004 | if (!FinalizeHists())
|
---|
1005 | return kFALSE;
|
---|
1006 |
|
---|
1007 |
|
---|
1008 | FinalizeBadPixels();
|
---|
1009 | CalcAverageSigma();
|
---|
1010 |
|
---|
1011 | return kTRUE;
|
---|
1012 | }
|
---|
1013 |
|
---|
1014 | // -------------------------------------------------------------
|
---|
1015 | //
|
---|
1016 | // If MBadPixelsPix::IsUnsuitable(MBadPixelsPix::kUnsuitableRun):
|
---|
1017 | // - calls MHCalibrationPix::SetExcluded()
|
---|
1018 | //
|
---|
1019 | // Calls:
|
---|
1020 | // - MHGausEvents::InitBins()
|
---|
1021 | // - MHCalibrationPix::ChangeHistId(i)
|
---|
1022 | // - MHCalibrationPix::SetEventFrequency(fPulserFrequency)
|
---|
1023 | //
|
---|
1024 | void MHCalibrationCam::InitHists(MHCalibrationPix &hist, MBadPixelsPix &bad, const Int_t i)
|
---|
1025 | {
|
---|
1026 |
|
---|
1027 | if (bad.IsUnsuitable(MBadPixelsPix::kUnsuitableRun))
|
---|
1028 | hist.SetExcluded();
|
---|
1029 |
|
---|
1030 | hist.InitBins();
|
---|
1031 | hist.SetEventFrequency(fPulserFrequency);
|
---|
1032 | }
|
---|
1033 |
|
---|
1034 | // -------------------------------------------------------------
|
---|
1035 | //
|
---|
1036 | // - Searches for the CalibrationIntensity*Cam corresponding to 'name'.
|
---|
1037 | // - In case, it does not exist in the parameter list, it searches
|
---|
1038 | // for the corresponding MCalibration*Cam.
|
---|
1039 | // - Initializes the MCalibration*Cam, if not yet done.
|
---|
1040 | //
|
---|
1041 | Bool_t MHCalibrationCam::InitCams( MParList *plist, const TString name )
|
---|
1042 | {
|
---|
1043 |
|
---|
1044 | TString intensname = "MCalibrationIntensity";
|
---|
1045 | intensname += name;
|
---|
1046 | intensname += "Cam";
|
---|
1047 |
|
---|
1048 | TString ordname = "MCalibration";
|
---|
1049 | ordname += name;
|
---|
1050 | ordname += "Cam";
|
---|
1051 |
|
---|
1052 | fIntensCam = (MCalibrationIntensityCam*)plist->FindObject(AddSerialNumber(intensname));
|
---|
1053 | if (fIntensCam)
|
---|
1054 | *fLog << inf << "Found " << intensname << "... " << endl;
|
---|
1055 | else
|
---|
1056 | {
|
---|
1057 | fCam = (MCalibrationCam*)plist->FindObject(AddSerialNumber(ordname));
|
---|
1058 | if (!fCam)
|
---|
1059 | {
|
---|
1060 | fCam = (MCalibrationCam*)plist->FindCreateObj(AddSerialNumber(ordname));
|
---|
1061 | if (!fCam)
|
---|
1062 | return kFALSE;
|
---|
1063 |
|
---|
1064 | fCam->Init(*fGeom);
|
---|
1065 | }
|
---|
1066 | }
|
---|
1067 | return kTRUE;
|
---|
1068 | }
|
---|
1069 |
|
---|
1070 | // --------------------------------------------------------------------------
|
---|
1071 | //
|
---|
1072 | // Calls FitHiGainHists for every entry in:
|
---|
1073 | // - fHiGainArray
|
---|
1074 | // - fAverageHiGainAreas
|
---|
1075 | // - fAverageHiGainSectors
|
---|
1076 | //
|
---|
1077 | void MHCalibrationCam::FitHiGainArrays(MCalibrationCam &calcam, MBadPixelsCam &badcam,
|
---|
1078 | MBadPixelsPix::UncalibratedType_t fittyp,
|
---|
1079 | MBadPixelsPix::UncalibratedType_t osctyp)
|
---|
1080 | {
|
---|
1081 | fIsHiGainFitRanges = TMath::Abs(fUpperFitLimitHiGain - fLowerFitLimitHiGain) > 1E-5;
|
---|
1082 |
|
---|
1083 | for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
|
---|
1084 | {
|
---|
1085 |
|
---|
1086 | MHCalibrationPix &hist = (*this)[i];
|
---|
1087 |
|
---|
1088 | if (hist.IsExcluded())
|
---|
1089 | continue;
|
---|
1090 |
|
---|
1091 | MCalibrationPix &pix = calcam[i];
|
---|
1092 | MBadPixelsPix &bad = badcam[i];
|
---|
1093 |
|
---|
1094 | FitHiGainHists(hist,pix,bad,fittyp,osctyp);
|
---|
1095 | }
|
---|
1096 |
|
---|
1097 | if (!IsAverageing())
|
---|
1098 | return;
|
---|
1099 |
|
---|
1100 | for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
|
---|
1101 | {
|
---|
1102 |
|
---|
1103 | MHCalibrationPix &hist = GetAverageHiGainArea(j);
|
---|
1104 | MCalibrationPix &pix = calcam.GetAverageArea(j);
|
---|
1105 | MBadPixelsPix &bad = calcam.GetAverageBadArea(j);
|
---|
1106 |
|
---|
1107 | FitHiGainHists(hist,pix,bad,fittyp,osctyp);
|
---|
1108 | }
|
---|
1109 |
|
---|
1110 | for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
|
---|
1111 | {
|
---|
1112 | MHCalibrationPix &hist = GetAverageHiGainSector(j);
|
---|
1113 | MCalibrationPix &pix = calcam.GetAverageSector(j);
|
---|
1114 | MBadPixelsPix &bad = calcam.GetAverageBadSector(j);
|
---|
1115 |
|
---|
1116 | FitHiGainHists(hist,pix,bad,fittyp,osctyp);
|
---|
1117 | }
|
---|
1118 | }
|
---|
1119 |
|
---|
1120 | // --------------------------------------------------------------------------
|
---|
1121 | //
|
---|
1122 | // Calls FitLoGainHists for every entry in:
|
---|
1123 | // - fLoGainArray
|
---|
1124 | // - fAverageLoGainAreas
|
---|
1125 | // - fAverageLoGainSectors
|
---|
1126 | //
|
---|
1127 | void MHCalibrationCam::FitLoGainArrays(MCalibrationCam &calcam, MBadPixelsCam &badcam,
|
---|
1128 | MBadPixelsPix::UncalibratedType_t fittyp,
|
---|
1129 | MBadPixelsPix::UncalibratedType_t osctyp)
|
---|
1130 | {
|
---|
1131 | fIsLoGainFitRanges = TMath::Abs(fUpperFitLimitLoGain - fLowerFitLimitLoGain) > 1E-5;
|
---|
1132 |
|
---|
1133 | if (!IsLoGain())
|
---|
1134 | return;
|
---|
1135 |
|
---|
1136 | for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
|
---|
1137 | {
|
---|
1138 |
|
---|
1139 | MHCalibrationPix &hist = (*this)(i);
|
---|
1140 |
|
---|
1141 | if (hist.IsExcluded())
|
---|
1142 | continue;
|
---|
1143 |
|
---|
1144 | MCalibrationPix &pix = calcam[i];
|
---|
1145 | MBadPixelsPix &bad = badcam[i];
|
---|
1146 |
|
---|
1147 | FitLoGainHists(hist,pix,bad,fittyp,osctyp);
|
---|
1148 |
|
---|
1149 | }
|
---|
1150 |
|
---|
1151 | if (!IsAverageing())
|
---|
1152 | return;
|
---|
1153 |
|
---|
1154 | for (Int_t j=0; j<fAverageLoGainAreas->GetSize(); j++)
|
---|
1155 | {
|
---|
1156 |
|
---|
1157 | MHCalibrationPix &hist = GetAverageLoGainArea(j);
|
---|
1158 | MCalibrationPix &pix = calcam.GetAverageArea(j);
|
---|
1159 | MBadPixelsPix &bad = calcam.GetAverageBadArea(j);
|
---|
1160 |
|
---|
1161 | FitLoGainHists(hist,pix,bad,fittyp,osctyp);
|
---|
1162 | }
|
---|
1163 |
|
---|
1164 | for (Int_t j=0; j<fAverageLoGainSectors->GetSize(); j++)
|
---|
1165 | {
|
---|
1166 |
|
---|
1167 | MHCalibrationPix &hist = GetAverageLoGainSector(j);
|
---|
1168 | MCalibrationPix &pix = calcam.GetAverageSector(j);
|
---|
1169 | MBadPixelsPix &bad = calcam.GetAverageBadSector(j);
|
---|
1170 |
|
---|
1171 | FitLoGainHists(hist,pix,bad,fittyp,osctyp);
|
---|
1172 | }
|
---|
1173 | }
|
---|
1174 |
|
---|
1175 | //------------------------------------------------------------
|
---|
1176 | //
|
---|
1177 | // For all averaged areas, the fitted sigma is multiplied with the square root of
|
---|
1178 | // the number involved pixels
|
---|
1179 | //
|
---|
1180 | void MHCalibrationCam::CalcAverageSigma()
|
---|
1181 | {
|
---|
1182 |
|
---|
1183 | if (!fCam)
|
---|
1184 | return;
|
---|
1185 |
|
---|
1186 | if (!IsAverageing())
|
---|
1187 | return;
|
---|
1188 |
|
---|
1189 | for (UInt_t j=0; j<fGeom->GetNumAreas(); j++)
|
---|
1190 | {
|
---|
1191 |
|
---|
1192 | MCalibrationPix &pix = fCam->GetAverageArea(j);
|
---|
1193 |
|
---|
1194 | const Float_t numsqr = TMath::Sqrt((Float_t)fAverageAreaNum[j]);
|
---|
1195 | fAverageAreaSigma[j] = pix.GetSigma () * numsqr;
|
---|
1196 | fAverageAreaSigmaVar[j] = pix.GetSigmaErr () * pix.GetSigmaErr() * numsqr;
|
---|
1197 |
|
---|
1198 | pix.SetSigma (fAverageAreaSigma[j]);
|
---|
1199 | pix.SetSigmaVar(fAverageAreaSigmaVar[j]);
|
---|
1200 |
|
---|
1201 | fAverageAreaRelSigma [j] = fAverageAreaSigma[j] / pix.GetMean();
|
---|
1202 | fAverageAreaRelSigmaVar[j] = fAverageAreaSigmaVar[j] / (fAverageAreaSigma[j]*fAverageAreaSigma[j]);
|
---|
1203 | fAverageAreaRelSigmaVar[j] += pix.GetMeanRelVar();
|
---|
1204 | fAverageAreaRelSigmaVar[j] *= fAverageAreaRelSigma[j];
|
---|
1205 | }
|
---|
1206 | }
|
---|
1207 |
|
---|
1208 | // ---------------------------------------------------------------------------
|
---|
1209 | //
|
---|
1210 | // Returns if the histogram is empty and sets the following flag:
|
---|
1211 | // - MBadPixelsPix::SetUnsuitable(MBadPixelsPix::kUnsuitableRun)
|
---|
1212 | //
|
---|
1213 | // Fits the histograms with a Gaussian, in case of failure
|
---|
1214 | // calls MHCalibrationPix::RepeatFit(), in case of repeated failure
|
---|
1215 | // calls MHCalibrationPix::BypassFit() and sets the following flags:
|
---|
1216 | // - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t fittyp )
|
---|
1217 | // - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
|
---|
1218 | //
|
---|
1219 | // Creates the fourier spectrum and tests MHGausEvents::IsFourierSpectrumOK().
|
---|
1220 | // In case no, sets the following flags:
|
---|
1221 | // - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t osctyp )
|
---|
1222 | // - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
|
---|
1223 | //
|
---|
1224 | // Retrieves the results and store them in MCalibrationPix
|
---|
1225 | //
|
---|
1226 | void MHCalibrationCam::FitHiGainHists(MHCalibrationPix &hist,
|
---|
1227 | MCalibrationPix &pix,
|
---|
1228 | MBadPixelsPix &bad,
|
---|
1229 | MBadPixelsPix::UncalibratedType_t fittyp,
|
---|
1230 | MBadPixelsPix::UncalibratedType_t osctyp)
|
---|
1231 | {
|
---|
1232 |
|
---|
1233 | if (hist.IsEmpty() || hist.IsOnlyOverflow() || hist.IsOnlyUnderflow())
|
---|
1234 | {
|
---|
1235 | *fLog << warn << GetDescriptor()
|
---|
1236 | << ": Only overflow or underflow in hi-gain pixel: " << pix.GetPixId() << endl;
|
---|
1237 | return;
|
---|
1238 | }
|
---|
1239 | //
|
---|
1240 | // 2) Fit the Hi Gain histograms with a Gaussian
|
---|
1241 | //
|
---|
1242 | if (fIsHiGainFitRanges)
|
---|
1243 | {
|
---|
1244 | if (!hist.FitGaus("R",fLowerFitLimitHiGain,fUpperFitLimitHiGain))
|
---|
1245 | bad.SetUncalibrated( fittyp );
|
---|
1246 | }
|
---|
1247 | else
|
---|
1248 | if (!hist.FitGaus())
|
---|
1249 | //
|
---|
1250 | // 3) In case of failure set the bit Fitted to false and take histogram means and RMS
|
---|
1251 | //
|
---|
1252 | if (!hist.RepeatFit())
|
---|
1253 | {
|
---|
1254 | hist.BypassFit();
|
---|
1255 | bad.SetUncalibrated( fittyp );
|
---|
1256 | }
|
---|
1257 |
|
---|
1258 |
|
---|
1259 | //
|
---|
1260 | // 4) Check for oscillations
|
---|
1261 | //
|
---|
1262 | if (IsOscillations())
|
---|
1263 | {
|
---|
1264 | hist.CreateFourierSpectrum();
|
---|
1265 |
|
---|
1266 | if (!hist.IsFourierSpectrumOK())
|
---|
1267 | bad.SetUncalibrated( osctyp );
|
---|
1268 | }
|
---|
1269 |
|
---|
1270 | //
|
---|
1271 | // 5) Retrieve the results and store them in this class
|
---|
1272 | //
|
---|
1273 | pix.SetHiGainMean ( hist.GetMean() );
|
---|
1274 | pix.SetHiGainMeanVar ( hist.GetMeanErr() * hist.GetMeanErr() );
|
---|
1275 | pix.SetHiGainRms ( hist.GetHistRms() );
|
---|
1276 | pix.SetHiGainSigma ( hist.GetSigma() );
|
---|
1277 | pix.SetHiGainSigmaVar ( hist.GetSigmaErr()* hist.GetSigmaErr() );
|
---|
1278 | pix.SetHiGainProb ( hist.GetProb() );
|
---|
1279 | pix.SetHiGainNumBlackout( hist.GetBlackout() );
|
---|
1280 | pix.SetHiGainNumPickup ( hist.GetPickup() );
|
---|
1281 |
|
---|
1282 | if (IsDebug())
|
---|
1283 | {
|
---|
1284 | *fLog << dbginf << GetDescriptor() << ": ID " << GetName()
|
---|
1285 | << " HiGainSaturation: " << pix.IsHiGainSaturation()
|
---|
1286 | << " HiGainMean: " << hist.GetMean ()
|
---|
1287 | << " HiGainMeanErr: " << hist.GetMeanErr ()
|
---|
1288 | << " HiGainMeanSigma: " << hist.GetSigma ()
|
---|
1289 | << " HiGainMeanSigmaErr: " << hist.GetSigmaErr()
|
---|
1290 | << " HiGainMeanProb: " << hist.GetProb ()
|
---|
1291 | << " HiGainNumBlackout: " << hist.GetBlackout()
|
---|
1292 | << " HiGainNumPickup : " << hist.GetPickup ()
|
---|
1293 | << endl;
|
---|
1294 | }
|
---|
1295 |
|
---|
1296 | }
|
---|
1297 |
|
---|
1298 |
|
---|
1299 | // ---------------------------------------------------------------------------
|
---|
1300 | //
|
---|
1301 | // Returns if the histogram is empty and sets the following flag:
|
---|
1302 | // - MBadPixelsPix::SetUnsuitable(MBadPixelsPix::kUnsuitableRun)
|
---|
1303 | //
|
---|
1304 | // Fits the histograms with a Gaussian, in case of failure
|
---|
1305 | // calls MHCalibrationPix::RepeatFit(), in case of repeated failure
|
---|
1306 | // calls MHCalibrationPix::BypassFit() and sets the following flags:
|
---|
1307 | // - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t fittyp )
|
---|
1308 | // - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
|
---|
1309 | //
|
---|
1310 | // Creates the fourier spectrum and tests MHGausEvents::IsFourierSpectrumOK().
|
---|
1311 | // In case no, sets the following flags:
|
---|
1312 | // - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::UncalibratedType_t osctyp )
|
---|
1313 | // - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
|
---|
1314 | //
|
---|
1315 | // Retrieves the results and store them in MCalibrationPix
|
---|
1316 | //
|
---|
1317 | void MHCalibrationCam::FitLoGainHists(MHCalibrationPix &hist,
|
---|
1318 | MCalibrationPix &pix,
|
---|
1319 | MBadPixelsPix &bad,
|
---|
1320 | MBadPixelsPix::UncalibratedType_t fittyp,
|
---|
1321 | MBadPixelsPix::UncalibratedType_t osctyp)
|
---|
1322 | {
|
---|
1323 |
|
---|
1324 | if (hist.IsEmpty() || hist.IsOnlyOverflow() || hist.IsOnlyUnderflow())
|
---|
1325 | {
|
---|
1326 | *fLog << warn << GetDescriptor()
|
---|
1327 | << ": Only overflow or underflow in lo-gain pixel: " << pix.GetPixId() << endl;
|
---|
1328 | return;
|
---|
1329 | }
|
---|
1330 | //
|
---|
1331 | // 2) Fit the Hi Gain histograms with a Gaussian
|
---|
1332 | //
|
---|
1333 | if (fIsLoGainFitRanges)
|
---|
1334 | {
|
---|
1335 | if (!hist.FitGaus("R",fLowerFitLimitLoGain,fUpperFitLimitLoGain))
|
---|
1336 | bad.SetUncalibrated( fittyp );
|
---|
1337 | }
|
---|
1338 | else
|
---|
1339 | if (!hist.FitGaus())
|
---|
1340 | //
|
---|
1341 | // 3) In case of failure set the bit Fitted to false and take histogram means and RMS
|
---|
1342 | //
|
---|
1343 | if (!hist.RepeatFit())
|
---|
1344 | {
|
---|
1345 | hist.BypassFit();
|
---|
1346 | if (pix.IsHiGainSaturation())
|
---|
1347 | bad.SetUncalibrated( fittyp );
|
---|
1348 | }
|
---|
1349 |
|
---|
1350 | //
|
---|
1351 | // 4) Check for oscillations
|
---|
1352 | //
|
---|
1353 | if (IsOscillations())
|
---|
1354 | {
|
---|
1355 | hist.CreateFourierSpectrum();
|
---|
1356 |
|
---|
1357 | if (!hist.IsFourierSpectrumOK())
|
---|
1358 | bad.SetUncalibrated( osctyp );
|
---|
1359 | }
|
---|
1360 |
|
---|
1361 | //
|
---|
1362 | // 5) Retrieve the results and store them in this class
|
---|
1363 | //
|
---|
1364 | pix.SetLoGainMean ( hist.GetMean() );
|
---|
1365 | pix.SetLoGainMeanVar ( hist.GetMeanErr() * hist.GetMeanErr() );
|
---|
1366 | pix.SetLoGainRms ( hist.GetHistRms() );
|
---|
1367 | pix.SetLoGainSigma ( hist.GetSigma() );
|
---|
1368 | pix.SetLoGainSigmaVar ( hist.GetSigmaErr() * hist.GetSigmaErr() );
|
---|
1369 | pix.SetLoGainProb ( hist.GetProb() );
|
---|
1370 | pix.SetLoGainNumBlackout( hist.GetBlackout() );
|
---|
1371 | pix.SetLoGainNumPickup ( hist.GetPickup() );
|
---|
1372 |
|
---|
1373 | if (IsDebug())
|
---|
1374 | {
|
---|
1375 | *fLog << dbginf << GetDescriptor() << "ID: " << hist.GetName()
|
---|
1376 | << " HiGainSaturation: " << pix.IsHiGainSaturation()
|
---|
1377 | << " LoGainMean: " << hist.GetMean ()
|
---|
1378 | << " LoGainMeanErr: " << hist.GetMeanErr ()
|
---|
1379 | << " LoGainMeanSigma: " << hist.GetSigma ()
|
---|
1380 | << " LoGainMeanSigmaErr: " << hist.GetSigmaErr()
|
---|
1381 | << " LoGainMeanProb: " << hist.GetProb ()
|
---|
1382 | << " LoGainNumBlackout: " << hist.GetBlackout()
|
---|
1383 | << " LoGainNumPickup : " << hist.GetPickup ()
|
---|
1384 | << endl;
|
---|
1385 | }
|
---|
1386 |
|
---|
1387 | }
|
---|
1388 |
|
---|
1389 |
|
---|
1390 |
|
---|
1391 | // -----------------------------------------------------------------------------
|
---|
1392 | //
|
---|
1393 | // Default draw:
|
---|
1394 | //
|
---|
1395 | // Displays the averaged areas, both High Gain and Low Gain
|
---|
1396 | //
|
---|
1397 | // Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
|
---|
1398 | //
|
---|
1399 | void MHCalibrationCam::Draw(const Option_t *opt)
|
---|
1400 | {
|
---|
1401 |
|
---|
1402 | if (!IsAverageing())
|
---|
1403 | return;
|
---|
1404 |
|
---|
1405 | const Int_t nareas = fAverageHiGainAreas->GetSize();
|
---|
1406 | if (nareas == 0)
|
---|
1407 | return;
|
---|
1408 |
|
---|
1409 | TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
|
---|
1410 | pad->SetBorderMode(0);
|
---|
1411 |
|
---|
1412 | pad->Divide(IsLoGain() ? 2 : 1,nareas);
|
---|
1413 |
|
---|
1414 | for (Int_t i=0; i<nareas;i++)
|
---|
1415 | {
|
---|
1416 |
|
---|
1417 | pad->cd(IsLoGain() ? 2*(i+1)-1 : i+1);
|
---|
1418 | GetAverageHiGainArea(i).Draw(opt);
|
---|
1419 |
|
---|
1420 | if (!fAverageAreaSat[i])
|
---|
1421 | DrawAverageSigma(fAverageAreaSat[i], i,
|
---|
1422 | fAverageAreaSigma[i], fAverageAreaSigmaVar[i],
|
---|
1423 | fAverageAreaRelSigma[i], fAverageAreaRelSigmaVar[i]);
|
---|
1424 |
|
---|
1425 | if (IsLoGain())
|
---|
1426 | {
|
---|
1427 | pad->cd(2*(i+1));
|
---|
1428 | GetAverageLoGainArea(i).Draw(opt);
|
---|
1429 | }
|
---|
1430 |
|
---|
1431 | if (fAverageAreaSat[i])
|
---|
1432 | DrawAverageSigma(fAverageAreaSat[i], i,
|
---|
1433 | fAverageAreaSigma[i], fAverageAreaSigmaVar[i],
|
---|
1434 | fAverageAreaRelSigma[i], fAverageAreaRelSigmaVar[i]);
|
---|
1435 | }
|
---|
1436 |
|
---|
1437 | }
|
---|
1438 |
|
---|
1439 | // -----------------------------------------------------------------------------
|
---|
1440 | //
|
---|
1441 | // Default draw:
|
---|
1442 | //
|
---|
1443 | // Displays a TPaveText with the re-normalized sigmas of the average area
|
---|
1444 | //
|
---|
1445 | void MHCalibrationCam::DrawAverageSigma(Bool_t sat, Bool_t inner,
|
---|
1446 | Float_t sigma, Float_t sigmavar,
|
---|
1447 | Float_t relsigma, Float_t relsigmavar) const
|
---|
1448 | {
|
---|
1449 |
|
---|
1450 | if (sigma != 0 && sigmavar >= 0 && relsigmavar >= 0.)
|
---|
1451 | {
|
---|
1452 |
|
---|
1453 | TPad *newpad = new TPad("newpad","transparent",0,0,1,1);
|
---|
1454 | newpad->SetFillStyle(4000);
|
---|
1455 | newpad->Draw();
|
---|
1456 | newpad->cd();
|
---|
1457 |
|
---|
1458 | TPaveText *text = new TPaveText(sat? 0.1 : 0.35,0.7,sat ? 0.4 : 0.7,1.0);
|
---|
1459 | text->SetTextSize(0.05);
|
---|
1460 | const TString line1 = Form("%s%s%s",inner ? "Outer" : "Inner",
|
---|
1461 | " Pixels ", sat ? "Low Gain" : "High Gain");
|
---|
1462 | TText *txt1 = text->AddText(line1.Data());
|
---|
1463 | const TString line2 = Form("#sigma per pix: %2.2f #pm %2.2f",sigma,TMath::Sqrt(sigmavar));
|
---|
1464 | TText *txt2 = text->AddText(line2.Data());
|
---|
1465 | const TString line3 = Form("Rel. #sigma per pix: %2.2f #pm %2.2f",relsigma,TMath::Sqrt(relsigmavar));
|
---|
1466 | TText *txt3 = text->AddText(line3.Data());
|
---|
1467 | text->Draw("");
|
---|
1468 |
|
---|
1469 | text->SetBit(kCanDelete);
|
---|
1470 | txt1->SetBit(kCanDelete);
|
---|
1471 | txt2->SetBit(kCanDelete);
|
---|
1472 | txt3->SetBit(kCanDelete);
|
---|
1473 | newpad->SetBit(kCanDelete);
|
---|
1474 | }
|
---|
1475 | }
|
---|
1476 |
|
---|
1477 | // -----------------------------------------------------------------------------
|
---|
1478 | //
|
---|
1479 | // Available options
|
---|
1480 | // Debug
|
---|
1481 | // LoGain
|
---|
1482 | // Oscillations
|
---|
1483 | // SizeCheck
|
---|
1484 | // Averageing
|
---|
1485 | // Nbins
|
---|
1486 | // First
|
---|
1487 | // Last
|
---|
1488 | // ProbLimit
|
---|
1489 | // OverflowLimit
|
---|
1490 | // PulserFrequency
|
---|
1491 | // LowerFitLimitHiGain
|
---|
1492 | // UpperFitLimitHiGain
|
---|
1493 | // LowerFitLimitLoGain
|
---|
1494 | // UpperFitLimitLoGain
|
---|
1495 | //
|
---|
1496 | Int_t MHCalibrationCam::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
|
---|
1497 | {
|
---|
1498 |
|
---|
1499 | Bool_t rc = kFALSE;
|
---|
1500 | if (IsEnvDefined(env, prefix, "Debug", print))
|
---|
1501 | {
|
---|
1502 | SetDebug(GetEnvValue(env, prefix, "Debug", IsDebug()));
|
---|
1503 | rc = kTRUE;
|
---|
1504 | }
|
---|
1505 | if (IsEnvDefined(env, prefix, "LoGain", print))
|
---|
1506 | {
|
---|
1507 | SetDebug(GetEnvValue(env, prefix, "LoGain", IsLoGain()));
|
---|
1508 | rc = kTRUE;
|
---|
1509 | }
|
---|
1510 | if (IsEnvDefined(env, prefix, "Oscillations", print))
|
---|
1511 | {
|
---|
1512 | SetOscillations(GetEnvValue(env, prefix, "Oscillations", IsOscillations()));
|
---|
1513 | rc = kTRUE;
|
---|
1514 | }
|
---|
1515 | if (IsEnvDefined(env, prefix, "SizeCheck", print))
|
---|
1516 | {
|
---|
1517 | SetSizeCheck(GetEnvValue(env, prefix, "SizeCheck", IsSizeCheck()));
|
---|
1518 | rc = kTRUE;
|
---|
1519 | }
|
---|
1520 | if (IsEnvDefined(env, prefix, "Averageing", print))
|
---|
1521 | {
|
---|
1522 | SetAverageing(GetEnvValue(env, prefix, "Averageing", IsAverageing()));
|
---|
1523 | rc = kTRUE;
|
---|
1524 | }
|
---|
1525 |
|
---|
1526 | if (IsEnvDefined(env, prefix, "Nbins", print))
|
---|
1527 | {
|
---|
1528 | SetNbins(GetEnvValue(env, prefix, "Nbins", fNbins));
|
---|
1529 | rc = kTRUE;
|
---|
1530 | }
|
---|
1531 | if (IsEnvDefined(env, prefix, "First", print))
|
---|
1532 | {
|
---|
1533 | SetFirst(GetEnvValue(env, prefix, "First", fFirst));
|
---|
1534 | rc = kTRUE;
|
---|
1535 | }
|
---|
1536 |
|
---|
1537 | if (IsEnvDefined(env, prefix, "Last", print))
|
---|
1538 | {
|
---|
1539 | SetLast(GetEnvValue(env, prefix, "Last", fLast));
|
---|
1540 | rc = kTRUE;
|
---|
1541 | }
|
---|
1542 |
|
---|
1543 | if (IsEnvDefined(env, prefix, "ProbLimit", print))
|
---|
1544 | {
|
---|
1545 | SetProbLimit(GetEnvValue(env, prefix, "ProbLimit", fProbLimit));
|
---|
1546 | rc = kTRUE;
|
---|
1547 | }
|
---|
1548 |
|
---|
1549 | if (IsEnvDefined(env, prefix, "OverflowLimit", print))
|
---|
1550 | {
|
---|
1551 | SetOverflowLimit(GetEnvValue(env, prefix, "OverflowLimit", fOverflowLimit));
|
---|
1552 | rc = kTRUE;
|
---|
1553 | }
|
---|
1554 |
|
---|
1555 | if (IsEnvDefined(env, prefix, "PulserFrequency", print))
|
---|
1556 | {
|
---|
1557 | SetPulserFrequency(GetEnvValue(env, prefix, "PulserFrequency", fPulserFrequency));
|
---|
1558 | rc = kTRUE;
|
---|
1559 | }
|
---|
1560 |
|
---|
1561 | if (IsEnvDefined(env, prefix, "LowerFitLimitHiGain", print))
|
---|
1562 | {
|
---|
1563 | SetLowerFitLimitHiGain(GetEnvValue(env, prefix, "LowerFitLimitHiGain", fLowerFitLimitHiGain));
|
---|
1564 | rc = kTRUE;
|
---|
1565 | }
|
---|
1566 |
|
---|
1567 | if (IsEnvDefined(env, prefix, "UpperFitLimitHiGain", print))
|
---|
1568 | {
|
---|
1569 | SetUpperFitLimitHiGain(GetEnvValue(env, prefix, "UpperFitLimitHiGain", fUpperFitLimitHiGain));
|
---|
1570 | rc = kTRUE;
|
---|
1571 | }
|
---|
1572 |
|
---|
1573 | if (IsEnvDefined(env, prefix, "LowerFitLimitLoGain", print))
|
---|
1574 | {
|
---|
1575 | SetLowerFitLimitLoGain(GetEnvValue(env, prefix, "LowerFitLimitLoGain", fLowerFitLimitLoGain));
|
---|
1576 | rc = kTRUE;
|
---|
1577 | }
|
---|
1578 |
|
---|
1579 | if (IsEnvDefined(env, prefix, "UpperFitLimitLoGain", print))
|
---|
1580 | {
|
---|
1581 | SetUpperFitLimitLoGain(GetEnvValue(env, prefix, "UpperFitLimitLoGain", fUpperFitLimitLoGain));
|
---|
1582 | rc = kTRUE;
|
---|
1583 | }
|
---|
1584 |
|
---|
1585 |
|
---|
1586 | return rc;
|
---|
1587 | }
|
---|
1588 |
|
---|