source: trunk/MagicSoft/Mars/mpedestal/MPedCalcFromLoGain.cc@ 5245

Last change on this file since 5245 was 5204, checked in by gaug, 20 years ago
*** empty log message ***
File size: 20.6 KB
Line 
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): Josep Flix 04/2001 <mailto:jflix@ifae.es>
19! Author(s): Thomas Bretz 05/2001 <mailto:tbretz@astro.uni-wuerzburg.de>
20! Author(s): Sebastian Commichau 12/2003
21! Author(s): Javier Rico 01/2004 <mailto:jrico@ifae.es>
22! Author(s): Markus Gaug 01/2004 <mailto:markus@ifae.es>
23! Author(s): Florian Goebel 06/2004 <mailto:fgoebel@mppmu.mpg.de>
24!
25! Copyright: MAGIC Software Development, 2000-2004
26!
27!
28\* ======================================================================== */
29
30/////////////////////////////////////////////////////////////////////////////
31//
32// MPedCalcLoGain
33//
34//
35// This task is derived form MPedCalcPedRun, described below. However, It
36// calculates the pedstals using the low gain slices, whenever the difference
37// between the highest and the lowest slice in the high gain
38// slices is below a given threshold (SetMaxHiGainVar). In this case the receiver
39// boards do not switch to lo gain and the so called lo gain slices are actually
40// high gain slices.
41//
42// MPedCalcLoGain also fills the ABoffset in MPedestalPix which allows to correct
43// the 150 MHz clock noise.
44//
45// This task takes a pedestal run file and fills MPedestalCam during
46// the Process() with the pedestal and rms computed in an event basis.
47// In the PostProcess() MPedestalCam is finally filled with the pedestal
48// mean and rms computed in a run basis.
49// More than one run (file) can be merged
50//
51// MPedCalcPedRun applies the following formula (1):
52//
53// Pedestal per slice = sum(x_i) / n / slices
54// PedRMS per slice = Sqrt( ( sum(x_i^2) - sum(x_i)^2/n ) / n-1 / slices )
55//
56// where x_i is the sum of "slices" FADC slices and sum means the sum over all
57// events. "n" is the number of events, "slices" is the number of summed FADC samples.
58//
59// Note that the slice-to-slice fluctuations are not Gaussian, but Poissonian, thus
60// asymmetric and they are correlated.
61//
62// It is important to know that the Pedestal per slice and PedRMS per slice depend
63// on the number of used FADC slices, as seen in the following plots:
64//
65//Begin_Html
66/*
67<img src="images/PedestalStudyInner.gif">
68*/
69//End_Html
70//
71//Begin_Html
72/*
73<img src="images/PedestalStudyOuter.gif">
74*/
75//
76// The plots show the inner and outer pixels, respectivly and have the following meaning:
77//
78// 1) The calculated mean pedestal per slice (from MPedCalcFromLoGain)
79// 2) The fitted mean pedestal per slice (from MHPedestalCam)
80// 3) The calculated pedestal RMS per slice (from MPedCalcFromLoGain)
81// 4) The fitted sigma of the pedestal distribution per slice
82// (from MHPedestalCam)
83// 5) The relative difference between calculation and histogram fit
84// for the mean
85// 6) The relative difference between calculation and histogram fit
86// for the sigma or RMS, respectively.
87//
88// The calculated means do not change significantly except for the case of 2 slices,
89// however the RMS changes from 5.7 per slice in the case of 2 extracted slices
90// to 8.3 per slice in the case of 26 extracted slices. This change is very significant.
91//
92// The plots have been produced on run 20123. You can reproduce them using
93// the macro pedestalstudies.C
94//
95// Usage of this class:
96// ====================
97//
98// Call: SetRange(higainfirst, higainlast, logainfirst, logainlast)
99// to modify the ranges in which the window is allowed to move.
100// Defaults are:
101//
102// fHiGainFirst = fgHiGainFirst = 0
103// fHiGainLast = fgHiGainLast = 29
104// fLoGainFirst = fgLoGainFirst = 0
105// fLoGainLast = fgLoGainLast = 14
106//
107// Call: SetWindowSize(windowhigain, windowlogain)
108// to modify the sliding window widths. Windows have to be an even number.
109// In case of odd numbers, the window will be modified.
110//
111// Defaults are:
112//
113// fHiGainWindowSize = fgHiGainWindowSize = 14
114// fLoGainWindowSize = fgLoGainWindowSize = 0
115//
116// Variables:
117// fgHiGainFirst; First FADC slice Hi-Gain (currently set to: 3)
118// fgHiGainLast: Last FADC slice Hi-Gain (currently set to: 14)
119// fgLoGainFirst: First FADC slice Lo-Gain (currently set to: 3)
120// fgLoGainLast: Last FADC slice Lo-Gain (currently set to: 14)
121// fgHiGainWindowSize: The extraction window Hi-Gain
122// fgLoGainWindowSize: The extraction window Lo-Gain
123// fgMaxHiGainVar: The maximum difference between the highest and lowest slice
124// in the high gain window allowed in order to use low gain
125//
126// Input Containers:
127// MRawEvtData
128// MRawRunHeader
129// MGeomCam
130//
131// Output Containers:
132// MPedestalCam
133//
134// See also: MPedestalCam, MPedestalPix, MHPedestalCam, MExtractor
135//
136/////////////////////////////////////////////////////////////////////////////
137#include "MPedCalcFromLoGain.h"
138#include "MExtractor.h"
139
140#include "MParList.h"
141
142#include "MLog.h"
143#include "MLogManip.h"
144
145#include "MRawRunHeader.h"
146#include "MRawEvtPixelIter.h"
147#include "MRawEvtData.h"
148
149#include "MPedestalPix.h"
150#include "MPedestalCam.h"
151
152#include "MGeomPix.h"
153#include "MGeomCam.h"
154
155ClassImp(MPedCalcFromLoGain);
156
157using namespace std;
158
159const Byte_t MPedCalcFromLoGain::fgHiGainFirst = 0;
160const Byte_t MPedCalcFromLoGain::fgHiGainLast = 11;
161const Byte_t MPedCalcFromLoGain::fgLoGainFirst = 1;
162const Byte_t MPedCalcFromLoGain::fgLoGainLast = 14;
163const Byte_t MPedCalcFromLoGain::fgHiGainWindowSize = 12;
164const Byte_t MPedCalcFromLoGain::fgLoGainWindowSize = 14;
165const Byte_t MPedCalcFromLoGain::fgMaxHiGainVar = 40;
166
167// --------------------------------------------------------------------------
168//
169// Default constructor:
170//
171// Sets:
172// - all pointers to NULL
173// - fWindowSizeHiGain to fgHiGainWindowSize
174// - fWindowSizeLoGain to fgLoGainWindowSize
175//
176// Calls:
177// - AddToBranchList("fHiGainPixId");
178// - AddToBranchList("fHiGainFadcSamples");
179// - SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast)
180// - Clear()
181//
182MPedCalcFromLoGain::MPedCalcFromLoGain(const char *name, const char *title)
183 : fWindowSizeHiGain(fgHiGainWindowSize),
184 fWindowSizeLoGain(fgLoGainWindowSize),
185 fGeom(NULL), fPedContainerName("MPedestalCam")
186{
187 fName = name ? name : "MPedCalcFromLoGain";
188 fTitle = title ? title : "Task to calculate pedestals from pedestal runs raw data";
189
190 AddToBranchList("fHiGainPixId");
191 AddToBranchList("fLoGainPixId");
192 AddToBranchList("fHiGainFadcSamples");
193 AddToBranchList("fLoGainFadcSamples");
194
195 SetRange();
196 SetMaxHiGainVar();
197 SetPedestalUpdate(kTRUE);
198 SetNumEventsDump(500);
199
200 Clear();
201}
202
203// --------------------------------------------------------------------------
204//
205// Sets:
206// - fNumSamplesTot to 0
207// - fRawEvt to NULL
208// - fRunHeader to NULL
209// - fPedestals to NULL
210//
211void MPedCalcFromLoGain::Clear(const Option_t *o)
212{
213 fRawEvt = NULL;
214 fRunHeader = NULL;
215 fPedestals = NULL;
216
217 // If the size is yet set, set the size
218 if (fSumx.GetSize()>0)
219 {
220 // Reset contents of arrays.
221 fSumx.Reset();
222 fSumx2.Reset();
223 fSumAB0.Reset();
224 fSumAB1.Reset();
225 fNumEventsUsed.Reset();
226 fTotalCounter.Reset();
227 }
228}
229
230// --------------------------------------------------------------------------
231//
232// SetRange:
233//
234// Calls:
235// - MExtractor::SetRange(hifirst,hilast,lofirst,lolast);
236// - SetWindowSize(fWindowSizeHiGain,fWindowSizeLoGain);
237//
238void MPedCalcFromLoGain::SetRange(Byte_t hifirst, Byte_t hilast, Byte_t lofirst, Byte_t lolast)
239{
240 MExtractor::SetRange(hifirst, hilast, lofirst, lolast);
241
242 //
243 // Redo the checks if the window is still inside the ranges
244 //
245 SetWindowSize(fWindowSizeHiGain,fWindowSizeLoGain);
246}
247
248// --------------------------------------------------------------------------
249//
250void MPedCalcFromLoGain::SetMaxHiGainVar(Byte_t maxvar)
251{
252 fMaxHiGainVar = maxvar;
253}
254
255// --------------------------------------------------------------------------
256//
257// Checks:
258// - if a window is odd, subtract one
259// - if a window is bigger than the one defined by the ranges, set it to the available range
260// - if a window is smaller than 2, set it to 2
261//
262// Sets:
263// - fNumHiGainSamples to: (Float_t)fWindowSizeHiGain
264// - fNumLoGainSamples to: (Float_t)fWindowSizeLoGain
265// - fSqrtHiGainSamples to: TMath::Sqrt(fNumHiGainSamples)
266// - fSqrtLoGainSamples to: TMath::Sqrt(fNumLoGainSamples)
267//
268void MPedCalcFromLoGain::SetWindowSize(Byte_t windowh, Byte_t windowl)
269{
270
271 fWindowSizeHiGain = windowh & ~1;
272 fWindowSizeLoGain = windowl & ~1;
273
274 if (fWindowSizeHiGain != windowh)
275 {
276 *fLog << warn;
277 *fLog << GetDescriptor() << ": HiGain window has to be even, set to: ";
278 *fLog << int(fWindowSizeHiGain) << " samples " << endl;
279 }
280
281 if (fWindowSizeLoGain != windowl)
282 {
283 *fLog << warn;
284 *fLog << GetDescriptor() << ": Lo Gain window has to be even, set to: ";
285 *fLog << int(fWindowSizeLoGain) << " samples " << endl;
286 }
287
288 if (fWindowSizeHiGain == 0)
289 {
290 *fLog << warn;
291 *fLog << GetDescriptor() << ": HiGain window currently set to 0, will set it to 2 samples ";
292 fWindowSizeHiGain = 2;
293 }
294
295 if (fWindowSizeLoGain == 0)
296 {
297 *fLog << warn;
298 *fLog << GetDescriptor() << ": LoGain window currently set to 0, will set it to 2 samples ";
299 fWindowSizeLoGain = 2;
300 }
301
302 const Byte_t availhirange = (fHiGainLast-fHiGainFirst+1) & ~1;
303 const Byte_t availlorange = (fLoGainLast-fLoGainFirst+1) & ~1;
304
305 if (fWindowSizeHiGain > availhirange)
306 {
307 *fLog << warn;
308 *fLog << GetDescriptor() << ": HiGain window " << (int)fWindowSizeHiGain;
309 *fLog << " out of range [" << (int)fHiGainFirst;
310 *fLog << "," << (int)fHiGainLast << "]" << endl;
311 *fLog << "Will set window size to " << (int)availhirange << endl;
312 fWindowSizeHiGain = availhirange;
313 }
314
315 if (fWindowSizeLoGain > availlorange)
316 {
317 *fLog << warn;
318 *fLog << GetDescriptor() << ": LoGain window " << (int)fWindowSizeLoGain;
319 *fLog << " out of range [" << (int)fLoGainFirst;
320 *fLog << "," << (int)fLoGainLast << "]" << endl;
321 *fLog << "Will set window size to " << (int)availlorange << endl;
322 fWindowSizeLoGain = availlorange;
323 }
324 /*
325 fNumHiGainSamples = (Float_t)fWindowSizeHiGain;
326 fNumLoGainSamples = (Float_t)fWindowSizeLoGain;
327
328 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
329 fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
330 */
331}
332
333// --------------------------------------------------------------------------
334//
335// Look for the following input containers:
336//
337// - MRawEvtData
338// - MRawRunHeader
339// - MGeomCam
340//
341// The following output containers are also searched and created if
342// they were not found:
343//
344// - MPedestalCam
345//
346Int_t MPedCalcFromLoGain::PreProcess(MParList *pList)
347{
348 Clear();
349
350 fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
351 if (!fRawEvt)
352 {
353 *fLog << err << "MRawEvtData not found... aborting." << endl;
354 return kFALSE;
355 }
356
357 fRunHeader = (MRawRunHeader*)pList->FindObject(AddSerialNumber("MRawRunHeader"));
358 if (!fRunHeader)
359 {
360 *fLog << err << AddSerialNumber("MRawRunHeader") << " not found... aborting." << endl;
361 return kFALSE;
362 }
363
364 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
365 if (!fGeom)
366 {
367 *fLog << err << "MGeomCam not found... aborting." << endl;
368 return kFALSE;
369 }
370
371 fPedestals = (MPedestalCam*)pList->FindCreateObj("MPedestalCam", AddSerialNumber(fPedContainerName));
372 if (!fPedestals)
373 return kFALSE;
374
375 if (fNumEventsDump<=0 && fPedestalUpdate)
376 {
377 *fLog << warn << "Pedestal Update switched on and Number of Events to dump <= 0... fNumEventsDump=1000" << endl;
378 fNumEventsDump=1000;
379 }
380
381 *fLog << inf << "Continous update switched " << (fPedestalUpdate?"on":"off");
382 if (fPedestalUpdate)
383 *fLog << " (dump each " << fNumEventsDump << " events)" << endl;
384 *fLog << endl;
385
386 return kTRUE;
387}
388
389// --------------------------------------------------------------------------
390//
391// The ReInit searches for:
392// - MRawRunHeader::GetNumSamplesHiGain()
393// - MRawRunHeader::GetNumSamplesLoGain()
394//
395// In case that the variables fHiGainLast and fLoGainLast are smaller than
396// the even part of the number of samples obtained from the run header, a
397// warning is given an the range is set back accordingly. A call to:
398// - SetRange(fHiGainFirst, fHiGainLast-diff, fLoGainFirst, fLoGainLast) or
399// - SetRange(fHiGainFirst, fHiGainLast, fLoGainFirst, fLoGainLast-diff)
400// is performed in that case. The variable diff means here the difference
401// between the requested range (fHiGainLast) and the available one. Note that
402// the functions SetRange() are mostly overloaded and perform more checks,
403// modifying the ranges again, if necessary.
404//
405Bool_t MPedCalcFromLoGain::ReInit(MParList *pList)
406{
407 if (fRunHeader->GetNumSamplesHiGain()<1)
408 {
409 *fLog << err << "ERROR - Number of available HiGainSamples<1... abort." << endl;
410 return kFALSE;
411 }
412 if (fRunHeader->GetNumSamplesLoGain()<1)
413 {
414 *fLog << err << "ERROR - Number of available LoGainSamples<1... abort." << endl;
415 return kFALSE;
416 }
417
418 fHiGainFirst = TMath::Max(fHiGainFirst, 0);
419 fHiGainLast = TMath::Min(fHiGainLast, fRunHeader->GetNumSamplesHiGain()-1);
420
421 fLoGainFirst = TMath::Max(fLoGainFirst, 0);
422 fLoGainLast = TMath::Min(fLoGainLast, fRunHeader->GetNumSamplesLoGain()-1);
423
424 const Double_t wh = fWindowSizeHiGain;
425 const Double_t wl = fWindowSizeLoGain;
426 const Double_t fh = fHiGainFirst;
427 const Double_t fl = fLoGainFirst;;
428
429 fWindowSizeHiGain = TMath::Min(fWindowSizeHiGain, fHiGainLast-fHiGainFirst+1);
430 fWindowSizeLoGain = TMath::Min(fWindowSizeLoGain, fLoGainLast-fLoGainFirst+1);
431
432 SetRange(fHiGainFirst, fHiGainLast, fLoGainFirst, fLoGainLast);
433
434 if (wh!=fWindowSizeHiGain || fh!=fHiGainFirst || wl!=fWindowSizeLoGain || fl!=fLoGainFirst)
435 {
436 *fLog << inf << endl;
437 *fLog << "Taking " << Form("%2d", (int)fWindowSizeHiGain) << " slices of Hi-Gain starting at slice " << (int)fHiGainFirst << endl;
438 *fLog << "Taking " << Form("%2d", (int)fWindowSizeLoGain) << " slices of Lo-Gain starting at slice " << (int)fLoGainFirst << endl;
439 }
440
441 if (fWindowSizeHiGain==0)
442 {
443 *fLog << err << "ERROR - HiGain windows size == 0... abort." << endl;
444 return kFALSE;
445 }
446 if (fWindowSizeLoGain==0)
447 {
448 *fLog << err << "ERROR - LoGain windows size == 0... abort." << endl;
449 return kFALSE;
450 }
451
452 // If the size is not yet set, set the size
453 if (fSumx.GetSize()==0)
454 {
455 const Int_t npixels = fPedestals->GetSize();
456
457 fSumx. Set(npixels);
458 fSumx2.Set(npixels);
459 fSumAB0.Set(npixels);
460 fSumAB1.Set(npixels);
461 fNumEventsUsed.Set(npixels);
462 fTotalCounter.Set(npixels);
463
464 // Reset contents of arrays.
465 fSumx.Reset();
466 fSumx2.Reset();
467 fSumAB0.Reset();
468 fSumAB1.Reset();
469 fNumEventsUsed.Reset();
470 fTotalCounter.Reset();
471 }
472
473 return kTRUE;
474}
475
476void MPedCalcFromLoGain::Calc(ULong_t n, UInt_t idx)
477{
478 const ULong_t nsamplestot = n*fWindowSizeLoGain;
479
480 const Float_t sum = fSumx.At(idx);
481 const Float_t sum2 = fSumx2.At(idx);
482 const Float_t ped = sum/nsamplestot;
483
484 // 1. Calculate the Variance of the sums:
485 Float_t var = (sum2-sum*sum/n)/(n-1.);
486
487 // 2. Scale the variance to the number of slices:
488 var /= (Float_t)(fWindowSizeLoGain);
489
490 // 3. Calculate the RMS from the Variance:
491 const Float_t rms = var<0 ? 0 : TMath::Sqrt(var);
492
493 // 4. Calculate the amplitude of the 150MHz "AB" noise
494 const Float_t abOffs = (fSumAB0[idx] - fSumAB1[idx]) / nsamplestot;
495
496 (*fPedestals)[idx].Set(ped, rms, abOffs, n);
497
498 fTotalCounter[idx]++;
499}
500
501// --------------------------------------------------------------------------
502//
503// Fill the MPedestalCam container with the signal mean and rms for the event.
504// Store the measured signal in arrays fSumx and fSumx2 so that we can
505// calculate the overall mean and rms in the PostProcess()
506//
507Int_t MPedCalcFromLoGain::Process()
508{
509 MRawEvtPixelIter pixel(fRawEvt);
510
511 while (pixel.Next())
512 {
513 const UInt_t idx = pixel.GetPixelId();
514
515 Byte_t *ptr = pixel.GetHiGainSamples() + fHiGainFirst;
516 Byte_t *end = ptr + fWindowSizeHiGain;
517
518 UInt_t sum = 0;
519 UInt_t sqr = 0;
520
521 UInt_t max = 0;
522 UInt_t min = 255;
523
524 // Find the maximum and minimum signal per slice in the high gain window
525 do {
526 if (*ptr > max)
527 max = *ptr;
528 if (*ptr < min)
529 min = *ptr;
530 } while (++ptr != end);
531
532 // If the maximum in the high gain window is smaller than
533 if (max-min>=fMaxHiGainVar || max>=255)
534 continue;
535
536 ptr = pixel.GetLoGainSamples() + fLoGainFirst;
537 end = ptr + fWindowSizeLoGain;
538
539 Byte_t *firstSlice = ptr;
540
541 do {
542 sum += *ptr;
543 sqr += *ptr * *ptr;
544 } while (++ptr != end);
545
546 const Float_t msum = (Float_t)sum;
547 const Float_t sqrsum = msum*msum;
548
549 fSumx[idx] += msum;
550 fSumx2[idx] += sqrsum;
551 fNumEventsUsed[idx]++;
552
553 // Calculate the amplitude of the 150MHz "AB" noise
554
555 if (pixel.IsABFlagValid())
556 {
557 const Int_t abFlag = (fRunHeader->GetNumSamplesHiGain()
558 + fLoGainFirst + pixel.HasABFlag()) & 0x1;
559 for (Int_t islice=0; islice<fWindowSizeLoGain; islice+=2)
560 {
561 const Int_t sliceAB0 = islice + abFlag;
562 const Int_t sliceAB1 = islice - abFlag + 1;
563 fSumAB0[idx] += firstSlice[sliceAB0];
564 fSumAB1[idx] += firstSlice[sliceAB1];
565 }
566 }
567
568 if (!fPedestalUpdate || fNumEventsUsed[idx]<fNumEventsDump)
569 continue;
570
571 Calc(fNumEventsDump, idx);
572
573 fNumEventsUsed[idx]=0;
574 fSumx[idx]=0;
575 fSumx2[idx]=0;
576 fSumAB0[idx]=0;
577 fSumAB1[idx]=0;
578 }
579
580 if (fPedestalUpdate)
581 fPedestals->SetReadyToSave();
582
583 return kTRUE;
584}
585
586// --------------------------------------------------------------------------
587//
588// Compute signal mean and rms in the whole run and store it in MPedestalCam
589//
590Int_t MPedCalcFromLoGain::PostProcess()
591{
592 // Compute pedestals and rms from the whole run
593 if (fPedestalUpdate || GetNumExecutions()<1)
594 return kTRUE;
595
596 *fLog << flush << inf << "Calculating pedestals..." << flush;
597
598 Double_t sum = 0;
599
600 const Int_t npix = fGeom->GetNumPixels();
601 for (Int_t idx=0; idx<npix; idx++)
602 {
603 const ULong_t n = fNumEventsUsed[idx];
604 if (n>1)
605 Calc(n, idx);
606 sum += n;
607 }
608
609 *fLog << flush << inf << "Calculating means..." << flush;
610
611 fPedestals->SetTotalEntries((UInt_t)(sum/npix*(fWindowSizeLoGain+fWindowSizeHiGain)));
612 fPedestals->SetReadyToSave();
613 fPedestals->ReCalc(*fGeom);
614 return kTRUE;
615}
616
617Int_t MPedCalcFromLoGain::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
618{
619 if (MExtractor::ReadEnv(env, prefix, print)==kERROR)
620 return kERROR;
621
622 Bool_t rc=kFALSE;
623
624 Byte_t hw = fWindowSizeHiGain;
625 Byte_t lw = fWindowSizeLoGain;
626 if (IsEnvDefined(env, prefix, "WindowSizeHiGain", print))
627 {
628 hw = GetEnvValue(env, prefix, "WindowSizeHiGain", hw);
629 rc = kTRUE;
630 }
631 if (IsEnvDefined(env, prefix, "WindowSizeLoGain", print))
632 {
633 lw = GetEnvValue(env, prefix, "WindowSizeLoGain", lw);
634 rc = kTRUE;
635 }
636
637 if (rc)
638 SetWindowSize(hw, lw);
639
640 if (IsEnvDefined(env, prefix, "NumEventsDump", print))
641 {
642 SetNumEventsDump(GetEnvValue(env, prefix, "NumEventsDump", fNumEventsDump));
643 rc = kTRUE;
644 }
645
646 if (IsEnvDefined(env, prefix, "MaxHiGainVar", print))
647 {
648 SetMaxHiGainVar(GetEnvValue(env, prefix, "MaxHiGainVar", fMaxHiGainVar));
649 rc = kTRUE;
650 }
651
652 if (IsEnvDefined(env, prefix, "PedestalUpdate", print))
653 {
654 SetPedestalUpdate(GetEnvValue(env, prefix, "PedestalUpdate", fPedestalUpdate));
655 rc = kTRUE;
656 }
657
658 return rc;
659}
Note: See TracBrowser for help on using the repository browser.