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

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