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

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