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@uni-sw.gwdg.de>
|
---|
20 | !
|
---|
21 | ! Copyright: MAGIC Software Development, 2000-2001
|
---|
22 | !
|
---|
23 | !
|
---|
24 | \* ======================================================================== */
|
---|
25 |
|
---|
26 | /////////////////////////////////////////////////////////////////////////////
|
---|
27 | // //
|
---|
28 | // MPedCalcPedRun //
|
---|
29 | // //
|
---|
30 | // Input Containers: //
|
---|
31 | // MRawEvtData //
|
---|
32 | // //
|
---|
33 | // Output Containers: //
|
---|
34 | // MPedestalCam //
|
---|
35 | // //
|
---|
36 | /////////////////////////////////////////////////////////////////////////////
|
---|
37 |
|
---|
38 | #include "MPedCalcPedRun.h"
|
---|
39 |
|
---|
40 | #include "MParList.h"
|
---|
41 |
|
---|
42 | #include "MLog.h"
|
---|
43 | #include "MLogManip.h"
|
---|
44 |
|
---|
45 | #include "MRawRunHeader.h"
|
---|
46 | #include "MRawEvtPixelIter.h"
|
---|
47 | #include "MRawEvtData.h"
|
---|
48 |
|
---|
49 | #include "MPedestalPix.h"
|
---|
50 | #include "MPedestalCam.h"
|
---|
51 |
|
---|
52 | #include "MGeomCamMagic.h"
|
---|
53 |
|
---|
54 | ClassImp(MPedCalcPedRun);
|
---|
55 |
|
---|
56 | using namespace std;
|
---|
57 |
|
---|
58 | MPedCalcPedRun::MPedCalcPedRun(const char *name, const char *title)
|
---|
59 | {
|
---|
60 | fName = name ? name : "MPedCalcPedRun";
|
---|
61 | fTitle = title ? title : "Task to calculate pedestals from pedestal runs raw data";
|
---|
62 |
|
---|
63 | AddToBranchList("fHiGainPixId");
|
---|
64 | AddToBranchList("fHiGainFadcSamples");
|
---|
65 | }
|
---|
66 |
|
---|
67 | Int_t MPedCalcPedRun::PreProcess( MParList *pList )
|
---|
68 | {
|
---|
69 | fRawEvt = (MRawEvtData*)pList->FindObject("MRawEvtData");
|
---|
70 | if (!fRawEvt)
|
---|
71 | {
|
---|
72 | *fLog << dbginf << "MRawEvtData not found... aborting." << endl;
|
---|
73 | return kFALSE;
|
---|
74 | }
|
---|
75 |
|
---|
76 | fPedestals = (MPedestalCam*)pList->FindCreateObj("MPedestalCam");
|
---|
77 | if (!fPedestals)
|
---|
78 | return kFALSE;
|
---|
79 |
|
---|
80 | MGeomCamMagic magiccam;
|
---|
81 |
|
---|
82 | fSumx.Set(magiccam.GetNumPixels());
|
---|
83 | fSumx2.Set(magiccam.GetNumPixels());
|
---|
84 |
|
---|
85 | for(UInt_t i=0;i<magiccam.GetNumPixels();i++)
|
---|
86 | {
|
---|
87 | fSumx.AddAt(0,i);
|
---|
88 | fSumx2.AddAt(0,i);
|
---|
89 | }
|
---|
90 |
|
---|
91 | return kTRUE;
|
---|
92 | }
|
---|
93 |
|
---|
94 | Bool_t MPedCalcPedRun::ReInit(MParList *pList )
|
---|
95 | {
|
---|
96 |
|
---|
97 | fRunheader = (MRawRunHeader*)pList->FindObject("MRawRunHeader");
|
---|
98 | if (!fRunheader)
|
---|
99 | {
|
---|
100 | *fLog << warn << dbginf <<
|
---|
101 | "Warning - cannot check file type, MRawRunHeader not found." << endl;
|
---|
102 | }
|
---|
103 | else
|
---|
104 | if (fRunheader->GetRunType() == kRTMonteCarlo)
|
---|
105 | {
|
---|
106 | return kTRUE;
|
---|
107 | }
|
---|
108 |
|
---|
109 | fNumHiGainSamples = fRunheader->GetNumSamplesHiGain();
|
---|
110 |
|
---|
111 | fPedestals->InitSize(fRunheader->GetNumPixel());
|
---|
112 |
|
---|
113 | return kTRUE;
|
---|
114 |
|
---|
115 | }
|
---|
116 |
|
---|
117 |
|
---|
118 | Int_t MPedCalcPedRun::Process()
|
---|
119 | {
|
---|
120 |
|
---|
121 | MRawEvtPixelIter pixel(fRawEvt);
|
---|
122 |
|
---|
123 | while (pixel.Next())
|
---|
124 | {
|
---|
125 | Byte_t shift=(fNumHiGainSamples/2*2==fNumHiGainSamples) ? 0:1;
|
---|
126 | Byte_t *ptr = pixel.GetHiGainSamples();
|
---|
127 | const Byte_t *end = ptr + fRawEvt->GetNumHiGainSamples()-shift;
|
---|
128 |
|
---|
129 | const Float_t higainped = CalcHiGainMean(ptr, end);
|
---|
130 | const Float_t higainrms = CalcHiGainRms(ptr, end, higainped);
|
---|
131 |
|
---|
132 | const UInt_t pixid = pixel.GetPixelId();
|
---|
133 | MPedestalPix &pix = (*fPedestals)[pixid];
|
---|
134 |
|
---|
135 | // cumulate the sum of pedestals and of pedestal squares
|
---|
136 | fSumx.AddAt(higainped+fSumx.At(pixid),pixid);
|
---|
137 | fSumx2.AddAt(GetSumx2(ptr, end)+fSumx2.At(pixid),pixid);
|
---|
138 |
|
---|
139 | // set the value of the pedestal and rms computed from the processed event
|
---|
140 | pix.Set(higainped, higainrms);
|
---|
141 | }
|
---|
142 |
|
---|
143 | fPedestals->SetReadyToSave();
|
---|
144 |
|
---|
145 | return kTRUE;
|
---|
146 | }
|
---|
147 |
|
---|
148 | Int_t MPedCalcPedRun::PostProcess()
|
---|
149 | {
|
---|
150 | // Compute pedestals and rms from the whole run
|
---|
151 |
|
---|
152 | MRawEvtPixelIter pixel(fRawEvt);
|
---|
153 |
|
---|
154 | while (pixel.Next())
|
---|
155 | {
|
---|
156 | const UInt_t pixid = pixel.GetPixelId();
|
---|
157 | MPedestalPix &pix = (*fPedestals)[pixid];
|
---|
158 |
|
---|
159 | const Int_t N = GetNumExecutions();
|
---|
160 | const Float_t sum = fSumx.At(pixid);
|
---|
161 | const Float_t sum2 = fSumx2.At(pixid);
|
---|
162 | const Float_t higainped = sum/N;
|
---|
163 | const Float_t higainrms = sqrt(1./(N-1.)*(sum2-sum*sum/N));
|
---|
164 | pix.Set(higainped,higainrms);
|
---|
165 |
|
---|
166 | }
|
---|
167 | return kTRUE;
|
---|
168 |
|
---|
169 | }
|
---|
170 |
|
---|
171 |
|
---|
172 | Float_t MPedCalcPedRun::CalcHiGainMean(Byte_t *ptr, const Byte_t *end) const
|
---|
173 | {
|
---|
174 | Int_t sum=0;
|
---|
175 | Byte_t EvenNumSamples=(fNumHiGainSamples/2*2==fNumHiGainSamples) ? fNumHiGainSamples:fNumHiGainSamples-1;
|
---|
176 |
|
---|
177 | do sum += *ptr;
|
---|
178 | while (++ptr != end);
|
---|
179 |
|
---|
180 | return (Float_t)sum/EvenNumSamples;
|
---|
181 | }
|
---|
182 |
|
---|
183 |
|
---|
184 |
|
---|
185 | Float_t MPedCalcPedRun::GetSumx2(Byte_t *ptr, const Byte_t *end) const
|
---|
186 | {
|
---|
187 | Float_t square = 0;
|
---|
188 |
|
---|
189 | // Take an even number of time slices to avoid biases due to A/B effect
|
---|
190 | Byte_t EvenNumSamples=(fNumHiGainSamples/2*2==fNumHiGainSamples) ? fNumHiGainSamples:fNumHiGainSamples-1;
|
---|
191 |
|
---|
192 | do
|
---|
193 | {
|
---|
194 | const Float_t val = (Float_t)(*ptr);
|
---|
195 |
|
---|
196 | square += val*val;
|
---|
197 | } while (++ptr != end);
|
---|
198 |
|
---|
199 | return square/EvenNumSamples;
|
---|
200 | }
|
---|
201 |
|
---|
202 |
|
---|
203 | Float_t MPedCalcPedRun::CalcHiGainRms(Byte_t *ptr, const Byte_t *end, Float_t higainped) const
|
---|
204 | {
|
---|
205 | Float_t rms = 0;
|
---|
206 | Byte_t EvenNumSamples=(fNumHiGainSamples/2*2==fNumHiGainSamples) ? fNumHiGainSamples:fNumHiGainSamples-1;
|
---|
207 |
|
---|
208 | do
|
---|
209 | {
|
---|
210 | const Float_t diff = (Float_t)(*ptr)-higainped;
|
---|
211 |
|
---|
212 | rms += diff*diff;
|
---|
213 | } while (++ptr != end);
|
---|
214 |
|
---|
215 | return sqrt(rms/(EvenNumSamples-1));
|
---|
216 | }
|
---|
217 |
|
---|
218 |
|
---|
219 |
|
---|
220 | /*
|
---|
221 | Float_t MPedCalcPedRun::CalcHiGainMeanErr(Float_t higainrms) const
|
---|
222 | {
|
---|
223 | return higainrms/sqrt((float)fNumHiGainSamples);
|
---|
224 | }
|
---|
225 |
|
---|
226 | Float_t MPedCalcPedRun::CalcHiGainRmsErr(Float_t higainrms) const
|
---|
227 | {
|
---|
228 | return higainrms/sqrt(2.*fNumHiGainSamples);
|
---|
229 | }
|
---|
230 | */
|
---|