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): Markus Gaug, 02/2004 <mailto:markus@ifae.es>
|
---|
19 | !
|
---|
20 | ! Copyright: MAGIC Software Development, 2000-2004
|
---|
21 | !
|
---|
22 | !
|
---|
23 | \* ======================================================================== */
|
---|
24 |
|
---|
25 | //////////////////////////////////////////////////////////////////////////////
|
---|
26 | //
|
---|
27 | // MExtractPINDiode
|
---|
28 | //
|
---|
29 | // Extracts the signal from a fixed window in a given range.
|
---|
30 | //
|
---|
31 | // Call: SetRange(fHiGainFirst, fHiGainLast, fLoGainFirst, fLoGainLast)
|
---|
32 | // to modify the ranges.
|
---|
33 | //
|
---|
34 | // Defaults are:
|
---|
35 | //
|
---|
36 | // fHiGainFirst = fgHiGainFirst = 0
|
---|
37 | // fHiGainLast = fgHiGainLast = 29
|
---|
38 | // fLoGainFirst = fgLoGainFirst = 0
|
---|
39 | // fLoGainLast = fgLoGainLast = 0
|
---|
40 | //
|
---|
41 | // The FADC slices are fit by a Gaussian around the pulse maximum.
|
---|
42 | // The following figures show two typical (high-intensity and low-intensity)
|
---|
43 | // pulses together with the applied fit:
|
---|
44 | //
|
---|
45 | //Begin_Html
|
---|
46 | /*
|
---|
47 | <img src="images/PINDiode_pulse_high.gif">
|
---|
48 | <img src="images/PINDiode_pulse_low.gif">
|
---|
49 | */
|
---|
50 | //End_Html
|
---|
51 | //
|
---|
52 | // The fit ranges can be modified with the functions:
|
---|
53 | // - SetLowerFitLimit()
|
---|
54 | // - SetUpperFitLimit()
|
---|
55 | //
|
---|
56 | // Defaults are:
|
---|
57 | // - fLowerFitLimit: 2
|
---|
58 | // - fUpperFitLimit: 5
|
---|
59 | //
|
---|
60 | //////////////////////////////////////////////////////////////////////////////
|
---|
61 | #include "MExtractPINDiode.h"
|
---|
62 |
|
---|
63 | #include <fstream>
|
---|
64 |
|
---|
65 | #include <TF1.h>
|
---|
66 | #include <TH1.h>
|
---|
67 | #include <TPad.h>
|
---|
68 |
|
---|
69 | #include "MLog.h"
|
---|
70 | #include "MLogManip.h"
|
---|
71 |
|
---|
72 | #include "MParList.h"
|
---|
73 |
|
---|
74 | #include "MRawEvtData.h"
|
---|
75 | #include "MRawRunHeader.h"
|
---|
76 | #include "MRawEvtPixelIter.h"
|
---|
77 |
|
---|
78 | #include "MPedestalCam.h"
|
---|
79 | #include "MPedestalPix.h"
|
---|
80 |
|
---|
81 | #include "MExtractedSignalPINDiode.h"
|
---|
82 |
|
---|
83 | ClassImp(MExtractPINDiode);
|
---|
84 |
|
---|
85 | using namespace std;
|
---|
86 |
|
---|
87 | const UInt_t MExtractPINDiode::fgPINDiodeIdx = 3;
|
---|
88 | const Byte_t MExtractPINDiode::fgHiGainFirst = 0;
|
---|
89 | const Byte_t MExtractPINDiode::fgHiGainLast = 29;
|
---|
90 | const Byte_t MExtractPINDiode::fgLoGainFirst = 0;
|
---|
91 | const Byte_t MExtractPINDiode::fgLoGainLast = 0;
|
---|
92 | const Byte_t MExtractPINDiode::fgLowerFitLimit = 2;
|
---|
93 | const Byte_t MExtractPINDiode::fgUpperFitLimit = 5;
|
---|
94 |
|
---|
95 | // --------------------------------------------------------------------------
|
---|
96 | //
|
---|
97 | // Default constructor.
|
---|
98 | //
|
---|
99 | // Calls:
|
---|
100 | // - SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast)
|
---|
101 | //
|
---|
102 | // - Set fPINDiodeIdx to fgPINDiodeIdx
|
---|
103 | // - Set fLowerFitLimit to fgLowerFitLimit
|
---|
104 | // - Set fUpperFitLimit to fgUpperFitLimit
|
---|
105 | //
|
---|
106 | MExtractPINDiode::MExtractPINDiode(const char *name, const char *title)
|
---|
107 | : fSlices(NULL)
|
---|
108 | {
|
---|
109 |
|
---|
110 | fName = name ? name : "MExtractPINDiode";
|
---|
111 | fTitle = title ? title : "Task to extract the signal from the FADC slices";
|
---|
112 |
|
---|
113 | SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast);
|
---|
114 | SetPINDiodeIdx();
|
---|
115 |
|
---|
116 | SetLowerFitLimit();
|
---|
117 | SetUpperFitLimit();
|
---|
118 |
|
---|
119 | fPedMean.Set(2);
|
---|
120 | }
|
---|
121 |
|
---|
122 | // --------------------------------------------------------------------------
|
---|
123 | //
|
---|
124 | // - delete the histogram fSlices, if it exists
|
---|
125 | //
|
---|
126 | MExtractPINDiode::~MExtractPINDiode()
|
---|
127 | {
|
---|
128 | if (fSlices)
|
---|
129 | delete fSlices;
|
---|
130 | }
|
---|
131 |
|
---|
132 | // --------------------------------------------------------------------------
|
---|
133 | //
|
---|
134 | // SetRange:
|
---|
135 | //
|
---|
136 | // Checks:
|
---|
137 | // - if the Hi Gain window is smaller than 4, set fHiGainLast to fHiGainFirst+3
|
---|
138 | //
|
---|
139 | // Calls:
|
---|
140 | // - MExtractor::SetRange(hifirst,hilast,lofirst,lolast);
|
---|
141 | //
|
---|
142 | // Sets:
|
---|
143 | // - fNumHiGainSamples to: (Float_t)(fHiGainLast-fHiGainFirst+1)
|
---|
144 | // - fNumLoGainSamples to: 0.
|
---|
145 | // - fSqrtHiGainSamples to: TMath::Sqrt(fNumHiGainSamples)
|
---|
146 | // - fSqrtLoGainSamples to: 0.
|
---|
147 | //
|
---|
148 | void MExtractPINDiode::SetRange(Byte_t hifirst, Byte_t hilast, Int_t lofirst, Byte_t lolast)
|
---|
149 | {
|
---|
150 |
|
---|
151 | lofirst = 0;
|
---|
152 | lolast = 0;
|
---|
153 |
|
---|
154 | const Byte_t window = hilast-hifirst+1;
|
---|
155 |
|
---|
156 | if (window<4)
|
---|
157 | {
|
---|
158 | *fLog << warn << GetDescriptor()
|
---|
159 | << Form("%s%2i%s%2i",": Total window is smaller than 4 FADC sampes, set last slice from"
|
---|
160 | ,(int)lolast," to ",(int)(lofirst+3)) << endl;
|
---|
161 | hilast = hifirst+3;
|
---|
162 | }
|
---|
163 |
|
---|
164 | MExtractor::SetRange(hifirst,hilast,lofirst,lolast);
|
---|
165 |
|
---|
166 | fNumHiGainSamples = (Float_t)(fHiGainLast-fHiGainFirst+1);
|
---|
167 | fNumLoGainSamples = 0.;
|
---|
168 |
|
---|
169 | fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
|
---|
170 | fSqrtLoGainSamples = 0.;
|
---|
171 |
|
---|
172 | }
|
---|
173 |
|
---|
174 | // --------------------------------------------------------------------------
|
---|
175 | //
|
---|
176 | // Calls:
|
---|
177 | // - MExtractor::PreProcess
|
---|
178 | //
|
---|
179 | // The following output containers are also searched and created if
|
---|
180 | // they were not found:
|
---|
181 | //
|
---|
182 | // - MRawEvtData2
|
---|
183 | // - MExtractedPINDiode
|
---|
184 | //
|
---|
185 | // Initializes fPedMean to:
|
---|
186 | // - fPedMean[0]: pedestal + AB-offset
|
---|
187 | // - fPedMean[1]: pedestal - AB-offset
|
---|
188 | //
|
---|
189 | // Initializes TH1F fSlices to [fHiGainFirst-0.5,fHiGainLast+0.5]
|
---|
190 | //
|
---|
191 | Int_t MExtractPINDiode::PreProcess(MParList *pList)
|
---|
192 | {
|
---|
193 |
|
---|
194 | if (!MExtractor::PreProcess(pList))
|
---|
195 | return kFALSE;
|
---|
196 |
|
---|
197 | fRawEvt = NULL;
|
---|
198 | fRawEvt = (MRawEvtData*)pList->FindObject(AddSerialNumber("MRawEvtData2"));
|
---|
199 | if (!fRawEvt)
|
---|
200 | {
|
---|
201 | *fLog << err << AddSerialNumber("MRawEvtData2") << " not found... aborting." << endl;
|
---|
202 | return kFALSE;
|
---|
203 | }
|
---|
204 |
|
---|
205 | fPINDiode = (MExtractedSignalPINDiode*)pList->FindCreateObj(AddSerialNumber("MExtractedSignalPINDiode"));
|
---|
206 | if (!fPINDiode)
|
---|
207 | return kFALSE;
|
---|
208 |
|
---|
209 | fPedMean.Reset();
|
---|
210 | /*
|
---|
211 | const MPedestalPix &ped = (*fPedestals)[fPINDiodeIdx];
|
---|
212 |
|
---|
213 | if (&ped)
|
---|
214 | {
|
---|
215 | fPedMean[0] = ped.GetPedestal() + ped.GetPedestalABoffset();
|
---|
216 | fPedMean[1] = ped.GetPedestal() - ped.GetPedestalABoffset();
|
---|
217 | }
|
---|
218 | else
|
---|
219 | {
|
---|
220 | *fLog << err << " Cannot find MPedestalPix of the PIN Diode (idx="
|
---|
221 | << fPINDiodeIdx << ")" << endl;
|
---|
222 | return kFALSE;
|
---|
223 | }
|
---|
224 | */
|
---|
225 | if (fSlices)
|
---|
226 | delete fSlices;
|
---|
227 |
|
---|
228 | fSlices = new TH1F("PINDiode","PIN Diode fitted slices",(Int_t)(fHiGainLast-fHiGainFirst+1),
|
---|
229 | fHiGainFirst-0.5,fHiGainLast+0.5);
|
---|
230 | fSlices->SetDirectory(NULL);
|
---|
231 |
|
---|
232 | return kTRUE;
|
---|
233 | }
|
---|
234 |
|
---|
235 | // --------------------------------------------------------------------------
|
---|
236 | //
|
---|
237 | // The ReInit calls:
|
---|
238 | // - MExtractor::ReInit()
|
---|
239 | // - fPINDiode->SetUsedFADCSlices(fHiGainFirst, fLoGainLast);
|
---|
240 | //
|
---|
241 | Bool_t MExtractPINDiode::ReInit(MParList *pList)
|
---|
242 | {
|
---|
243 |
|
---|
244 | MExtractor::ReInit(pList);
|
---|
245 |
|
---|
246 | fPINDiode->SetUsedFADCSlices(fHiGainFirst, fLoGainLast);
|
---|
247 |
|
---|
248 | *fLog << endl;
|
---|
249 | *fLog << inf << "Taking " << fNumHiGainSamples
|
---|
250 | << " HiGain samples from slice " << (Int_t)fHiGainFirst
|
---|
251 | << " to " << (Int_t)(fHiGainLast/*+fHiLoLast*/) << " incl" << endl;
|
---|
252 |
|
---|
253 | return kTRUE;
|
---|
254 | }
|
---|
255 |
|
---|
256 |
|
---|
257 |
|
---|
258 |
|
---|
259 | // ----------------------------------------------------------------------------------
|
---|
260 | //
|
---|
261 | // Extracts the (pedestal-subtracted) FADC slices between fHiGainFirst and fHiGainLast
|
---|
262 | // and fills them into the histogram fSlices. Memorizes the position of maximum at
|
---|
263 | // maxpos.
|
---|
264 | //
|
---|
265 | // Checks for saturation
|
---|
266 | //
|
---|
267 | // Fits fSlices to a Gaussian in the ranges: maxpos-fLowerFitLimit, maxpos+fUpperFitLimit
|
---|
268 | //
|
---|
269 | // Writes fit results into MExtractedSignalPINDiode
|
---|
270 | //
|
---|
271 | Int_t MExtractPINDiode::Process()
|
---|
272 | {
|
---|
273 |
|
---|
274 |
|
---|
275 | MRawEvtPixelIter pixel(fRawEvt);
|
---|
276 |
|
---|
277 | fPINDiode->Clear();
|
---|
278 | fSlices->Reset();
|
---|
279 |
|
---|
280 | pixel.Jump(fPINDiodeIdx);
|
---|
281 |
|
---|
282 | Byte_t sat = 0;
|
---|
283 |
|
---|
284 | const Int_t higainsamples = fRunHeader->GetNumSamplesHiGain();
|
---|
285 | const Int_t logainsamples = fRunHeader->GetNumSamplesLoGain();
|
---|
286 |
|
---|
287 | const Bool_t higainabflag = pixel.HasABFlag();
|
---|
288 | Byte_t *ptr = pixel.GetHiGainSamples()+fHiGainFirst;
|
---|
289 | Byte_t *end = ptr+higainsamples;
|
---|
290 |
|
---|
291 | Int_t cnt=0;
|
---|
292 |
|
---|
293 | Float_t max = 0.;
|
---|
294 | Int_t maxpos = 0;
|
---|
295 |
|
---|
296 | while (ptr<end)
|
---|
297 | {
|
---|
298 |
|
---|
299 | if (*ptr >= fSaturationLimit)
|
---|
300 | {
|
---|
301 | sat++;
|
---|
302 | break;
|
---|
303 | }
|
---|
304 |
|
---|
305 | const Float_t cont = (Float_t)*ptr - fPedMean[(cnt + higainabflag) & 0x1];
|
---|
306 | fSlices->Fill(cnt,cont);
|
---|
307 |
|
---|
308 | if (cont > max)
|
---|
309 | {
|
---|
310 | max = cont;
|
---|
311 | maxpos = cnt;
|
---|
312 | }
|
---|
313 |
|
---|
314 | ptr++;
|
---|
315 | cnt++;
|
---|
316 | }
|
---|
317 |
|
---|
318 | cnt = 0;
|
---|
319 |
|
---|
320 | if (logainsamples>0 && !sat)
|
---|
321 | {
|
---|
322 |
|
---|
323 | ptr = pixel.GetLoGainSamples();
|
---|
324 | end = ptr+logainsamples;
|
---|
325 |
|
---|
326 | const Bool_t logainabflag = (higainabflag + higainsamples) & 0x1;
|
---|
327 |
|
---|
328 | while (ptr<end)
|
---|
329 | {
|
---|
330 |
|
---|
331 | if (*ptr >= fSaturationLimit)
|
---|
332 | {
|
---|
333 | sat++;
|
---|
334 | break;
|
---|
335 | }
|
---|
336 |
|
---|
337 | const Float_t cont = (Float_t)*ptr - fPedMean[(cnt + logainabflag) & 0x1];
|
---|
338 |
|
---|
339 | fSlices->Fill(cnt+ higainsamples,cont);
|
---|
340 |
|
---|
341 | if (cont > max)
|
---|
342 | {
|
---|
343 | max = cont;
|
---|
344 | maxpos = cnt+higainsamples;
|
---|
345 | }
|
---|
346 | ptr++;
|
---|
347 | cnt++;
|
---|
348 | }
|
---|
349 | }
|
---|
350 |
|
---|
351 | if (sat)
|
---|
352 | {
|
---|
353 | fPINDiode->SetSaturation(1);
|
---|
354 | return kTRUE;
|
---|
355 | }
|
---|
356 |
|
---|
357 | fSlices->Fit("gaus", "Q0", "", maxpos-fLowerFitLimit,maxpos+fUpperFitLimit);
|
---|
358 |
|
---|
359 | TF1 &gaus = *fSlices->GetFunction("gaus");
|
---|
360 |
|
---|
361 | fPINDiode->SetExtractedSignal(gaus.GetParameter(0), gaus.GetParError(0));
|
---|
362 | fPINDiode->SetExtractedTime (gaus.GetParameter(1), gaus.GetParError(1));
|
---|
363 | fPINDiode->SetExtractedSigma (gaus.GetParameter(2), gaus.GetParError(2));
|
---|
364 | fPINDiode->SetExtractedChi2 (gaus.GetChisquare());
|
---|
365 | fPINDiode->SetReadyToSave();
|
---|
366 |
|
---|
367 | return kTRUE;
|
---|
368 | }
|
---|
369 |
|
---|
370 | // ----------------------------------------------------------------------------------
|
---|
371 | //
|
---|
372 | // deletes fSlices and sets pointer to NULL
|
---|
373 | //
|
---|
374 | Int_t MExtractPINDiode::PostProcess()
|
---|
375 | {
|
---|
376 |
|
---|
377 | delete fSlices;
|
---|
378 | fSlices = NULL;
|
---|
379 |
|
---|
380 | return kTRUE;
|
---|
381 | }
|
---|