source: trunk/MagicSoft/Mars/manalysis/MPedCalcPedRun.cc@ 3005

Last change on this file since 3005 was 3005, checked in by gaug, 22 years ago
*** empty log message ***
File size: 8.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
22!
23! Copyright: MAGIC Software Development, 2000-2004
24!
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//
72// Input Containers:
73// MRawEvtData
74//
75// Output Containers:
76// MPedestalCam
77//
78//
79/////////////////////////////////////////////////////////////////////////////
80#include "MPedCalcPedRun.h"
81
82#include "MParList.h"
83
84#include "MLog.h"
85#include "MLogManip.h"
86
87#include "MRawRunHeader.h"
88#include "MRawEvtPixelIter.h"
89#include "MRawEvtData.h"
90
91#include "MPedestalPix.h"
92#include "MPedestalCam.h"
93
94#include "MExtractedSignalPix.h"
95#include "MExtractedSignalCam.h"
96
97
98#include "MGeomCamMagic.h"
99
100ClassImp(MPedCalcPedRun);
101
102using namespace std;
103
104// --------------------------------------------------------------------------
105//
106// default constructor
107//
108MPedCalcPedRun::MPedCalcPedRun(const char *name, const char *title)
109 : fRawEvt(NULL), fPedestals(NULL), fSignals(NULL)
110{
111 fName = name ? name : "MPedCalcPedRun";
112 fTitle = title ? title : "Task to calculate pedestals from pedestal runs raw data";
113
114 AddToBranchList("fHiGainPixId");
115 AddToBranchList("fHiGainFadcSamples");
116
117 Clear();
118}
119
120void MPedCalcPedRun::Clear(const Option_t *o)
121{
122
123 fNumHiGainSamples = 0;
124 fNumPixels = 0;
125 fNumSamplesTot = 0;
126 fUseHists = kFALSE;
127
128 fRawEvt = NULL;
129 fPedestals = NULL;
130 fSignals = NULL;
131
132 return;
133
134}
135
136
137// --------------------------------------------------------------------------
138//
139// Look for the following input containers:
140//
141// - MRawEvtData
142//
143// The following output containers are also searched and created if
144// they were not found:
145//
146// - MPedestalCam
147//
148Int_t MPedCalcPedRun::PreProcess( MParList *pList )
149{
150 fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
151 if (!fRawEvt)
152 {
153 *fLog << err << "MRawEvtData not found... aborting." << endl;
154 return kFALSE;
155 }
156
157 fPedestals = (MPedestalCam*)pList->FindCreateObj("MPedestalCam");
158 if (!fPedestals)
159 return kFALSE;
160
161 fSignals = (MExtractedSignalCam*)pList->FindObject("MExtractedSignalCam");
162
163 if (!fSignals && fUseHists)
164 {
165 *fLog << warn << "Cannot find MExtractedSignalCam... will not use histograms!" << endl;
166 fUseHists = kFALSE;
167 }
168
169 fNumSamplesTot=0;
170
171 return kTRUE;
172}
173
174// --------------------------------------------------------------------------
175//
176// The ReInit searches for the following input containers:
177// - MRawRunHeader
178//
179// It also initializes the data arrays fSumx and fSumx2
180// (only for the first read file)
181//
182Bool_t MPedCalcPedRun::ReInit(MParList *pList)
183{
184 const MRawRunHeader *runheader = (MRawRunHeader*)pList->FindObject("MRawRunHeader");
185 if (!runheader)
186 {
187 *fLog << warn << dbginf;
188 *fLog << "Warning - cannot check file type, MRawRunHeader not found." << endl;
189 }
190 else
191 if (runheader->GetRunType() == kRTMonteCarlo)
192 return kTRUE;
193
194 fNumPixels = fPedestals->GetSize();
195
196 if(fSumx.GetSize()==0)
197 {
198 fSumx.Set(fNumPixels);
199 fSumx2.Set(fNumPixels);
200
201 fSumx.Reset();
202 fSumx2.Reset();
203 }
204
205 // Calculate an even number for the hi gain samples to avoid
206 // biases due to the fluctuation in pedestal from one slice to
207 // the other one
208 fNumHiGainSamples = runheader->GetNumSamplesHiGain() & ~1;
209
210 if (fUseHists)
211 fPedestals->InitUseHists();
212
213 return kTRUE;
214}
215// --------------------------------------------------------------------------
216//
217// Fill the MPedestalCam container with the signal mean and rms for the event.
218// Store the measured signal in arrays fSumx and fSumx2 so that we can
219// calculate the overall mean and rms in the PostProcess()
220//
221Int_t MPedCalcPedRun::Process()
222{
223 MRawEvtPixelIter pixel(fRawEvt);
224
225 while (pixel.Next())
226 {
227 Byte_t *ptr = pixel.GetHiGainSamples();
228 const Byte_t *end = ptr + fNumHiGainSamples;
229
230 UInt_t sum = 0;
231 UInt_t sqr = 0;
232
233 do
234 {
235 sum += *ptr;
236 sqr += *ptr * *ptr;
237 }
238 while (++ptr != end);
239
240 const Float_t msum = (Float_t)sum;
241
242 const UInt_t idx = pixel.GetPixelId();
243 //
244 // These three lines have been uncommented by Markus Gaug
245 // If anybody needs them, please contact me!!
246 //
247 // const Float_t higainped = msum/fNumHiGainSamples;
248 // const Float_t higainrms = TMath::Sqrt((msqr-msum*msum/fNumHiGainSamples)/(fNumHiGainSamples-1.));
249 // (*fPedestals)[idx].Set(higainped, higainrms);
250
251 fSumx[idx] += msum;
252 //
253 // The old version:
254 //
255 // const Float_t msqr = (Float_t)sqr;
256 // fSumx2[idx] += msqr;
257 //
258 // The new version:
259 //
260 fSumx2[idx] += msum*msum;
261
262 if (fUseHists)
263 {
264 MExtractedSignalPix &sig = (*fSignals)[idx];
265 const Float_t signal = sig.GetExtractedSignalHiGain();
266 (*fPedestals)[idx].FillHists(signal);
267 }
268 }
269
270 fPedestals->SetReadyToSave();
271 fNumSamplesTot += fNumHiGainSamples;
272
273
274 return kTRUE;
275}
276
277// --------------------------------------------------------------------------
278//
279// Compute signal mean and rms in the whole run and store it in MPedestalCam
280//
281Int_t MPedCalcPedRun::PostProcess()
282{
283
284 // Compute pedestals and rms from the whole run
285 const ULong_t n = fNumSamplesTot;
286 const ULong_t nevts = GetNumExecutions();
287
288 MRawEvtPixelIter pixel(fRawEvt);
289
290 while (pixel.Next())
291 {
292 const UInt_t pixid = pixel.GetPixelId();
293
294 const Float_t sum = fSumx.At(pixid);
295 const Float_t sum2 = fSumx2.At(pixid);
296
297 const Float_t higainped = sum/n;
298 //
299 // The old version:
300 //
301 // const Float_t higainrms = TMath::Sqrt((sum2-sum*sum/n)/(n-1.));
302 //
303 // The new version:
304 //
305 const Float_t higainrms = TMath::Sqrt((sum2-sum*sum/nevts)/(nevts-1.)/(Float_t)fNumHiGainSamples);
306
307 (*fPedestals)[pixid].Set(higainped, higainrms);
308
309 if (fUseHists)
310 (*fPedestals)[pixid].FitCharge();
311
312 }
313
314 fPedestals->SetNumTotSlices(fNumSamplesTot);
315
316 if (fUseHists)
317 fPedestals->SetNumExtractSlices(fSignals->GetNumUsedHiGainFADCSlices());
318
319 return kTRUE;
320}
321
322
Note: See TracBrowser for help on using the repository browser.