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

Last change on this file since 5338 was 5338, checked in by otte, 20 years ago
*** empty log message ***
File size: 16.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//
34//
35// MPedCalcLoGain
36//
37//
38// This task is devide form MPedCalcPedRun, described below. However, It
39// calculates the pedstals using the low gain slices, whenever the difference
40// between the highest and the lowest slice in the high gain
41// slices is below a given threshold (SetMaxHiGainVar). In this case the receiver
42// boards do not switch to lo gain and the so called lo gain slices are actually
43// high gain slices.
44//
45// MPedCalcLoGain also fills the ABoffset in MPedestalPix which allows to correct
46// the 150 MHz clock noise.
47//
48// This task takes a pedestal run file and fills MPedestalCam during
49// the Process() with the pedestal and rms computed in an event basis.
50// In the PostProcess() MPedestalCam is finally filled with the pedestal
51// mean and rms computed in a run basis.
52// More than one run (file) can be merged
53//
54// MPedCalcPedRun applies the following formula (1):
55//
56// Pedestal per slice = sum(x_i) / n / slices
57// PedRMS per slice = Sqrt( ( sum(x_i^2) - sum(x_i)^2/n ) / n-1 / slices )
58//
59// where x_i is the sum of "slices" FADC slices and sum means the sum over all
60// events. "n" is the number of events, "slices" is the number of summed FADC samples.
61//
62// Note that the slice-to-slice fluctuations are not Gaussian, but Poissonian, thus
63// asymmetric and they are correlated.
64//
65// It is important to know that the Pedestal per slice and PedRMS per slice depend
66// on the number of used FADC slices, as seen in the following plots:
67//
68//Begin_Html
69/*
70<img src="images/PedestalStudyInner.gif">
71*/
72//End_Html
73//
74//Begin_Html
75/*
76<img src="images/PedestalStudyOuter.gif">
77*/
78//
79// The plots show the inner and outer pixels, respectivly and have the following meaning:
80//
81// 1) The calculated mean pedestal per slice (from MPedCalcFromLoGain)
82// 2) The fitted mean pedestal per slice (from MHPedestalCam)
83// 3) The calculated pedestal RMS per slice (from MPedCalcFromLoGain)
84// 4) The fitted sigma of the pedestal distribution per slice
85// (from MHPedestalCam)
86// 5) The relative difference between calculation and histogram fit
87// for the mean
88// 6) The relative difference between calculation and histogram fit
89// for the sigma or RMS, respectively.
90//
91// The calculated means do not change significantly except for the case of 2 slices,
92// however the RMS changes from 5.7 per slice in the case of 2 extracted slices
93// to 8.3 per slice in the case of 26 extracted slices. This change is very significant.
94//
95// The plots have been produced on run 20123. You can reproduce them using
96// the macro pedestalstudies.C
97//
98// Usage of this class:
99// ====================
100//
101//
102// fCheckWinFirst = fgCheckWinFirst = 0
103// fHiGainLast = fgCheckWinLast = 29
104// fExtractWinFirst = fgExtractWinFirst = 0
105// fExtractWinSize = fgExtractWinSize = 6
106// fMaxSignalVar = fgMaxSignalVar = 40;
107//
108// Call:
109// SetCheckRange(fCheckWinFirst,fCheckWinLast);
110// to set the Window in which a signal is searched
111//
112// SetExtractWindow(fExtractWindFirst,fExtractWinSize);
113// to set the Window from which a signal is extracted
114//
115// SetMaxSignalVar(fMaxSignalVar);
116// set the maximum allowed difference between maximum and minimal signal in CheckWindow
117//
118// Variables:
119// fgCheckWinFirst; First FADC slice to check for signal (currently set to: 0)
120// fgCheckWinLast: Last FADC slice to check for signal (currently set to: 29)
121// fgExtractWinFirst: First FADC slice to be used for pedestal extraction (currently set to: 15)
122// fgExtractWinSize: Window size in slices used for pedestal extraction (currently set to: 6)
123// fgMaxSignalVar: The maximum difference between the highest and lowest slice
124// in the check window allowed in order to use event
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 Int_t MPedCalcFromLoGain::fgCheckWinFirst = 0;
160const Int_t MPedCalcFromLoGain::fgCheckWinLast = 29;
161const Int_t MPedCalcFromLoGain::fgExtractWinFirst = 15;
162const Int_t MPedCalcFromLoGain::fgExtractWinSize = 6;
163const Int_t MPedCalcFromLoGain::fgMaxSignalVar = 40;
164
165// --------------------------------------------------------------------------
166//
167// Default constructor:
168//
169// Sets:
170// - all pointers to NULL
171//
172// Calls:
173// - AddToBranchList("fHiGainPixId");
174// - AddToBranchList("fHiGainFadcSamples");
175// - SetCheckRange(fgCheckWinFirst, fgCheckWinLast, fgExtractWinFirst, fgExtractWinSize)
176// - Clear()
177//
178MPedCalcFromLoGain::MPedCalcFromLoGain(const char *name, const char *title)
179 : fGeom(NULL), fPedContainerName("MPedestalCam")
180{
181 fName = name ? name : "MPedCalcFromLoGain";
182 fTitle = title ? title : "Task to calculate pedestals from pedestal runs raw data";
183
184 AddToBranchList("fHiGainPixId");
185 AddToBranchList("fHiGainFadcSamples");
186
187 SetCheckRange(fgCheckWinFirst, fgCheckWinLast);
188 SetExtractWindow(fgExtractWinFirst, fgExtractWinSize);
189
190 SetMaxSignalVar(fgMaxSignalVar);
191 SetPedestalUpdate(kTRUE);
192 Clear();
193}
194
195// --------------------------------------------------------------------------
196//
197// Sets:
198// - fNumSamplesTot to 0
199// - fRawEvt to NULL
200// - fRunHeader to NULL
201// - fPedestals to NULL
202//
203void MPedCalcFromLoGain::Clear(const Option_t *o)
204{
205 fRawEvt = NULL;
206 fRunHeader = NULL;
207 fPedestals = NULL;
208
209 // If the size is yet set, set the size
210 if (fSumx.GetSize()>0)
211 {
212 // Reset contents of arrays.
213 fSumx.Reset();
214 fSumx2.Reset();
215 fSumAB0.Reset();
216 fSumAB1.Reset();
217 fNumEventsUsed.Reset();
218 fTotalCounter.Reset();
219 }
220}
221
222// --------------------------------------------------------------------------
223//
224// SetCheckRange:
225//
226// Calls:
227// - MExtractor::SetCheckRange(hifirst,hilast,lofirst,lolast);
228// - SetExtractWindow(fWindowSizeHiGain,fExtractWinSize);
229//
230void MPedCalcFromLoGain::SetCheckRange(Int_t chfirst, Int_t chlast)
231{
232
233
234 if(chfirst<0){
235 *fLog << warn << GetDescriptor()
236 << Form(": First slice in window to check for Signal <0, adjust:")<< endl;
237 exit(-1);
238 }
239
240 if(chlast<=chfirst){
241 *fLog << warn << GetDescriptor()
242 << Form(": Last slice in Check window smaller than first slice in window, adjust:")<< endl;
243 exit(-1);
244 }
245
246 fCheckWinFirst = chfirst;
247 fCheckWinLast = chlast;
248
249
250}
251
252// --------------------------------------------------------------------------
253//
254void MPedCalcFromLoGain::SetMaxSignalVar(Int_t maxvar)
255{
256 fMaxSignalVar = maxvar;
257}
258
259// --------------------------------------------------------------------------
260//
261// Checks:
262// - if a window is odd
263//
264
265void MPedCalcFromLoGain::SetExtractWindow(Int_t windowf, Int_t windows)
266{
267
268 if(windowf<0){
269 *fLog << warn << GetDescriptor()
270 << Form(": First slice in Extract window has to be >0, adjust:")<< endl;
271 exit(-1);
272 }
273
274 Int_t odd = windows & 0x1;
275
276
277 if (odd||(windows==0)){
278 *fLog << warn << GetDescriptor() << ": Extract window size has to be even and larger 0, adjust!"<< endl;
279 exit(-1);
280 }
281
282 fExtractWinSize = windows;
283 fExtractWinFirst = windowf;
284 fExtractWinLast = fExtractWinFirst+fExtractWinSize-1;
285
286
287}
288
289// --------------------------------------------------------------------------
290//
291// Look for the following input containers:
292//
293// - MRawEvtData
294// - MRawRunHeader
295// - MGeomCam
296//
297// The following output containers are also searched and created if
298// they were not found:
299//
300// - MPedestalCam
301//
302Int_t MPedCalcFromLoGain::PreProcess(MParList *pList)
303{
304 Clear();
305
306 fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
307 if (!fRawEvt)
308 {
309 *fLog << err << "MRawEvtData not found... aborting." << endl;
310 return kFALSE;
311 }
312
313 fRunHeader = (MRawRunHeader*)pList->FindObject(AddSerialNumber("MRawRunHeader"));
314 if (!fRunHeader)
315 {
316 *fLog << err << AddSerialNumber("MRawRunHeader") << " not found... aborting." << endl;
317 return kFALSE;
318 }
319
320 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
321 if (!fGeom)
322 {
323 *fLog << err << "MGeomCam not found... aborting." << endl;
324 return kFALSE;
325 }
326
327
328 fPedestals = (MPedestalCam*)pList->FindCreateObj("MPedestalCam", AddSerialNumber(fPedContainerName));
329 if (!fPedestals)
330 return kFALSE;
331
332 return kTRUE;
333}
334
335// --------------------------------------------------------------------------
336//
337// The ReInit searches for:
338// - MRawRunHeader::GetNumSamplesHiGain()
339// - MRawRunHeader::GetNumSamplesLoGain()
340//
341
342Bool_t MPedCalcFromLoGain::ReInit(MParList *pList)
343{
344
345 Int_t lastavailableslice = (Int_t)fRunHeader->GetNumSamplesHiGain()+(Int_t)fRunHeader->GetNumSamplesLoGain()-1;
346
347 Int_t lastextractslice = fExtractWinSize+ fExtractWinFirst - 1;
348
349 if ( lastextractslice > lastavailableslice)//changed to override check
350 {
351 *fLog << endl;
352 *fLog << warn << GetDescriptor()
353 << Form(": Selected Extract Window ranges out of the available limits adjust. Last available slice is %2i",
354 lastavailableslice) << endl;
355 exit(-1);
356 }
357
358 if ( fCheckWinLast > lastavailableslice)//changed to override check
359 {
360 *fLog << endl;
361 *fLog << warn << GetDescriptor()
362 << Form(": Last Check Window slice is out of the available limits adjust. Last available slice is %2i",
363 lastavailableslice) << endl;
364 exit(-1);
365 }
366
367
368
369 // If the size is not yet set, set the size
370 if (fSumx.GetSize()==0)
371 {
372 const Int_t npixels = fPedestals->GetSize();
373
374 fSumx. Set(npixels);
375 fSumx2.Set(npixels);
376 fSumAB0.Set(npixels);
377 fSumAB1.Set(npixels);
378 fNumEventsUsed.Set(npixels);
379 fTotalCounter.Set(npixels);
380
381 // Reset contents of arrays.
382 fSumx.Reset();
383 fSumx2.Reset();
384 fSumAB0.Reset();
385 fSumAB1.Reset();
386 fNumEventsUsed.Reset();
387 fTotalCounter.Reset();
388 }
389
390 return kTRUE;
391}
392
393void MPedCalcFromLoGain::Calc(ULong_t n, UInt_t idx)
394{
395 const ULong_t nsamplestot = n*fExtractWinSize;
396
397 const Float_t sum = fSumx.At(idx);
398 const Float_t sum2 = fSumx2.At(idx);
399 const Float_t ped = sum/nsamplestot;
400
401 // 1. Calculate the Variance of the sums:
402 Float_t var = (sum2-sum*sum/n)/(n-1.);
403
404 // 2. Scale the variance to the number of slices:
405 var /= (Float_t)(fExtractWinSize);
406
407 // 3. Calculate the RMS from the Variance:
408 const Float_t rms = var<0 ? 0 : TMath::Sqrt(var);
409
410 // 4. Calculate the amplitude of the 150MHz "AB" noise
411 const Float_t abOffs = (fSumAB0[idx] - fSumAB1[idx]) / nsamplestot;
412
413 (*fPedestals)[idx].Set(ped, rms, abOffs, n);
414
415 fTotalCounter[idx]++;
416}
417
418UInt_t MPedCalcFromLoGain::GetSlice(MRawEvtPixelIter *pixel, UInt_t slice)
419{
420
421 const UInt_t nh = (Int_t)fRunHeader->GetNumSamplesHiGain();
422
423 Byte_t *ptr;
424
425 if(slice>=nh)
426 {
427 ptr = pixel->GetLoGainSamples();
428 ptr += slice - nh;
429 }
430 else
431 {
432 ptr = pixel->GetHiGainSamples();
433 ptr += slice;
434 }
435
436 return *ptr;
437}
438
439// --------------------------------------------------------------------------
440//
441// Fill the MPedestalCam container with the signal mean and rms for the event.
442// Store the measured signal in arrays fSumx and fSumx2 so that we can
443// calculate the overall mean and rms in the PostProcess()
444//
445Int_t MPedCalcFromLoGain::Process()
446{
447 MRawEvtPixelIter pixel(fRawEvt);
448
449 while (pixel.Next())
450 {
451
452 const UInt_t idx = pixel.GetPixelId();
453
454 Int_t max = 0;
455 Int_t min = 1025;
456
457 // Find the maximum and minimum signal per slice in the high gain window
458
459 for(Int_t slice=fCheckWinFirst; slice<=fCheckWinLast; slice++)
460 {
461
462 Int_t svalue = GetSlice(&pixel,slice);
463 if (svalue > max) {
464 max = svalue;
465 }
466 if (svalue < min) {
467 min = svalue;
468 }
469 }
470
471 // If the maximum in the high gain window is smaller than
472 if (max-min>=fMaxSignalVar || max>=250)
473 continue;
474
475 UInt_t sum = 0;
476 UInt_t sqr = 0;
477
478 //extract pedestal
479 for(Int_t slice=fExtractWinFirst; slice<=fExtractWinLast; slice++)
480 {
481 UInt_t svalue = GetSlice(&pixel,slice);
482 sum += svalue;
483 sqr += svalue*svalue;
484 }
485
486 const Float_t msum = (Float_t)sum;
487 const Float_t sqrsum = msum*msum;
488
489 fSumx[idx] += msum;
490 fSumx2[idx] += sqrsum;
491 fNumEventsUsed[idx]++;
492
493 // Calculate the amplitude of the 150MHz "AB" noise
494
495 Int_t abFlag = (fExtractWinFirst + pixel.HasABFlag()) & 0x1;
496
497 //cout << " MPedCalcFromLoGain: idx: " << idx << " abFlag: " << abFlag << endl;
498
499 for (Int_t islice=fExtractWinFirst; islice<=fExtractWinLast; islice+=2)
500 {
501 Int_t sliceAB0 = islice + abFlag;
502 Int_t sliceAB1 = islice - abFlag + 1;
503 fSumAB0[idx] += GetSlice(&pixel, sliceAB0);
504 fSumAB1[idx] += GetSlice(&pixel, sliceAB1);
505 }
506
507 if (!fPedestalUpdate || fNumEventsUsed[idx]<fNumEventsDump)
508 continue;
509
510 Calc(fNumEventsDump, idx);
511
512 fNumEventsUsed[idx]=0;
513 fSumx[idx]=0;
514 fSumx2[idx]=0;
515 fSumAB0[idx]=0;
516 fSumAB1[idx]=0;
517 }
518
519 if (fPedestalUpdate)
520 fPedestals->SetReadyToSave();
521
522 return kTRUE;
523}
524
525// --------------------------------------------------------------------------
526//
527// Compute signal mean and rms in the whole run and store it in MPedestalCam
528//
529Int_t MPedCalcFromLoGain::PostProcess()
530{
531 // Compute pedestals and rms from the whole run
532 if (fPedestalUpdate)
533 return kTRUE;
534
535 *fLog << flush << inf << "Calculating Pedestals..." << flush;
536
537 const Int_t npix = fGeom->GetNumPixels();
538 for (Int_t idx=0; idx<npix; idx++)
539 {
540 const ULong_t n = fNumEventsUsed[idx];
541 if (n>1)
542 Calc(n, idx);
543 }
544
545 fPedestals->SetReadyToSave();
546 return kTRUE;
547}
548
549Int_t MPedCalcFromLoGain::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
550{
551 if (MExtractor::ReadEnv(env, prefix, print)==kERROR)
552 return kERROR;
553
554 Bool_t rc=kFALSE;
555
556 Int_t lw = fExtractWinSize;
557 Int_t wf = fExtractWinFirst;
558
559 if (IsEnvDefined(env, prefix, "ExtractWindowSize", print))
560 {
561 lw = GetEnvValue(env, prefix, "ExtractWindowSize", lw);
562 wf = GetEnvValue(env, prefix, "ExtractWindowFirst", wf);
563 rc = kTRUE;
564
565 }
566 SetExtractWindow(wf,lw);
567
568 Int_t num = fNumEventsDump;
569 if (IsEnvDefined(env, prefix, "NumEventsDump", print))
570 {
571 num = GetEnvValue(env, prefix, "NumEventsDump", num);
572 rc = kTRUE;
573 }
574 SetNumEventsDump(num);
575
576 Byte_t max = fMaxSignalVar;
577 if (IsEnvDefined(env, prefix, "MaxSignalVar", print))
578 {
579 max = GetEnvValue(env, prefix, "MaxSignalVar", max);
580 rc = kTRUE;
581 }
582 SetMaxSignalVar(max);
583
584 Bool_t upd = fPedestalUpdate;
585 if (IsEnvDefined(env, prefix, "PedestalUpdate", print))
586 {
587 upd = GetEnvValue(env, prefix, "PedestalUpdate", upd);
588 rc = kTRUE;
589 }
590 SetPedestalUpdate(upd);
591
592 return rc;
593}
Note: See TracBrowser for help on using the repository browser.