source: trunk/MagicSoft/Mars/msignal/MExtractPINDiode.cc@ 7344

Last change on this file since 7344 was 7059, checked in by tbretz, 19 years ago
*** empty log message ***
File size: 9.5 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): 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 "MRawEvtPixelIter.h"
76
77#include "MPedestalCam.h"
78#include "MPedestalPix.h"
79
80#include "MExtractedSignalPINDiode.h"
81
82ClassImp(MExtractPINDiode);
83
84using namespace std;
85
86const UInt_t MExtractPINDiode::fgPINDiodeIdx = 3;
87const Byte_t MExtractPINDiode::fgHiGainFirst = 0;
88const Byte_t MExtractPINDiode::fgHiGainLast = 29;
89const Byte_t MExtractPINDiode::fgLoGainFirst = 0;
90const Byte_t MExtractPINDiode::fgLoGainLast = 0;
91const Byte_t MExtractPINDiode::fgLowerFitLimit = 2;
92const Byte_t MExtractPINDiode::fgUpperFitLimit = 5;
93
94// --------------------------------------------------------------------------
95//
96// Default constructor.
97//
98// Calls:
99// - SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast)
100//
101// - Set fPINDiodeIdx to fgPINDiodeIdx
102// - Set fLowerFitLimit to fgLowerFitLimit
103// - Set fUpperFitLimit to fgUpperFitLimit
104//
105MExtractPINDiode::MExtractPINDiode(const char *name, const char *title)
106 : fSlices(NULL)
107{
108
109 fName = name ? name : "MExtractPINDiode";
110 fTitle = title ? title : "Task to extract the signal from the FADC slices";
111
112 SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast);
113 SetPINDiodeIdx();
114
115 SetLowerFitLimit();
116 SetUpperFitLimit();
117
118 fPedMean.Set(2);
119}
120
121// --------------------------------------------------------------------------
122//
123// - delete the histogram fSlices, if it exists
124//
125MExtractPINDiode::~MExtractPINDiode()
126{
127 if (fSlices)
128 delete fSlices;
129}
130
131// --------------------------------------------------------------------------
132//
133// SetRange:
134//
135// Checks:
136// - if the Hi Gain window is smaller than 4, set fHiGainLast to fHiGainFirst+3
137//
138// Calls:
139// - MExtractor::SetRange(hifirst,hilast,lofirst,lolast);
140//
141// Sets:
142// - fNumHiGainSamples to: (Float_t)(fHiGainLast-fHiGainFirst+1)
143// - fNumLoGainSamples to: 0.
144// - fSqrtHiGainSamples to: TMath::Sqrt(fNumHiGainSamples)
145// - fSqrtLoGainSamples to: 0.
146//
147void MExtractPINDiode::SetRange(Byte_t hifirst, Byte_t hilast, Byte_t lofirst, Byte_t lolast)
148{
149
150 lofirst = 0;
151 lolast = 0;
152
153 const Byte_t window = hilast-hifirst+1;
154
155 if (window<4)
156 {
157 *fLog << warn << GetDescriptor()
158 << Form("%s%2i%s%2i",": Total window is smaller than 4 FADC sampes, set last slice from"
159 ,(int)lolast," to ",(int)(lofirst+3)) << endl;
160 hilast = hifirst+3;
161 }
162
163 MExtractor::SetRange(hifirst,hilast,lofirst,lolast);
164
165 fNumHiGainSamples = (Float_t)(fHiGainLast-fHiGainFirst+1);
166 fNumLoGainSamples = 0.;
167
168 fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
169 fSqrtLoGainSamples = 0.;
170
171}
172
173// --------------------------------------------------------------------------
174//
175// Calls:
176// - MExtractor::PreProcess
177//
178// The following output containers are also searched and created if
179// they were not found:
180//
181// - MRawEvtData2
182// - MExtractedPINDiode
183//
184// Initializes fPedMean to:
185// - fPedMean[0]: pedestal + AB-offset
186// - fPedMean[1]: pedestal - AB-offset
187//
188// Initializes TH1F fSlices to [fHiGainFirst-0.5,fHiGainLast+0.5]
189//
190Int_t MExtractPINDiode::PreProcess(MParList *pList)
191{
192
193 if (!MExtractor::PreProcess(pList))
194 return kFALSE;
195
196 fRawEvt = NULL;
197 fRawEvt = (MRawEvtData*)pList->FindObject(AddSerialNumber("MRawEvtData2"));
198 if (!fRawEvt)
199 {
200 *fLog << err << AddSerialNumber("MRawEvtData2") << " not found... aborting." << endl;
201 return kFALSE;
202 }
203
204 fPINDiode = (MExtractedSignalPINDiode*)pList->FindCreateObj(AddSerialNumber("MExtractedSignalPINDiode"));
205 if (!fPINDiode)
206 return kFALSE;
207
208 fPedMean.Reset();
209
210 const MPedestalPix &ped = (*fPedestals)[fPINDiodeIdx];
211
212 if (&ped)
213 {
214 fPedMean[0] = ped.GetPedestal() + ped.GetPedestalABoffset();
215 fPedMean[1] = ped.GetPedestal() - ped.GetPedestalABoffset();
216 }
217 else
218 {
219 *fLog << err << " Cannot find MPedestalPix of the PIN Diode (idx="
220 << fPINDiodeIdx << ")" << endl;
221 return kFALSE;
222 }
223
224 if (fSlices)
225 delete fSlices;
226
227 fSlices = new TH1F("PINDiode","PIN Diode fitted slices",(Int_t)(fHiGainLast-fHiGainFirst+1),
228 fHiGainFirst-0.5,fHiGainLast+0.5);
229 fSlices->SetDirectory(NULL);
230
231 return kTRUE;
232}
233
234// --------------------------------------------------------------------------
235//
236// The ReInit calls:
237// - MExtractor::ReInit()
238// - fPINDiode->SetUsedFADCSlices(fHiGainFirst, fLoGainLast);
239//
240Bool_t MExtractPINDiode::ReInit(MParList *pList)
241{
242
243 MExtractor::ReInit(pList);
244
245 fPINDiode->SetUsedFADCSlices(fHiGainFirst, fLoGainLast);
246
247 *fLog << endl;
248 *fLog << inf << "Taking " << fNumHiGainSamples
249 << " HiGain samples from slice " << (Int_t)fHiGainFirst
250 << " to " << (Int_t)(fHiGainLast+fHiLoLast) << " incl" << endl;
251
252 return kTRUE;
253}
254
255
256
257
258// ----------------------------------------------------------------------------------
259//
260// Extracts the (pedestal-subtracted) FADC slices between fHiGainFirst and fHiGainLast
261// and fills them into the histogram fSlices. Memorizes the position of maximum at
262// maxpos.
263//
264// Checks for saturation
265//
266// Fits fSlices to a Gaussian in the ranges: maxpos-fLowerFitLimit, maxpos+fUpperFitLimit
267//
268// Writes fit results into MExtractedSignalPINDiode
269//
270Int_t MExtractPINDiode::Process()
271{
272
273
274 MRawEvtPixelIter pixel(fRawEvt);
275
276 fPINDiode->Clear();
277 fSlices->Reset();
278
279 pixel.Jump(fPINDiodeIdx);
280
281 Byte_t sat = 0;
282
283 Int_t higainsamples = pixel.GetNumHiGainSamples();
284 Int_t logainsamples = pixel.GetNumLoGainSamples();
285
286 const Bool_t higainabflag = pixel.HasABFlag();
287 Byte_t *ptr = pixel.GetHiGainSamples()+fHiGainFirst;
288 Byte_t *end = ptr+higainsamples;
289
290 Int_t cnt=0;
291
292 Float_t max = 0.;
293 Int_t maxpos = 0;
294
295 while (ptr<end)
296 {
297
298 if (*ptr >= fSaturationLimit)
299 {
300 sat++;
301 break;
302 }
303
304 const Float_t cont = (Float_t)*ptr - fPedMean[(cnt + higainabflag) & 0x1];
305 fSlices->Fill(cnt,cont);
306
307 if (cont > max)
308 {
309 max = cont;
310 maxpos = cnt;
311 }
312
313 ptr++;
314 cnt++;
315 }
316
317 cnt = 0;
318
319 if (pixel.HasLoGain() && !sat)
320 {
321
322 ptr = pixel.GetLoGainSamples();
323 end = ptr+logainsamples;
324
325 const Bool_t logainabflag = (higainabflag + pixel.GetNumHiGainSamples()) & 0x1;
326
327 while (ptr<end)
328 {
329
330 if (*ptr >= fSaturationLimit)
331 {
332 sat++;
333 break;
334 }
335
336 const Float_t cont = (Float_t)*ptr - fPedMean[(cnt + logainabflag) & 0x1];
337
338 fSlices->Fill(cnt+ higainsamples,cont);
339
340 if (cont > max)
341 {
342 max = cont;
343 maxpos = cnt+higainsamples;
344 }
345 ptr++;
346 cnt++;
347 }
348 }
349
350 if (sat)
351 {
352 fPINDiode->SetSaturation(1);
353 return kTRUE;
354 }
355
356 fSlices->Fit("gaus", "Q0", "", maxpos-fLowerFitLimit,maxpos+fUpperFitLimit);
357
358 TF1 &gaus = *fSlices->GetFunction("gaus");
359
360 fPINDiode->SetExtractedSignal(gaus.GetParameter(0), gaus.GetParError(0));
361 fPINDiode->SetExtractedTime (gaus.GetParameter(1), gaus.GetParError(1));
362 fPINDiode->SetExtractedSigma (gaus.GetParameter(2), gaus.GetParError(2));
363 fPINDiode->SetExtractedChi2 (gaus.GetChisquare());
364 fPINDiode->SetReadyToSave();
365
366 return kTRUE;
367}
368
369// ----------------------------------------------------------------------------------
370//
371// deletes fSlices and sets pointer to NULL
372//
373Int_t MExtractPINDiode::PostProcess()
374{
375
376 delete fSlices;
377 fSlices = NULL;
378
379 return kTRUE;
380}
Note: See TracBrowser for help on using the repository browser.