source: trunk/MagicSoft/Mars/mhcalib/MHCalibrationHiLoCam.cc@ 7454

Last change on this file since 7454 was 7189, checked in by tbretz, 19 years ago
*** empty log message ***
File size: 32.7 KB
Line 
1
2/* ======================================================================== *\
3!
4! *
5! * This file is part of MARS, the MAGIC Analysis and Reconstruction
6! * Software. It is distributed to you in the hope that it can be a useful
7! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
8! * It is distributed WITHOUT ANY WARRANTY.
9! *
10! * Permission to use, copy, modify and distribute this software and its
11! * documentation for any purpose is hereby granted without fee,
12! * provided that the above copyright notice appear in all copies and
13! * that both that copyright notice and this permission notice appear
14! * in supporting documentation. It is provided "as is" without express
15! * or implied warranty.
16! *
17!
18!
19! Author(s): Markus Gaug 02/2004 <mailto:markus@ifae.es>
20!
21! Copyright: MAGIC Software Development, 2000-2004
22!
23!
24\* ======================================================================== */
25/////////////////////////////////////////////////////////////////////////////
26//
27// MHCalibrationHiLoCam
28//
29// Fills the extracted high-gain low-gain charge ratios of MArrivalTimeCam into
30// the MHCalibrationPix-classes MHCalibrationPix for every:
31//
32// - Pixel, stored in the TObjArray's MHCalibrationCam::fHiGainArray
33// or MHCalibrationCam::fHiGainArray, respectively, depending if
34// MArrivalTimePix::IsLoGainUsed() is set.
35//
36// - Average pixel per AREA index (e.g. inner and outer for the MAGIC camera),
37// stored in the TObjArray's MHCalibrationCam::fAverageHiGainAreas and
38// MHCalibrationCam::fAverageHiGainAreas
39//
40// - Average pixel per camera SECTOR (e.g. sectors 1-6 for the MAGIC camera),
41// stored in the TObjArray's MHCalibrationCam::fAverageHiGainSectors
42// and MHCalibrationCam::fAverageHiGainSectors
43//
44// The histograms are fitted to a Gaussian, mean and sigma with its errors
45// and the fit probability are extracted. If none of these values are NaN's and
46// if the probability is bigger than MHGausEvents::fProbLimit (default: 0.5%),
47// the fit is declared valid.
48// Otherwise, the fit is repeated within ranges of the previous mean
49// +- MHCalibrationPix::fPickupLimit (default: 5) sigma (see MHCalibrationPix::RepeatFit())
50// In case this does not make the fit valid, the histogram means and RMS's are
51// taken directly (see MHCalibrationPix::BypassFit()) and the following flags are set:
52// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiLoNotFitted ) and
53// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
54//
55// Outliers of more than MHCalibrationPix::fPickupLimit (default: 5) sigmas
56// from the mean are counted as Pickup events (stored in MHCalibrationPix::fPickup)
57//
58// The class also fills arrays with the signal vs. event number, creates a fourier
59// spectrum (see MHGausEvents::CreateFourierSpectrum()) and investigates if the
60// projected fourier components follow an exponential distribution.
61// In case that the probability of the exponential fit is less than
62// MHGausEvents::fProbLimit (default: 0.5%), the following flags are set:
63// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiLoOscillating ) and
64// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
65//
66// This same procedure is performed for the average pixels.
67//
68// The following results are written into MCalibrationHiLoCam:
69//
70// - MCalibrationPix::SetMean()
71// - MCalibrationPix::SetMeanErr()
72// - MCalibrationPix::SetSigma()
73// - MCalibrationPix::SetSigmaErr()
74// - MCalibrationPix::SetProb()
75// - MCalibrationPix::SetNumPickup()
76//
77// For all averaged areas, the fitted sigma is multiplied with the square root of
78// the number involved pixels in order to be able to compare it to the average of
79// sigmas in the camera.
80//
81/////////////////////////////////////////////////////////////////////////////
82#include "MHCalibrationHiLoCam.h"
83
84#include <TOrdCollection.h>
85#include <TPad.h>
86#include <TVirtualPad.h>
87#include <TCanvas.h>
88#include <TStyle.h>
89#include <TF1.h>
90#include <TLine.h>
91#include <TLatex.h>
92#include <TLegend.h>
93#include <TGraph.h>
94#include <TProfile.h>
95
96#include "MHCalibrationHiLoPix.h"
97
98#include "MLog.h"
99#include "MLogManip.h"
100
101#include "MParList.h"
102
103#include "MCalibrationHiLoCam.h"
104#include "MCalibrationHiLoPix.h"
105#include "MCalibrationCam.h"
106#include "MCalibrationIntensityCam.h"
107#include "MCalibrationPix.h"
108
109#include "MExtractedSignalCam.h"
110#include "MExtractedSignalPix.h"
111#include "MArrivalTimeCam.h"
112#include "MArrivalTimePix.h"
113
114#include "MGeomCam.h"
115#include "MGeomPix.h"
116
117#include "MBadPixelsIntensityCam.h"
118#include "MBadPixelsCam.h"
119#include "MBadPixelsPix.h"
120
121ClassImp(MHCalibrationHiLoCam);
122
123using namespace std;
124
125const Int_t MHCalibrationHiLoCam::fgNbins = 175;
126const Axis_t MHCalibrationHiLoCam::fgFirst = -5.1;
127const Axis_t MHCalibrationHiLoCam::fgLast = 29.9;
128const Float_t MHCalibrationHiLoCam::fgProbLimit = 0.;
129const Int_t MHCalibrationHiLoCam::fgHivsLoNbins = 90;
130const Axis_t MHCalibrationHiLoCam::fgHivsLoFirst = 95.;
131const Axis_t MHCalibrationHiLoCam::fgHivsLoLast = 995.;
132const Axis_t MHCalibrationHiLoCam::fgLowerFitLimitProfile = 480.;
133const Axis_t MHCalibrationHiLoCam::fgUpperFitLimitProfile = 680.;
134const TString MHCalibrationHiLoCam::gsHistName = "HiLo";
135const TString MHCalibrationHiLoCam::gsHistTitle = "HiGain vs. LoGain";
136const TString MHCalibrationHiLoCam::gsHistXTitle = "Amplification Ratio [1]";
137const TString MHCalibrationHiLoCam::gsHistYTitle = "Nr. events";
138const TString MHCalibrationHiLoCam::gsHivsLoHistName = "HivsLo";
139const TString MHCalibrationHiLoCam::gsHivsLoHistTitle = "High-gain vs. Low-gain Charge";
140const TString MHCalibrationHiLoCam::gsHivsLoHistXTitle = "Q High-Gain [FADC counts]";
141const TString MHCalibrationHiLoCam::gsHivsLoHistYTitle = "Q Low-Gain [FADC counts]";
142
143// --------------------------------------------------------------------------
144//
145// Default Constructor.
146//
147// Sets:
148// - fNbins to fgNbins
149// - fFirst to fgFirst
150// - fLast to fgLast
151//
152// - fHistName to gsHistName
153// - fHistTitle to gsHistTitle
154// - fHistXTitle to gsHistXTitle
155// - fHistYTitle to gsHistYTitle
156//
157// - fLowerLimt to fgLowerLim
158// - fUpperLimt to fgUpperLim
159//
160MHCalibrationHiLoCam::MHCalibrationHiLoCam(const char *name, const char *title)
161 : fArrTimes(NULL), fHivsLoResults("Results","Fit Results high-gain vs. low-gain",
162 200,-10.,10.,200,0.,20.),
163 fUsedLoGainSlices(0)
164{
165
166 fName = name ? name : "MHCalibrationHiLoCam";
167 fTitle = title ? title : "Histogram class for the high-gain vs. low-gain amplification ratio calibration";
168
169 SetNbins(fgNbins);
170 SetFirst(fgFirst);
171 SetLast (fgLast );
172
173 SetProbLimit(fgProbLimit);
174
175 SetHistName (gsHistName .Data());
176 SetHistTitle (gsHistTitle .Data());
177 SetHistXTitle(gsHistXTitle.Data());
178 SetHistYTitle(gsHistYTitle.Data());
179
180 SetHivsLoNbins(fgHivsLoNbins);
181 SetHivsLoFirst(fgHivsLoFirst);
182 SetHivsLoLast (fgHivsLoLast );
183
184 SetLowerFitLimitProfile();
185 SetUpperFitLimitProfile();
186
187 SetHivsLoHistName (gsHivsLoHistName .Data());
188 SetHivsLoHistTitle (gsHivsLoHistTitle .Data());
189 SetHivsLoHistXTitle(gsHivsLoHistXTitle.Data());
190 SetHivsLoHistYTitle(gsHivsLoHistYTitle.Data());
191
192 SetOscillations(kFALSE);
193
194 fHivsLoResults.GetXaxis()->SetTitle("Offset per FADC slices [FADC cnts]");
195 fHivsLoResults.GetYaxis()->SetTitle("Gains ratio [1]");
196 fHivsLoResults.SetDirectory(0);
197
198}
199
200// --------------------------------------------------------------------------
201//
202// Creates new MHCalibrationHiLoCam only with the averaged areas:
203// the rest has to be retrieved directly, e.g. via:
204// MHCalibrationHiLoCam *cam = MParList::FindObject("MHCalibrationHiLoCam");
205// - cam->GetAverageSector(5).DrawClone();
206// - (*cam)[100].DrawClone()
207//
208TObject *MHCalibrationHiLoCam::Clone(const char *) const
209{
210
211 MHCalibrationHiLoCam *cam = new MHCalibrationHiLoCam();
212
213 //
214 // Copy the data members
215 //
216 cam->fColor = fColor;
217 cam->fRunNumbers = fRunNumbers;
218 cam->fPulserFrequency = fPulserFrequency;
219 cam->fFlags = fFlags;
220 cam->fNbins = fNbins;
221 cam->fFirst = fFirst;
222 cam->fLast = fLast;
223
224 //
225 // Copy the MArrays
226 //
227 cam->fAverageAreaRelSigma = fAverageAreaRelSigma;
228 cam->fAverageAreaRelSigmaVar = fAverageAreaRelSigmaVar;
229 cam->fAverageAreaSat = fAverageAreaSat;
230 cam->fAverageAreaSigma = fAverageAreaSigma;
231 cam->fAverageAreaSigmaVar = fAverageAreaSigmaVar;
232 cam->fAverageAreaNum = fAverageAreaNum;
233 cam->fAverageSectorNum = fAverageSectorNum;
234
235 if (!IsAverageing())
236 return cam;
237
238 const Int_t navhi = fAverageHiGainAreas->GetSize();
239 const Int_t navlo = fAverageLoGainAreas->GetSize();
240
241 for (int i=0; i<navhi; i++)
242 cam->fAverageHiGainAreas->AddAt(GetAverageHiGainArea(i).Clone(),i);
243
244 for (int i=0; i<navlo; i++)
245 cam->fAverageLoGainAreas->AddAt(GetAverageLoGainArea(i).Clone(),i);
246
247 return cam;
248}
249
250// --------------------------------------------------------------------------
251//
252// Gets or creates the pointers to:
253// - MCalibrationHiLoCam
254//
255// Searches pointer to:
256// - MExtractedSignalCam
257// - MArrivalTimeCam
258//
259// Calls:
260// - MHCalibrationCam::InitHiGainArrays()
261// - MHCalibrationCam::InitLoGainArrays()
262//
263// Sets:
264// - fSumarea to nareas
265// - fSumsector to nareas
266// - fNumarea to nareas
267// - fNumsector to nareas
268//
269Bool_t MHCalibrationHiLoCam::ReInitHists(MParList *pList)
270{
271
272 fCam = (MCalibrationCam*)pList->FindObject(AddSerialNumber("MCalibrationHiLoCam"));
273 if (!fCam)
274 {
275 fCam = (MCalibrationCam*)pList->FindCreateObj(AddSerialNumber("MCalibrationHiLoCam"));
276 if (!fCam)
277 return kFALSE;
278 fCam->Init(*fGeom);
279 }
280
281 MExtractedSignalCam *signal = (MExtractedSignalCam*)pList->FindObject("MExtractedSignalCam");
282 if (!signal)
283 {
284 *fLog << err << "MExtractedSignalCam not found... abort." << endl;
285 return kFALSE;
286 }
287
288 fUsedLoGainSlices = signal->GetNumUsedLoGainFADCSlices();
289
290 fArrTimes = (MArrivalTimeCam*)pList->FindObject("MArrivalTimeCam");
291 if (!fArrTimes)
292 {
293 *fLog << warn << "MArrivalTimeCam not found... cannot calibrated arrival times between "
294 <<"high and low-gain" << endl;
295 SetLoGain(kFALSE);
296 }
297
298 const Int_t npixels = fGeom->GetNumPixels();
299 const Int_t nsectors = fGeom->GetNumSectors();
300 const Int_t nareas = fGeom->GetNumAreas();
301
302 InitHiGainArrays(npixels,nareas,nsectors);
303 InitLoGainArrays(npixels,nareas,nsectors);
304
305 fSumareahi .Set(nareas);
306 fSumsectorhi.Set(nsectors);
307 fNumareahi .Set(nareas);
308 fNumsectorhi.Set(nsectors);
309 if (IsLoGain())
310 {
311 fSumarealo .Set(nareas);
312 fSumsectorlo.Set(nsectors);
313 fNumarealo .Set(nareas);
314 fNumsectorlo.Set(nsectors);
315 }
316 return kTRUE;
317}
318
319// --------------------------------------------------------------------------
320//
321// Retrieve:
322// - fRunHeader->GetNumSamplesHiGain();
323//
324// Initializes the High Gain Arrays:
325//
326// - For every entry in the expanded arrays:
327// * Initialize an MHCalibrationHiLoPix
328// * Set Binning from fNbins, fFirst and fLast
329// * Set Binning of Abs Times histogram from fAbsNbins, fAbsFirst and fAbsLast
330// * Set Histgram names and titles from fHistName and fHistTitle
331// * Set Abs Times Histgram names and titles from fAbsHistName and fAbsHistTitle
332// * Set X-axis and Y-axis titles from fHistXTitle and fHistYTitle
333// * Set X-axis and Y-axis titles of Abs Times Histogram from fAbsHistXTitle and fAbsHistYTitle
334// * Call InitHists
335//
336//
337void MHCalibrationHiLoCam::InitHiGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
338{
339
340 TProfile *h;
341
342 if (fHiGainArray->GetSize()==0)
343 {
344 for (Int_t i=0; i<npixels; i++)
345 {
346 fHiGainArray->AddAt(new MHCalibrationHiLoPix(Form("%sHiGainPix%04d",fHistName.Data(),i),
347 Form("%s High Gain Pixel%04d",fHistTitle.Data(),i)),i);
348
349 MHCalibrationHiLoPix &pix = (MHCalibrationHiLoPix&)(*this)[i];
350
351 pix.SetNbins(fNbins);
352 pix.SetFirst(fFirst);
353 pix.SetLast (fLast);
354
355 pix.SetProbLimit(fProbLimit);
356
357 pix.SetHivsLoNbins(fHivsLoNbins);
358 pix.SetHivsLoFirst(fHivsLoFirst);
359 pix.SetHivsLoLast (fHivsLoLast);
360
361 InitHists(pix,(*fBadPixels)[i],i);
362
363 if (fCam)
364 (*fCam)[i].SetPixId(i);
365
366 h = pix.GetHivsLo();
367
368 h->SetName (Form("H%sHiGainPix%04d",fHivsLoHistName.Data(),i));
369 h->SetTitle(Form("%s High Gain Pixel %04d",fHivsLoHistTitle.Data(),i));
370 h->SetXTitle(fHivsLoHistXTitle.Data());
371 h->SetYTitle(fHivsLoHistYTitle.Data());
372 h->SetDirectory(0);
373 }
374 }
375
376
377 if (fAverageHiGainAreas->GetSize()==0)
378 {
379 for (Int_t j=0; j<nareas; j++)
380 {
381 fAverageHiGainAreas->AddAt(new MHCalibrationHiLoPix(Form("%sHiGainArea%d",fHistName.Data(),j),
382 Form("%s High Gain Area Idx %d",fHistTitle.Data(),j)),j);
383
384 MHCalibrationHiLoPix &pix = (MHCalibrationHiLoPix&)GetAverageHiGainArea(j);
385
386 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
387 pix.SetFirst(fFirst);
388 pix.SetLast (fLast);
389
390 pix.SetHivsLoNbins(fHivsLoNbins);
391 pix.SetHivsLoFirst(fHivsLoFirst);
392 pix.SetHivsLoLast (fHivsLoLast);
393
394 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadArea(j) : fCam->GetAverageBadArea(j),j);
395
396 if (fCam)
397 fCam->GetAverageArea(j).SetPixId(j);
398
399 h = pix.GetHivsLo();
400
401 h->SetName (Form("H%sHiGainArea%d",fHivsLoHistName.Data(),j));
402 h->SetTitle(Form("%s averaged on event-by-event basis High Gain Area Idx %d",
403 fHivsLoHistTitle.Data(), j));
404 h->SetXTitle(fHivsLoHistXTitle.Data());
405 h->SetYTitle(fHivsLoHistYTitle.Data());
406 h->SetDirectory(0);
407 }
408 }
409
410 if (fAverageHiGainSectors->GetSize()==0)
411 {
412 for (Int_t j=0; j<nsectors; j++)
413 {
414 fAverageHiGainSectors->AddAt(new MHCalibrationHiLoPix(Form("%sHiGainSector%02d",fHistName.Data(),j),
415 Form("%s High Gain Sector %02d",fHistTitle.Data(),j)),j);
416
417 MHCalibrationHiLoPix &pix = (MHCalibrationHiLoPix&)GetAverageHiGainSector(j);
418
419 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
420 pix.SetFirst(fFirst);
421 pix.SetLast (fLast);
422
423 pix.SetHivsLoNbins(fHivsLoNbins);
424 pix.SetHivsLoFirst(fHivsLoFirst);
425 pix.SetHivsLoLast (fHivsLoLast);
426
427 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadSector(j) : fCam->GetAverageBadSector(j),j);
428
429 if (fCam)
430 fCam->GetAverageSector(j).SetPixId(j);
431
432 h = pix.GetHivsLo();
433
434 h->SetName (Form("H%sHiGainSector%02d",fHivsLoHistName.Data(),j));
435 h->SetTitle(Form("%s averaged on event-by-event basis High Gain Area Sector %02d",
436 fHivsLoHistTitle.Data(),j));
437 h->SetXTitle(fHivsLoHistXTitle.Data());
438 h->SetYTitle(fHivsLoHistYTitle.Data());
439 h->SetDirectory(0);
440 }
441 }
442}
443
444//--------------------------------------------------------------------------------------
445//
446// Return, if IsLoGain() is kFALSE
447//
448// Retrieve:
449// - fRunHeader->GetNumSamplesHiGain();
450//
451// Initializes the Low Gain Arrays:
452//
453// - For every entry in the expanded arrays:
454// * Initialize an MHCalibrationHiLoPix
455// * Set Binning from fNbins, fFirst and fLast
456// * Set Binning of HivsLo Times histogram from fHivsLoNbins, fHivsLoFirst and fHivsLoLast
457// * Set Histgram names and titles from fHistName and fHistTitle
458// * Set HivsLo Times Histgram names and titles from fHivsLoHistName and fHivsLoHistTitle
459// * Set X-axis and Y-axis titles from fHistXTitle and fHistYTitle
460// * Set X-axis and Y-axis titles of HivsLo Times Histogram from fHivsLoHistXTitle and fHivsLoHistYTitle
461// * Call InitHists
462//
463void MHCalibrationHiLoCam::InitLoGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
464{
465 if (!IsLoGain())
466 return;
467
468 TProfile *h;
469
470 if (fLoGainArray->GetSize()==0 )
471 {
472 for (Int_t i=0; i<npixels; i++)
473 {
474 fLoGainArray->AddAt(new MHCalibrationHiLoPix(Form("%sLoGainPix%04d",fHistName.Data(),i),
475 Form("%s Low Gain Pixel %04d",fHistTitle.Data(),i)),i);
476
477 MHCalibrationHiLoPix &pix = (MHCalibrationHiLoPix&)(*this)(i);
478
479 pix.SetNbins(fNbins);
480 pix.SetFirst(fFirst);
481 pix.SetLast (fLast);
482
483 pix.SetProbLimit(fProbLimit);
484
485 pix.SetHivsLoNbins(fHivsLoNbins);
486 pix.SetHivsLoFirst(fHivsLoFirst);
487 pix.SetHivsLoLast (fHivsLoLast );
488
489 InitHists(pix,(*fBadPixels)[i],i);
490
491 h = pix.GetHivsLo();
492
493 h->SetName (Form("H%sLoGainPix%04d",fHivsLoHistName.Data(),i));
494 h->SetTitle(Form("%s Low Gain Pixel %04d",fHivsLoHistTitle.Data(),i));
495 h->SetXTitle(fHivsLoHistXTitle.Data());
496 h->SetYTitle(fHivsLoHistYTitle.Data());
497 h->SetDirectory(0);
498 }
499 }
500
501 if (fAverageLoGainAreas->GetSize()==0)
502 {
503 for (Int_t j=0; j<nareas; j++)
504 {
505 fAverageLoGainAreas->AddAt(new MHCalibrationHiLoPix(Form("%sLoGainArea%d",fHistName.Data(),j),
506 Form("%s Low Gain Area Idx %d",fHistTitle.Data(),j)),j);
507
508 MHCalibrationHiLoPix &pix = (MHCalibrationHiLoPix&)GetAverageLoGainArea(j);
509
510 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
511 pix.SetFirst(fFirst);
512 pix.SetLast (fLast);
513
514 pix.SetHivsLoNbins(fHivsLoNbins);
515 pix.SetHivsLoFirst(fHivsLoFirst);
516 pix.SetHivsLoLast (fHivsLoLast );
517
518 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadArea(j) : fCam->GetAverageBadArea(j),j);
519
520 h = pix.GetHivsLo();
521
522 h->SetName (Form("H%sLoGainArea%02d",fHivsLoHistName.Data(),j));
523 h->SetTitle(Form("%s%s%02d",fHivsLoHistTitle.Data(),
524 " averaged on event-by-event basis Low Gain Area Idx ",j));
525 h->SetXTitle(fHivsLoHistXTitle.Data());
526 h->SetYTitle(fHivsLoHistYTitle.Data());
527 h->SetDirectory(0);
528 }
529 }
530
531
532 if (fAverageLoGainSectors->GetSize()==0 && IsLoGain())
533 {
534 for (Int_t j=0; j<nsectors; j++)
535 {
536 fAverageLoGainSectors->AddAt(new MHCalibrationHiLoPix(Form("%sLoGainSector%02d",fHistName.Data(),j),
537 Form("%s Low Gain Sector %02d",fHistTitle.Data(),j)),j);
538
539 MHCalibrationHiLoPix &pix = (MHCalibrationHiLoPix&)GetAverageLoGainSector(j);
540
541 pix.SetNbins(fNbins*(Int_t)TMath::Sqrt((Float_t)npixels/nareas));
542 pix.SetFirst(fFirst);
543 pix.SetLast (fLast);
544
545 pix.SetHivsLoNbins(fHivsLoNbins);
546 pix.SetHivsLoFirst(fHivsLoFirst);
547 pix.SetHivsLoLast (fHivsLoLast);
548
549 InitHists(pix,fIntensCam ? fIntensCam->GetAverageBadSector(j) : fCam->GetAverageBadSector(j),j);
550
551 h = pix.GetHivsLo();
552
553 h->SetName (Form("H%sLoGainSector%02d",fHivsLoHistName.Data(),j));
554 h->SetTitle(Form("%s%s%02d",fHivsLoHistTitle.Data(),
555 " averaged on event-by-event basis Low Gain Area Sector ",j));
556 h->SetXTitle(fHivsLoHistXTitle.Data());
557 h->SetYTitle(fHivsLoHistYTitle.Data());
558 h->SetDirectory(0);
559 }
560 }
561}
562
563
564
565// -------------------------------------------------------------------------------
566//
567// Retrieves pointer to MExtractedSignalCam:
568//
569// Retrieves from MGeomCam:
570// - number of pixels
571// - number of pixel areas
572// - number of sectors
573//
574// Fills histograms (MHGausEvents::FillHistAndArray()) with:
575// - MExtractedSignalPix::GetExtractedSignalHiGain(pixid) / MExtractedSignalPix::GetExtractedSignalLoGain;
576// if the high-gain signal does not show high-gain saturation, but the low-gain
577// has been extracted.
578// - MArrivalTimePix::GetArrivalTimeHiGain(pixid) / MArrivalTimePix::GetArrivalTimeLoGain;
579// if the high-gain signal does not show high-gain saturation, but the low-gain
580// has been extracted.
581//
582Bool_t MHCalibrationHiLoCam::FillHists(const MParContainer *par, const Stat_t w)
583{
584
585 MExtractedSignalCam *signal = (MExtractedSignalCam*)par;
586 if (!signal)
587 {
588 gLog << err << "No argument in MExtractedSignal::Fill... abort." << endl;
589 return kFALSE;
590 }
591
592 const Int_t npixels = fGeom->GetNumPixels();
593 const Int_t nareas = fGeom->GetNumAreas();
594 const Int_t nsectors = fGeom->GetNumSectors();
595
596 fSumareahi .Reset();
597 fSumsectorhi.Reset();
598 fNumareahi .Reset();
599 fNumsectorhi.Reset();
600 fSumarealo .Reset();
601 fSumsectorlo.Reset();
602 fNumarealo .Reset();
603 fNumsectorlo.Reset();
604
605 for (Int_t i=0; i<npixels; i++)
606 {
607 const MExtractedSignalPix &pix = (*signal)[i];
608 const Int_t aidx = (*fGeom)[i].GetAidx();
609 const Int_t sector = (*fGeom)[i].GetSector();
610
611 const Float_t siglo = pix.GetExtractedSignalLoGain();
612
613 //
614 // Skip all pixels with:
615 // - Saturated high-gain
616 // - Not extracted low-gain
617 // (see MExtractTimeAndCharge::fLoGainSwitch for setting the criteria)
618 //
619 if (siglo < 0.5 || pix.GetNumHiGainSaturated() > 0)
620 continue;
621
622 const Float_t sighi = pix.GetExtractedSignalHiGain();
623
624 // *fLog << err << sighi << " " << siglo << endl;
625 const Float_t ratio = sighi / siglo;
626
627 MHCalibrationHiLoPix &histhi = (MHCalibrationHiLoPix&)(*this)[i];
628
629 histhi.FillHist(ratio);
630 histhi.FillHivsLo(sighi,siglo);
631
632 if (IsAverageing())
633 {
634 MHCalibrationHiLoPix &histhi = (MHCalibrationHiLoPix&)GetAverageHiGainArea(aidx);
635 histhi.FillHivsLo(sighi,siglo);
636 }
637
638 fSumareahi [aidx] += ratio;
639 fNumareahi [aidx] ++;
640 fSumsectorhi[sector] += ratio;
641 fNumsectorhi[sector] ++;
642
643 if (IsLoGain())
644 {
645 const MArrivalTimePix &tix = (*fArrTimes)[i];
646 MHCalibrationPix &histlo = (*this)(i);
647
648 const Float_t diff = tix.GetArrivalTimeLoGain() - tix.GetArrivalTimeHiGain();
649
650 histlo.FillHist(diff);
651 fSumarealo [aidx] += diff;
652 fNumarealo [aidx] ++;
653 fSumsectorlo[sector] += diff;
654 fNumsectorlo[sector] ++;
655 }
656 }
657
658 if (!IsAverageing())
659 return kTRUE;
660
661 for (Int_t j=0; j<nareas; j++)
662 {
663
664 MHCalibrationHiLoPix &histhi = (MHCalibrationHiLoPix&)GetAverageHiGainArea(j);
665
666 if (IsOscillations())
667 histhi.FillHistAndArray(fNumareahi[j] == 0 ? 0. : fSumareahi[j]/fNumareahi[j]);
668 else
669 histhi.FillHist(fNumareahi[j] == 0 ? 0. : fSumareahi[j]/fNumareahi[j]);
670
671 if (IsLoGain())
672 {
673 MHCalibrationPix &histlo = GetAverageLoGainArea(j);
674 if (IsOscillations())
675 histlo.FillHistAndArray(fNumarealo[j] == 0 ? 0. : fSumarealo[j]/fNumarealo[j]);
676 else
677 histlo.FillHist(fNumarealo[j] == 0 ? 0. : fSumarealo[j]/fNumarealo[j]);
678 }
679 }
680
681 for (Int_t j=0; j<nsectors; j++)
682 {
683 MHCalibrationPix &hist = GetAverageHiGainSector(j);
684
685 if (IsOscillations())
686 hist.FillHistAndArray(fNumsectorhi[j] == 0 ? 0. : fSumsectorhi[j]/fNumsectorhi[j]);
687 else
688 hist.FillHist(fNumsectorhi[j] == 0 ? 0. : fSumsectorhi[j]/fNumsectorhi[j]);
689
690 if (IsLoGain())
691 {
692
693 MHCalibrationPix &histlo = GetAverageLoGainSector(j);
694
695 if (IsOscillations())
696 histlo.FillHistAndArray(fNumsectorlo[j] == 0 ? 0. : fSumsectorlo[j]/fNumsectorlo[j]);
697 else
698 histlo.FillHist(fNumsectorlo[j] == 0 ? 0. : fSumsectorlo[j]/fNumsectorlo[j]);
699 }
700 }
701
702 return kTRUE;
703}
704
705// --------------------------------------------------------------------------
706//
707// Calls:
708// - MHCalibrationCam::FitHiGainArrays() with flags:
709// MBadPixelsPix::kHiLoNotFitted and MBadPixelsPix::kHiLoOscillating
710// - MHCalibrationCam::FitLoGainArrays() with flags:
711// MBadPixelsPix::kHiLoNotFitted and MBadPixelsPix::kHiLoOscillating
712//
713Bool_t MHCalibrationHiLoCam::FinalizeHists()
714{
715
716 *fLog << endl;
717
718 MCalibrationCam *hilocam = fCam;
719 MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
720
721 const Int_t nareas = fAverageHiGainAreas->GetSize();
722
723 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
724 {
725
726 MHCalibrationHiLoPix &hist = (MHCalibrationHiLoPix&)(*this)[i];
727
728 if (hist.IsExcluded())
729 continue;
730
731 CheckOverflow(hist);
732
733 TProfile *h = hist.GetHivsLo();
734 h->Fit("pol1","RQ","",fLowerFitLimitProfile,fUpperFitLimitProfile);
735
736 TF1 *fit = h->GetFunction("pol1");
737
738 const Float_t gainr = fit->GetParameter(1) > 0.001
739 ? 1./fit->GetParameter(1)
740 : 0.;
741 const Float_t offset = fit->GetParameter(0)/fUsedLoGainSlices;
742
743 fHivsLoResults.Fill(offset,gainr);
744
745 MCalibrationHiLoPix &pix = (MCalibrationHiLoPix&)(*fCam)[i];
746 pix.SetOffsetPerSlice(offset);
747 pix.SetGainRatio (gainr );
748
749 }
750
751 //
752 // Check histogram overflow
753 //
754 if (IsAverageing())
755 {
756 for (Int_t j=0; j<nareas; j++)
757 CheckOverflow(GetAverageHiGainArea(j));
758
759 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
760 CheckOverflow(GetAverageHiGainSector(j));
761 }
762
763
764 FitHiGainArrays(*hilocam,*badcam,
765 MBadPixelsPix::kHiLoNotFitted,
766 MBadPixelsPix::kHiLoOscillating);
767
768 if (!IsLoGain())
769 return kTRUE;
770
771 for (Int_t i=0; i<fLoGainArray->GetSize(); i++)
772 {
773
774 MHCalibrationPix &hist = (*this)(i);
775
776 if (hist.IsExcluded())
777 continue;
778
779 CheckOverflow(hist);
780 }
781
782 if (IsAverageing())
783 {
784 for (Int_t j=0; j<nareas; j++)
785 {
786
787 MHCalibrationHiLoPix &hist = (MHCalibrationHiLoPix&)GetAverageHiGainArea(j);
788 //
789 // Check histogram overflow
790 //
791 CheckOverflow(hist);
792
793 TProfile *h = hist.GetHivsLo();
794 h->Fit("pol1","RQ","",fLowerFitLimitProfile,fUpperFitLimitProfile);
795
796 TF1 *fit = h->GetFunction("pol1");
797
798 const Float_t gainr = fit->GetParameter(1) > 0.001
799 ? 1./fit->GetParameter(1)
800 : 0.;
801 const Float_t offset = fit->GetParameter(0)/fUsedLoGainSlices;
802
803 MCalibrationHiLoPix &pix = (MCalibrationHiLoPix&)fCam->GetAverageArea(0);
804 pix.SetOffsetPerSlice(offset);
805 pix.SetGainRatio (gainr );
806
807 }
808
809 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
810 {
811
812 MHCalibrationHiLoPix &hist = (MHCalibrationHiLoPix&)GetAverageHiGainSector(j);
813 //
814 // Check histogram overflow
815 //
816 CheckOverflow(hist);
817
818 TProfile *h = hist.GetHivsLo();
819 h->Fit("pol1","RQ","",fLowerFitLimitProfile,fUpperFitLimitProfile);
820
821 TF1 *fit = h->GetFunction("pol1");
822
823 const Float_t gainr = fit->GetParameter(1) > 0.001
824 ? 1./fit->GetParameter(1)
825 : 0.;
826 const Float_t offset = fit->GetParameter(0)/fUsedLoGainSlices;
827
828 MCalibrationHiLoPix &pix = (MCalibrationHiLoPix&)fCam->GetAverageSector(0);
829 pix.SetOffsetPerSlice(offset);
830 pix.SetGainRatio (gainr );
831
832 }
833 }
834
835 FitLoGainArrays(*hilocam,*badcam,
836 MBadPixelsPix::kHiLoNotFitted,
837 MBadPixelsPix::kHiLoOscillating);
838
839 return kTRUE;
840}
841
842// --------------------------------------------------------------------------
843//
844// Sets all pixels to MBadPixelsPix::kUnreliableRun, if following flags are set:
845// - MBadPixelsPix::kHiLoNotFitted
846// - MBadPixelsPix::kHiLoOscillating
847//
848void MHCalibrationHiLoCam::FinalizeBadPixels()
849{
850
851 MBadPixelsCam *badcam = fIntensBad ? fIntensBad->GetCam() : fBadPixels;
852
853 for (Int_t i=0; i<badcam->GetSize(); i++)
854 {
855 MBadPixelsPix &bad = (*badcam)[i];
856
857 if (bad.IsUncalibrated( MBadPixelsPix::kHiLoNotFitted ))
858 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
859
860 if (bad.IsUncalibrated( MBadPixelsPix::kHiLoOscillating))
861 bad.SetUnsuitable( MBadPixelsPix::kUnreliableRun );
862
863 }
864}
865
866// --------------------------------------------------------------------------
867//
868// The types are as follows:
869//
870// Fitted values:
871// ==============
872//
873// 0: Fitted Mean High-Gain Low-Gain Charge Ratio in FADC slices (MHGausEvents::GetMean()
874// 1: Error Mean High-Gain Low-Gain Charge Ratio in FADC slices (MHGausEvents::GetMeanErr()
875// 2: Sigma fitted High-Gain Low-Gain Charge Ratio in FADC slices (MHGausEvents::GetSigma()
876// 3: Error Sigma High-Gain Low-Gain Charge Ratio in FADC slices (MHGausEvents::GetSigmaErr()
877//
878// Useful variables derived from the fit results:
879// =============================================
880//
881// 4: Returned probability of Gauss fit (calls: MHGausEvents::GetProb())
882//
883// Localized defects:
884// ==================
885//
886// 5: Gaus fit not OK (calls: MHGausEvents::IsGausFitOK())
887// 6: Fourier spectrum not OK (calls: MHGausEvents::IsFourierSpectrumOK())
888//
889Bool_t MHCalibrationHiLoCam::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const
890{
891
892 if (fHiGainArray->GetSize() <= idx)
893 return kFALSE;
894
895 const MHCalibrationPix &pixhi = (*this)[idx];
896 const MHCalibrationPix &pixlo = (*this)(idx);
897
898 switch (type)
899 {
900 case 0:
901 val = pixhi.GetMean();
902 break;
903 case 1:
904 val = pixhi.GetMeanErr();
905 break;
906 case 2:
907 val = pixhi.GetSigma();
908 break;
909 case 3:
910 val = pixhi.GetSigmaErr();
911 break;
912 case 4:
913 val = pixhi.GetProb();
914 break;
915 case 5:
916 if (!pixhi.IsGausFitOK())
917 val = 1.;
918 break;
919 case 6:
920 if (!pixhi.IsFourierSpectrumOK())
921 val = 1.;
922 break;
923 case 7:
924 if (!IsLoGain())
925 break;
926 val = pixlo.GetMean();
927 break;
928 case 8:
929 if (!IsLoGain())
930 break;
931 val = pixlo.GetMeanErr();
932 break;
933 case 9:
934 if (!IsLoGain())
935 break;
936 val = pixlo.GetSigma();
937 break;
938 case 10:
939 if (!IsLoGain())
940 break;
941 val = pixlo.GetSigmaErr();
942 break;
943 case 11:
944 if (!IsLoGain())
945 break;
946 val = pixlo.GetProb();
947 break;
948 case 12:
949 if (!IsLoGain())
950 break;
951 if (!pixlo.IsGausFitOK())
952 val = 1.;
953 break;
954 case 13:
955 if (!IsLoGain())
956 break;
957 if (!pixlo.IsFourierSpectrumOK())
958 val = 1.;
959 break;
960 default:
961 return kFALSE;
962 }
963 return kTRUE;
964}
965
966// --------------------------------------------------------------------------
967//
968// Calls MHCalibrationPix::DrawClone() for pixel idx
969//
970void MHCalibrationHiLoCam::DrawPixelContent(Int_t idx) const
971{
972 (*this)[idx].DrawClone();
973}
974
975void MHCalibrationHiLoCam::CheckOverflow( MHCalibrationPix &pix )
976{
977
978 if (pix.IsExcluded())
979 return;
980
981 TH1F *hist = pix.GetHGausHist();
982
983 Stat_t overflow = hist->GetBinContent(hist->GetNbinsX()+1);
984 if (overflow > fOverflowLimit*hist->GetEntries())
985 {
986 *fLog << warn << "Hist-overflow " << overflow
987 << " times in " << pix.GetName() << endl;
988 }
989
990 overflow = hist->GetBinContent(0);
991 if (overflow > fOverflowLimit*hist->GetEntries())
992 {
993 *fLog << warn << "Hist-underflow " << overflow
994 << " times in " << pix.GetName() << endl;
995 }
996}
997
998// -----------------------------------------------------------------------------
999//
1000// Default draw:
1001//
1002// Displays the averaged areas, both amplification ratio as time difference
1003//
1004void MHCalibrationHiLoCam::Draw(const Option_t *opt)
1005{
1006
1007 if (!IsAverageing())
1008 return;
1009
1010 const Int_t nareas = fAverageHiGainAreas->GetSize();
1011 if (nareas == 0)
1012 return;
1013
1014 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
1015 pad->SetBorderMode(0);
1016
1017 pad->Divide(IsLoGain() ? 2 : 1,nareas);
1018
1019 for (Int_t i=0; i<nareas;i++)
1020 {
1021
1022 pad->cd(IsLoGain() ? 2*(i+1)-1 : i+1);
1023
1024 GetAverageHiGainArea(i).Draw(opt);
1025
1026 if (IsLoGain())
1027 {
1028 pad->cd(2*(i+1));
1029
1030 TH1F *hist = GetAverageLoGainArea(i).GetHGausHist();
1031 hist->SetXTitle("Extracted Time Difference [FADC sl.]");
1032 GetAverageLoGainArea(i).Draw(opt);
1033 }
1034 }
1035}
1036
1037Int_t MHCalibrationHiLoCam::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
1038{
1039
1040 Bool_t rc = kFALSE;
1041
1042 if (MHCalibrationCam::ReadEnv(env,prefix,print))
1043 rc = kTRUE;
1044
1045 if (IsEnvDefined(env, prefix, "LowerFitLimitProfile", print))
1046 {
1047 SetLowerFitLimitProfile(GetEnvValue(env, prefix, "LowerFitLimitProfile", fLowerFitLimitProfile));
1048 rc = kTRUE;
1049 }
1050
1051 if (IsEnvDefined(env, prefix, "UpperFitLimitProfile", print))
1052 {
1053 SetUpperFitLimitProfile(GetEnvValue(env, prefix, "UpperFitLimitProfile", fUpperFitLimitProfile));
1054 rc = kTRUE;
1055 }
1056
1057 if (IsEnvDefined(env, prefix, "HivsLoNbins", print))
1058 {
1059 SetHivsLoNbins(GetEnvValue(env, prefix, "HivsLoNbins", fHivsLoNbins));
1060 rc = kTRUE;
1061 }
1062
1063 if (IsEnvDefined(env, prefix, "HivsLoFirst", print))
1064 {
1065 SetHivsLoFirst(GetEnvValue(env, prefix, "HivsLoFirst", fHivsLoFirst));
1066 rc = kTRUE;
1067 }
1068
1069 if (IsEnvDefined(env, prefix, "HivsLoLast", print))
1070 {
1071 SetHivsLoLast(GetEnvValue(env, prefix, "HivsLoLast", fHivsLoLast));
1072 rc = kTRUE;
1073 }
1074
1075 return rc;
1076}
Note: See TracBrowser for help on using the repository browser.