source: trunk/MagicSoft/Mars/mhcalib/MHCalibrationPulseTimeCam.cc@ 8550

Last change on this file since 8550 was 8519, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 21.1 KB
Line 
1/* ======================================================================== *\
2! $Name: not supported by cvs2svn $:$Id: MHCalibrationPulseTimeCam.cc,v 1.39 2007-05-16 13:56:17 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! Author(s): Thomas Bretz <mailto:tbretz@astro.uni-wuerzburg.de>
22!
23! Copyright: MAGIC Software Development, 2000-2007
24!
25!
26\* ======================================================================== */
27/////////////////////////////////////////////////////////////////////////////
28//
29// MHCalibrationPulseTimeCam
30//
31// Fills the extracted signals of MExtractedSignalCam into the MHCalibrationPix-classes
32// MHCalibrationPulseTimeHiGainPix and MHCalibrationPulseTimeLoGainPix for every:
33//
34// - Pixel, stored in the TOrdCollection's MHCalibrationCam::fHiGainArray and
35// MHCalibrationCam::fLoGainArray
36//
37// - Average pixel per AREA index (e.g. inner and outer for the MAGIC camera),
38// stored in the TOrdCollection's MHCalibrationCam::fAverageHiGainAreas and
39// MHCalibrationCam::fAverageLoGainAreas
40//
41// - Average pixel per camera SECTOR (e.g. sectors 1-6 for the MAGIC camera),
42// stored in the TOrdCollection's MHCalibrationCam::fAverageHiGainSectors and
43// MHCalibrationCam::fAverageLoGainSectors
44//
45// Every signal is taken from MExtractedSignalCam and filled into a histogram and
46// an array, in order to perform a Fourier analysis (see MHGausEvents).
47// The signals are moreover averaged on an event-by-event basis and written into
48// the corresponding average pixels.
49//
50// Additionally, the (FADC slice) position of the maximum is stored in an Absolute
51// Arrival Time histogram. This histogram serves for a rough cross-check if the
52// signal does not lie at or outside the edges of the extraction window.
53//
54// The PulseTime histograms are fitted to a Gaussian, mean and sigma with its errors
55// and the fit probability are extracted. If none of these values are NaN's and
56// if the probability is bigger than MHGausEvents::fProbLimit (default: 0.5%),
57// the fit is declared valid.
58// Otherwise, the fit is repeated within ranges of the previous mean
59// +- MHCalibrationPix::fPickupLimit (default: 5) sigma (see MHCalibrationPix::RepeatFit())
60// In case this does not make the fit valid, the histogram means and RMS's are
61// taken directly (see MHCalibrationPix::BypassFit()) and the following flags are set:
62// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainNotFitted ) or
63// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainNotFitted ) and
64// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
65//
66// Outliers of more than MHCalibrationPix::fPickupLimit (default: 5) sigmas
67// from the mean are counted as Pickup events (stored in MHCalibrationPix::fPickup)
68//
69// Unless more than fNumHiGainSaturationLimit (default: 1%) of the overall FADC
70// slices show saturation, the following flag is set:
71// - MCalibrationPulseTimePix::SetHiGainSaturation();
72// In that case, the calibration constants are derived from the low-gain results.
73//
74// If more than fNumLoGainSaturationLimit (default: 1%) of the overall
75// low-gain FADC slices saturate, the following flags are set:
76// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainSaturation ) and
77// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnsuitableRun )
78//
79// The class also fills arrays with the signal vs. event number, creates a fourier
80// spectrum and investigates if the projected fourier components follow an exponential
81// distribution. In case that the probability of the exponential fit is less than
82// MHGausEvents::fProbLimit (default: 0.5%), the following flags are set:
83// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainOscillating ) or
84// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainOscillating ) and
85// - MBadPixelsPix::SetUnsuitable( MBadPixelsPix::kUnreliableRun )
86//
87// This same procedure is performed for the average pixels.
88//
89// The following results are written into MCalibrationPulseTimeCam:
90//
91// - MCalibrationPix::SetHiGainSaturation()
92// - MCalibrationPix::SetHiGainMean()
93// - MCalibrationPix::SetHiGainMeanErr()
94// - MCalibrationPix::SetHiGainSigma()
95// - MCalibrationPix::SetHiGainSigmaErr()
96// - MCalibrationPix::SetHiGainProb()
97// - MCalibrationPix::SetHiGainNumPickup()
98//
99// For all averaged areas, the fitted sigma is multiplied with the square root of
100// the number involved pixels in order to be able to compare it to the average of
101// sigmas in the camera.
102//
103/////////////////////////////////////////////////////////////////////////////
104#include "MHCalibrationPulseTimeCam.h"
105
106#include <TEnv.h>
107#include <TLine.h>
108#include <TGraph.h>
109#include <TLegend.h>
110#include <TCanvas.h>
111#include <TOrdCollection.h>
112
113#include "MLog.h"
114#include "MLogManip.h"
115
116#include "MParList.h"
117
118#include "MArrayD.h"
119#include "MRawRunHeader.h"
120#include "MExtractedSignalCam.h"
121#include "MPedestalSubtractedEvt.h"
122
123#include "MCalibrationPix.h"
124#include "MHCalibrationPix.h"
125
126#include "MGeomCam.h"
127#include "MGeomPix.h"
128
129#include "MBadPixelsCam.h"
130
131ClassImp(MHCalibrationPulseTimeCam);
132
133using namespace std;
134
135const Byte_t MHCalibrationPulseTimeCam::fgSaturationLimit = 245;
136const Byte_t MHCalibrationPulseTimeCam::fgLowerSignalLimit = 85;
137const Int_t MHCalibrationPulseTimeCam::fgNumPixelsRequired = 2;
138const Int_t MHCalibrationPulseTimeCam::fgHiGainNbins = 20;
139const Axis_t MHCalibrationPulseTimeCam::fgHiGainFirst = -0.5;
140const Axis_t MHCalibrationPulseTimeCam::fgHiGainLast = 19.5;
141const Float_t MHCalibrationPulseTimeCam::fgProbLimit = 0.001;
142const TString MHCalibrationPulseTimeCam::gsHistName = "PulseTime";
143const TString MHCalibrationPulseTimeCam::gsHistTitle = "Extracted Times";
144const TString MHCalibrationPulseTimeCam::gsHistXTitle = "Time [FADC slices]";
145const TString MHCalibrationPulseTimeCam::gsHistYTitle = "Nr. events";
146const TString MHCalibrationPulseTimeCam::fgReferenceFile = "mjobs/signalref.rc";
147
148// --------------------------------------------------------------------------
149//
150// Default Constructor.
151//
152// Sets:
153// - all pointers to NULL
154//
155// - fNbins to fgHiGainNbins
156// - fFirst to fgHiGainFirst
157// - fLast to fgHiGainLast
158//
159// - fHistName to gsHistName
160// - fHistTitle to gsHistTitle
161// - fHistXTitle to gsHistXTitle
162// - fHistYTitle to gsHistYTitle
163//
164// - fSaturationLimit to fgSaturationLimit
165// - fLowerSignalLimit to fgLowerSignalLimit
166// - fNumPixelsRequired to fgNumPixelsRequired
167//
168MHCalibrationPulseTimeCam::MHCalibrationPulseTimeCam(const char *name, const char *title)
169 : fBadPixels(NULL)
170{
171
172 fName = name ? name : "MHCalibrationPulseTimeCam";
173 fTitle = title ? title : "Class to fill the extracted pulse times for cosmics ";
174
175 SetBinning(fgHiGainNbins, fgHiGainFirst, fgHiGainLast);
176
177 SetProbLimit(fgProbLimit);
178
179 SetHistName (gsHistName .Data());
180 SetHistTitle (gsHistTitle .Data());
181 SetHistXTitle(gsHistXTitle.Data());
182 SetHistYTitle(gsHistYTitle.Data());
183
184 SetReferenceFile();
185 SetLoGain(kFALSE);
186 SetOscillations(kFALSE);
187
188 SetSaturationLimit();
189 SetLowerSignalLimit();
190 SetNumPixelsRequired();
191
192 fInnerRefTime = 5.;
193 fOuterRefTime = 5.;
194}
195// --------------------------------------------------------------------------
196//
197// Creates new MHCalibrationPulseTimeCam only with the averaged areas:
198// the rest has to be retrieved directly, e.g. via:
199// MHCalibrationPulseTimeCam *cam = MParList::FindObject("MHCalibrationPulseTimeCam");
200// - cam->GetAverageSector(5).DrawClone();
201// - (*cam)[100].DrawClone()
202//
203TObject *MHCalibrationPulseTimeCam::Clone(const char *) const
204{
205
206 MHCalibrationPulseTimeCam *cam = new MHCalibrationPulseTimeCam();
207
208 //
209 // Copy the data members
210 //
211 cam->fColor = fColor;
212 cam->fRunNumbers = fRunNumbers;
213 cam->fPulserFrequency = fPulserFrequency;
214 cam->fFlags = fFlags;
215 cam->fNbins = fNbins;
216 cam->fFirst = fFirst;
217 cam->fLast = fLast;
218 cam->fReferenceFile = fReferenceFile;
219 cam->fInnerRefTime = fInnerRefTime;
220 cam->fOuterRefTime = fOuterRefTime;
221
222 if (!IsAverageing())
223 return cam;
224
225 const Int_t navhi = fAverageHiGainAreas->GetSize();
226
227 for (int i=0; i<navhi; i++)
228 cam->fAverageHiGainAreas->AddAt(GetAverageHiGainArea(i).Clone(),i);
229
230 return cam;
231}
232
233// --------------------------------------------------------------------------
234//
235// Gets the pointers to:
236// - MRawEvtData
237//
238Bool_t MHCalibrationPulseTimeCam::SetupHists(const MParList *pList)
239{
240
241 fBadPixels = (MBadPixelsCam*)pList->FindObject("MBadPixelsCam");
242 if (!fBadPixels)
243 {
244 *fLog << warn << GetDescriptor() << "MBadPixelsCam not found... " << endl;
245 }
246
247 return kTRUE;
248}
249
250// --------------------------------------------------------------------------
251//
252// Gets or creates the pointers to:
253// - MExtractedSignalCam
254// - MCalibrationPulseTimeCam
255// - MBadPixelsCam
256//
257// Initializes the number of used FADC slices from MExtractedSignalCam
258// into MCalibrationPulseTimeCam and test for changes in that variable
259//
260// Calls:
261// - InitHiGainArrays()
262//
263// Sets:
264// - fSumhiarea to nareas
265// - fSumloarea to nareas
266// - fSumhisector to nsectors
267// - fSumlosector to nsectors
268//
269Bool_t MHCalibrationPulseTimeCam::ReInitHists(MParList *pList)
270{
271
272 MExtractedSignalCam *signal =
273 (MExtractedSignalCam*)pList->FindObject(AddSerialNumber("MExtractedSignalCam"));
274 if (!signal)
275 {
276 *fLog << err << "MExtractedSignalCam not found... abort." << endl;
277 return kFALSE;
278 }
279
280 if (!InitCams(pList,"PulseTime"))
281 return kFALSE;
282
283 const Int_t npixels = fGeom->GetNumPixels();
284 const Int_t nsectors = fGeom->GetNumSectors();
285 const Int_t nareas = fGeom->GetNumAreas();
286
287 InitHiGainArrays(npixels,nareas,nsectors);
288
289 return kTRUE;
290}
291
292void MHCalibrationPulseTimeCam::InitHiGainArrays(const Int_t npixels, const Int_t nareas, const Int_t nsectors)
293{
294
295 const Int_t samples = fRunHeader->GetNumSamples();
296 SetBinning(samples, -0.5, samples-0.5);
297
298 if (fHiGainArray->GetSize()==0)
299 {
300 for (Int_t i=0; i<npixels; i++)
301 {
302 fHiGainArray->AddAt(new MHCalibrationPix(Form("%sHiGainPix%04d",fHistName.Data(),i),
303 Form("%s High Gain Pixel %4d",fHistTitle.Data(),i)),i);
304
305 MHCalibrationPix &pix = (*this)[i];
306 pix.SetBinning(fNbins, fFirst, fLast);
307
308 InitHists(pix, (*fBadPixels)[i], i);
309
310 if (fCam)
311 (*fCam)[i].SetPixId(i);
312 }
313 }
314
315 if (!IsAverageing())
316 return;
317
318 if (fAverageHiGainAreas->GetSize()==0)
319 {
320 for (Int_t j=0; j<nareas; j++)
321 {
322 fAverageHiGainAreas->AddAt(new MHCalibrationPix(Form("%sHiGainArea%d",fHistName.Data(),j),
323 Form("%s High Gain Area Idx %d",fHistTitle.Data(),j)),j);
324
325 MHCalibrationPix &pix = GetAverageHiGainArea(j);
326 pix.SetBinning(fNbins, fFirst, fLast);
327
328 InitHists(pix, fCam->GetAverageBadArea(j), j);
329 }
330 }
331
332 if (fAverageHiGainSectors->GetSize()==0)
333 {
334 for (Int_t j=0; j<nsectors; j++)
335 {
336 fAverageHiGainSectors->AddAt(new MHCalibrationPix(Form("%sHiGainSector%02d",fHistName.Data(),j),
337 Form("%s High Gain Sector %02d",fHistTitle.Data(),j)),j);
338 MHCalibrationPix &pix = GetAverageHiGainSector(j);
339 pix.SetBinning(fNbins, fFirst, fLast);
340
341 InitHists(pix, fCam->GetAverageBadSector(j), j);
342 }
343 }
344}
345
346// --------------------------------------------------------------------------
347//
348// Retrieves from MExtractedSignalCam:
349// - first used LoGain FADC slice
350//
351// Retrieves from MGeomCam:
352// - number of pixels
353// - number of pixel areas
354// - number of sectors
355//
356// For all TOrdCollection's (including the averaged ones), the following steps are performed:
357//
358// 1) Fill PulseTimes histograms (MHGausEvents::FillHistAndArray()) with:
359// - MExtractedSignalPix::GetExtractedSignalHiGain();
360// - MExtractedSignalPix::GetExtractedSignalLoGain();
361//
362Bool_t MHCalibrationPulseTimeCam::FillHists(const MParContainer *par, const Stat_t w)
363{
364 MPedestalSubtractedEvt *evt = (MPedestalSubtractedEvt*)par;
365 if (!evt)
366 {
367 *fLog << err << "No argument in MHCalibrationPulseTimeCam::Fill... abort." << endl;
368 return kFALSE;
369 }
370
371 const UInt_t nareas = fGeom->GetNumAreas();
372 const UInt_t nsectors = fGeom->GetNumSectors();
373
374 MArrayD sumarea(nareas);
375 MArrayD sumsector(nsectors);
376
377 fAverageAreaNum.Reset();
378 fAverageSectorNum.Reset();
379
380 const Int_t npix = evt->GetNumPixels();
381 for (int idx=0; idx<npix; idx++)
382 {
383 if (fBadPixels)
384 if ((*fBadPixels)[idx].IsUnsuitable())
385 continue;
386
387 // Check for saturation
388 if (evt->GetSaturation(idx, fSaturationLimit)>0)
389 continue;
390
391 // Get position of maximum
392 Float_t max;
393 const Int_t pos = evt->GetMax(idx, max);
394
395 // check if maximum is high enough
396 if (max<fLowerSignalLimit)
397 continue;
398
399 const Int_t maxpos = pos-1;
400
401 (*this)[idx].FillHist(maxpos);
402
403 const Int_t aidx = (*fGeom)[idx].GetAidx();
404 const Int_t sector = (*fGeom)[idx].GetSector();
405
406 sumarea[aidx] += maxpos;
407 sumsector[sector] += maxpos;
408
409 fAverageAreaNum[aidx]++;
410 fAverageSectorNum[sector]++;
411 }
412
413 for (UInt_t j=0; j<nareas; j++)
414 {
415 if (fAverageAreaNum[j] > fNumPixelsRequired)
416 {
417 sumarea[j] /= fAverageAreaNum[j];
418
419 if (IsOscillations())
420 GetAverageHiGainArea(j).FillHistAndArray(sumarea[j]);
421 else
422 GetAverageHiGainArea(j).FillHist(sumarea[j]);
423
424 }
425 }
426
427 for (UInt_t j=0; j<nsectors; j++)
428 {
429 if (fAverageSectorNum[j] > 0)
430 {
431 sumsector[j] /= fAverageSectorNum[j];
432
433 if (IsOscillations())
434 GetAverageHiGainSector(j).FillHistAndArray(sumsector[j]);
435 else
436 GetAverageHiGainSector(j).FillHist(sumsector[j]);
437 }
438 }
439
440 return kTRUE;
441}
442
443// --------------------------------------------------------------------------
444//
445// For all TOrdCollection's (including the averaged ones), the following steps are performed:
446//
447// 1) Returns if the pixel is excluded.
448// 2) Tests saturation. In case yes, set the flag: MCalibrationPix::SetHiGainSaturation()
449// or the flag: MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainSaturated )
450// 3) Store the absolute arrival times in the MCalibrationPulseTimePix's. If flag
451// MCalibrationPix::IsHiGainSaturation() is set, the Low-Gain arrival times are stored,
452// otherwise the Hi-Gain ones.
453// 4) Calls to MHCalibrationCam::FitHiGainArrays() and MCalibrationCam::FitLoGainArrays()
454// with the flags:
455// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainNotFitted )
456// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainNotFitted )
457// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kHiGainOscillating )
458// - MBadPixelsPix::SetUncalibrated( MBadPixelsPix::kLoGainOscillating )
459//
460Bool_t MHCalibrationPulseTimeCam::FinalizeHists()
461{
462
463 *fLog << endl;
464
465 //
466 // Perform the fitting for the High Gain (done in MHCalibrationCam)
467 //
468 for (Int_t i=0; i<fHiGainArray->GetSize(); i++)
469 {
470
471 MHCalibrationPix &hist = (*this)[i];
472 if (hist.IsExcluded())
473 continue;
474
475 CalcHists(hist, (*fCam)[i]);
476 }
477
478 if (!IsAverageing())
479 return kTRUE;
480
481 for (Int_t j=0; j<fAverageHiGainAreas->GetSize(); j++)
482 CalcHists(GetAverageHiGainArea(j), fCam->GetAverageArea(j));
483
484 for (Int_t j=0; j<fAverageHiGainSectors->GetSize(); j++)
485 CalcHists(GetAverageHiGainSector(j), fCam->GetAverageSector(j));
486
487 return kTRUE;
488}
489
490void MHCalibrationPulseTimeCam::CalcHists(MHCalibrationPix &hist, MCalibrationPix &pix) const
491{
492
493 if (hist.IsEmpty())
494 {
495 *fLog << warn << hist.GetName() << ": Histogram empty." << endl;
496 return;
497 }
498 if (hist.IsOnlyOverflow())
499 {
500 *fLog << warn << hist.GetName() << ": Histogram contains only overflows." << endl;
501 return;
502 }
503 if (hist.IsOnlyUnderflow())
504 {
505 *fLog << warn << hist.GetName() << ": Histogram contains only underflows." << endl;
506 return;
507 }
508
509 hist.BypassFit();
510
511 pix.SetHiGainMean ( hist.GetMean() );
512 pix.SetHiGainMeanVar ( hist.GetMeanErr() * hist.GetMeanErr() );
513 pix.SetHiGainRms ( hist.GetHistRms() );
514 pix.SetHiGainSigma ( hist.GetSigma() );
515 pix.SetHiGainSigmaVar ( hist.GetSigmaErr()* hist.GetSigmaErr() );
516
517 if (IsDebug())
518 {
519 *fLog << dbginf << GetDescriptor() << ": ID " << GetName()
520 << " "<<pix.GetPixId()
521 << " Mean: " << hist.GetMean ()
522 << " MeanErr: " << hist.GetMeanErr ()
523 << " MeanSigma: " << hist.GetSigma ()
524 << " MeanSigmaErr: " << hist.GetSigmaErr()
525 << " Prob: " << hist.GetProb ()
526 << endl;
527 }
528}
529
530
531// --------------------------------------------------------------------------
532//
533// Calls MHCalibrationPix::DrawClone() for pixel idx
534//
535void MHCalibrationPulseTimeCam::DrawPixelContent(Int_t idx) const
536{
537 (*this)[idx].DrawClone();
538}
539
540
541// -----------------------------------------------------------------------------
542//
543// Default draw:
544//
545// Displays the averaged areas, both High Gain and Low Gain
546//
547// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
548//
549void MHCalibrationPulseTimeCam::Draw(const Option_t *opt)
550{
551
552 const Int_t nareas = fAverageHiGainAreas->GetSize();
553 if (nareas == 0)
554 return;
555
556 TString option(opt);
557 option.ToLower();
558
559 TVirtualPad *pad = gPad ? gPad : MH::MakeDefCanvas(this);
560 pad->SetBorderMode(0);
561 pad->Divide(1,nareas, 1e-10, 1e-10);
562
563 //
564 // Loop over inner and outer pixels
565 //
566 for (Int_t i=0; i<nareas;i++)
567 {
568 pad->cd(i+1);
569
570 MHCalibrationPix &hipix = GetAverageHiGainArea(i);
571 DrawDataCheckPixel(hipix,i ? fOuterRefTime : fInnerRefTime);
572 }
573}
574
575// -----------------------------------------------------------------------------
576//
577// Draw the average pixel for the datacheck:
578//
579// Displays the averaged areas, both High Gain and Low Gain
580//
581// Calls the Draw of the fAverageHiGainAreas and fAverageLoGainAreas objects with options
582//
583void MHCalibrationPulseTimeCam::DrawDataCheckPixel(MHCalibrationPix &pix, const Float_t refline)
584{
585 gPad->SetBorderMode(0);
586 gPad->SetTicks();
587
588 TH1F *hist = pix.GetHGausHist();
589
590 //
591 // set the labels bigger
592 //
593 TAxis *xaxe = hist->GetXaxis();
594 TAxis *yaxe = hist->GetYaxis();
595 xaxe->CenterTitle();
596 yaxe->CenterTitle();
597 xaxe->SetTitleSize(0.06);
598 yaxe->SetTitleSize(0.076);
599 xaxe->SetTitleOffset(0.6);
600 yaxe->SetTitleOffset(0.65);
601 xaxe->SetLabelSize(0.06);
602 yaxe->SetLabelSize(0.06);
603 xaxe->SetTitle(hist->GetXaxis()->GetTitle());
604 yaxe->SetTitle(hist->GetYaxis()->GetTitle());
605 xaxe->SetRange(hist->GetMaximumBin()-30, hist->GetMaximumBin()+30);
606
607 hist->Draw();
608
609 DisplayRefLines(hist, refline);
610}
611
612void MHCalibrationPulseTimeCam::DisplayRefLines(const TH1F *hist, const Float_t refline) const
613{
614 TLine *line = new TLine(refline, 0, refline, hist->GetMaximum());
615 line->SetLineColor(106);
616 line->SetLineStyle(2);
617 line->SetLineWidth(3);
618 line->SetBit(kCanDelete);
619 line->Draw();
620
621 TLegend *leg = new TLegend(0.8,0.35,0.99,0.65);
622 leg->AddEntry(line, "Reference", "l");
623 leg->SetBit(kCanDelete);
624 leg->Draw();
625}
626
627Int_t MHCalibrationPulseTimeCam::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
628{
629
630 Int_t rc = MHCalibrationCam::ReadEnv(env,prefix,print);
631 if (rc==kERROR)
632 return kERROR;
633
634 if (IsEnvDefined(env, prefix, "SaturationLimit", print))
635 {
636 SetSaturationLimit(GetEnvValue(env, prefix, "SaturationLimit", fSaturationLimit));
637 rc = kTRUE;
638 }
639
640 if (IsEnvDefined(env, prefix, "LowerSignalLimit", print))
641 {
642 SetLowerSignalLimit(GetEnvValue(env,prefix,"LowerSignalLimit",fLowerSignalLimit));
643 rc = kTRUE;
644 }
645
646 if (IsEnvDefined(env, prefix, "NumPixelsRequired", print))
647 {
648 SetNumPixelsRequired(GetEnvValue(env,prefix,"NumPixelsRequired",fNumPixelsRequired));
649 rc = kTRUE;
650 }
651
652 if (IsEnvDefined(env, prefix, "ReferenceFile", print))
653 {
654 SetReferenceFile(GetEnvValue(env,prefix,"ReferenceFile",fReferenceFile.Data()));
655 rc = kTRUE;
656 }
657
658 TEnv refenv(fReferenceFile);
659
660 fInnerRefTime = refenv.GetValue("InnerRefTime",fInnerRefTime);
661 fOuterRefTime = refenv.GetValue("OuterRefTime",fOuterRefTime);
662
663 return rc;
664}
Note: See TracBrowser for help on using the repository browser.