source: trunk/MagicSoft/Mars/mpedestal/MPedCalcFromData.cc@ 8154

Last change on this file since 8154 was 8106, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 5.3 KB
Line 
1/* ======================================================================== *\
2! $Name: not supported by cvs2svn $:$Id: MPedCalcFromData.cc,v 1.3 2006-10-17 17:16:01 tbretz Exp $
3! --------------------------------------------------------------------------
4!
5! *
6! * This file is part of MARS, the MAGIC Analysis and Reconstruction
7! * Software. It is distributed to you in the hope that it can be a useful
8! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
9! * It is distributed WITHOUT ANY WARRANTY.
10! *
11! * Permission to use, copy, modify and distribute this software and its
12! * documentation for any purpose is hereby granted without fee,
13! * provided that the above copyright notice appear in all copies and
14! * that both that copyright notice and this permission notice appear
15! * in supporting documentation. It is provided "as is" without expressed
16! * or implied warranty.
17! *
18!
19!
20! Author(s): Josep Flix 06/2004 <mailto:jflix@ifae.es>
21!
22! Copyright: MAGIC Software Development, 2000-2004
23!
24!
25\* ======================================================================== */
26
27/////////////////////////////////////////////////////////////////////////////
28// //
29// MPedCalcFromData //
30// //
31// Input Containers: //
32// MRawEvtData //
33// //
34// Output Containers: //
35// MPedestalCam //
36// //
37/////////////////////////////////////////////////////////////////////////////
38
39#include "MPedCalcFromData.h"
40#include "MExtractor.h"
41
42#include "MParList.h"
43
44#include "MLog.h"
45#include "MLogManip.h"
46
47#include "MRawRunHeader.h"
48#include "MRawEvtPixelIter.h"
49#include "MRawEvtData.h"
50
51#include "MPedestalPix.h"
52#include "MPedestalCam.h"
53
54ClassImp(MPedCalcFromData);
55
56using namespace std;
57
58MPedCalcFromData::MPedCalcFromData(const char *name, const char *title)
59{
60 fName = name ? name : "MPedCalcFromData";
61 fTitle = title ? title : "Task to calculate pedestals from data runs";
62
63 AddToBranchList("fHiGainPixId");
64 AddToBranchList("fHiGainFadcSamples");
65
66 SetRange(fLoGainFirst, fLoGainLast);
67 Clear();
68}
69
70void MPedCalcFromData::Clear(const Option_t *o)
71{
72
73 fRawEvt = NULL;
74 fRunHeader = NULL;
75
76}
77
78void MPedCalcFromData::SetLoRange(Byte_t lofirst, Byte_t lolast)
79{
80 fLoGainFirst = lofirst;
81 fLoGainLast = lolast;
82
83 fWindowSizeLoGain = lolast - lofirst;
84
85}
86
87Int_t MPedCalcFromData::PreProcess( MParList *pList )
88{
89 Clear();
90
91 fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
92 if (!fRawEvt)
93 {
94 *fLog << err << "MRawEvtData not found... aborting." << endl;
95 return kFALSE;
96 }
97
98 fRunHeader = (MRawRunHeader*)pList->FindObject(AddSerialNumber("MRawRunHeader"));
99 if (!fRunHeader)
100 {
101 *fLog << err << AddSerialNumber("MRawRunHeader") << " not found... aborting." << endl;
102 return kFALSE;
103 }
104
105 fPedestals = (MPedestalCam*)pList->FindCreateObj("MPedestalCam");
106 if (!fPedestals)
107 return kFALSE;
108
109 return kTRUE;
110}
111
112
113Bool_t MPedCalcFromData::ReInit(MParList *pList)
114{
115
116 Int_t npixels = fPedestals->GetSize();
117
118 if (fSumx.GetSize()==0)
119 {
120 fSumx. Set(npixels);
121 fSumx2.Set(npixels);
122 fEvtCounter.Set(npixels);
123 fTotalCounter.Set(npixels);
124
125 fEvtCounter.Reset();
126 fTotalCounter.Reset();
127 fSumx.Reset();
128 fSumx2.Reset();
129 }
130
131 return kTRUE;
132
133}
134
135Int_t MPedCalcFromData::Process()
136{
137
138 MRawEvtPixelIter pixel(fRawEvt);
139
140 while (pixel.Next())
141 {
142
143 const UInt_t idx = pixel.GetPixelId();
144
145 if ( (UInt_t)pixel.GetMaxHiGainSample() < fHiGainThreshold ) {
146
147 fEvtCounter[idx]++;
148
149 Byte_t *ptr = pixel.GetLoGainSamples() + fLoGainFirst;
150 Byte_t *end = ptr + fWindowSizeLoGain;
151
152 UInt_t sum = 0;
153 UInt_t sqr = 0;
154
155 if (fWindowSizeLoGain != 0)
156 {
157 do
158 {
159 sum += *ptr;
160 sqr += *ptr * *ptr;
161 }
162 while (++ptr != end);
163 }
164
165 const Float_t msum = (Float_t)sum;
166 fSumx[idx] += msum;
167
168 const Float_t sqrsum = msum*msum;
169 fSumx2[idx] += sqrsum;
170
171 if ((UInt_t)fEvtCounter[idx] == fDump){
172
173 // Compute pedestals and rms from the sample
174 const ULong_t n = fWindowSizeLoGain*fDump;
175
176 const Float_t sum1 = fSumx.At(idx);
177 const Float_t sum2 = fSumx2.At(idx);
178 const Float_t higainped = sum1/n;
179
180 // 1. Calculate the Variance of the sums:
181 Float_t higainVar = (sum2-sum1*sum1/fDump)/(fDump-1.);
182 // 2. Scale the variance to the number of slices:
183 higainVar /= (Float_t)(fWindowSizeLoGain);
184 // 3. Calculate the RMS from the Variance:
185 (*fPedestals)[idx].Set(higainped, higainVar < 0 ? 0. : TMath::Sqrt(higainVar));
186
187 fTotalCounter[idx]++;
188 fEvtCounter[idx]=0;
189 fSumx[idx]=0;
190 fSumx2[idx]=0;
191 };
192 }
193
194 };
195
196 fPedestals->SetReadyToSave();
197 return kTRUE;
198}
199
Note: See TracBrowser for help on using the repository browser.