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

Last change on this file since 5489 was 5395, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 19.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! Author(s): Nepomuk Otte 10/2004 <mailto:otte@mppmu.mpg.de>
25!
26! Copyright: MAGIC Software Development, 2000-2004
27!
28!
29\* ======================================================================== */
30
31/////////////////////////////////////////////////////////////////////////////
32//
33// MPedCalcLoGain
34//
35// This task is devide 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//
99// fCheckWinFirst = fgCheckWinFirst = 0
100// fHiGainLast = fgCheckWinLast = 29
101// fExtractWinFirst = fgExtractWinFirst = 0
102// fExtractWinSize = fgExtractWinSize = 6
103// fMaxSignalVar = fgMaxSignalVar = 40;
104//
105// Call:
106// SetCheckRange(fCheckWinFirst,fCheckWinLast);
107// to set the Window in which a signal is searched
108//
109// SetExtractWindow(fExtractWindFirst,fExtractWinSize);
110// to set the Window from which a signal is extracted
111//
112// SetMaxSignalVar(fMaxSignalVar);
113// set the maximum allowed difference between maximum and minimal signal in CheckWindow
114//
115// Variables:
116// fgCheckWinFirst; First FADC slice to check for signal (currently set to: 0)
117// fgCheckWinLast: Last FADC slice to check for signal (currently set to: 29)
118// fgExtractWinFirst: First FADC slice to be used for pedestal extraction (currently set to: 15)
119// fgExtractWinSize: Window size in slices used for pedestal extraction (currently set to: 6)
120// fgMaxSignalVar: The maximum difference between the highest and lowest slice
121// in the check window allowed in order to use event
122//
123// Input Containers:
124// MRawEvtData
125// MRawRunHeader
126// MGeomCam
127//
128// Output Containers:
129// MPedestalCam
130//
131// See also: MPedestalCam, MPedestalPix, MHPedestalCam, MExtractor
132//
133/////////////////////////////////////////////////////////////////////////////
134#include "MPedCalcFromLoGain.h"
135#include "MExtractor.h"
136
137#include "MParList.h"
138
139#include "MLog.h"
140#include "MLogManip.h"
141
142#include "MRawRunHeader.h"
143#include "MRawEvtPixelIter.h"
144#include "MRawEvtData.h"
145
146#include "MPedestalPix.h"
147#include "MPedestalCam.h"
148
149#include "MGeomPix.h"
150#include "MGeomCam.h"
151
152ClassImp(MPedCalcFromLoGain);
153
154using namespace std;
155
156const UShort_t MPedCalcFromLoGain::fgCheckWinFirst = 0;
157const UShort_t MPedCalcFromLoGain::fgCheckWinLast = 29;
158const UShort_t MPedCalcFromLoGain::fgExtractWinFirst = 15;
159const UShort_t MPedCalcFromLoGain::fgExtractWinSize = 6;
160const UShort_t MPedCalcFromLoGain::fgMaxSignalVar = 40;
161
162// --------------------------------------------------------------------------
163//
164// Default constructor:
165//
166// Sets:
167// - all pointers to NULL
168//
169// Calls:
170// - AddToBranchList("fHiGainPixId");
171// - AddToBranchList("fHiGainFadcSamples");
172// - SetCheckRange(fgCheckWinFirst, fgCheckWinLast, fgExtractWinFirst, fgExtractWinSize)
173// - Clear()
174//
175MPedCalcFromLoGain::MPedCalcFromLoGain(const char *name, const char *title)
176 : fGeom(NULL), fNamePedestalCam("MPedestalCam")
177{
178 fName = name ? name : "MPedCalcFromLoGain";
179 fTitle = title ? title : "Task to calculate pedestals from lo-gains";
180
181 AddToBranchList("fHiGainPixId");
182 AddToBranchList("fLoGainPixId");
183 AddToBranchList("fHiGainFadcSamples");
184 AddToBranchList("fLoGainFadcSamples");
185
186 SetCheckRange(fgCheckWinFirst, fgCheckWinLast);
187 SetExtractWindow(fgExtractWinFirst, fgExtractWinSize);
188
189 SetMaxSignalVar(fgMaxSignalVar);
190 SetPedestalUpdate(kTRUE);
191
192 Clear();
193}
194
195void MPedCalcFromLoGain::ResetArrays()
196{
197 // Reset contents of arrays.
198 fSumx.Reset();
199 fSumx2.Reset();
200 fSumAB0.Reset();
201 fSumAB1.Reset();
202 fNumEventsUsed.Reset();
203 fTotalCounter.Reset();
204}
205
206// --------------------------------------------------------------------------
207//
208// Sets:
209// - fRawEvt to NULL
210// - fRunHeader to NULL
211// - fPedestals to NULL
212//
213// Resets:
214// - fSumx
215// - fSumx2
216// - fSumAB0
217// - fSumAB1
218// - fNumEventsUsed
219// - fTotalCounter
220//
221void MPedCalcFromLoGain::Clear(const Option_t *o)
222{
223
224 fRawEvt = NULL;
225 fRunHeader = NULL;
226 fPedestals = NULL;
227
228 // If the size is yet set, set the size
229 if (fSumx.GetSize()>0)
230 ResetArrays();
231}
232
233// --------------------------------------------------------------------------
234//
235// SetCheckRange:
236//
237// Exits, if the first argument is smaller than 0
238// Exits, if the the last argument is smaller than the first
239//
240Bool_t MPedCalcFromLoGain::SetCheckRange(UShort_t chfirst, UShort_t chlast)
241{
242 Bool_t rc = kTRUE;
243
244 if (chlast<=chfirst)
245 {
246 *fLog << warn << GetDescriptor();
247 *fLog << " - WARNING: Last slice in SetCheckRange smaller than first slice... set to first+2" << endl;
248 chlast = chfirst+1;
249 rc = kFALSE;
250 }
251
252 fCheckWinFirst = chfirst;
253 fCheckWinLast = chlast;
254
255 return rc;
256}
257
258// --------------------------------------------------------------------------
259//
260// Checks:
261// - if a window is odd
262//
263Bool_t MPedCalcFromLoGain::SetExtractWindow(UShort_t windowf, UShort_t windows)
264{
265 Bool_t rc = kTRUE;
266
267 const Int_t odd = windows & 0x1;
268
269 if (odd || windows==0)
270 {
271 *fLog << warn << GetDescriptor();
272 *fLog << " - WARNING: Window size in SetExtraxtWindow has to be even and > 0... adjusting!" << endl;
273 windows += 1;
274 rc = kFALSE;
275 }
276
277 fExtractWinSize = windows;
278 fExtractWinFirst = windowf;
279 fExtractWinLast = fExtractWinFirst+fExtractWinSize-1;
280
281//
282// NO RANGE CHECK IMPLEMENTED, YET
283//
284/*
285 const Byte_t availhirange = (fHiGainLast-fHiGainFirst+1) & ~1;
286 const Byte_t availlorange = (fLoGainLast-fLoGainFirst+1) & ~1;
287
288 if (fWindowSizeHiGain > availhirange)
289 {
290 *fLog << warn;
291 *fLog << GetDescriptor() << ": HiGain window " << (int)fWindowSizeHiGain;
292 *fLog << " out of range [" << (int)fHiGainFirst;
293 *fLog << "," << (int)fHiGainLast << "]" << endl;
294 *fLog << "Will set window size to " << (int)availhirange << endl;
295 fWindowSizeHiGain = availhirange;
296 }
297
298 if (fWindowSizeLoGain > availlorange)
299 {
300 *fLog << warn;
301 *fLog << GetDescriptor() << ": LoGain window " << (int)fWindowSizeLoGain;
302 *fLog << " out of range [" << (int)fLoGainFirst;
303 *fLog << "," << (int)fLoGainLast << "]" << endl;
304 *fLog << "Will set window size to " << (int)availlorange << endl;
305 fWindowSizeLoGain = availlorange;
306 }
307 */
308
309 return rc;
310}
311
312// --------------------------------------------------------------------------
313//
314// Look for the following input containers:
315//
316// - MRawEvtData
317// - MRawRunHeader
318// - MGeomCam
319//
320// The following output containers are also searched and created if
321// they were not found:
322//
323// - MPedestalCam with the name fPedContainerName
324//
325Int_t MPedCalcFromLoGain::PreProcess(MParList *pList)
326{
327 Clear();
328
329 fRawEvt = (MRawEvtData*)pList->FindObject(AddSerialNumber("MRawEvtData"));
330 if (!fRawEvt)
331 {
332 *fLog << err << AddSerialNumber("MRawEvtData") << " not found... aborting." << endl;
333 return kFALSE;
334 }
335
336 fRunHeader = (MRawRunHeader*)pList->FindObject(AddSerialNumber("MRawRunHeader"));
337 if (!fRunHeader)
338 {
339 *fLog << err << AddSerialNumber("MRawRunHeader") << " not found... aborting." << endl;
340 return kFALSE;
341 }
342
343 fGeom = (MGeomCam*)pList->FindObject(AddSerialNumber("MGeomCam"));
344 if (!fGeom)
345 {
346 *fLog << err << AddSerialNumber("MGeomCam") << " not found... aborting." << endl;
347 return kFALSE;
348 }
349
350 fPedestals = (MPedestalCam*)pList->FindCreateObj("MPedestalCam", AddSerialNumber(fNamePedestalCam));
351 if (!fPedestals)
352 return kFALSE;
353
354 *fLog << inf;
355 Print();
356
357 return kTRUE;
358}
359
360// ---------------------------------------------------------------------------------
361//
362// Checks:
363// - if the number of available slices
364// (fRunHeader->GetNumSamplesHiGain()+(Int_t)fRunHeader->GetNumSamplesLoGain()-1)
365// is smaller than the number of used slices
366// (fExtractWinSize+ fExtractWinFirst-1) or
367// fCheckWinLast
368//
369// Sets the size (from MPedestalCam::GetSize() ) and resets the following arrays:
370// - fSumx
371// - fSumx2
372// - fSumAB0
373// - fSumAB1
374// - fNumEventsUsed
375// - fTotalCounter
376//
377Bool_t MPedCalcFromLoGain::ReInit(MParList *pList)
378{
379 Int_t lastavailableslice = fRunHeader->GetNumSamplesHiGain()+fRunHeader->GetNumSamplesLoGain()-1;
380 Int_t lastextractslice = fExtractWinSize + fExtractWinFirst - 1;
381
382 if (lastextractslice>lastavailableslice)//changed to override check
383 {
384 *fLog << warn << GetDescriptor();
385 *fLog << " - WARNING: Selected Extract Window ranges out of range...adjusting to last available slice ";
386 *fLog << lastavailableslice << endl;
387
388 lastextractslice=lastavailableslice;
389 }
390
391 if (fCheckWinLast>lastavailableslice)//changed to override check
392 {
393 *fLog << warn << GetDescriptor();
394 *fLog << " - WARNING: Last Check Window slice out of range...adjusting to last available slice ";
395 *fLog << lastavailableslice << endl;
396
397 fCheckWinLast=lastavailableslice;
398 }
399
400
401 // If the size is not yet set, set the size
402 if (fSumx.GetSize()==0)
403 {
404 const Int_t npixels = fPedestals->GetSize();
405
406 fSumx. Set(npixels);
407 fSumx2.Set(npixels);
408 fSumAB0.Set(npixels);
409 fSumAB1.Set(npixels);
410 fNumEventsUsed.Set(npixels);
411 fTotalCounter.Set(npixels);
412
413 ResetArrays();
414 }
415
416 return kTRUE;
417}
418
419// ---------------------------------------------------------------------------------
420//
421// Calculates for pixel "idx":
422//
423// Ped per slice = sum / n / fExtractWinSize;
424// RMS per slice = sqrt { (sum2 - sum*sum/n) / (n-1) / fExtractWinSize }
425// ABOffset per slice = (fSumAB0[idx] - fSumAB1[idx]) / n / fExtractWinSize;
426//
427// Sets fTotalCounter up by 1.
428//
429void MPedCalcFromLoGain::Calc(ULong_t n, UInt_t idx)
430{
431
432 const ULong_t nsamplestot = n*fExtractWinSize;
433
434 const Float_t sum = fSumx.At(idx);
435 const Float_t sum2 = fSumx2.At(idx);
436 const Float_t ped = sum/nsamplestot;
437
438 // 1. Calculate the Variance of the sums:
439 Float_t var = (sum2-sum*sum/n)/(n-1.);
440
441 // 2. Scale the variance to the number of slices:
442 var /= fExtractWinSize;
443
444 // 3. Calculate the RMS from the Variance:
445 const Float_t rms = var<0 ? 0 : TMath::Sqrt(var);
446
447 // 4. Calculate the amplitude of the 150MHz "AB" noise
448 const Float_t abOffs = (fSumAB0[idx] - fSumAB1[idx]) / nsamplestot;
449
450 (*fPedestals)[idx].Set(ped, rms, abOffs, n);
451
452 fTotalCounter[idx]++;
453}
454
455// ---------------------------------------------------------------------------------
456//
457// Returns the pointer to slice "slice".
458//
459UShort_t MPedCalcFromLoGain::GetSlice(MRawEvtPixelIter *pixel, UInt_t slice)
460{
461 const UShort_t nh = (Int_t)fRunHeader->GetNumSamplesHiGain();
462
463 Byte_t *ptr;
464
465 if(slice<nh)
466 ptr = pixel->GetHiGainSamples() + slice;
467 else
468 ptr = pixel->GetLoGainSamples() + slice - nh;
469
470 return *ptr;
471}
472
473// --------------------------------------------------------------------------
474//
475// Fill the MPedestalCam container with the signal mean and rms for the event.
476// Store the measured signal in arrays fSumx and fSumx2 so that we can
477// calculate the overall mean and rms in the PostProcess()
478//
479Int_t MPedCalcFromLoGain::Process()
480{
481 MRawEvtPixelIter pixel(fRawEvt);
482
483 while (pixel.Next())
484 {
485 const UInt_t idx = pixel.GetPixelId();
486
487 UShort_t max = 0;
488 UShort_t min = (UShort_t)-1;
489
490 // Find the maximum and minimum signal per slice in the high gain window
491 for (Int_t slice=fCheckWinFirst; slice<=fCheckWinLast; slice++)
492 {
493 const UShort_t svalue = GetSlice(&pixel,slice);
494 if (svalue > max)
495 max = svalue;
496 if (svalue < min)
497 min = svalue;
498 }
499
500 // If the maximum in the high gain window is smaller than
501 if (max-min>=fMaxSignalVar || max>=250)
502 continue;
503
504 UInt_t sum = 0;
505 UInt_t sqr = 0;
506
507 //extract pedestal
508 for(Int_t slice=fExtractWinFirst; slice<=fExtractWinLast; slice++)
509 {
510 const UInt_t svalue = GetSlice(&pixel,slice);
511 sum += svalue;
512 sqr += svalue*svalue;
513 }
514
515 fSumx[idx] += sum;
516 fSumx2[idx] += sum*sum;
517
518 fNumEventsUsed[idx]++;
519
520 // Calculate the amplitude of the 150MHz "AB" noise
521 UShort_t abFlag = (fExtractWinFirst + pixel.HasABFlag()) & 0x1;
522
523 for (UShort_t islice=fExtractWinFirst; islice<=fExtractWinLast; islice+=2)
524 {
525 const UShort_t sliceAB0 = islice + abFlag;
526 const UShort_t sliceAB1 = islice - abFlag + 1;
527 fSumAB0[idx] += GetSlice(&pixel, sliceAB0);
528 fSumAB1[idx] += GetSlice(&pixel, sliceAB1);
529 }
530
531 if (!fPedestalUpdate || (UInt_t)fNumEventsUsed[idx]<fNumEventsDump)
532 continue;
533
534 Calc(fNumEventsDump, idx);
535
536 fNumEventsUsed[idx]=0;
537 fSumx[idx]=0;
538 fSumx2[idx]=0;
539 fSumAB0[idx]=0;
540 fSumAB1[idx]=0;
541 }
542
543 if (fPedestalUpdate)
544 fPedestals->SetReadyToSave();
545
546 return kTRUE;
547}
548
549// --------------------------------------------------------------------------
550//
551// Compute signal mean and rms in the whole run and store it in MPedestalCam
552//
553Int_t MPedCalcFromLoGain::PostProcess()
554{
555 // Compute pedestals and rms from the whole run
556 if (fPedestalUpdate)
557 return kTRUE;
558
559 *fLog << flush << inf << "Calculating Pedestals..." << flush;
560
561 const Int_t npix = fGeom->GetNumPixels();
562 for (Int_t idx=0; idx<npix; idx++)
563 {
564 const ULong_t n = fNumEventsUsed[idx];
565 if (n>1)
566 Calc(n, idx);
567 }
568
569 fPedestals->SetReadyToSave();
570 return kTRUE;
571}
572
573
574// --------------------------------------------------------------------------
575//
576// The following resources are available:
577// FirstCheckWindowSlice: 0
578// LastCheckWindowSlice: 29
579// ExtractWindowFirst: 15
580// ExtractWindowSize: 6
581// NumEventsDump: 500
582// MaxSignalVar: 40
583// PedestalUpdate: yes
584//
585Int_t MPedCalcFromLoGain::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
586{
587 Bool_t rc=kFALSE;
588
589 // Find resources for CheckWindow
590 Int_t fs = fCheckWinFirst;
591 Int_t ls = fCheckWinLast;
592 if (IsEnvDefined(env, prefix, "CheckWinFirst", print))
593 {
594 fs = GetEnvValue(env, prefix, "CheckWinFirst", fs);
595 rc = kTRUE;
596 }
597 if (IsEnvDefined(env, prefix, "CheckWinLast", print))
598 {
599 ls = GetEnvValue(env, prefix, "CheckWinLast", ls);
600 rc = kTRUE;
601 }
602
603 SetCheckRange(fs,ls);
604
605 // Find resources for ExtractWindow
606 Int_t lw = fExtractWinSize;
607 Int_t wf = fExtractWinFirst;
608 if (IsEnvDefined(env, prefix, "ExtractWinSize", print))
609 {
610 lw = GetEnvValue(env, prefix, "ExtractWinSize", lw);
611 rc = kTRUE;
612 }
613
614 if (IsEnvDefined(env, prefix, "ExtractWinFirst", print))
615 {
616 wf = GetEnvValue(env, prefix, "ExtractWinFirst", wf);
617 rc = kTRUE;
618 }
619 SetExtractWindow(wf,lw);
620
621 // find resource for numeventsdump
622 if (IsEnvDefined(env, prefix, "NumEventsDump", print))
623 {
624 SetNumEventsDump(GetEnvValue(env, prefix, "NumEventsDump", (Int_t)fNumEventsDump));
625 rc = kTRUE;
626 }
627
628 // find resource for maximum signal variation
629 if (IsEnvDefined(env, prefix, "MaxSignalVar", print))
630 {
631 SetMaxSignalVar(GetEnvValue(env, prefix, "MaxSignalVar", fMaxSignalVar));
632 rc = kTRUE;
633 }
634
635 // find resource for pedestal update
636 if (IsEnvDefined(env, prefix, "PedestalUpdate", print))
637 {
638 SetPedestalUpdate(GetEnvValue(env, prefix, "PedestalUpdate", fPedestalUpdate));
639 rc = kTRUE;
640 }
641
642 return rc;
643}
644
645void MPedCalcFromLoGain::Print(Option_t *o) const
646{
647 *fLog << GetDescriptor() << ":" << endl;
648 *fLog << "Parameters used for pedestal calculation from " << fNamePedestalCam << ":"<<endl;
649 *fLog << "CheckWindow from slice " << fCheckWinFirst << " to " << fCheckWinLast << " included." << endl;
650 *fLog << "ExtractWindow from slice " << fExtractWinFirst << " to " << fExtractWinLast << " included." << endl;
651 *fLog << "Max allowed signal variation: " << fMaxSignalVar << endl;
652 *fLog << "Number of events to calculate pedestal from: " << fNumEventsDump << endl;
653 *fLog << "Pedestal Update is " << (fPedestalUpdate?"on":"off") << endl;
654}
Note: See TracBrowser for help on using the repository browser.