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 | // MHCalibrationChargeCam
|
---|
27 | //
|
---|
28 | // Fills the extracted signals of MExtractedSignalCam into the MHCalibrationPix-classes
|
---|
29 | // MHCalibrationChargeHiGainPix and MHCalibrationChargeLoGainPix for every:
|
---|
30 | //
|
---|
31 | // - Pixel, stored in the TOrdCollection's MHCalibrationCam::fHiGainArray and
|
---|
32 | // MHCalibrationCam::fLoGainArray
|
---|
33 | //
|
---|
34 | // - Average pixel per AREA index (e.g. inner and outer for the MAGIC camera),
|
---|
35 | // stored in the TOrdCollection's MHCalibrationCam::fAverageHiGainAreas and
|
---|
36 | // MHCalibrationCam::fAverageLoGainAreas
|
---|
37 | //
|
---|
38 | // - Average pixel per camera SECTOR (e.g. sectors 1-6 for the MAGIC camera),
|
---|
39 | // stored in the TOrdCollection's MHCalibrationCam::fAverageHiGainSectors and
|
---|
40 | // MHCalibrationCam::fAverageLoGainSectors
|
---|
41 | //
|
---|
42 | // Every signal is taken from MExtractedSignalCam and filled into a histogram and
|
---|
43 | // an array, in order to perform a Fourier analysis (see MHGausEvents).
|
---|
44 | // The signals are moreover averaged on an event-by-event basis and written into
|
---|
45 | // the corresponding average pixels.
|
---|
46 | //
|
---|
47 | // Additionally, the (FADC slice) position of the maximum is stored in an Absolute
|
---|
48 | // Arrival Time histogram. This histogram serves for a rough cross-check if the
|
---|
49 | // signal does not lie at or outside the edges of the extraction window.
|
---|
50 | //
|
---|
51 | // The Charge histograms are fitted to a Gaussian, mean and sigma with its errors
|
---|
52 | // and the fit probability are extracted. If none of these values are NaN's and
|
---|
53 | // if the probability is bigger than MHGausEvents::fProbLimit (default: 0.5%),
|
---|
54 | // the fit is declared valid.
|
---|
55 | // Otherwise, the fit is repeated within ranges of the previous mean
|
---|
56 | // +- MHCalibrationPix::fPickupLimit (default: 5) sigma (see MHCalibrationPix::RepeatFit())
|
---|
57 | // In case this does not make the fit valid, the histogram means and RMS's are
|
---|
58 | // taken directly (see MHCalibrationPix::BypassFit()) and the following flags are set:
|
---|
59 | // - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainNotFitted ) or
|
---|
60 | // - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainNotFitted ) and
|
---|
61 | // - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
|
---|
62 | //
|
---|
63 | // Outliers of more than MHCalibrationPix::fPickupLimit (default: 5) sigmas
|
---|
64 | // from the mean are counted as Pickup events (stored in MHCalibrationPix::fPickup)
|
---|
65 | //
|
---|
66 | // Unless more than fNumHiGainSaturationLimit (default: 1%) of the overall FADC
|
---|
67 | // slices show saturation, the following flag is set:
|
---|
68 | // - MCalibrationChargePix::SetHiGainSaturation();
|
---|
69 | // In that case, the calibration constants are derived from the low-gain results.
|
---|
70 | //
|
---|
71 | // If more than fNumLoGainSaturationLimit (default: 1%) of the overall
|
---|
72 | // low-gain FADC slices saturate, the following flags are set:
|
---|
73 | // - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainSaturation ) and
|
---|
74 | // - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnsuitableRun )
|
---|
75 | //
|
---|
76 | // The class also fills arrays with the signal vs. event number, creates a fourier
|
---|
77 | // spectrum and investigates if the projected fourier components follow an exponential
|
---|
78 | // distribution. In case that the probability of the exponential fit is less than
|
---|
79 | // MHGausEvents::fProbLimit (default: 0.5%), the following flags are set:
|
---|
80 | // - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainOscillating ) or
|
---|
81 | // - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainOscillating ) and
|
---|
82 | // - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
|
---|
83 | //
|
---|
84 | // This same procedure is performed for the average pixels.
|
---|
85 | //
|
---|
86 | // The following results are written into MCalibrationChargeCam:
|
---|
87 | //
|
---|
88 | // - MCalibrationPix::SetHiGainSaturation()
|
---|
89 | // - MCalibrationPix::SetHiGainMean()
|
---|
90 | // - MCalibrationPix::SetHiGainMeanErr()
|
---|
91 | // - MCalibrationPix::SetHiGainSigma()
|
---|
92 | // - MCalibrationPix::SetHiGainSigmaErr()
|
---|
93 | // - MCalibrationPix::SetHiGainProb()
|
---|
94 | // - MCalibrationPix::SetHiGainNumPickup()
|
---|
95 | //
|
---|
96 | // - MCalibrationPix::SetLoGainMean()
|
---|
97 | // - MCalibrationPix::SetLoGainMeanErr()
|
---|
98 | // - MCalibrationPix::SetLoGainSigma()
|
---|
99 | // - MCalibrationPix::SetLoGainSigmaErr()
|
---|
100 | // - MCalibrationPix::SetLoGainProb()
|
---|
101 | // - MCalibrationPix::SetLoGainNumPickup()
|
---|
102 | //
|
---|
103 | // - MCalibrationChargePix::SetAbsTimeMean()
|
---|
104 | // - MCalibrationChargePix::SetAbsTimeRms()
|
---|
105 | //
|
---|
106 | // For all averaged areas, the fitted sigma is multiplied with the square root of
|
---|
107 | // the number involved pixels in order to be able to compare it to the average of
|
---|
108 | // sigmas in the camera.
|
---|
109 | //
|
---|
110 | /////////////////////////////////////////////////////////////////////////////
|
---|
111 | #include "MHCalibrationChargeCam.h"
|
---|
112 | #include "MHCalibrationCam.h"
|
---|
113 |
|
---|
114 | #include "MLog.h"
|
---|
115 | #include "MLogManip.h"
|
---|
116 |
|
---|
117 | #include "MParList.h"
|
---|
118 |
|
---|
119 | #include "MHCalibrationChargePix.h"
|
---|
120 | #include "MHCalibrationPix.h"
|
---|
121 |
|
---|
122 | #include "MCalibrationIntensityCam.h"
|
---|
123 | #include "MCalibrationChargeCam.h"
|
---|
124 | #include "MCalibrationChargePix.h"
|
---|
125 |
|
---|
126 | #include "MGeomCam.h"
|
---|
127 | #include "MGeomPix.h"
|
---|
128 |
|
---|
129 | #include "MBadPixelsIntensityCam.h"
|
---|
130 | #include "MBadPixelsCam.h"
|
---|
131 | #include "MBadPixelsPix.h"
|
---|
132 |
|
---|
133 | #include "MRawEvtData.h"
|
---|
134 | #include "MRawRunHeader.h"
|
---|
135 | #include "MRawEvtPixelIter.h"
|
---|
136 |
|
---|
137 | #include "MExtractedSignalCam.h"
|
---|
138 | #include "MExtractedSignalPix.h"
|
---|
139 |
|
---|
140 | #include "MArrayI.h"
|
---|
141 | #include "MArrayD.h"
|
---|
142 |
|
---|
143 | #include <TOrdCollection.h>
|
---|
144 | #include <TPad.h>
|
---|
145 | #include <TVirtualPad.h>
|
---|
146 | #include <TCanvas.h>
|
---|
147 | #include <TStyle.h>
|
---|
148 | #include <TF1.h>
|
---|
149 | #include <TLatex.h>
|
---|
150 | #include <TLegend.h>
|
---|
151 | #include <TGraph.h>
|
---|
152 |
|
---|
153 | ClassImp(MHCalibrationChargeCam);
|
---|
154 |
|
---|
155 | using namespace std;
|
---|
156 |
|
---|
157 | const Int_t MHCalibrationChargeCam::fgChargeHiGainNbins = 500;
|
---|
158 | const Axis_t MHCalibrationChargeCam::fgChargeHiGainFirst = -100.5;
|
---|
159 | const Axis_t MHCalibrationChargeCam::fgChargeHiGainLast = 1899.5;
|
---|
160 | const Int_t MHCalibrationChargeCam::fgChargeLoGainNbins = 620;
|
---|
161 | const Axis_t MHCalibrationChargeCam::fgChargeLoGainFirst = -350.5;
|
---|
162 | const Axis_t MHCalibrationChargeCam::fgChargeLoGainLast = 2049.5;
|
---|
163 | const TString MHCalibrationChargeCam::gsHistName = "Charge";
|
---|
164 | const TString MHCalibrationChargeCam::gsHistTitle = "Signals";
|
---|
165 | const TString MHCalibrationChargeCam::gsHistXTitle = "Signal [FADC counts]";
|
---|
166 | const TString MHCalibrationChargeCam::gsHistYTitle = "Nr. events";
|
---|
167 | const TString MHCalibrationChargeCam::gsAbsHistName = "AbsTime";
|
---|
168 | const TString MHCalibrationChargeCam::gsAbsHistTitle = "Abs. Arr. Times";
|
---|
169 | const TString MHCalibrationChargeCam::gsAbsHistXTitle = "Time [FADC slices]";
|
---|
170 | const TString MHCalibrationChargeCam::gsAbsHistYTitle = "Nr. events";
|
---|
171 | const Float_t MHCalibrationChargeCam::fgNumHiGainSaturationLimit = 0.01;
|
---|
172 | const Float_t MHCalibrationChargeCam::fgNumLoGainSaturationLimit = 0.005;
|
---|
173 | const Float_t MHCalibrationChargeCam::fgTimeLowerLimit = 1.;
|
---|
174 | const Float_t MHCalibrationChargeCam::fgTimeUpperLimit = 2.;
|
---|
175 | // 1Led Green, 1 LED blue, 5 LEDs blue, 10 LEDs blue, 10 LEDs UV, CT1, 5Leds Green
|
---|
176 | const Float_t MHCalibrationChargeCam::gkHiGainInnerRefLines[7] = { 245., 323. , 1065., 1467., 180., 211. , 533.5};
|
---|
177 | const Float_t MHCalibrationChargeCam::gkHiGainOuterRefLines[7] = { 217., 307.5, 932. , 1405., 167., 183.5, 405.5};
|
---|
178 | const Float_t MHCalibrationChargeCam::gkLoGainInnerRefLines[7] = { 20.8, 28.0 , 121. , 200.2, 16.5, 13.5 , 41.7 };
|
---|
179 | const Float_t MHCalibrationChargeCam::gkLoGainOuterRefLines[7] = { 18.9, 26.0 , 108.3, 198. , 14.0, 11. , 42. };
|
---|
180 | // --------------------------------------------------------------------------
|
---|
181 | //
|
---|
182 | // Default Constructor.
|
---|
183 | //
|
---|
184 | // Sets:
|
---|
185 | // - all pointers to NULL
|
---|
186 | //
|
---|
187 | // Initializes:
|
---|
188 | // - fNumHiGainSaturationLimit to fgNumHiGainSaturationLimit
|
---|
189 | // - fNumLoGainSaturationLimit to fgNumLoGainSaturationLimit
|
---|
190 | // - fTimeLowerLimit to fgTimeLowerLimit
|
---|
191 | // - fTimeUpperLimit to fgTimeUpperLimit
|
---|
192 | //
|
---|
193 | // - fNbins to fgChargeHiGainNbins
|
---|
194 | // - fFirst to fgChargeHiGainFirst
|
---|
195 | // - fLast to fgChargeHiGainLast
|
---|
196 | //
|
---|
197 | // - fLoGainNbins to fgChargeLoGainNbins
|
---|
198 | // - fLoGainFirst to fgChargeLoGainFirst
|
---|
199 | // - fLoGainLast to fgChargeLoGainLast
|
---|
200 | //
|
---|
201 | // - fHistName to gsHistName
|
---|
202 | // - fHistTitle to gsHistTitle
|
---|
203 | // - fHistXTitle to gsHistXTitle
|
---|
204 | // - fHistYTitle to gsHistYTitle
|
---|
205 | //
|
---|
206 | // - fAbsHistName to gsAbsHistName
|
---|
207 | // - fAbsHistTitle to gsAbsHistTitle
|
---|
208 | // - fAbsHistXTitle to gsAbsHistXTitle
|
---|
209 | // - fAbsHistYTitle to gsAbsHistYTitle
|
---|
210 | //
|
---|
211 | MHCalibrationChargeCam::MHCalibrationChargeCam(const char *name, const char *title)
|
---|
212 | : fRawEvt(NULL)
|
---|
213 | {
|
---|
214 |
|
---|
215 | fName = name ? name : "MHCalibrationChargeCam";
|
---|
216 | fTitle = title ? title : "Class to fill the calibration histograms ";
|
---|
217 |
|
---|
218 | SetNumHiGainSaturationLimit(fgNumHiGainSaturationLimit);
|
---|
219 | SetNumLoGainSaturationLimit(fgNumLoGainSaturationLimit);
|
---|
220 | SetTimeLowerLimit();
|
---|
221 | SetTimeUpperLimit();
|
---|
222 |
|
---|
223 | SetNbins(fgChargeHiGainNbins);
|
---|
224 | SetFirst(fgChargeHiGainFirst);
|
---|
225 | SetLast (fgChargeHiGainLast );
|
---|
226 |
|
---|
227 | SetLoGainNbins(fgChargeLoGainNbins);
|
---|
228 | SetLoGainFirst(fgChargeLoGainFirst);
|
---|
229 | SetLoGainLast (fgChargeLoGainLast );
|
---|
230 |
|
---|
231 | SetHistName (gsHistName .Data());
|
---|
232 | SetHistTitle (gsHistTitle .Data());
|
---|
233 | SetHistXTitle(gsHistXTitle.Data());
|
---|
234 | SetHistYTitle(gsHistYTitle.Data());
|
---|
235 |
|
---|
236 | SetAbsHistName (gsAbsHistName .Data());
|
---|
237 | SetAbsHistTitle (gsAbsHistTitle .Data());
|
---|
238 | SetAbsHistXTitle(gsAbsHistXTitle.Data());
|
---|
239 | SetAbsHistYTitle(gsAbsHistYTitle.Data());
|
---|
240 | }
|
---|
241 |
|
---|
242 | // --------------------------------------------------------------------------
|
---|
243 | //
|
---|
244 | // Creates new MHCalibrationChargeCam only with the averaged areas:
|
---|
245 | // the rest has to be retrieved directly, e.g. via:
|
---|
246 | // MHCalibrationChargeCam *cam = MParList::FindObject("MHCalibrationChargeCam");
|
---|
247 | // - cam->GetAverageSector(5).DrawClone();
|
---|
248 | // - (*cam)[100].DrawClone()
|
---|
249 | //
|
---|
250 | TObject *MHCalibrationChargeCam::Clone(const char *) const
|
---|
251 | {
|
---|
252 |
|
---|
253 | MHCalibrationChargeCam *cam = new MHCalibrationChargeCam();
|
---|
254 |
|
---|
255 | //
|
---|
256 | // Copy the data members
|
---|
257 | //
|
---|
258 | cam->fColor = fColor;
|
---|
259 | cam->fRunNumbers = fRunNumbers;
|
---|
260 | cam->fPulserFrequency = fPulserFrequency;
|
---|
261 | cam->fFlags = fFlags;
|
---|
262 | cam->fNbins = fNbins;
|
---|
263 | cam->fFirst = fFirst;
|
---|
264 | cam->fLast = fLast;
|
---|
265 | cam->fLoGainNbins = fLoGainNbins;
|
---|
266 | cam->fLoGainFirst = fLoGainFirst;
|
---|
267 | cam->fLoGainLast = fLoGainLast;
|
---|
268 |
|
---|
269 | //
|
---|
270 | // Copy the MArrays
|
---|
271 | //
|
---|
272 | cam->fAverageAreaRelSigma = fAverageAreaRelSigma;
|
---|
273 | cam->fAverageAreaRelSigmaVar = fAverageAreaRelSigmaVar;
|
---|
274 | cam->fAverageAreaSat = fAverageAreaSat;
|
---|
275 | cam->fAverageAreaSigma = fAverageAreaSigma;
|
---|
276 | cam->fAverageAreaSigmaVar = fAverageAreaSigmaVar;
|
---|
277 | cam->fAverageAreaNum = fAverageAreaNum;
|
---|
278 | cam->fAverageSectorNum = fAverageSectorNum;
|
---|
279 |
|
---|
280 | if (!IsAverageing())
|
---|
281 | return cam;
|
---|
282 |
|
---|
283 | const Int_t navhi = fAverageHiGainAreas->GetSize();
|
---|
284 |
|
---|
285 | for (int i=0; i<navhi; i++)
|
---|
286 | cam->fAverageHiGainAreas->AddAt(GetAverageHiGainArea(i).Clone(),i);
|
---|
287 |
|
---|
288 | if (IsLoGain())
|
---|
289 | {
|
---|
290 |
|
---|
291 | const Int_t navlo = fAverageLoGainAreas->GetSize();
|
---|
292 | for (int i=0; i<navlo; i++)
|
---|
293 | cam->fAverageLoGainAreas->AddAt(GetAverageLoGainArea(i).Clone(),i);
|
---|
294 |
|
---|
295 | }
|
---|
296 |
|
---|
297 | return cam;
|
---|
298 | }
|
---|
299 |
|
---|
300 | // --------------------------------------------------------------------------
|
---|
301 | //
|
---|
302 | // Gets the pointers to:
|
---|
303 | // - MRawEvtData
|
---|
304 | //
|
---|
305 | Bool_t MHCalibrationChargeCam::SetupHists(const MParList *pList)
|
---|
306 | {
|
---|
307 |
|
---|
308 | fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
|
---|
309 | if (!fRawEvt)
|
---|
310 | {
|
---|
311 | *fLog << err << dbginf << "MRawEvtData not found... aborting." << endl;
|
---|
312 | return kFALSE;
|
---|
313 | }
|
---|
314 |
|
---|
315 | return kTRUE;
|
---|
316 | }
|
---|
317 |
|
---|
318 | // --------------------------------------------------------------------------
|
---|
319 | //
|
---|
320 | // Gets or creates the pointers to:
|
---|
321 | // - MExtractedSignalCam
|
---|
322 | // - MCalibrationChargeCam or MCalibrationIntensityChargeCam
|
---|
323 | // - MBadPixelsCam
|
---|
324 | //
|
---|
325 | // Initializes the number of used FADC slices from MExtractedSignalCam
|
---|
326 | // into MCalibrationChargeCam and test for changes in that variable
|
---|
327 | //
|
---|
328 | // Calls:
|
---|
329 | // - InitHiGainArrays()
|
---|
330 | // - InitLoGainArrays()
|
---|
331 | //
|
---|
332 | // Sets:
|
---|
333 | // - fSumhiarea to nareas
|
---|
334 | // - fSumloarea to nareas
|
---|
335 | // - fTimehiarea to nareas
|
---|
336 | // - fTimeloarea to nareas
|
---|
337 | // - fSumhisector to nsectors
|
---|
338 | // - fSumlosector to nsectors
|
---|
339 | // - fTimehisector to nsectors
|
---|
340 | // - fTimelosector to nsectors
|
---|
341 | // - fSathiarea to nareas
|
---|
342 | // - fSatloarea to nareas
|
---|
343 | // - fSathisector to nsectors
|
---|
344 | // - fSatlosector to nsectors
|
---|
345 | //
|
---|
346 | Bool_t MHCalibrationChargeCam::ReInitHists(MParList *pList)
|
---|
347 | {
|
---|
348 |
|
---|
349 | MExtractedSignalCam *signal =
|
---|
350 | (MExtractedSignalCam*)pList->FindObject(AddSerialNumber("MExtractedSignalCam"));
|
---|
351 | if (!signal)
|
---|
352 | {
|
---|
353 | *fLog << err << "MExtractedSignalCam not found... abort." << endl;
|
---|
354 | return kFALSE;
|
---|
355 | }
|
---|
356 |
|
---|
357 | fIntensCam = (MCalibrationIntensityCam*)pList->FindObject(AddSerialNumber("MCalibrationIntensityChargeCam"));
|
---|
358 | if (fIntensCam)
|
---|
359 | *fLog << inf << "Found MCalibrationIntensityChargeCam ... " << endl;
|
---|
360 | else
|
---|
361 | {
|
---|
362 | fCam = (MCalibrationCam*)pList->FindObject(AddSerialNumber("MCalibrationChargeCam"));
|
---|
363 | if (!fCam)
|
---|
364 | {
|
---|
365 | fCam = (MCalibrationCam*)pList->FindCreateObj(AddSerialNumber("MCalibrationChargeCam"));
|
---|
366 | if (!fCam)
|
---|
367 | {
|
---|
368 | *fLog << err << "Cannot find nor create MCalibrationChargeCam ... abort." << endl;
|
---|
369 | return kFALSE;
|
---|
370 | }
|
---|
371 | fCam->Init(*fGeom);
|
---|
372 | }
|
---|
373 | }
|
---|
374 |
|
---|
375 | fFirstHiGain = signal->GetFirstUsedSliceHiGain();
|
---|
376 | fLastHiGain = signal->GetLastUsedSliceHiGain();
|
---|
377 | fFirstLoGain = signal->GetFirstUsedSliceLoGain();
|
---|
378 | fLastLoGain = signal->GetLastUsedSliceLoGain();
|
---|
379 |
|
---|
380 | /*
|
---|
381 | const Float_t numhigain = signal->GetNumUsedHiGainFADCSlices();
|
---|
382 | const Float_t numlogain = signal->GetNumUsedLoGainFADCSlices();
|
---|
383 |
|
---|
384 | if (fCam)
|
---|
385 | {
|
---|
386 | if (fCam->GetNumHiGainFADCSlices() == 0.)
|
---|
387 | fCam->SetNumHiGainFADCSlices ( numhigain );
|
---|
388 | else if (fCam->GetNumHiGainFADCSlices() != numhigain)
|
---|
389 | {
|
---|
390 | *fLog << err << GetDescriptor()
|
---|
391 | << ": Number of High Gain FADC extraction slices has changed, abort..." << endl;
|
---|
392 | return kFALSE;
|
---|
393 | }
|
---|
394 |
|
---|
395 | if (fCam->GetNumLoGainFADCSlices() == 0.)
|
---|
396 | fCam->SetNumLoGainFADCSlices ( numlogain );
|
---|
397 | else if (fCam->GetNumLoGainFADCSlices() != numlogain)
|
---|
398 | {
|
---|
399 | *fLog << err << GetDescriptor()
|
---|
400 | << ": Number of Low Gain FADC extraction slices has changes, abort..." << endl;
|
---|
401 | return kFALSE;
|
---|
402 | }
|
---|
403 | }
|
---|
404 | */
|
---|
405 |
|
---|
406 | const Int_t npixels = fGeom->GetNumPixels();
|
---|
407 | const Int_t nsectors = fGeom->GetNumSectors();
|
---|
408 | const Int_t nareas = fGeom->GetNumAreas();
|
---|
409 |
|
---|
410 | //
|
---|
411 | // In case of the intense blue, double the range
|
---|
412 | //
|
---|
413 | if (fGeom->InheritsFrom("MGeomCamMagic"))
|
---|
414 | if ( fColor == MCalibrationCam::kBLUE)
|
---|
415 | SetLoGainLast(2.*fLoGainLast - fLoGainFirst);
|
---|
416 |
|
---|
417 | InitHiGainArrays(npixels,nareas,nsectors);
|
---|
418 | InitLoGainArrays(npixels,nareas,nsectors);
|
---|
419 |
|
---|
420 | fSumhiarea .Set(nareas);
|
---|
421 | fSumloarea .Set(nareas);
|
---|
422 | fTimehiarea .Set(nareas);
|
---|
423 | fTimeloarea .Set(nareas);
|
---|
424 | fSumhisector.Set(nsectors);
|
---|
425 | fSumlosector.Set(nsectors);
|
---|
426 | fTimehisector.Set(nsectors);
|
---|
427 | fTimelosector.Set(nsectors);
|
---|
428 |
|
---|
429 | fSathiarea .Set(nareas);
|
---|
430 | fSatloarea .Set(nareas);
|
---|
431 | fSathisector.Set(nsectors);
|
---|
432 | fSatlosector.Set(nsectors);
|
---|
433 |
|
---|
434 | return kTRUE;
|
---|
435 | }
|
---|
436 |
|
---|
437 | // --------------------------------------------------------------------------
|
---|
438 | //
|
---|
439 | // Retrieve:
|
---|
440 | // - fRunHeader->GetNumSamplesHiGain();
|
---|
441 | //
|
---|
442 | // Initializes the High Gain Arrays:
|
---|
443 | //
|
---|
444 | // - For every entry in the expanded arrays:
|
---|
445 | // * Initialize an MHCalibrationChargePix
|
---|
446 | // * Set Binning from fNbins, fFirst and fLast
|
---|
447 | // * Set Binning of Abs Times histogram from fAbsNbins, fAbsFirst and fAbsLast
|
---|
448 | // * Set Histgram names and titles from fHistName and fHistTitle
|
---|
449 | // * Set Abs Times Histgram names and titles from fAbsHistName and fAbsHistTitle
|
---|
450 | // * Set X-axis and Y-axis titles from fHistXTitle and fHistYTitle
|
---|
451 | // * Set X-axis and Y-axis titles of Abs Times Histogram from fAbsHistXTitle and fAbsHistYTitle
|
---|
452 | // * Call InitHists
|
---|
453 | //
|
---|
454 | //
|
---|
455 | void MHCalibrationChargeCam::InitHiGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
|
---|
456 | {
|
---|
457 |
|
---|
458 | TH1F *h;
|
---|
459 |
|
---|
460 | const Int_t higainsamples = fRunHeader->GetNumSamplesHiGain();
|
---|
461 |
|
---|
462 | if (fHiGainArray->GetSize()==0)
|
---|
463 | {
|
---|
464 | for (Int_t i=0; i<npixels; i++)
|
---|
465 | {
|
---|
466 | fHiGainArray->AddAt(new MHCalibrationChargePix(Form("%s%s%4i",fHistName.Data(),"HiGainPix",i),
|
---|
467 | Form("%s%s%4i",fHistTitle.Data()," High Gain Pixel",i)),i);
|
---|
468 |
|
---|
469 | MHCalibrationChargePix &pix = (MHCalibrationChargePix&)(*this)[i];
|
---|
470 |
|
---|
471 | pix.SetNbins(fNbins);
|
---|
472 | pix.SetFirst(fFirst);
|
---|
473 | pix.SetLast (fLast);
|
---|
474 |
|
---|
475 | pix.SetAbsTimeNbins(higainsamples);
|
---|
476 | pix.SetAbsTimeFirst(-0.5);
|
---|
477 | pix.SetAbsTimeLast(higainsamples-0.5);
|
---|
478 |
|
---|
479 | h = pix.GetHAbsTime();
|
---|
480 |
|
---|
481 | h->SetName (Form("%s%s%s%4i","H",fAbsHistName.Data(),"HiGainPix",i));
|
---|
482 | h->SetTitle(Form("%s%s%4i",fAbsHistTitle.Data()," High Gain Pixel ",i));
|
---|
483 | h->SetXTitle(fAbsHistXTitle.Data());
|
---|
484 | h->SetYTitle(fAbsHistYTitle.Data());
|
---|
485 |
|
---|
486 | MBadPixelsPix &bad = fIntensBad ? (*fIntensBad)[i] : (*fBadPixels)[i];
|
---|
487 | InitHists((*this)[i],bad,i);
|
---|
488 | }
|
---|
489 | }
|
---|
490 |
|
---|
491 |
|
---|
492 | if (fAverageHiGainAreas->GetSize()==0)
|
---|
493 | {
|
---|
494 | for (Int_t j=0; j<nareas; j++)
|
---|
495 | {
|
---|
496 | fAverageHiGainAreas->AddAt(new MHCalibrationChargePix(Form("%s%s%d",fHistName.Data(),"HiGainArea",j),
|
---|
497 | Form("%s%s%d",fHistTitle.Data()," High Gain Area Idx ",j)),j);
|
---|
498 |
|
---|
499 | MHCalibrationChargePix &pix = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
|
---|
500 |
|
---|
501 | pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
|
---|
502 | pix.SetFirst(fFirst);
|
---|
503 | pix.SetLast (fLast);
|
---|
504 |
|
---|
505 | pix.SetAbsTimeNbins(higainsamples);
|
---|
506 | pix.SetAbsTimeFirst(-0.5);
|
---|
507 | pix.SetAbsTimeLast(higainsamples-0.5);
|
---|
508 |
|
---|
509 | InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadArea(j) : fCam->GetAverageBadArea(j),j);
|
---|
510 |
|
---|
511 | h = pix.GetHAbsTime();
|
---|
512 |
|
---|
513 | h->SetName (Form("%s%s%s%d","H",fAbsHistName.Data(),"HiGainArea",j));
|
---|
514 | h->SetTitle(Form("%s%s%d",fAbsHistTitle.Data(),
|
---|
515 | " averaged on event-by-event basis High Gain Area Idx ",j));
|
---|
516 | h->SetXTitle(fAbsHistXTitle.Data());
|
---|
517 | h->SetYTitle(fAbsHistYTitle.Data());
|
---|
518 | }
|
---|
519 | }
|
---|
520 |
|
---|
521 | if (fAverageHiGainSectors->GetSize()==0)
|
---|
522 | {
|
---|
523 | for (Int_t j=0; j<nsectors; j++)
|
---|
524 | {
|
---|
525 | fAverageHiGainSectors->AddAt(new MHCalibrationChargePix(Form("%s%s%2i",fHistName.Data(),"HiGainSector",j),
|
---|
526 | Form("%s%s%2i",fHistTitle.Data()," High Gain Sector ",j)),j);
|
---|
527 |
|
---|
528 | MHCalibrationChargePix &pix = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
|
---|
529 |
|
---|
530 | pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
|
---|
531 | pix.SetFirst(fFirst);
|
---|
532 | pix.SetLast (fLast);
|
---|
533 |
|
---|
534 | pix.SetAbsTimeNbins(higainsamples);
|
---|
535 | pix.SetAbsTimeFirst(-0.5);
|
---|
536 | pix.SetAbsTimeLast(higainsamples-0.5);
|
---|
537 |
|
---|
538 | h = pix.GetHAbsTime();
|
---|
539 |
|
---|
540 | h->SetName (Form("%s%s%s%2i","H",fAbsHistName.Data(),"HiGainSector",j));
|
---|
541 | h->SetTitle(Form("%s%s%2i",fAbsHistTitle.Data(),
|
---|
542 | " averaged on event-by-event basis High Gain Area Sector ",j));
|
---|
543 | h->SetXTitle(fAbsHistXTitle.Data());
|
---|
544 | h->SetYTitle(fAbsHistYTitle.Data());
|
---|
545 |
|
---|
546 | InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadSector(j) : fCam->GetAverageBadSector(j),j);
|
---|
547 | }
|
---|
548 | }
|
---|
549 | }
|
---|
550 |
|
---|
551 | //--------------------------------------------------------------------------------------
|
---|
552 | //
|
---|
553 | // Return, if IsLoGain() is kFALSE
|
---|
554 | //
|
---|
555 | // Retrieve:
|
---|
556 | // - fRunHeader->GetNumSamplesHiGain();
|
---|
557 | //
|
---|
558 | // Initializes the Low Gain Arrays:
|
---|
559 | //
|
---|
560 | // - For every entry in the expanded arrays:
|
---|
561 | // * Initialize an MHCalibrationChargePix
|
---|
562 | // * Set Binning from fNbins, fFirst and fLast
|
---|
563 | // * Set Binning of Abs Times histogram from fAbsNbins, fAbsFirst and fAbsLast
|
---|
564 | // * Set Histgram names and titles from fHistName and fHistTitle
|
---|
565 | // * Set Abs Times Histgram names and titles from fAbsHistName and fAbsHistTitle
|
---|
566 | // * Set X-axis and Y-axis titles from fHistXTitle and fHistYTitle
|
---|
567 | // * Set X-axis and Y-axis titles of Abs Times Histogram from fAbsHistXTitle and fAbsHistYTitle
|
---|
568 | // * Call InitHists
|
---|
569 | //
|
---|
570 | void MHCalibrationChargeCam::InitLoGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
|
---|
571 | {
|
---|
572 |
|
---|
573 | if (!IsLoGain())
|
---|
574 | return;
|
---|
575 |
|
---|
576 | const Int_t logainsamples = fRunHeader->GetNumSamplesLoGain();
|
---|
577 |
|
---|
578 | TH1F *h;
|
---|
579 |
|
---|
580 | if (fLoGainArray->GetSize()==0 )
|
---|
581 | {
|
---|
582 | for (Int_t i=0; i<npixels; i++)
|
---|
583 | {
|
---|
584 | fLoGainArray->AddAt(new MHCalibrationChargePix(Form("%s%s%4i",fHistName.Data(),"LoGainPix",i),
|
---|
585 | Form("%s%s%4i",fHistTitle.Data()," Low Gain Pixel",i)),i);
|
---|
586 |
|
---|
587 | MHCalibrationChargePix &pix = (MHCalibrationChargePix&)(*this)(i);
|
---|
588 |
|
---|
589 | pix.SetNbins(fLoGainNbins);
|
---|
590 | pix.SetFirst(fLoGainFirst);
|
---|
591 | pix.SetLast (fLoGainLast);
|
---|
592 |
|
---|
593 | pix.SetAbsTimeNbins(logainsamples);
|
---|
594 | pix.SetAbsTimeFirst(-0.5);
|
---|
595 | pix.SetAbsTimeLast(logainsamples-0.5);
|
---|
596 |
|
---|
597 | h = pix.GetHAbsTime();
|
---|
598 |
|
---|
599 | h->SetName (Form("%s%s%s%4i","H",fAbsHistName.Data(),"HiGainPix",i));
|
---|
600 | h->SetTitle(Form("%s%s%4i",fAbsHistTitle.Data()," High Gain Pixel ",i));
|
---|
601 | h->SetXTitle(fAbsHistXTitle.Data());
|
---|
602 | h->SetYTitle(fAbsHistYTitle.Data());
|
---|
603 |
|
---|
604 | MBadPixelsPix &bad = fIntensBad ? (*fIntensBad)[i] : (*fBadPixels)[i];
|
---|
605 | InitHists(pix,bad,i);
|
---|
606 | }
|
---|
607 | }
|
---|
608 |
|
---|
609 | if (fAverageLoGainAreas->GetSize()==0)
|
---|
610 | {
|
---|
611 | for (Int_t j=0; j<nareas; j++)
|
---|
612 | {
|
---|
613 | fAverageLoGainAreas->AddAt(new MHCalibrationChargePix(Form("%s%s%d",fHistName.Data(),"LoGainArea",j),
|
---|
614 | Form("%s%s%d",fHistTitle.Data()," Low Gain Area Idx ",j)),j);
|
---|
615 |
|
---|
616 | MHCalibrationChargePix &pix = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
|
---|
617 |
|
---|
618 | pix.SetNbins(fLoGainNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
|
---|
619 | pix.SetFirst(fLoGainFirst);
|
---|
620 | pix.SetLast (fLoGainLast);
|
---|
621 |
|
---|
622 | pix.SetAbsTimeNbins(logainsamples);
|
---|
623 | pix.SetAbsTimeFirst(-0.5);
|
---|
624 | pix.SetAbsTimeLast(logainsamples-0.5);
|
---|
625 |
|
---|
626 | h = pix.GetHAbsTime();
|
---|
627 |
|
---|
628 | h->SetName (Form("%s%s%s%2i","H",fAbsHistName.Data(),"LoGainArea",j));
|
---|
629 | h->SetTitle(Form("%s%s%2i",fAbsHistTitle.Data(),
|
---|
630 | " averaged on event-by-event basis Low Gain Area Idx ",j));
|
---|
631 | h->SetXTitle(fAbsHistXTitle.Data());
|
---|
632 | h->SetYTitle(fAbsHistYTitle.Data());
|
---|
633 |
|
---|
634 | InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadArea(j) : fCam->GetAverageBadArea(j),j);
|
---|
635 | }
|
---|
636 | }
|
---|
637 |
|
---|
638 |
|
---|
639 | if (fAverageLoGainSectors->GetSize()==0 && IsLoGain())
|
---|
640 | {
|
---|
641 | for (Int_t j=0; j<nsectors; j++)
|
---|
642 | {
|
---|
643 | fAverageLoGainSectors->AddAt(new MHCalibrationChargePix(Form("%s%s%2i",fHistName.Data(),"LoGainSector",j),
|
---|
644 | Form("%s%s%2i",fHistTitle.Data()," Low Gain Sector ",j)),j);
|
---|
645 |
|
---|
646 | MHCalibrationChargePix &pix = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
|
---|
647 |
|
---|
648 | pix.SetNbins(fLoGainNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
|
---|
649 | pix.SetFirst(fLoGainFirst);
|
---|
650 | pix.SetLast (fLoGainLast);
|
---|
651 |
|
---|
652 | pix.SetAbsTimeNbins(logainsamples);
|
---|
653 | pix.SetAbsTimeFirst(-0.5);
|
---|
654 | pix.SetAbsTimeLast(logainsamples-0.5);
|
---|
655 |
|
---|
656 | h = pix.GetHAbsTime();
|
---|
657 |
|
---|
658 | h->SetName (Form("%s%s%s%2i","H",fAbsHistName.Data(),"LoGainSector",j));
|
---|
659 | h->SetTitle(Form("%s%s%2i",fAbsHistTitle.Data(),
|
---|
660 | " averaged on event-by-event basis Low Gain Area Sector ",j));
|
---|
661 | h->SetXTitle(fAbsHistXTitle.Data());
|
---|
662 | h->SetYTitle(fAbsHistYTitle.Data());
|
---|
663 |
|
---|
664 | InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadSector(j) : fCam->GetAverageBadSector(j),j);
|
---|
665 | }
|
---|
666 | }
|
---|
667 | }
|
---|
668 |
|
---|
669 |
|
---|
670 | // --------------------------------------------------------------------------
|
---|
671 | //
|
---|
672 | // Retrieves from MExtractedSignalCam:
|
---|
673 | // - first used LoGain FADC slice
|
---|
674 | //
|
---|
675 | // Retrieves from MGeomCam:
|
---|
676 | // - number of pixels
|
---|
677 | // - number of pixel areas
|
---|
678 | // - number of sectors
|
---|
679 | //
|
---|
680 | // For all TOrdCollection's (including the averaged ones), the following steps are performed:
|
---|
681 | //
|
---|
682 | // 1) Fill Charges histograms (MHGausEvents::FillHistAndArray()) with:
|
---|
683 | // - MExtractedSignalPix::GetExtractedSignalHiGain();
|
---|
684 | // - MExtractedSignalPix::GetExtractedSignalLoGain();
|
---|
685 | //
|
---|
686 | // 2) Set number of saturated slices (MHCalibrationChargePix::AddSaturated()) with:
|
---|
687 | // - MExtractedSignalPix::GetNumHiGainSaturated();
|
---|
688 | // - MExtractedSignalPix::GetNumLoGainSaturated();
|
---|
689 | //
|
---|
690 | // 3) Fill AbsTime histograms (MHCalibrationChargePix::FillAbsTime()) with:
|
---|
691 | // - MRawEvtPixelIter::GetIdxMaxHiGainSample();
|
---|
692 | // - MRawEvtPixelIter::GetIdxMaxLoGainSample(first slice);
|
---|
693 | //
|
---|
694 | Bool_t MHCalibrationChargeCam::FillHists(const MParContainer *par, const Stat_t w)
|
---|
695 | {
|
---|
696 |
|
---|
697 | MExtractedSignalCam *signal = (MExtractedSignalCam*)par;
|
---|
698 | if (!signal)
|
---|
699 | {
|
---|
700 | *fLog << err << "No argument in MExtractedSignalCam::Fill... abort." << endl;
|
---|
701 | return kFALSE;
|
---|
702 | }
|
---|
703 |
|
---|
704 | const UInt_t npixels = fGeom->GetNumPixels();
|
---|
705 | const UInt_t nareas = fGeom->GetNumAreas();
|
---|
706 | const UInt_t nsectors = fGeom->GetNumSectors();
|
---|
707 | const UInt_t lofirst = signal->GetFirstUsedSliceLoGain();
|
---|
708 |
|
---|
709 | fSumhiarea .Reset();
|
---|
710 | fSumloarea .Reset();
|
---|
711 | fTimehiarea .Reset();
|
---|
712 | fTimeloarea .Reset();
|
---|
713 | fSumhisector.Reset();
|
---|
714 | fSumlosector.Reset();
|
---|
715 | fTimehisector.Reset();
|
---|
716 | fTimelosector.Reset();
|
---|
717 |
|
---|
718 | fSathiarea .Reset();
|
---|
719 | fSatloarea .Reset();
|
---|
720 | fSathisector.Reset();
|
---|
721 | fSatlosector.Reset();
|
---|
722 |
|
---|
723 | for (UInt_t i=0; i<npixels; i++)
|
---|
724 | {
|
---|
725 |
|
---|
726 | MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
|
---|
727 |
|
---|
728 | if (histhi.IsExcluded())
|
---|
729 | continue;
|
---|
730 |
|
---|
731 | const MExtractedSignalPix &pix = (*signal)[i];
|
---|
732 |
|
---|
733 | const Float_t sumhi = pix.GetExtractedSignalHiGain();
|
---|
734 | const Int_t sathi = (Int_t)pix.GetNumHiGainSaturated();
|
---|
735 |
|
---|
736 | if (IsOscillations())
|
---|
737 | histhi.FillHistAndArray(sumhi);
|
---|
738 | else
|
---|
739 | histhi.FillHist(sumhi);
|
---|
740 |
|
---|
741 | histhi.AddSaturated(sathi);
|
---|
742 |
|
---|
743 | const Int_t aidx = (*fGeom)[i].GetAidx();
|
---|
744 | const Int_t sector = (*fGeom)[i].GetSector();
|
---|
745 |
|
---|
746 | fSumhiarea[aidx] += sumhi;
|
---|
747 | fSathiarea[aidx] += sathi;
|
---|
748 |
|
---|
749 | fSumhisector[sector] += sumhi;
|
---|
750 | fSathisector[sector] += sathi;
|
---|
751 |
|
---|
752 | if (IsLoGain())
|
---|
753 | {
|
---|
754 | MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
|
---|
755 | const Float_t sumlo = pix.GetExtractedSignalLoGain();
|
---|
756 | const Int_t satlo = (Int_t)pix.GetNumLoGainSaturated();
|
---|
757 |
|
---|
758 | if (IsOscillations())
|
---|
759 | histlo.FillHistAndArray(sumlo);
|
---|
760 | else
|
---|
761 | histlo.FillHist(sumlo);
|
---|
762 |
|
---|
763 | histlo.AddSaturated(satlo);
|
---|
764 |
|
---|
765 | fSumloarea[aidx] += sumlo;
|
---|
766 | fSatloarea[aidx] += satlo;
|
---|
767 | fSumlosector[sector] += sumlo;
|
---|
768 | fSatlosector[sector] += satlo;
|
---|
769 | }
|
---|
770 |
|
---|
771 | }
|
---|
772 |
|
---|
773 | MRawEvtPixelIter pixel(fRawEvt);
|
---|
774 | while (pixel.Next())
|
---|
775 | {
|
---|
776 |
|
---|
777 | const UInt_t pixid = pixel.GetPixelId();
|
---|
778 |
|
---|
779 | MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[pixid];
|
---|
780 |
|
---|
781 | if (histhi.IsExcluded())
|
---|
782 | continue;
|
---|
783 |
|
---|
784 | const Float_t timehi = (Float_t)pixel.GetIdxMaxHiGainSample();
|
---|
785 |
|
---|
786 | histhi.FillAbsTime(timehi);
|
---|
787 |
|
---|
788 | const Int_t aidx = (*fGeom)[pixid].GetAidx();
|
---|
789 | const Int_t sector = (*fGeom)[pixid].GetSector();
|
---|
790 |
|
---|
791 | fTimehiarea [aidx] += timehi;
|
---|
792 | fTimehisector[sector] += timehi;
|
---|
793 |
|
---|
794 | if (IsLoGain())
|
---|
795 | {
|
---|
796 | MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(pixid);
|
---|
797 |
|
---|
798 | const Float_t timelo = (Float_t)pixel.GetIdxMaxLoGainSample(lofirst);
|
---|
799 | histlo.FillAbsTime(timelo);
|
---|
800 |
|
---|
801 | fTimeloarea[aidx] += timelo;
|
---|
802 | fTimelosector[sector] += timelo;
|
---|
803 | }
|
---|
804 | }
|
---|
805 |
|
---|
806 | for (UInt_t j=0; j<nareas; j++)
|
---|
807 | {
|
---|
808 |
|
---|
809 | const Int_t npix = fAverageAreaNum[j];
|
---|
810 |
|
---|
811 | if (npix == 0)
|
---|
812 | continue;
|
---|
813 |
|
---|
814 | MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
|
---|
815 |
|
---|
816 |
|
---|
817 | if (IsOscillations())
|
---|
818 | hipix.FillHistAndArray(fSumhiarea [j]/npix);
|
---|
819 | else
|
---|
820 | hipix.FillHist(fSumhiarea[j]/npix);
|
---|
821 |
|
---|
822 | hipix.AddSaturated ((Float_t)fSathiarea [j]/npix > 0.5 ? 1 : 0);
|
---|
823 | hipix.FillAbsTime (fTimehiarea[j]/npix);
|
---|
824 |
|
---|
825 | if (IsLoGain())
|
---|
826 | {
|
---|
827 | MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
|
---|
828 | if (IsOscillations())
|
---|
829 | lopix.FillHistAndArray(fSumloarea [j]/npix);
|
---|
830 | else
|
---|
831 | lopix.FillHist(fSumloarea [j]/npix);
|
---|
832 | lopix.AddSaturated ((Float_t)fSatloarea [j]/npix > 0.5 ? 1 : 0);
|
---|
833 | lopix.FillAbsTime (fTimeloarea[j]/npix);
|
---|
834 | }
|
---|
835 | }
|
---|
836 |
|
---|
837 | for (UInt_t j=0; j<nsectors; j++)
|
---|
838 | {
|
---|
839 |
|
---|
840 | const Int_t npix = fAverageSectorNum[j];
|
---|
841 |
|
---|
842 | if (npix == 0)
|
---|
843 | continue;
|
---|
844 |
|
---|
845 | MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
|
---|
846 |
|
---|
847 | if (IsOscillations())
|
---|
848 | hipix.FillHistAndArray(fSumhisector [j]/npix);
|
---|
849 | else
|
---|
850 | hipix.FillHist(fSumhisector [j]/npix);
|
---|
851 |
|
---|
852 | hipix.AddSaturated ((Float_t)fSathisector[j]/npix > 0.5 ? 1 : 0);
|
---|
853 | hipix.FillAbsTime (fTimehisector[j]/npix);
|
---|
854 |
|
---|
855 | if (IsLoGain())
|
---|
856 | {
|
---|
857 | MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
|
---|
858 |
|
---|
859 | if (IsOscillations())
|
---|
860 | lopix.FillHistAndArray(fSumlosector [j]/npix);
|
---|
861 | else
|
---|
862 | lopix.FillHist(fSumlosector [j]/npix);
|
---|
863 |
|
---|
864 | lopix.AddSaturated ((Float_t)fSatlosector[j]/npix > 0.5 ? 1 : 0);
|
---|
865 | lopix.FillAbsTime (fTimelosector[j]/npix);
|
---|
866 | }
|
---|
867 | }
|
---|
868 |
|
---|
869 | return kTRUE;
|
---|
870 | }
|
---|
871 |
|
---|
872 | // --------------------------------------------------------------------------
|
---|
873 | //
|
---|
874 | // For all TOrdCollection's (including the averaged ones), the following steps are performed:
|
---|
875 | //
|
---|
876 | // 1) Returns if the pixel is excluded.
|
---|
877 | // 2) Tests saturation. In case yes, set the flag: MCalibrationPix::SetHiGainSaturation()
|
---|
878 | // or the flag: MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainSaturated )
|
---|
879 | // 3) Store the absolute arrival times in the MCalibrationChargePix's. If flag
|
---|
880 | // MCalibrationPix::IsHiGainSaturation() is set, the Low-Gain arrival times are stored,
|
---|
881 | // otherwise the Hi-Gain ones.
|
---|
882 | // 4) Calls to MHCalibrationCam::FitHiGainArrays() and MCalibrationCam::FitLoGainArrays()
|
---|
883 | // with the flags:
|
---|
884 | // - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainNotFitted )
|
---|
885 | // - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainNotFitted )
|
---|
886 | // - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainOscillating )
|
---|
887 | // - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainOscillating )
|
---|
888 | //
|
---|
889 | Bool_t MHCalibrationChargeCam::FinalizeHists()
|
---|
890 | {
|
---|
891 |
|
---|
892 | *fLog << endl;
|
---|
893 |
|
---|
894 | TH1F *h = NULL;
|
---|
895 |
|
---|
896 | MCalibrationCam *chargecam = fIntensCam ? fIntensCam->GetCam() : fCam;
|
---|
897 | MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
|
---|
898 |
|
---|
899 |
|
---|
900 | for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
|
---|
901 | {
|
---|
902 |
|
---|
903 | MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)(*this)[i];
|
---|
904 |
|
---|
905 | if (histhi.IsExcluded())
|
---|
906 | continue;
|
---|
907 |
|
---|
908 | MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i] ;
|
---|
909 |
|
---|
910 | if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
|
---|
911 | {
|
---|
912 | pix.SetHiGainSaturation();
|
---|
913 | if (IsOscillations())
|
---|
914 | histhi.CreateFourierSpectrum();
|
---|
915 | continue;
|
---|
916 | }
|
---|
917 |
|
---|
918 | MBadPixelsPix &bad = (*badcam)[i];
|
---|
919 |
|
---|
920 | h = histhi.GetHGausHist();
|
---|
921 |
|
---|
922 | Stat_t overflow = h->GetBinContent(h->GetNbinsX()+1);
|
---|
923 | if (overflow > 0.1)
|
---|
924 | {
|
---|
925 | *fLog << warn
|
---|
926 | << "HiGain Hist-overflow " << overflow
|
---|
927 | << " times in pix: " << i << " (w/o saturation!) " << endl;
|
---|
928 | bad.SetUncalibrated( MBadPixelsPix::kHiGainOverFlow );
|
---|
929 | }
|
---|
930 |
|
---|
931 | overflow = h->GetBinContent(0);
|
---|
932 | if (overflow > 0.1)
|
---|
933 | {
|
---|
934 | *fLog << warn
|
---|
935 | << "HiGain Hist-underflow " << overflow
|
---|
936 | << " times in pix: " << i << " (w/o saturation!) " << endl;
|
---|
937 | bad.SetUncalibrated( MBadPixelsPix::kHiGainOverFlow );
|
---|
938 | }
|
---|
939 |
|
---|
940 | FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
|
---|
941 | }
|
---|
942 |
|
---|
943 | if (IsLoGain())
|
---|
944 | for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
|
---|
945 | {
|
---|
946 |
|
---|
947 | MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)(*this)(i);
|
---|
948 | MBadPixelsPix &bad = (*badcam)[i];
|
---|
949 |
|
---|
950 | if (histlo.IsExcluded())
|
---|
951 | continue;
|
---|
952 |
|
---|
953 | if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
|
---|
954 | {
|
---|
955 | *fLog << warn << "Saturated Lo Gain histogram in pixel: " << i << endl;
|
---|
956 | bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
|
---|
957 | if (IsOscillations())
|
---|
958 | histlo.CreateFourierSpectrum();
|
---|
959 | continue;
|
---|
960 | }
|
---|
961 |
|
---|
962 | h = histlo.GetHGausHist();
|
---|
963 |
|
---|
964 | Stat_t overflow = h->GetBinContent(h->GetNbinsX()+1);
|
---|
965 | if (overflow > 0.1)
|
---|
966 | {
|
---|
967 | *fLog << warn
|
---|
968 | << "LoGain Hist-overflow " << overflow
|
---|
969 | << " times in pix: " << i << " (w/o saturation!) " << endl;
|
---|
970 | bad.SetUncalibrated( MBadPixelsPix::kLoGainOverFlow );
|
---|
971 | }
|
---|
972 |
|
---|
973 | overflow = h->GetBinContent(0);
|
---|
974 | if (overflow > 0.1)
|
---|
975 | {
|
---|
976 | *fLog << warn
|
---|
977 | << "LoGain Hist-underflow " << overflow
|
---|
978 | << " times in pix: " << i << " (w/o saturation!) " << endl;
|
---|
979 | bad.SetUncalibrated( MBadPixelsPix::kLoGainOverFlow );
|
---|
980 | }
|
---|
981 |
|
---|
982 | MCalibrationChargePix &pix = (MCalibrationChargePix&)(*chargecam)[i] ;
|
---|
983 |
|
---|
984 | if (pix.IsHiGainSaturation())
|
---|
985 | FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
|
---|
986 | }
|
---|
987 |
|
---|
988 | for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
|
---|
989 | {
|
---|
990 |
|
---|
991 | MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainArea(j);
|
---|
992 | MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageArea(j);
|
---|
993 |
|
---|
994 | if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
|
---|
995 | {
|
---|
996 | pix.SetHiGainSaturation();
|
---|
997 | if (IsOscillations())
|
---|
998 | histhi.CreateFourierSpectrum();
|
---|
999 | continue;
|
---|
1000 | }
|
---|
1001 |
|
---|
1002 | MBadPixelsPix &bad = chargecam->GetAverageBadArea(j);
|
---|
1003 | FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
|
---|
1004 | }
|
---|
1005 |
|
---|
1006 | if (IsLoGain())
|
---|
1007 | for (Int_t j=0; j<fAverageLoGainAreas->GetSize(); j++)
|
---|
1008 | {
|
---|
1009 |
|
---|
1010 | MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainArea(j);
|
---|
1011 |
|
---|
1012 | if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
|
---|
1013 | {
|
---|
1014 | *fLog << warn << "Saturated Lo Gain histogram in area idx: " << j << endl;
|
---|
1015 | histlo.CreateFourierSpectrum();
|
---|
1016 | continue;
|
---|
1017 | }
|
---|
1018 |
|
---|
1019 | MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageArea(j);
|
---|
1020 |
|
---|
1021 | if (pix.IsHiGainSaturation())
|
---|
1022 | {
|
---|
1023 | MBadPixelsPix &bad = chargecam->GetAverageBadArea(j);
|
---|
1024 | FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
|
---|
1025 | }
|
---|
1026 |
|
---|
1027 | }
|
---|
1028 |
|
---|
1029 | for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
|
---|
1030 | {
|
---|
1031 |
|
---|
1032 | MHCalibrationChargePix &histhi = (MHCalibrationChargePix&)GetAverageHiGainSector(j);
|
---|
1033 | MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageSector(j);
|
---|
1034 |
|
---|
1035 | if (histhi.GetSaturated() > fNumHiGainSaturationLimit*histhi.GetHGausHist()->GetEntries())
|
---|
1036 | {
|
---|
1037 | pix.SetHiGainSaturation();
|
---|
1038 | if (IsOscillations())
|
---|
1039 | histhi.CreateFourierSpectrum();
|
---|
1040 | continue;
|
---|
1041 | }
|
---|
1042 |
|
---|
1043 | MBadPixelsPix &bad = chargecam->GetAverageBadSector(j);
|
---|
1044 | FinalizeAbsTimes(histhi, pix, bad, fFirstHiGain, fLastHiGain);
|
---|
1045 | }
|
---|
1046 |
|
---|
1047 | if (IsLoGain())
|
---|
1048 | for (Int_t j=0; j<fAverageLoGainSectors->GetSize(); j++)
|
---|
1049 | {
|
---|
1050 |
|
---|
1051 | MHCalibrationChargePix &histlo = (MHCalibrationChargePix&)GetAverageLoGainSector(j);
|
---|
1052 | MBadPixelsPix &bad = chargecam->GetAverageBadSector(j);
|
---|
1053 |
|
---|
1054 | if (histlo.GetSaturated() > fNumLoGainSaturationLimit*histlo.GetHGausHist()->GetEntries())
|
---|
1055 | {
|
---|
1056 | *fLog << warn << "Saturated Lo Gain histogram in sector: " << j << endl;
|
---|
1057 | bad.SetUncalibrated( MBadPixelsPix::kLoGainSaturation );
|
---|
1058 | if (IsOscillations())
|
---|
1059 | histlo.CreateFourierSpectrum();
|
---|
1060 | continue;
|
---|
1061 | }
|
---|
1062 |
|
---|
1063 | MCalibrationChargePix &pix = (MCalibrationChargePix&)chargecam->GetAverageSector(j);
|
---|
1064 |
|
---|
1065 | if (pix.IsHiGainSaturation())
|
---|
1066 | FinalizeAbsTimes(histlo, pix, bad, fFirstLoGain, fLastLoGain);
|
---|
1067 | }
|
---|
1068 |
|
---|
1069 | //
|
---|
1070 | // Perform the fitting for the High Gain (done in MHCalibrationCam)
|
---|
1071 | //
|
---|
1072 | FitHiGainArrays(*chargecam, *badcam,
|
---|
1073 | MBadPixelsPix::kHiGainNotFitted,
|
---|
1074 | MBadPixelsPix::kHiGainOscillating);
|
---|
1075 |
|
---|
1076 | //
|
---|
1077 | // Perform the fitting for the Low Gain (done in MHCalibrationCam)
|
---|
1078 | //
|
---|
1079 | if (IsLoGain())
|
---|
1080 | FitLoGainArrays(*chargecam, *badcam,
|
---|
1081 | MBadPixelsPix::kLoGainNotFitted,
|
---|
1082 | MBadPixelsPix::kLoGainOscillating);
|
---|
1083 |
|
---|
1084 |
|
---|
1085 | return kTRUE;
|
---|
1086 | }
|
---|
1087 |
|
---|
1088 | // --------------------------------------------------------------------------------
|
---|
1089 | //
|
---|
1090 | // Fill the absolute time results into MCalibrationChargePix
|
---|
1091 | //
|
---|
1092 | // Check absolute time validity:
|
---|
1093 | // - Mean arrival time is at least fTimeLowerLimit slices from the lower edge
|
---|
1094 | // - Mean arrival time is at least fUpperLimit slices from the upper edge
|
---|
1095 | //
|
---|
1096 | void MHCalibrationChargeCam::FinalizeAbsTimes(MHCalibrationChargePix &hist, MCalibrationChargePix &pix, MBadPixelsPix &bad,
|
---|
1097 | Byte_t first, Byte_t last)
|
---|
1098 | {
|
---|
1099 |
|
---|
1100 | const Float_t mean = hist.GetAbsTimeMean();
|
---|
1101 | const Float_t rms = hist.GetAbsTimeRms();
|
---|
1102 |
|
---|
1103 | pix.SetAbsTimeMean ( mean );
|
---|
1104 | pix.SetAbsTimeRms ( rms );
|
---|
1105 |
|
---|
1106 | const Float_t lowerlimit = (Float_t)first + fTimeLowerLimit;
|
---|
1107 | const Float_t upperlimit = (Float_t)last - fTimeUpperLimit;
|
---|
1108 |
|
---|
1109 | if ( mean < lowerlimit)
|
---|
1110 | {
|
---|
1111 | *fLog << warn
|
---|
1112 | << Form("%s%3.1f%s%2.1f%s%2.0f%s%s","Mean ArrivalTime: ",mean," < ",fTimeLowerLimit,
|
---|
1113 | " slices from lower edge: ",lowerlimit," in pixel ",hist.GetName()) << endl;
|
---|
1114 | bad.SetUncalibrated( MBadPixelsPix::kMeanTimeInFirstBin );
|
---|
1115 | }
|
---|
1116 |
|
---|
1117 | if ( mean > upperlimit )
|
---|
1118 | {
|
---|
1119 | *fLog << warn
|
---|
1120 | << Form("%s%3.1f%s%2.1f%s%2.0f%s%s","Mean ArrivalTime: ",mean," > ",fTimeUpperLimit,
|
---|
1121 | " slices from upper edge: ",upperlimit," in pixel ",hist.GetName()) << endl;
|
---|
1122 | bad.SetUncalibrated( MBadPixelsPix::kMeanTimeInLast2Bins );
|
---|
1123 | }
|
---|
1124 | }
|
---|
1125 |
|
---|
1126 | // --------------------------------------------------------------------------
|
---|
1127 | //
|
---|
1128 | // Sets all pixels to MBadPixelsPix::kUnsuitableRun, if following flags are set:
|
---|
1129 | // - MBadPixelsPix::kLoGainSaturation
|
---|
1130 | //
|
---|
1131 | // Sets all pixels to MBadPixelsPix::kUnreliableRun, if following flags are set:
|
---|
1132 | // - if MBadPixelsPix::kHiGainNotFitted and !MCalibrationPix::IsHiGainSaturation()
|
---|
1133 | // - if MBadPixelsPix::kHiGainOscillating and !MCalibrationPix::IsHiGainSaturation()
|
---|
1134 | // - if MBadPixelsPix::kLoGainNotFitted and MCalibrationPix::IsLoGainSaturation()
|
---|
1135 | // - if MBadPixelsPix::kLoGainOscillating and MCalibrationPix::IsLoGainSaturation()
|
---|
1136 | //
|
---|
1137 | void MHCalibrationChargeCam::FinalizeBadPixels()
|
---|
1138 | {
|
---|
1139 |
|
---|
1140 | MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
|
---|
1141 | MCalibrationCam *chargecam = fIntensCam ? fIntensCam->GetCam() : fCam;
|
---|
1142 |
|
---|
1143 | for (Int_t i=0; i<badcam->GetSize(); i++)
|
---|
1144 | {
|
---|
1145 |
|
---|
1146 | MBadPixelsPix &bad = (*badcam)[i];
|
---|
1147 | MCalibrationPix &pix = (*chargecam)[i];
|
---|
1148 |
|
---|
1149 | if (bad.IsUncalibrated( MBadPixelsPix::kHiGainNotFitted ))
|
---|
1150 | if (!pix.IsHiGainSaturation())
|
---|
1151 | bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
|
---|
1152 |
|
---|
1153 | if (bad.IsUncalibrated( MBadPixelsPix::kLoGainNotFitted ))
|
---|
1154 | if (pix.IsHiGainSaturation())
|
---|
1155 | bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
|
---|
1156 |
|
---|
1157 | if (bad.IsUncalibrated( MBadPixelsPix::kLoGainSaturation ))
|
---|
1158 | bad.SetUnsuitable( MBadPixelsPix::kUnsuitableRun );
|
---|
1159 |
|
---|
1160 | if (IsOscillations())
|
---|
1161 | {
|
---|
1162 | if (bad.IsUncalibrated( MBadPixelsPix::kHiGainOscillating ))
|
---|
1163 | bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
|
---|
1164 |
|
---|
1165 | if (bad.IsUncalibrated( MBadPixelsPix::kLoGainOscillating ))
|
---|
1166 | if (pix.IsHiGainSaturation())
|
---|
1167 | bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
|
---|
1168 | }
|
---|
1169 | }
|
---|
1170 |
|
---|
1171 |
|
---|
1172 |
|
---|
1173 | }
|
---|
1174 |
|
---|
1175 | // --------------------------------------------------------------------------
|
---|
1176 | //
|
---|
1177 | // Calls MHCalibrationPix::DrawClone() for pixel idx
|
---|
1178 | //
|
---|
1179 | void MHCalibrationChargeCam::DrawPixelContent(Int_t idx) const
|
---|
1180 | {
|
---|
1181 | (*this)[idx].DrawClone();
|
---|
1182 | }
|
---|
1183 |
|
---|
1184 |
|
---|
1185 | // -----------------------------------------------------------------------------
|
---|
1186 | //
|
---|
1187 | // Default draw:
|
---|
1188 | //
|
---|
1189 | // Displays the averaged areas, both High Gain and Low Gain
|
---|
1190 | //
|
---|
1191 | // Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
|
---|
1192 | //
|
---|
1193 | void MHCalibrationChargeCam::Draw(const Option_t *opt)
|
---|
1194 | {
|
---|
1195 |
|
---|
1196 | const Int_t nareas = fAverageHiGainAreas->GetSize();
|
---|
1197 | if (nareas == 0)
|
---|
1198 | return;
|
---|
1199 |
|
---|
1200 | TString option(opt);
|
---|
1201 | option.ToLower();
|
---|
1202 |
|
---|
1203 | if (!option.Contains("datacheck"))
|
---|
1204 | {
|
---|
1205 | MHCalibrationCam::Draw(opt);
|
---|
1206 | return;
|
---|
1207 | }
|
---|
1208 |
|
---|
1209 | //
|
---|
1210 | // From here on , the datacheck - Draw
|
---|
1211 | //
|
---|
1212 | TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
|
---|
1213 | pad->SetBorderMode(0);
|
---|
1214 | pad->Divide(1,nareas);
|
---|
1215 |
|
---|
1216 | //
|
---|
1217 | // Loop over inner and outer pixels
|
---|
1218 | //
|
---|
1219 | for (Int_t i=0; i<nareas;i++)
|
---|
1220 | {
|
---|
1221 |
|
---|
1222 | pad->cd(i+1);
|
---|
1223 |
|
---|
1224 | MHCalibrationChargePix &hipix = (MHCalibrationChargePix&)GetAverageHiGainArea(i);
|
---|
1225 | //
|
---|
1226 | // Ask for Hi-Gain saturation
|
---|
1227 | //
|
---|
1228 | *fLog << err << hipix.GetSaturated() << " " << fNumHiGainSaturationLimit*hipix.GetHGausHist()->GetEntries() << endl;
|
---|
1229 | if (hipix.GetSaturated() > fNumHiGainSaturationLimit*hipix.GetHGausHist()->GetEntries() && IsLoGain())
|
---|
1230 | {
|
---|
1231 | MHCalibrationChargePix &lopix = (MHCalibrationChargePix&)GetAverageLoGainArea(i);
|
---|
1232 | DrawDataCheckPixel(lopix,i ? gkLoGainOuterRefLines : gkLoGainInnerRefLines);
|
---|
1233 | }
|
---|
1234 | else
|
---|
1235 | DrawDataCheckPixel(hipix,i ? gkHiGainOuterRefLines : gkHiGainInnerRefLines);
|
---|
1236 |
|
---|
1237 | }
|
---|
1238 | }
|
---|
1239 |
|
---|
1240 | // -----------------------------------------------------------------------------
|
---|
1241 | //
|
---|
1242 | // Draw the average pixel for the datacheck:
|
---|
1243 | //
|
---|
1244 | // Displays the averaged areas, both High Gain and Low Gain
|
---|
1245 | //
|
---|
1246 | // Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
|
---|
1247 | //
|
---|
1248 | void MHCalibrationChargeCam::DrawDataCheckPixel(MHCalibrationChargePix &pix, const Float_t refline[])
|
---|
1249 | {
|
---|
1250 |
|
---|
1251 | TVirtualPad *newpad = gPad;
|
---|
1252 | newpad->Divide(1,2);
|
---|
1253 | newpad->cd(1);
|
---|
1254 |
|
---|
1255 | gPad->SetTicks();
|
---|
1256 | if (!pix.IsEmpty() && !pix.IsOnlyOverflow() && !pix.IsOnlyUnderflow())
|
---|
1257 | gPad->SetLogy();
|
---|
1258 |
|
---|
1259 | TH1F *hist = pix.GetHGausHist();
|
---|
1260 |
|
---|
1261 | TH1F *null = new TH1F("Null",hist->GetTitle(),100,
|
---|
1262 | pix.GetFirst() > 0. ? pix.GetFirst() : 0.,
|
---|
1263 | pix.GetLast() > pix.GetFirst()
|
---|
1264 | ? ( pix.GetLast() > 450.
|
---|
1265 | ? ( fColor == MCalibrationCam::kBLUE ? 800. : 450. )
|
---|
1266 | : pix.GetLast() )
|
---|
1267 | : pix.GetFirst()*2.);
|
---|
1268 |
|
---|
1269 | null->SetMaximum(1.1*hist->GetMaximum());
|
---|
1270 | null->SetDirectory(NULL);
|
---|
1271 | null->SetBit(kCanDelete);
|
---|
1272 | null->SetStats(kFALSE);
|
---|
1273 | //
|
---|
1274 | // set the labels bigger
|
---|
1275 | //
|
---|
1276 | TAxis *xaxe = null->GetXaxis();
|
---|
1277 | TAxis *yaxe = null->GetYaxis();
|
---|
1278 | xaxe->CenterTitle();
|
---|
1279 | yaxe->CenterTitle();
|
---|
1280 | xaxe->SetTitleSize(0.07);
|
---|
1281 | yaxe->SetTitleSize(0.07);
|
---|
1282 | xaxe->SetTitleOffset(0.7);
|
---|
1283 | yaxe->SetTitleOffset(0.55);
|
---|
1284 | xaxe->SetLabelSize(0.06);
|
---|
1285 | yaxe->SetLabelSize(0.06);
|
---|
1286 | xaxe->SetTitle(hist->GetXaxis()->GetTitle());
|
---|
1287 | yaxe->SetTitle(hist->GetYaxis()->GetTitle());
|
---|
1288 |
|
---|
1289 | null->Draw();
|
---|
1290 | hist->Draw("same");
|
---|
1291 |
|
---|
1292 | gStyle->SetOptFit();
|
---|
1293 |
|
---|
1294 | TF1 *fit = pix.GetFGausFit();
|
---|
1295 |
|
---|
1296 | if (fit)
|
---|
1297 | {
|
---|
1298 | switch ( fColor )
|
---|
1299 | {
|
---|
1300 | case MCalibrationCam::kGREEN:
|
---|
1301 | fit->SetLineColor(kGreen);
|
---|
1302 | break;
|
---|
1303 | case MCalibrationCam::kBLUE:
|
---|
1304 | fit->SetLineColor(kBlue);
|
---|
1305 | break;
|
---|
1306 | case MCalibrationCam::kUV:
|
---|
1307 | fit->SetLineColor(106);
|
---|
1308 | break;
|
---|
1309 | case MCalibrationCam::kCT1:
|
---|
1310 | fit->SetLineColor(006);
|
---|
1311 | break;
|
---|
1312 | default:
|
---|
1313 | fit->SetLineColor(kRed);
|
---|
1314 | }
|
---|
1315 | fit->Draw("same");
|
---|
1316 | }
|
---|
1317 |
|
---|
1318 | DisplayRefLines(null,refline);
|
---|
1319 |
|
---|
1320 | newpad->cd(2);
|
---|
1321 | gPad->SetTicks();
|
---|
1322 |
|
---|
1323 | TH1F *null2 = new TH1F("Null2",hist->GetTitle(),100,0.,pix.GetEvents()->GetSize()/pix.GetEventFrequency());
|
---|
1324 |
|
---|
1325 | null2->SetMinimum(pix.GetMean()-10.*pix.GetSigma());
|
---|
1326 | null2->SetMaximum(pix.GetMean()+10.*pix.GetSigma());
|
---|
1327 | null2->SetDirectory(NULL);
|
---|
1328 | null2->SetBit(kCanDelete);
|
---|
1329 | null2->SetStats(kFALSE);
|
---|
1330 | //
|
---|
1331 | // set the labels bigger
|
---|
1332 | //
|
---|
1333 | TAxis *xaxe2 = null2->GetXaxis();
|
---|
1334 | TAxis *yaxe2 = null2->GetYaxis();
|
---|
1335 | xaxe2->CenterTitle();
|
---|
1336 | yaxe2->CenterTitle();
|
---|
1337 | xaxe2->SetTitleSize(0.07);
|
---|
1338 | yaxe2->SetTitleSize(0.07);
|
---|
1339 | xaxe2->SetTitleOffset(0.7);
|
---|
1340 | yaxe2->SetTitleOffset(0.55);
|
---|
1341 | xaxe2->SetLabelSize(0.06);
|
---|
1342 | yaxe2->SetLabelSize(0.06);
|
---|
1343 |
|
---|
1344 | pix.CreateGraphEvents();
|
---|
1345 | TGraph *gr = pix.GetGraphEvents();
|
---|
1346 |
|
---|
1347 | xaxe2->SetTitle(gr->GetXaxis()->GetTitle());
|
---|
1348 | yaxe2->SetTitle(gr->GetYaxis()->GetTitle());
|
---|
1349 |
|
---|
1350 | null2->Draw();
|
---|
1351 |
|
---|
1352 | pix.DrawEvents("same");
|
---|
1353 |
|
---|
1354 | // newpad->cd(3);
|
---|
1355 | // pix.DrawPowerSpectrum(*newpad,4);
|
---|
1356 |
|
---|
1357 | return;
|
---|
1358 |
|
---|
1359 | }
|
---|
1360 |
|
---|
1361 |
|
---|
1362 | void MHCalibrationChargeCam::DisplayRefLines(const TH1F *hist, const Float_t refline[]) const
|
---|
1363 | {
|
---|
1364 |
|
---|
1365 | TGraph *green1 = new TGraph(2);
|
---|
1366 | green1->SetPoint(0,refline[0],0.1);
|
---|
1367 | green1->SetPoint(1,refline[0],hist->GetMaximum());
|
---|
1368 | green1->SetBit(kCanDelete);
|
---|
1369 | green1->SetLineColor(kGreen);
|
---|
1370 | green1->SetLineStyle(2);
|
---|
1371 | green1->SetLineWidth(3);
|
---|
1372 | green1->Draw("L");
|
---|
1373 |
|
---|
1374 | TGraph *green5 = new TGraph(2);
|
---|
1375 | green5->SetPoint(0,refline[6],0.1);
|
---|
1376 | green5->SetPoint(1,refline[6],hist->GetMaximum());
|
---|
1377 | green5->SetBit(kCanDelete);
|
---|
1378 | green5->SetLineColor(8);
|
---|
1379 | green5->SetLineStyle(2);
|
---|
1380 | green5->SetLineWidth(3);
|
---|
1381 | green5->Draw("L");
|
---|
1382 |
|
---|
1383 | TGraph *blue1 = new TGraph(2);
|
---|
1384 | blue1->SetPoint(0,refline[1],0.1);
|
---|
1385 | blue1->SetPoint(1,refline[1],hist->GetMaximum());
|
---|
1386 | blue1->SetBit(kCanDelete);
|
---|
1387 | blue1->SetLineColor(227);
|
---|
1388 | blue1->SetLineStyle(2);
|
---|
1389 | blue1->SetLineWidth(3);
|
---|
1390 | blue1->Draw("L");
|
---|
1391 |
|
---|
1392 | TGraph *blue5 = new TGraph(2);
|
---|
1393 | blue5->SetPoint(0,refline[2],0.1);
|
---|
1394 | blue5->SetPoint(1,refline[2],hist->GetMaximum());
|
---|
1395 | blue5->SetBit(kCanDelete);
|
---|
1396 | blue5->SetLineColor(68);
|
---|
1397 | blue5->SetLineStyle(2);
|
---|
1398 | blue5->SetLineWidth(3);
|
---|
1399 | blue5->Draw("L");
|
---|
1400 |
|
---|
1401 | TGraph *blue10 = new TGraph(2);
|
---|
1402 | blue10->SetPoint(0,refline[3],0.1);
|
---|
1403 | blue10->SetPoint(1,refline[3],hist->GetMaximum());
|
---|
1404 | blue10->SetBit(kCanDelete);
|
---|
1405 | blue10->SetLineColor(4);
|
---|
1406 | blue10->SetLineStyle(2);
|
---|
1407 | blue10->SetLineWidth(3);
|
---|
1408 | blue10->Draw("L");
|
---|
1409 |
|
---|
1410 | TGraph *uv10 = new TGraph(2);
|
---|
1411 | uv10->SetPoint(0,refline[4],0.1);
|
---|
1412 | uv10->SetPoint(1,refline[4],hist->GetMaximum());
|
---|
1413 | uv10->SetBit(kCanDelete);
|
---|
1414 | uv10->SetLineColor(106);
|
---|
1415 | uv10->SetLineStyle(2);
|
---|
1416 | uv10->SetLineWidth(3);
|
---|
1417 | uv10->Draw("L");
|
---|
1418 |
|
---|
1419 | TGraph *ct1 = new TGraph(2);
|
---|
1420 | ct1->SetPoint(0,refline[5],0.1);
|
---|
1421 | ct1->SetPoint(1,refline[5],hist->GetMaximum());
|
---|
1422 | ct1->SetBit(kCanDelete);
|
---|
1423 | ct1->SetLineColor(6);
|
---|
1424 | ct1->SetLineStyle(2);
|
---|
1425 | ct1->SetLineWidth(3);
|
---|
1426 | ct1->Draw("L");
|
---|
1427 |
|
---|
1428 | TLegend *leg = new TLegend(0.8,0.35,0.99,0.99);
|
---|
1429 | leg->SetBit(kCanDelete);
|
---|
1430 | leg->AddEntry(green1,"1 Led GREEN","l");
|
---|
1431 | leg->AddEntry(green5,"5 Leds GREEN","l");
|
---|
1432 | leg->AddEntry(blue1,"1 Led BLUE","l");
|
---|
1433 | leg->AddEntry(blue5,"5 Leds BLUE","l");
|
---|
1434 | leg->AddEntry(blue10,"10 Leds BLUE","l");
|
---|
1435 | leg->AddEntry(uv10,"10 Leds UV","l");
|
---|
1436 | leg->AddEntry(ct1,"CT1-Pulser","l");
|
---|
1437 |
|
---|
1438 | leg->Draw();
|
---|
1439 | }
|
---|
1440 |
|
---|
1441 | Int_t MHCalibrationChargeCam::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
|
---|
1442 | {
|
---|
1443 |
|
---|
1444 | Bool_t rc = kFALSE;
|
---|
1445 |
|
---|
1446 | if (MHCalibrationCam::ReadEnv(env,prefix,print))
|
---|
1447 | rc = kTRUE;
|
---|
1448 |
|
---|
1449 | if (IsEnvDefined(env, prefix, "LoGainNbins", print))
|
---|
1450 | {
|
---|
1451 | SetLoGainNbins(GetEnvValue(env, prefix, "LoGainNbins", fLoGainNbins));
|
---|
1452 | rc = kTRUE;
|
---|
1453 | }
|
---|
1454 |
|
---|
1455 | if (IsEnvDefined(env, prefix, "LoGainFirst", print))
|
---|
1456 | {
|
---|
1457 | SetLoGainFirst(GetEnvValue(env, prefix, "LoGainFirst", fLoGainFirst));
|
---|
1458 | rc = kTRUE;
|
---|
1459 | }
|
---|
1460 |
|
---|
1461 | if (IsEnvDefined(env, prefix, "LoGainLast", print))
|
---|
1462 | {
|
---|
1463 | SetLoGainLast(GetEnvValue(env, prefix, "LoGainLast", fLoGainLast));
|
---|
1464 | rc = kTRUE;
|
---|
1465 | }
|
---|
1466 |
|
---|
1467 | if (IsEnvDefined(env, prefix, "TimeLowerLimit", print))
|
---|
1468 | {
|
---|
1469 | SetTimeLowerLimit(GetEnvValue(env, prefix, "TimeLowerLimit", fTimeLowerLimit));
|
---|
1470 | rc = kTRUE;
|
---|
1471 | }
|
---|
1472 |
|
---|
1473 | if (IsEnvDefined(env, prefix, "TimeUpperLimit", print))
|
---|
1474 | {
|
---|
1475 | SetTimeUpperLimit(GetEnvValue(env, prefix, "TimeUpperLimit", fTimeUpperLimit));
|
---|
1476 | rc = kTRUE;
|
---|
1477 | }
|
---|
1478 |
|
---|
1479 |
|
---|
1480 | return rc;
|
---|
1481 | }
|
---|