source: trunk/MagicSoft/Mars/mpedestal/MPedCalcPedRun.cc@ 3967

Last change on this file since 3967 was 3931, checked in by gaug, 21 years ago
*** empty log message ***
File size: 18.4 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!
24! Copyright: MAGIC Software Development, 2000-2004
25!
26!
27\* ======================================================================== */
28/////////////////////////////////////////////////////////////////////////////
29//
30// MPedCalcPedRun
31//
32// This task takes a pedestal run file and fills MPedestalCam during
33// the Process() with the pedestal and rms computed in an event basis.
34// In the PostProcess() MPedestalCam is finally filled with the pedestal
35// mean and rms computed in a run basis.
36// More than one run (file) can be merged
37//
38//
39// Actually, MPedCalcPedRun applies the following formula (1):
40//
41// PedRMS = Sqrt( (sum(x_i^2) - sum(x_i)^2/n) / n-1 / 14 )
42//
43// where x_i is the sum of 14 FADC slices and sum means the sum over all
44// events, n is the number of events.
45//
46// For a high number of events, this formula is equivalent to formula (2):
47//
48// PedRMS = Sqrt( (<x_i*x_i> - <x_i>*<x_i>) / 14 )
49//
50// where <> is the mean over all events and x_i again the sum over the 14
51// slices.
52//
53// If you assume statistical equivalence of all slices (say, all have equal
54// offset and are not correlated and fluctuate Gaussian), it should also be
55// equivalent to (old formula) (3):
56//
57// PedRMS = Sqrt( (<p_i*p_i> - <p_i>*<p_i>) )
58//
59// which is the RMS per slice of a single slice (p_i) and
60// <> the mean over the total number of measurements, i.e. n*14.
61//
62// If we assume that at least our pairs fluctuate independently and Gaussian,
63// then we can use the actual formula (1) in order to get
64// fluctuations of pairs by the transformation:
65//
66// PedRMS/pair = PedRMS (form. (3)) * Sqrt(2)
67//
68// (However, we know that our slice-to-slice fluctuations are not Gaussian
69// (and moreover asymmetric) and that they are also correlated.)
70//
71// Call: SetRange(higainfirst, higainlast, logainfirst, logainlast)
72// to modify the ranges in which the window is allowed to move.
73// Defaults are:
74//
75// fHiGainFirst = fgHiGainFirst = 0
76// fHiGainLast = fgHiGainLast = 14
77// fLoGainFirst = fgLoGainFirst = 0
78// fLoGainLast = fgLoGainLast = 14
79//
80// Call: SetWindowSize(windowhigain, windowlogain)
81// to modify the sliding window widths. Windows have to be an even number.
82// In case of odd numbers, the window will be modified.
83//
84// Defaults are:
85//
86// fHiGainWindowSize = fgHiGainWindowSize = 14
87// fLoGainWindowSize = fgLoGainWindowSize = 0
88//
89//
90// Input Containers:
91// MRawEvtData
92// MRawRunHeader
93// MGeomCam
94//
95// Output Containers:
96// MPedestalCam
97//
98//
99/////////////////////////////////////////////////////////////////////////////
100#include "MPedCalcPedRun.h"
101#include "MExtractor.h"
102
103#include "MParList.h"
104
105#include "MLog.h"
106#include "MLogManip.h"
107
108#include "MRawRunHeader.h"
109#include "MRawEvtPixelIter.h"
110#include "MRawEvtData.h"
111
112#include "MPedestalPix.h"
113#include "MPedestalCam.h"
114
115#include "MGeomPix.h"
116#include "MGeomCam.h"
117
118#include "MBadPixelsPix.h"
119#include "MBadPixelsCam.h"
120
121#include "MGeomCamMagic.h"
122
123ClassImp(MPedCalcPedRun);
124
125using namespace std;
126
127const Byte_t MPedCalcPedRun::fgHiGainFirst = 0;
128const Byte_t MPedCalcPedRun::fgHiGainLast = 14;
129const Byte_t MPedCalcPedRun::fgLoGainFirst = 0;
130const Byte_t MPedCalcPedRun::fgLoGainLast = 14;
131const Byte_t MPedCalcPedRun::fgHiGainWindowSize = 14;
132const Byte_t MPedCalcPedRun::fgLoGainWindowSize = 0;
133// --------------------------------------------------------------------------
134//
135// Default constructor:
136//
137// Sets:
138// - all pointers to NULL
139// - fWindowSizeHiGain to fgHiGainWindowSize
140// - fWindowSizeLoGain to fgLoGainWindowSize
141//
142// Calls:
143// - AddToBranchList("fHiGainPixId");
144// - AddToBranchList("fHiGainFadcSamples");
145// - SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast)
146// - Clear()
147//
148MPedCalcPedRun::MPedCalcPedRun(const char *name, const char *title)
149 : fWindowSizeHiGain(fgHiGainWindowSize),
150 fWindowSizeLoGain(fgLoGainWindowSize),
151 fGeom(NULL),fBad(NULL)
152{
153 fName = name ? name : "MPedCalcPedRun";
154 fTitle = title ? title : "Task to calculate pedestals from pedestal runs raw data";
155
156 AddToBranchList("fHiGainPixId");
157 AddToBranchList("fHiGainFadcSamples");
158
159 SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast);
160 Clear();
161}
162
163// --------------------------------------------------------------------------
164//
165// Sets:
166// - fNumSamplesTot to 0
167// - fRawEvt to NULL
168// - fRunHeader to NULL
169// - fPedestals to NULL
170//
171void MPedCalcPedRun::Clear(const Option_t *o)
172{
173
174 fNumSamplesTot = 0;
175
176 fRawEvt = NULL;
177 fRunHeader = NULL;
178 fPedestals = NULL;
179}
180
181// --------------------------------------------------------------------------
182//
183// SetRange:
184//
185// Calls:
186// - MExtractor::SetRange(hifirst,hilast,lofirst,lolast);
187// - SetWindowSize(fWindowSizeHiGain,fWindowSizeLoGain);
188//
189void MPedCalcPedRun::SetRange(Byte_t hifirst, Byte_t hilast, Byte_t lofirst, Byte_t lolast)
190{
191
192 MExtractor::SetRange(hifirst,hilast,lofirst,lolast);
193
194 //
195 // Redo the checks if the window is still inside the ranges
196 //
197 SetWindowSize(fWindowSizeHiGain,fWindowSizeLoGain);
198
199}
200
201
202// --------------------------------------------------------------------------
203//
204// Checks:
205// - if a window is odd, subtract one
206// - if a window is bigger than the one defined by the ranges, set it to the available range
207// - if a window is smaller than 2, set it to 2
208//
209// Sets:
210// - fNumHiGainSamples to: (Float_t)fWindowSizeHiGain
211// - fNumLoGainSamples to: (Float_t)fWindowSizeLoGain
212// - fSqrtHiGainSamples to: TMath::Sqrt(fNumHiGainSamples)
213// - fSqrtLoGainSamples to: TMath::Sqrt(fNumLoGainSamples)
214//
215void MPedCalcPedRun::SetWindowSize(Byte_t windowh, Byte_t windowl)
216{
217
218 fWindowSizeHiGain = windowh & ~1;
219 fWindowSizeLoGain = windowl & ~1;
220
221 if (fWindowSizeHiGain != windowh)
222 *fLog << warn << GetDescriptor() << ": Hi Gain window size has to be even, set to: "
223 << int(fWindowSizeHiGain) << " samples " << endl;
224
225 if (fWindowSizeLoGain != windowl)
226 *fLog << warn << GetDescriptor() << ": Lo Gain window size has to be even, set to: "
227 << int(fWindowSizeLoGain) << " samples " << endl;
228
229 const Byte_t availhirange = (fHiGainLast-fHiGainFirst+1) & ~1;
230 const Byte_t availlorange = (fLoGainLast-fLoGainFirst+1) & ~1;
231
232 if (fWindowSizeHiGain > availhirange)
233 {
234 *fLog << warn << GetDescriptor()
235 << Form("%s%2i%s%2i%s%2i%s",": Hi Gain window size: ",(int)fWindowSizeHiGain,
236 " is bigger than available range: [",(int)fHiGainFirst,",",(int)fHiGainLast,"]") << endl;
237 *fLog << warn << GetDescriptor()
238 << ": Will set window size to: " << (int)availhirange << endl;
239 fWindowSizeHiGain = availhirange;
240 }
241
242 if (fWindowSizeLoGain > availlorange)
243 {
244 *fLog << warn << GetDescriptor()
245 << Form("%s%2i%s%2i%s%2i%s",": Lo Gain window size: ",(int)fWindowSizeLoGain,
246 " is bigger than available range: [",(int)fLoGainFirst,",",(int)fLoGainLast,"]") << endl;
247 *fLog << warn << GetDescriptor()
248 << ": Will set window size to: " << (int)availlorange << endl;
249 fWindowSizeLoGain = availlorange;
250 }
251
252
253 fNumHiGainSamples = (Float_t)fWindowSizeHiGain;
254 fNumLoGainSamples = (Float_t)fWindowSizeLoGain;
255
256 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
257 fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
258
259}
260
261
262
263// --------------------------------------------------------------------------
264//
265// Look for the following input containers:
266//
267// - MRawEvtData
268// - MRawRunHeader
269// - MGeomCam
270// - MBadPixelsCam
271//
272// The following output containers are also searched and created if
273// they were not found:
274//
275// - MPedestalCam
276//
277Int_t MPedCalcPedRun::PreProcess( MParList *pList )
278{
279
280 Clear();
281
282 fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
283 if (!fRawEvt)
284 {
285 *fLog << err << "MRawEvtData not found... aborting." << endl;
286 return kFALSE;
287 }
288
289 fRunHeader = (MRawRunHeader*)pList->FindObject(AddSerialNumber("MRawRunHeader"));
290 if (!fRunHeader)
291 {
292 *fLog << err << AddSerialNumber("MRawRunHeader") << " not found... aborting." << endl;
293 return kFALSE;
294 }
295
296 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
297 if (!fGeom)
298 {
299 *fLog << err << "MGeomCam not found... aborting." << endl;
300 return kFALSE;
301 }
302
303 fPedestals = (MPedestalCam*)pList->FindCreateObj("MPedestalCam");
304 if (!fPedestals)
305 return kFALSE;
306
307 fBad = (MBadPixelsCam*)pList->FindObject("MBadPixelsCam");
308
309 return kTRUE;
310}
311
312// --------------------------------------------------------------------------
313//
314// The ReInit searches for:
315// - MRawRunHeader::GetNumSamplesHiGain()
316// - MRawRunHeader::GetNumSamplesLoGain()
317//
318// In case that the variables fHiGainLast and fLoGainLast are smaller than
319// the even part of the number of samples obtained from the run header, a
320// warning is given an the range is set back accordingly. A call to:
321// - SetRange(fHiGainFirst, fHiGainLast-diff, fLoGainFirst, fLoGainLast) or
322// - SetRange(fHiGainFirst, fHiGainLast, fLoGainFirst, fLoGainLast-diff)
323// is performed in that case. The variable diff means here the difference
324// between the requested range (fHiGainLast) and the available one. Note that
325// the functions SetRange() are mostly overloaded and perform more checks,
326// modifying the ranges again, if necessary.
327//
328// A loop over the MBadPixelsCam is performed and bad pixels are set
329// to MPedestalPix::SetValid(kFALSE);
330//
331Bool_t MPedCalcPedRun::ReInit(MParList *pList)
332{
333
334 Int_t lastdesired = (Int_t)fHiGainLast;
335 Int_t lastavailable = (Int_t)fRunHeader->GetNumSamplesHiGain()-1;
336
337 if (lastdesired > lastavailable)
338 {
339 const Int_t diff = lastdesired - lastavailable;
340 *fLog << endl;
341 *fLog << warn << GetDescriptor()
342 << Form("%s%2i%s%2i%s%2i%s",": Selected Hi Gain FADC Window [",
343 (int)fHiGainFirst,",",lastdesired,
344 "] ranges out of the available limits: [0,",lastavailable,"].") << endl;
345 *fLog << GetDescriptor() << ": Will reduce the upper edge to " << (int)(fHiGainLast - diff) << endl;
346 SetRange(fHiGainFirst, fHiGainLast-diff, fLoGainFirst, fLoGainLast);
347 }
348
349 lastdesired = (Int_t)(fLoGainLast);
350 lastavailable = (Int_t)fRunHeader->GetNumSamplesLoGain()-1;
351
352 if (lastdesired > lastavailable)
353 {
354 const Int_t diff = lastdesired - lastavailable;
355 *fLog << endl;
356 *fLog << warn << GetDescriptor()
357 << Form("%s%2i%s%2i%s%2i%s",": Selected Lo Gain FADC Window [",
358 (int)fLoGainFirst,",",lastdesired,
359 "] ranges out of the available limits: [0,",lastavailable,"].") << endl;
360 *fLog << GetDescriptor() << ": Will reduce the upper edge to " << (int)(fLoGainLast - diff) << endl;
361 SetRange(fHiGainFirst, fHiGainLast, fLoGainFirst, fLoGainLast-diff);
362 }
363
364 Int_t npixels = fPedestals->GetSize();
365 Int_t areas = fPedestals->GetAverageAreas();
366 Int_t sectors = fPedestals->GetAverageSectors();
367
368 if (fSumx.GetSize()==0)
369 {
370 fSumx. Set(npixels);
371 fSumx2.Set(npixels);
372
373 fAreaSumx. Set(areas);
374 fAreaSumx2.Set(areas);
375 fAreaValid.Set(areas);
376
377 fSectorSumx. Set(sectors);
378 fSectorSumx2.Set(sectors);
379 fSectorValid.Set(sectors);
380
381 fSumx.Reset();
382 fSumx2.Reset();
383 }
384
385 if (fWindowSizeHiGain == 0 && fWindowSizeLoGain == 0)
386 {
387 *fLog << err << GetDescriptor()
388 << ": Number of extracted Slices is 0, abort ... " << endl;
389 return kFALSE;
390 }
391
392
393 *fLog << endl;
394 *fLog << inf << GetDescriptor() << ": Taking " << (int)fWindowSizeHiGain
395 << " HiGain FADC samples starting with slice: " << (int)fHiGainFirst << endl;
396 *fLog << inf << GetDescriptor() << ": Taking " << (int)fWindowSizeLoGain
397 << " LoGain FADC samples starting with slice: " << (int)fLoGainFirst << endl;
398
399
400 if (fBad)
401 {
402 const Int_t nbads = fBad->GetSize();
403 for (Int_t i=0; i<(nbads>npixels?npixels:nbads);i++)
404 if ((*fBad)[i].IsBad())
405 (*fPedestals)[i].SetValid(kFALSE);
406 }
407
408 return kTRUE;
409
410}
411// --------------------------------------------------------------------------
412//
413// Fill the MPedestalCam container with the signal mean and rms for the event.
414// Store the measured signal in arrays fSumx and fSumx2 so that we can
415// calculate the overall mean and rms in the PostProcess()
416//
417Int_t MPedCalcPedRun::Process()
418{
419
420 MRawEvtPixelIter pixel(fRawEvt);
421
422 while (pixel.Next())
423 {
424
425 const UInt_t idx = pixel.GetPixelId();
426 const UInt_t aidx = (*fGeom)[idx].GetAidx();
427 const UInt_t sector = (*fGeom)[idx].GetSector();
428
429 Byte_t *ptr = pixel.GetHiGainSamples() + fHiGainFirst;
430 Byte_t *end = ptr + fWindowSizeHiGain;
431
432 UInt_t sum = 0;
433 UInt_t sqr = 0;
434
435 do
436 {
437 sum += *ptr;
438 sqr += *ptr * *ptr;
439 }
440 while (++ptr != end);
441
442 if (fWindowSizeLoGain != 0)
443 {
444
445 ptr = pixel.GetLoGainSamples() + fLoGainFirst;
446 end = ptr + fWindowSizeLoGain;
447
448 do
449 {
450 sum += *ptr;
451 sqr += *ptr * *ptr;
452 }
453 while (++ptr != end);
454
455 }
456
457 const Float_t msum = (Float_t)sum;
458
459 //
460 // These three lines have been uncommented by Markus Gaug
461 // If anybody needs them, please contact me!!
462 //
463 // const Float_t higainped = msum/fNumHiGainSlices;
464 // const Float_t higainrms = TMath::Sqrt((msqr-msum*msum/fNumHiGainSlices)/(fNumHiGainSlices-1.));
465 // (*fPedestals)[idx].Set(higainped, higainrms);
466
467 fSumx[idx] += msum;
468 fAreaSumx[aidx] += msum;
469 fSectorSumx[sector] += msum;
470 //
471 // The old version:
472 //
473 // const Float_t msqr = (Float_t)sqr;
474 // fSumx2[idx] += msqr;
475 //
476 // The new version:
477 //
478 const Float_t sqrsum = msum*msum;
479 fSumx2[idx] += sqrsum;
480 fAreaSumx2[aidx] += sqrsum;
481 fSectorSumx2[sector] += sqrsum;
482 }
483
484 fPedestals->SetReadyToSave();
485 fNumSamplesTot += fWindowSizeHiGain + fWindowSizeLoGain;
486
487 return kTRUE;
488}
489
490// --------------------------------------------------------------------------
491//
492// Compute signal mean and rms in the whole run and store it in MPedestalCam
493//
494Int_t MPedCalcPedRun::PostProcess()
495{
496
497 // Compute pedestals and rms from the whole run
498 const ULong_t n = fNumSamplesTot;
499 const ULong_t nevts = GetNumExecutions();
500
501 MRawEvtPixelIter pixel(fRawEvt);
502
503 while (pixel.Next())
504 {
505
506 const Int_t pixid = pixel.GetPixelId();
507 const UInt_t aidx = (*fGeom)[pixid].GetAidx();
508 const UInt_t sector = (*fGeom)[pixid].GetSector();
509
510 fAreaValid [aidx]++;
511 fSectorValid[sector]++;
512
513 const Float_t sum = fSumx.At(pixid);
514 const Float_t sum2 = fSumx2.At(pixid);
515
516 const Float_t higainped = sum/n;
517 //
518 // The old version:
519 //
520 // const Float_t higainrms = TMath::Sqrt((sum2-sum*sum/n)/(n-1.));
521 //
522 // The new version:
523 //
524 // 1. Calculate the Variance of the sums:
525 Float_t higainVar = (sum2-sum*sum/nevts)/(nevts-1.);
526 // 2. Scale the variance to the number of slices:
527 higainVar /= (Float_t)(fWindowSizeHiGain+fWindowSizeLoGain);
528 // 3. Calculate the RMS from the Variance:
529 const Float_t higainrms = TMath::Sqrt(higainVar);
530
531 (*fPedestals)[pixid].Set(higainped, higainrms);
532
533 }
534
535 //
536 // Loop over the (two) area indices to get the averaged pedestal per aidx
537 //
538 for (Int_t aidx=0; aidx<fAreaValid.GetSize(); aidx++)
539 {
540
541 const Int_t napix = fAreaValid.At(aidx);
542 const Float_t sum = fAreaSumx.At(aidx);
543 const Float_t sum2 = fAreaSumx2.At(aidx);
544 const ULong_t an = napix * n;
545 const ULong_t aevts = napix * nevts;
546
547 const Float_t higainped = sum/an;
548
549 // 1. Calculate the Variance of the sums:
550 Float_t higainVar = (sum2-sum*sum/aevts)/(aevts-1.);
551 // 2. Scale the variance to the number of slices:
552 higainVar /= fWindowSizeHiGain+fWindowSizeLoGain;
553 // 3. Calculate the RMS from the Variance:
554 Float_t higainrms = TMath::Sqrt(higainVar);
555 // 4. Re-scale it with the square root of the number of involved pixels
556 // in order to be comparable to the mean of pedRMS of that area
557 higainrms *= TMath::Sqrt((Float_t)napix);
558
559 fPedestals->GetAverageArea(aidx).Set(higainped, higainrms);
560 }
561
562 //
563 // Loop over the (six) sector indices to get the averaged pedestal per sector
564 //
565 for (Int_t sector=0; sector<fSectorValid.GetSize(); sector++)
566 {
567
568 const Int_t nspix = fSectorValid.At(sector);
569 const Float_t sum = fSectorSumx.At(sector);
570 const Float_t sum2 = fSectorSumx2.At(sector);
571 const ULong_t sn = nspix * n;
572 const ULong_t sevts = nspix * nevts;
573
574 const Float_t higainped = sum/sn;
575
576 // 1. Calculate the Variance of the sums:
577 Float_t higainVar = (sum2-sum*sum/sevts)/(sevts-1.);
578 // 2. Scale the variance to the number of slices:
579 higainVar /= fWindowSizeHiGain+fWindowSizeLoGain;
580 // 3. Calculate the RMS from the Variance:
581 Float_t higainrms = TMath::Sqrt(higainVar);
582 // 4. Re-scale it with the square root of the number of involved pixels
583 // in order to be comparable to the mean of pedRMS of that sector
584 higainrms *= TMath::Sqrt((Float_t)nspix);
585
586 fPedestals->GetAverageSector(sector).Set(higainped, higainrms);
587 }
588
589 fPedestals->SetTotalEntries(fNumSamplesTot);
590 fPedestals->SetReadyToSave();
591
592 return kTRUE;
593}
Note: See TracBrowser for help on using the repository browser.