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

Last change on this file since 3930 was 3889, checked in by gaug, 21 years ago
*** empty log message ***
File size: 17.9 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 "MGeomCamMagic.h"
119
120ClassImp(MPedCalcPedRun);
121
122using namespace std;
123
124const Byte_t MPedCalcPedRun::fgHiGainFirst = 0;
125const Byte_t MPedCalcPedRun::fgHiGainLast = 14;
126const Byte_t MPedCalcPedRun::fgLoGainFirst = 0;
127const Byte_t MPedCalcPedRun::fgLoGainLast = 14;
128const Byte_t MPedCalcPedRun::fgHiGainWindowSize = 14;
129const Byte_t MPedCalcPedRun::fgLoGainWindowSize = 0;
130// --------------------------------------------------------------------------
131//
132// Default constructor:
133//
134// Sets:
135// - fWindowSizeHiGain to fgHiGainWindowSize
136// - fWindowSizeLoGain to fgLoGainWindowSize
137//
138// Calls:
139// - AddToBranchList("fHiGainPixId");
140// - AddToBranchList("fHiGainFadcSamples");
141// - SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast)
142// - Clear()
143//
144MPedCalcPedRun::MPedCalcPedRun(const char *name, const char *title)
145 : fWindowSizeHiGain(fgHiGainWindowSize),
146 fWindowSizeLoGain(fgLoGainWindowSize)
147{
148 fName = name ? name : "MPedCalcPedRun";
149 fTitle = title ? title : "Task to calculate pedestals from pedestal runs raw data";
150
151 AddToBranchList("fHiGainPixId");
152 AddToBranchList("fHiGainFadcSamples");
153
154 SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast);
155 Clear();
156}
157
158// --------------------------------------------------------------------------
159//
160// Sets:
161// - fNumSamplesTot to 0
162// - fRawEvt to NULL
163// - fRunHeader to NULL
164// - fPedestals to NULL
165//
166void MPedCalcPedRun::Clear(const Option_t *o)
167{
168
169 fNumSamplesTot = 0;
170
171 fRawEvt = NULL;
172 fRunHeader = NULL;
173 fPedestals = NULL;
174}
175
176// --------------------------------------------------------------------------
177//
178// SetRange:
179//
180// Calls:
181// - MExtractor::SetRange(hifirst,hilast,lofirst,lolast);
182// - SetWindowSize(fWindowSizeHiGain,fWindowSizeLoGain);
183//
184void MPedCalcPedRun::SetRange(Byte_t hifirst, Byte_t hilast, Byte_t lofirst, Byte_t lolast)
185{
186
187 MExtractor::SetRange(hifirst,hilast,lofirst,lolast);
188
189 //
190 // Redo the checks if the window is still inside the ranges
191 //
192 SetWindowSize(fWindowSizeHiGain,fWindowSizeLoGain);
193
194}
195
196
197// --------------------------------------------------------------------------
198//
199// Checks:
200// - if a window is odd, subtract one
201// - if a window is bigger than the one defined by the ranges, set it to the available range
202// - if a window is smaller than 2, set it to 2
203//
204// Sets:
205// - fNumHiGainSamples to: (Float_t)fWindowSizeHiGain
206// - fNumLoGainSamples to: (Float_t)fWindowSizeLoGain
207// - fSqrtHiGainSamples to: TMath::Sqrt(fNumHiGainSamples)
208// - fSqrtLoGainSamples to: TMath::Sqrt(fNumLoGainSamples)
209//
210void MPedCalcPedRun::SetWindowSize(Byte_t windowh, Byte_t windowl)
211{
212
213 fWindowSizeHiGain = windowh & ~1;
214 fWindowSizeLoGain = windowl & ~1;
215
216 if (fWindowSizeHiGain != windowh)
217 *fLog << warn << GetDescriptor() << ": Hi Gain window size has to be even, set to: "
218 << int(fWindowSizeHiGain) << " samples " << endl;
219
220 if (fWindowSizeLoGain != windowl)
221 *fLog << warn << GetDescriptor() << ": Lo Gain window size has to be even, set to: "
222 << int(fWindowSizeLoGain) << " samples " << endl;
223
224 const Byte_t availhirange = (fHiGainLast-fHiGainFirst+1) & ~1;
225 const Byte_t availlorange = (fLoGainLast-fLoGainFirst+1) & ~1;
226
227 if (fWindowSizeHiGain > availhirange)
228 {
229 *fLog << warn << GetDescriptor()
230 << Form("%s%2i%s%2i%s%2i%s",": Hi Gain window size: ",(int)fWindowSizeHiGain,
231 " is bigger than available range: [",(int)fHiGainFirst,",",(int)fHiGainLast,"]") << endl;
232 *fLog << warn << GetDescriptor()
233 << ": Will set window size to: " << (int)availhirange << endl;
234 fWindowSizeHiGain = availhirange;
235 }
236
237 if (fWindowSizeLoGain > availlorange)
238 {
239 *fLog << warn << GetDescriptor()
240 << Form("%s%2i%s%2i%s%2i%s",": Lo Gain window size: ",(int)fWindowSizeLoGain,
241 " is bigger than available range: [",(int)fLoGainFirst,",",(int)fLoGainLast,"]") << endl;
242 *fLog << warn << GetDescriptor()
243 << ": Will set window size to: " << (int)availlorange << endl;
244 fWindowSizeLoGain = availlorange;
245 }
246
247
248 fNumHiGainSamples = (Float_t)fWindowSizeHiGain;
249 fNumLoGainSamples = (Float_t)fWindowSizeLoGain;
250
251 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
252 fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
253
254}
255
256
257
258// --------------------------------------------------------------------------
259//
260// Look for the following input containers:
261//
262// - MRawEvtData
263// - MRawRunHeader
264// - MGeomCam
265//
266// The following output containers are also searched and created if
267// they were not found:
268//
269// - MPedestalCam
270//
271Int_t MPedCalcPedRun::PreProcess( MParList *pList )
272{
273
274 Clear();
275
276 fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
277 if (!fRawEvt)
278 {
279 *fLog << err << "MRawEvtData not found... aborting." << endl;
280 return kFALSE;
281 }
282
283 fRunHeader = (MRawRunHeader*)pList->FindObject(AddSerialNumber("MRawRunHeader"));
284 if (!fRunHeader)
285 {
286 *fLog << err << AddSerialNumber("MRawRunHeader") << " not found... aborting." << endl;
287 return kFALSE;
288 }
289
290 fGeom = (MGeomCam*)pList->FindObject("MGeomCam");
291 if (!fGeom)
292 {
293 *fLog << err << "MGeomCam not found... aborting." << endl;
294 return kFALSE;
295 }
296
297 fPedestals = (MPedestalCam*)pList->FindCreateObj("MPedestalCam");
298 if (!fPedestals)
299 return kFALSE;
300
301 return kTRUE;
302}
303
304// --------------------------------------------------------------------------
305//
306// The ReInit searches for:
307// - MRawRunHeader::GetNumSamplesHiGain()
308// - MRawRunHeader::GetNumSamplesLoGain()
309//
310// In case that the variables fHiGainLast and fLoGainLast are smaller than
311// the even part of the number of samples obtained from the run header, a
312// warning is given an the range is set back accordingly. A call to:
313// - SetRange(fHiGainFirst, fHiGainLast-diff, fLoGainFirst, fLoGainLast) or
314// - SetRange(fHiGainFirst, fHiGainLast, fLoGainFirst, fLoGainLast-diff)
315// is performed in that case. The variable diff means here the difference
316// between the requested range (fHiGainLast) and the available one. Note that
317// the functions SetRange() are mostly overloaded and perform more checks,
318// modifying the ranges again, if necessary.
319//
320Bool_t MPedCalcPedRun::ReInit(MParList *pList)
321{
322
323 Int_t lastdesired = (Int_t)fHiGainLast;
324 Int_t lastavailable = (Int_t)fRunHeader->GetNumSamplesHiGain()-1;
325
326 if (lastdesired > lastavailable)
327 {
328 const Int_t diff = lastdesired - lastavailable;
329 *fLog << endl;
330 *fLog << warn << GetDescriptor()
331 << Form("%s%2i%s%2i%s%2i%s",": Selected Hi Gain FADC Window [",
332 (int)fHiGainFirst,",",lastdesired,
333 "] ranges out of the available limits: [0,",lastavailable,"].") << endl;
334 *fLog << GetDescriptor() << ": Will reduce the upper edge to " << (int)(fHiGainLast - diff) << endl;
335 SetRange(fHiGainFirst, fHiGainLast-diff, fLoGainFirst, fLoGainLast);
336 }
337
338 lastdesired = (Int_t)(fLoGainLast);
339 lastavailable = (Int_t)fRunHeader->GetNumSamplesLoGain()-1;
340
341 if (lastdesired > lastavailable)
342 {
343 const Int_t diff = lastdesired - lastavailable;
344 *fLog << endl;
345 *fLog << warn << GetDescriptor()
346 << Form("%s%2i%s%2i%s%2i%s",": Selected Lo Gain FADC Window [",
347 (int)fLoGainFirst,",",lastdesired,
348 "] ranges out of the available limits: [0,",lastavailable,"].") << endl;
349 *fLog << GetDescriptor() << ": Will reduce the upper edge to " << (int)(fLoGainLast - diff) << endl;
350 SetRange(fHiGainFirst, fHiGainLast, fLoGainFirst, fLoGainLast-diff);
351 }
352
353 Int_t npixels = fPedestals->GetSize();
354 Int_t areas = fPedestals->GetAverageAreas();
355 Int_t sectors = fPedestals->GetAverageSectors();
356
357 if (fSumx.GetSize()==0)
358 {
359 fSumx. Set(npixels);
360 fSumx2.Set(npixels);
361
362 fAreaSumx. Set(areas);
363 fAreaSumx2.Set(areas);
364 fAreaValid.Set(areas);
365
366 fSectorSumx. Set(sectors);
367 fSectorSumx2.Set(sectors);
368 fSectorValid.Set(sectors);
369
370 fSumx.Reset();
371 fSumx2.Reset();
372 }
373
374 if (fWindowSizeHiGain == 0 && fWindowSizeLoGain == 0)
375 {
376 *fLog << err << GetDescriptor()
377 << ": Number of extracted Slices is 0, abort ... " << endl;
378 return kFALSE;
379 }
380
381
382 *fLog << endl;
383 *fLog << inf << GetDescriptor() << ": Taking " << (int)fWindowSizeHiGain
384 << " HiGain FADC samples starting with slice: " << (int)fHiGainFirst << endl;
385 *fLog << inf << GetDescriptor() << ": Taking " << (int)fWindowSizeLoGain
386 << " LoGain FADC samples starting with slice: " << (int)fLoGainFirst << endl;
387
388 return kTRUE;
389
390}
391// --------------------------------------------------------------------------
392//
393// Fill the MPedestalCam container with the signal mean and rms for the event.
394// Store the measured signal in arrays fSumx and fSumx2 so that we can
395// calculate the overall mean and rms in the PostProcess()
396//
397Int_t MPedCalcPedRun::Process()
398{
399
400 MRawEvtPixelIter pixel(fRawEvt);
401
402 while (pixel.Next())
403 {
404
405 const UInt_t idx = pixel.GetPixelId();
406 const UInt_t aidx = (*fGeom)[idx].GetAidx();
407 const UInt_t sector = (*fGeom)[idx].GetSector();
408
409 Byte_t *ptr = pixel.GetHiGainSamples() + fHiGainFirst;
410 Byte_t *end = ptr + fWindowSizeHiGain;
411
412 UInt_t sum = 0;
413 UInt_t sqr = 0;
414
415 do
416 {
417 sum += *ptr;
418 sqr += *ptr * *ptr;
419 }
420 while (++ptr != end);
421
422 if (fWindowSizeLoGain != 0)
423 {
424
425 ptr = pixel.GetLoGainSamples() + fLoGainFirst;
426 end = ptr + fWindowSizeLoGain;
427
428 do
429 {
430 sum += *ptr;
431 sqr += *ptr * *ptr;
432 }
433 while (++ptr != end);
434
435 }
436
437 const Float_t msum = (Float_t)sum;
438
439 //
440 // These three lines have been uncommented by Markus Gaug
441 // If anybody needs them, please contact me!!
442 //
443 // const Float_t higainped = msum/fNumHiGainSlices;
444 // const Float_t higainrms = TMath::Sqrt((msqr-msum*msum/fNumHiGainSlices)/(fNumHiGainSlices-1.));
445 // (*fPedestals)[idx].Set(higainped, higainrms);
446
447 fSumx[idx] += msum;
448 fAreaSumx[aidx] += msum;
449 fSectorSumx[sector] += msum;
450 //
451 // The old version:
452 //
453 // const Float_t msqr = (Float_t)sqr;
454 // fSumx2[idx] += msqr;
455 //
456 // The new version:
457 //
458 const Float_t sqrsum = msum*msum;
459 fSumx2[idx] += sqrsum;
460 fAreaSumx2[aidx] += sqrsum;
461 fSectorSumx2[sector] += sqrsum;
462 }
463
464 fPedestals->SetReadyToSave();
465 fNumSamplesTot += fWindowSizeHiGain + fWindowSizeLoGain;
466
467 return kTRUE;
468}
469
470// --------------------------------------------------------------------------
471//
472// Compute signal mean and rms in the whole run and store it in MPedestalCam
473//
474Int_t MPedCalcPedRun::PostProcess()
475{
476
477 // Compute pedestals and rms from the whole run
478 const ULong_t n = fNumSamplesTot;
479 const ULong_t nevts = GetNumExecutions();
480
481 MRawEvtPixelIter pixel(fRawEvt);
482
483 while (pixel.Next())
484 {
485
486 const Int_t pixid = pixel.GetPixelId();
487 const UInt_t aidx = (*fGeom)[pixid].GetAidx();
488 const UInt_t sector = (*fGeom)[pixid].GetSector();
489
490 fAreaValid [aidx]++;
491 fSectorValid[sector]++;
492
493 const Float_t sum = fSumx.At(pixid);
494 const Float_t sum2 = fSumx2.At(pixid);
495
496 const Float_t higainped = sum/n;
497 //
498 // The old version:
499 //
500 // const Float_t higainrms = TMath::Sqrt((sum2-sum*sum/n)/(n-1.));
501 //
502 // The new version:
503 //
504 // 1. Calculate the Variance of the sums:
505 Float_t higainVar = (sum2-sum*sum/nevts)/(nevts-1.);
506 // 2. Scale the variance to the number of slices:
507 higainVar /= (Float_t)(fWindowSizeHiGain+fWindowSizeLoGain);
508 // 3. Calculate the RMS from the Variance:
509 const Float_t higainrms = TMath::Sqrt(higainVar);
510
511 (*fPedestals)[pixid].Set(higainped, higainrms);
512
513 }
514
515 //
516 // Loop over the (two) area indices to get the averaged pedestal per aidx
517 //
518 for (Int_t aidx=0; aidx<fAreaValid.GetSize(); aidx++)
519 {
520
521 const Int_t napix = fAreaValid.At(aidx);
522 const Float_t sum = fAreaSumx.At(aidx);
523 const Float_t sum2 = fAreaSumx2.At(aidx);
524 const ULong_t an = napix * n;
525 const ULong_t aevts = napix * nevts;
526
527 const Float_t higainped = sum/an;
528
529 // 1. Calculate the Variance of the sums:
530 Float_t higainVar = (sum2-sum*sum/aevts)/(aevts-1.);
531 // 2. Scale the variance to the number of slices:
532 higainVar /= fWindowSizeHiGain+fWindowSizeLoGain;
533 // 3. Calculate the RMS from the Variance:
534 Float_t higainrms = TMath::Sqrt(higainVar);
535 // 4. Re-scale it with the square root of the number of involved pixels
536 // in order to be comparable to the mean of pedRMS of that area
537 higainrms *= TMath::Sqrt((Float_t)napix);
538
539 fPedestals->GetAverageArea(aidx).Set(higainped, higainrms);
540 }
541
542 //
543 // Loop over the (six) sector indices to get the averaged pedestal per sector
544 //
545 for (Int_t sector=0; sector<fSectorValid.GetSize(); sector++)
546 {
547
548 const Int_t nspix = fSectorValid.At(sector);
549 const Float_t sum = fSectorSumx.At(sector);
550 const Float_t sum2 = fSectorSumx2.At(sector);
551 const ULong_t sn = nspix * n;
552 const ULong_t sevts = nspix * nevts;
553
554 const Float_t higainped = sum/sn;
555
556 // 1. Calculate the Variance of the sums:
557 Float_t higainVar = (sum2-sum*sum/sevts)/(sevts-1.);
558 // 2. Scale the variance to the number of slices:
559 higainVar /= fWindowSizeHiGain+fWindowSizeLoGain;
560 // 3. Calculate the RMS from the Variance:
561 Float_t higainrms = TMath::Sqrt(higainVar);
562 // 4. Re-scale it with the square root of the number of involved pixels
563 // in order to be comparable to the mean of pedRMS of that sector
564 higainrms *= TMath::Sqrt((Float_t)nspix);
565
566 fPedestals->GetAverageSector(sector).Set(higainped, higainrms);
567 }
568
569 fPedestals->SetTotalEntries(fNumSamplesTot);
570 fPedestals->SetReadyToSave();
571
572 return kTRUE;
573}
Note: See TracBrowser for help on using the repository browser.