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 analyzing 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 | ! Author(s): Markus Gaug 09/2004 <mailto:markus@ifae.es>
|
---|
18 | !
|
---|
19 | ! Copyright: MAGIC Software Development, 2002-2004
|
---|
20 | !
|
---|
21 | !
|
---|
22 | \* ======================================================================== */
|
---|
23 | //////////////////////////////////////////////////////////////////////////////
|
---|
24 | //
|
---|
25 | // MExtractTimeAndChargeSpline
|
---|
26 | //
|
---|
27 | // Fast Spline extractor using a cubic spline algorithm, adapted from
|
---|
28 | // Numerical Recipes in C++, 2nd edition, pp. 116-119.
|
---|
29 | //
|
---|
30 | // The coefficients "ya" are here denoted as "fHiGainSignal" and "fLoGainSignal"
|
---|
31 | // which means the FADC value subtracted by the clock-noise corrected pedestal.
|
---|
32 | //
|
---|
33 | // The coefficients "y2a" get immediately divided 6. and are called here
|
---|
34 | // "fHiGainSecondDeriv" and "fLoGainSecondDeriv" although they are now not exactly
|
---|
35 | // the second derivative coefficients any more.
|
---|
36 | //
|
---|
37 | // The calculation of the cubic-spline interpolated value "y" on a point
|
---|
38 | // "x" along the FADC-slices axis becomes:
|
---|
39 | //
|
---|
40 | // y = a*fHiGainSignal[klo] + b*fHiGainSignal[khi]
|
---|
41 | // + (a*a*a-a)*fHiGainSecondDeriv[klo] + (b*b*b-b)*fHiGainSecondDeriv[khi]
|
---|
42 | //
|
---|
43 | // with:
|
---|
44 | // a = (khi - x)
|
---|
45 | // b = (x - klo)
|
---|
46 | //
|
---|
47 | // and "klo" being the lower bin edge FADC index and "khi" the upper bin edge FADC index.
|
---|
48 | // fHiGainSignal[klo] and fHiGainSignal[khi] are the FADC values at "klo" and "khi".
|
---|
49 | //
|
---|
50 | // An analogues formula is used for the low-gain values.
|
---|
51 | //
|
---|
52 | // The coefficients fHiGainSecondDeriv and fLoGainSecondDeriv are calculated with the
|
---|
53 | // following simplified algorithm:
|
---|
54 | //
|
---|
55 | // for (Int_t i=1;i<range-1;i++) {
|
---|
56 | // pp = fHiGainSecondDeriv[i-1] + 4.;
|
---|
57 | // fHiGainFirstDeriv[i] = fHiGainSignal[i+1] - 2.*fHiGainSignal[i] + fHiGainSignal[i-1]
|
---|
58 | // fHiGainFirstDeriv[i] = (6.0*fHiGainFirstDeriv[i]-fHiGainFirstDeriv[i-1])/pp;
|
---|
59 | // }
|
---|
60 | //
|
---|
61 | // for (Int_t k=range-2;k>=0;k--)
|
---|
62 | // fHiGainSecondDeriv[k] = (fHiGainSecondDeriv[k]*fHiGainSecondDeriv[k+1] + fHiGainFirstDeriv[k])/6.;
|
---|
63 | //
|
---|
64 | //
|
---|
65 | // This algorithm takes advantage of the fact that the x-values are all separated by exactly 1
|
---|
66 | // which simplifies the Numerical Recipes algorithm.
|
---|
67 | // (Note that the variables "fHiGainFirstDeriv" are not real first derivative coefficients.)
|
---|
68 | //
|
---|
69 | // The algorithm to search the time proceeds as follows:
|
---|
70 | //
|
---|
71 | // 1) Calculate all fHiGainSignal from fHiGainFirst to fHiGainLast
|
---|
72 | // (note that an "overlap" to the low-gain arrays is possible: i.e. fHiGainLast>14 in the case of
|
---|
73 | // the MAGIC FADCs).
|
---|
74 | // 2) Remember the position of the slice with the highest content "fAbMax" at "fAbMaxPos".
|
---|
75 | // 3) If one or more slices are saturated or fAbMaxPos is less than 2 slices from fHiGainFirst,
|
---|
76 | // return fAbMaxPos as time and fAbMax as charge (note that the pedestal is subtracted here).
|
---|
77 | // 4) Calculate all fHiGainSecondDeriv from the fHiGainSignal array
|
---|
78 | // 5) Search for the maximum, starting in interval fAbMaxPos-1 in steps of 0.2 till fAbMaxPos-0.2.
|
---|
79 | // If no maximum is found, go to interval fAbMaxPos+1.
|
---|
80 | // --> 4 function evaluations
|
---|
81 | // 6) Search for the absolute maximum from fAbMaxPos to fAbMaxPos+1 in steps of 0.2
|
---|
82 | // --> 4 function evaluations
|
---|
83 | // 7) Try a better precision searching from new max. position fAbMaxPos-0.2 to fAbMaxPos+0.2
|
---|
84 | // in steps of 0.025 (83 psec. in the case of the MAGIC FADCs).
|
---|
85 | // --> 14 function evaluations
|
---|
86 | // 8) If Time Extraction Type kMaximum has been chosen, the position of the found maximum is
|
---|
87 | // returned, else:
|
---|
88 | // 9) The Half Maximum is calculated.
|
---|
89 | // 10) fHiGainSignal is called beginning from fAbMaxPos-1 backwards until a value smaller than fHalfMax
|
---|
90 | // is found at "klo".
|
---|
91 | // 11) Then, the spline value between "klo" and "klo"+1 is halfed by means of bisection as long as
|
---|
92 | // the difference between fHalfMax and spline evaluation is less than fResolution (default: 0.01).
|
---|
93 | // --> maximum 12 interations.
|
---|
94 | //
|
---|
95 | // The algorithm to search the charge proceeds as follows:
|
---|
96 | //
|
---|
97 | // 1) If Charge Type: kAmplitude was chosen, return the Maximum of the spline, found during the
|
---|
98 | // time search.
|
---|
99 | // 2) If Charge Type: kIntegral was chosen, sum the fHiGainSignal between:
|
---|
100 | // (Int_t)(fAbMaxPos - fRiseTimeHiGain) and
|
---|
101 | // (Int_t)(fAbMaxPos + fFallTimeHiGain)
|
---|
102 | // (default: fRiseTime: 1.5, fFallTime: 4.5)
|
---|
103 | // sum the fLoGainSignal between:
|
---|
104 | // (Int_t)(fAbMaxPos - fRiseTimeHiGain*fLoGainStretch) and
|
---|
105 | // (Int_t)(fAbMaxPos + fFallTimeHiGain*fLoGainStretch)
|
---|
106 | // (default: fLoGainStretch: 1.5)
|
---|
107 | //
|
---|
108 | // The values: fNumHiGainSamples and fNumLoGainSamples are set to:
|
---|
109 | // 1) If Charge Type: kAmplitude was chosen: 1.
|
---|
110 | // 2) If Charge Type: kIntegral was chosen: fRiseTimeHiGain + fFallTimeHiGain
|
---|
111 | // or: fNumHiGainSamples*fLoGainStretch in the case of the low-gain
|
---|
112 | //
|
---|
113 | // Call: SetRange(fHiGainFirst, fHiGainLast, fLoGainFirst, fLoGainLast)
|
---|
114 | // to modify the ranges.
|
---|
115 | //
|
---|
116 | // Defaults:
|
---|
117 | // fHiGainFirst = 2
|
---|
118 | // fHiGainLast = 14
|
---|
119 | // fLoGainFirst = 2
|
---|
120 | // fLoGainLast = 14
|
---|
121 | //
|
---|
122 | // Call: SetResolution() to define the resolution of the half-maximum search.
|
---|
123 | // Default: 0.01
|
---|
124 | //
|
---|
125 | // Call: SetRiseTime() and SetFallTime() to define the integration ranges
|
---|
126 | // for the case, the extraction type kIntegral has been chosen.
|
---|
127 | //
|
---|
128 | // Call: - SetChargeType(MExtractTimeAndChargeSpline::kAmplitude) for the
|
---|
129 | // computation of the amplitude at the maximum (default) and extraction
|
---|
130 | // the position of the maximum (default)
|
---|
131 | // --> no further function evaluation needed
|
---|
132 | // - SetChargeType(MExtractTimeAndChargeSpline::kIntegral) for the
|
---|
133 | // computation of the integral beneith the spline between fRiseTimeHiGain
|
---|
134 | // from the position of the maximum to fFallTimeHiGain after the position of
|
---|
135 | // the maximum. The Low Gain is computed with half a slice more at the rising
|
---|
136 | // edge and half a slice more at the falling edge.
|
---|
137 | // The time of the half maximum is returned.
|
---|
138 | // --> needs one function evaluations but is more precise
|
---|
139 | //
|
---|
140 | //////////////////////////////////////////////////////////////////////////////
|
---|
141 | #include "MExtractTimeAndChargeSpline.h"
|
---|
142 |
|
---|
143 | #include "MPedestalPix.h"
|
---|
144 |
|
---|
145 | #include "MLog.h"
|
---|
146 | #include "MLogManip.h"
|
---|
147 |
|
---|
148 | ClassImp(MExtractTimeAndChargeSpline);
|
---|
149 |
|
---|
150 | using namespace std;
|
---|
151 |
|
---|
152 | const Byte_t MExtractTimeAndChargeSpline::fgHiGainFirst = 2;
|
---|
153 | const Byte_t MExtractTimeAndChargeSpline::fgHiGainLast = 14;
|
---|
154 | const Byte_t MExtractTimeAndChargeSpline::fgLoGainFirst = 2;
|
---|
155 | const Byte_t MExtractTimeAndChargeSpline::fgLoGainLast = 14;
|
---|
156 | const Float_t MExtractTimeAndChargeSpline::fgResolution = 0.05;
|
---|
157 | const Float_t MExtractTimeAndChargeSpline::fgRiseTimeHiGain = 0.5;
|
---|
158 | const Float_t MExtractTimeAndChargeSpline::fgFallTimeHiGain = 1.5;
|
---|
159 | const Float_t MExtractTimeAndChargeSpline::fgLoGainStretch = 1.5;
|
---|
160 | const Float_t MExtractTimeAndChargeSpline::fgOffsetLoGain = 1.7; // 5 ns
|
---|
161 | // --------------------------------------------------------------------------
|
---|
162 | //
|
---|
163 | // Default constructor.
|
---|
164 | //
|
---|
165 | // Calls:
|
---|
166 | // - SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast)
|
---|
167 | //
|
---|
168 | // Initializes:
|
---|
169 | // - fResolution to fgResolution
|
---|
170 | // - fRiseTimeHiGain to fgRiseTimeHiGain
|
---|
171 | // - fFallTimeHiGain to fgFallTimeHiGain
|
---|
172 | // - Charge Extraction Type to kAmplitude
|
---|
173 | // - fLoGainStretch to fgLoGainStretch
|
---|
174 | //
|
---|
175 | MExtractTimeAndChargeSpline::MExtractTimeAndChargeSpline(const char *name, const char *title)
|
---|
176 | : fAbMax(0.), fAbMaxPos(0.), fHalfMax(0.),
|
---|
177 | fRiseTimeHiGain(fgRiseTimeHiGain), fFallTimeHiGain(fgFallTimeHiGain),
|
---|
178 | fRandomIter(0), fExtractionType(kIntegral)
|
---|
179 | {
|
---|
180 |
|
---|
181 | fName = name ? name : "MExtractTimeAndChargeSpline";
|
---|
182 | fTitle = title ? title : "Calculate photons arrival time using a fast spline";
|
---|
183 |
|
---|
184 | SetResolution();
|
---|
185 | SetLoGainStretch();
|
---|
186 | SetOffsetLoGain(fgOffsetLoGain);
|
---|
187 |
|
---|
188 | SetRange(fgHiGainFirst, fgHiGainLast, fgLoGainFirst, fgLoGainLast);
|
---|
189 | }
|
---|
190 |
|
---|
191 |
|
---|
192 | //-------------------------------------------------------------------
|
---|
193 | //
|
---|
194 | // Set the ranges
|
---|
195 | // In order to set the fNum...Samples variables correctly for the case,
|
---|
196 | // the integral is computed, have to overwrite this function and make an
|
---|
197 | // explicit call to SetChargeType().
|
---|
198 | //
|
---|
199 | void MExtractTimeAndChargeSpline::SetRange(Byte_t hifirst, Byte_t hilast, Byte_t lofirst, Byte_t lolast)
|
---|
200 | {
|
---|
201 |
|
---|
202 | MExtractor::SetRange(hifirst, hilast, lofirst, lolast);
|
---|
203 |
|
---|
204 | SetChargeType(fExtractionType);
|
---|
205 | }
|
---|
206 |
|
---|
207 | //-------------------------------------------------------------------
|
---|
208 | //
|
---|
209 | // Set the Charge Extraction type. Possible are:
|
---|
210 | // - kAmplitude: Search the value of the spline at the maximum
|
---|
211 | // - kIntegral: Integral the spline from fHiGainFirst to fHiGainLast,
|
---|
212 | // by counting the edge bins only half and setting the
|
---|
213 | // second derivative to zero, there.
|
---|
214 | //
|
---|
215 | void MExtractTimeAndChargeSpline::SetChargeType( ExtractionType_t typ )
|
---|
216 | {
|
---|
217 |
|
---|
218 | fExtractionType = typ;
|
---|
219 |
|
---|
220 | if (fExtractionType == kAmplitude)
|
---|
221 | {
|
---|
222 | fNumHiGainSamples = 1.;
|
---|
223 | fNumLoGainSamples = fLoGainLast ? 1. : 0.;
|
---|
224 | fSqrtHiGainSamples = 1.;
|
---|
225 | fSqrtLoGainSamples = 1.;
|
---|
226 | fWindowSizeHiGain = 1;
|
---|
227 | fWindowSizeLoGain = 1;
|
---|
228 | fRiseTimeHiGain = 0.5;
|
---|
229 |
|
---|
230 | return;
|
---|
231 | }
|
---|
232 |
|
---|
233 | if (fExtractionType == kIntegral)
|
---|
234 | {
|
---|
235 |
|
---|
236 | fNumHiGainSamples = fRiseTimeHiGain + fFallTimeHiGain;
|
---|
237 | fNumLoGainSamples = fLoGainLast ? fRiseTimeLoGain + fFallTimeLoGain : 0.;
|
---|
238 | // fNumLoGainSamples *= 0.75;
|
---|
239 |
|
---|
240 | fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
|
---|
241 | fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
|
---|
242 | fWindowSizeHiGain = (Int_t)(fRiseTimeHiGain + fFallTimeHiGain);
|
---|
243 | fWindowSizeLoGain = (Int_t)(fRiseTimeLoGain + fFallTimeLoGain);
|
---|
244 | // fNumLoGainSamples *= 0.75;
|
---|
245 | }
|
---|
246 | }
|
---|
247 |
|
---|
248 | // --------------------------------------------------------------------------
|
---|
249 | //
|
---|
250 | // InitArrays
|
---|
251 | //
|
---|
252 | // Gets called in the ReInit() and initialized the arrays
|
---|
253 | //
|
---|
254 | Bool_t MExtractTimeAndChargeSpline::InitArrays()
|
---|
255 | {
|
---|
256 |
|
---|
257 | Int_t range = fHiGainLast - fHiGainFirst + 1 + fHiLoLast;
|
---|
258 |
|
---|
259 | fHiGainSignal .Set(range);
|
---|
260 | fHiGainFirstDeriv .Set(range);
|
---|
261 | fHiGainSecondDeriv.Set(range);
|
---|
262 |
|
---|
263 | range = fLoGainLast - fLoGainFirst + 1;
|
---|
264 |
|
---|
265 | fLoGainSignal .Set(range);
|
---|
266 | fLoGainFirstDeriv .Set(range);
|
---|
267 | fLoGainSecondDeriv.Set(range);
|
---|
268 |
|
---|
269 | fHiGainSignal .Reset();
|
---|
270 | fHiGainFirstDeriv .Reset();
|
---|
271 | fHiGainSecondDeriv.Reset();
|
---|
272 |
|
---|
273 | fLoGainSignal .Reset();
|
---|
274 | fLoGainFirstDeriv .Reset();
|
---|
275 | fLoGainSecondDeriv.Reset();
|
---|
276 |
|
---|
277 | if (fExtractionType == kAmplitude)
|
---|
278 | {
|
---|
279 | fNumHiGainSamples = 1.;
|
---|
280 | fNumLoGainSamples = fLoGainLast ? 1. : 0.;
|
---|
281 | fSqrtHiGainSamples = 1.;
|
---|
282 | fSqrtLoGainSamples = 1.;
|
---|
283 | fWindowSizeHiGain = 1;
|
---|
284 | fWindowSizeLoGain = 1;
|
---|
285 | fRiseTimeHiGain = 0.5;
|
---|
286 | }
|
---|
287 |
|
---|
288 | fRiseTimeLoGain = fRiseTimeHiGain * fLoGainStretch;
|
---|
289 | fFallTimeLoGain = fFallTimeHiGain * fLoGainStretch;
|
---|
290 |
|
---|
291 | if (fExtractionType == kIntegral)
|
---|
292 | {
|
---|
293 |
|
---|
294 | fNumHiGainSamples = fRiseTimeHiGain + fFallTimeHiGain;
|
---|
295 | fNumLoGainSamples = fLoGainLast ? fRiseTimeLoGain + fFallTimeLoGain : 0.;
|
---|
296 | // fNumLoGainSamples *= 0.75;
|
---|
297 |
|
---|
298 | fSqrtHiGainSamples = TMath::Sqrt(fNumHiGainSamples);
|
---|
299 | fSqrtLoGainSamples = TMath::Sqrt(fNumLoGainSamples);
|
---|
300 | fWindowSizeHiGain = (Int_t)(fRiseTimeHiGain + fFallTimeHiGain);
|
---|
301 | fWindowSizeLoGain = (Int_t)(fRiseTimeLoGain + fFallTimeLoGain);
|
---|
302 | }
|
---|
303 |
|
---|
304 | return kTRUE;
|
---|
305 |
|
---|
306 | }
|
---|
307 |
|
---|
308 | // --------------------------------------------------------------------------
|
---|
309 | //
|
---|
310 | // Calculates the arrival time and charge for each pixel
|
---|
311 | //
|
---|
312 | void MExtractTimeAndChargeSpline::FindTimeAndChargeHiGain(Byte_t *first, Byte_t *logain, Float_t &sum, Float_t &dsum,
|
---|
313 | Float_t &time, Float_t &dtime,
|
---|
314 | Byte_t &sat, const MPedestalPix &ped, const Bool_t abflag)
|
---|
315 | {
|
---|
316 |
|
---|
317 | Int_t range = fHiGainLast - fHiGainFirst + 1;
|
---|
318 | const Byte_t *end = first + range;
|
---|
319 | Byte_t *p = first;
|
---|
320 |
|
---|
321 | sat = 0;
|
---|
322 |
|
---|
323 | const Float_t pedes = ped.GetPedestal();
|
---|
324 | const Float_t ABoffs = ped.GetPedestalABoffset();
|
---|
325 |
|
---|
326 | const Float_t pedmean[2] = { pedes + ABoffs, pedes - ABoffs };
|
---|
327 |
|
---|
328 | fAbMax = 0.;
|
---|
329 | fAbMaxPos = 0.;
|
---|
330 | fHalfMax = 0.;
|
---|
331 | fMaxBinContent = 0;
|
---|
332 | Int_t maxpos = 0;
|
---|
333 |
|
---|
334 | //
|
---|
335 | // Check for saturation in all other slices
|
---|
336 | //
|
---|
337 | Int_t ids = fHiGainFirst;
|
---|
338 | Float_t *sample = fHiGainSignal.GetArray();
|
---|
339 | while (p<end)
|
---|
340 | {
|
---|
341 |
|
---|
342 | *sample++ = (Float_t)*p - pedmean[(ids++ + abflag) & 0x1];
|
---|
343 |
|
---|
344 | if (*p > fMaxBinContent)
|
---|
345 | {
|
---|
346 | maxpos = ids-fHiGainFirst-1;
|
---|
347 | fMaxBinContent = *p;
|
---|
348 | }
|
---|
349 |
|
---|
350 | if (*p++ >= fSaturationLimit)
|
---|
351 | if (!sat)
|
---|
352 | sat = ids-2;
|
---|
353 |
|
---|
354 | }
|
---|
355 |
|
---|
356 | if (fHiLoLast != 0)
|
---|
357 | {
|
---|
358 |
|
---|
359 | end = logain + fHiLoLast;
|
---|
360 |
|
---|
361 | while (logain<end)
|
---|
362 | {
|
---|
363 |
|
---|
364 | *sample++ = (Float_t)*logain - pedmean[(ids++ + abflag) & 0x1];
|
---|
365 |
|
---|
366 | if (*logain > fMaxBinContent)
|
---|
367 | {
|
---|
368 | maxpos = ids-fHiGainFirst-1;
|
---|
369 | fMaxBinContent = *logain;
|
---|
370 | }
|
---|
371 |
|
---|
372 | if (*logain++ >= fSaturationLimit)
|
---|
373 | if (!sat)
|
---|
374 | sat = ids-2;
|
---|
375 |
|
---|
376 | range++;
|
---|
377 | }
|
---|
378 | }
|
---|
379 |
|
---|
380 | fAbMax = fHiGainSignal[maxpos];
|
---|
381 |
|
---|
382 | Float_t pp;
|
---|
383 |
|
---|
384 | fHiGainSecondDeriv[0] = 0.;
|
---|
385 | fHiGainFirstDeriv[0] = 0.;
|
---|
386 |
|
---|
387 | for (Int_t i=1;i<range-1;i++)
|
---|
388 | {
|
---|
389 | pp = fHiGainSecondDeriv[i-1] + 4.;
|
---|
390 | fHiGainSecondDeriv[i] = -1.0/pp;
|
---|
391 | fHiGainFirstDeriv [i] = fHiGainSignal[i+1] - fHiGainSignal[i] - fHiGainSignal[i] + fHiGainSignal[i-1];
|
---|
392 | fHiGainFirstDeriv [i] = (6.0*fHiGainFirstDeriv[i]-fHiGainFirstDeriv[i-1])/pp;
|
---|
393 | }
|
---|
394 |
|
---|
395 | fHiGainSecondDeriv[range-1] = 0.;
|
---|
396 |
|
---|
397 | for (Int_t k=range-2;k>=0;k--)
|
---|
398 | fHiGainSecondDeriv[k] = fHiGainSecondDeriv[k]*fHiGainSecondDeriv[k+1] + fHiGainFirstDeriv[k];
|
---|
399 | for (Int_t k=range-2;k>=0;k--)
|
---|
400 | fHiGainSecondDeriv[k] /= 6.;
|
---|
401 |
|
---|
402 | if (IsNoiseCalculation())
|
---|
403 | {
|
---|
404 |
|
---|
405 | if (fRandomIter == int(1./fResolution))
|
---|
406 | fRandomIter = 0;
|
---|
407 |
|
---|
408 | const Float_t nsx = fRandomIter * fResolution;
|
---|
409 |
|
---|
410 | if (fExtractionType == kAmplitude)
|
---|
411 | {
|
---|
412 | const Float_t b = nsx;
|
---|
413 | const Float_t a = 1. - nsx;
|
---|
414 |
|
---|
415 | sum = a*fHiGainSignal[1]
|
---|
416 | + b*fHiGainSignal[2]
|
---|
417 | + (a*a*a-a)*fHiGainSecondDeriv[1]
|
---|
418 | + (b*b*b-b)*fHiGainSecondDeriv[2];
|
---|
419 | }
|
---|
420 | else
|
---|
421 | {
|
---|
422 | Float_t start = 2. + nsx;
|
---|
423 | Float_t last = start + fRiseTimeHiGain + fFallTimeHiGain;
|
---|
424 |
|
---|
425 | if (int(last) > range)
|
---|
426 | {
|
---|
427 | const Int_t diff = range - int(last);
|
---|
428 | last -= diff;
|
---|
429 | start -= diff;
|
---|
430 | }
|
---|
431 |
|
---|
432 | CalcIntegralHiGain(sum, start, last);
|
---|
433 | }
|
---|
434 | fRandomIter++;
|
---|
435 | return;
|
---|
436 | }
|
---|
437 |
|
---|
438 | //
|
---|
439 | // Allow no saturated slice
|
---|
440 | // and
|
---|
441 | // Don't start if the maxpos is too close to the limits.
|
---|
442 | //
|
---|
443 | if (sat || maxpos < TMath::Ceil(fRiseTimeHiGain) || maxpos > range-2)
|
---|
444 | {
|
---|
445 |
|
---|
446 | dtime = 1.0;
|
---|
447 | if (fExtractionType == kAmplitude)
|
---|
448 | {
|
---|
449 | sum = fAbMax;
|
---|
450 | time = (Float_t)(fHiGainFirst + maxpos);
|
---|
451 | return;
|
---|
452 | }
|
---|
453 |
|
---|
454 | if (maxpos > range - 2)
|
---|
455 | CalcIntegralHiGain(sum, (Float_t)range - fRiseTimeHiGain - fFallTimeHiGain, (Float_t)range - 0.001);
|
---|
456 | else
|
---|
457 | CalcIntegralHiGain(sum, 0.001, fRiseTimeHiGain + fFallTimeHiGain);
|
---|
458 |
|
---|
459 | time = (Float_t)(fHiGainFirst + maxpos - 1);
|
---|
460 | return;
|
---|
461 | }
|
---|
462 |
|
---|
463 | dtime = fResolution;
|
---|
464 |
|
---|
465 | //
|
---|
466 | // Now find the maximum
|
---|
467 | //
|
---|
468 | Float_t step = 0.2; // start with step size of 1ns and loop again with the smaller one
|
---|
469 | Float_t lower = -1. + maxpos;
|
---|
470 | Float_t upper = (Float_t)maxpos;
|
---|
471 | fAbMaxPos = upper;
|
---|
472 | Float_t x = lower;
|
---|
473 | Float_t y = 0.;
|
---|
474 | Float_t a = 1.;
|
---|
475 | Float_t b = 0.;
|
---|
476 | Int_t klo = maxpos-1;
|
---|
477 | Int_t khi = maxpos;
|
---|
478 |
|
---|
479 | //
|
---|
480 | // Search for the maximum, starting in interval maxpos-1 in steps of 0.2 till maxpos-0.2.
|
---|
481 | // If no maximum is found, go to interval maxpos+1.
|
---|
482 | //
|
---|
483 | while ( x < upper - 0.3 )
|
---|
484 | {
|
---|
485 |
|
---|
486 | x += step;
|
---|
487 | a -= step;
|
---|
488 | b += step;
|
---|
489 |
|
---|
490 | y = a*fHiGainSignal[klo]
|
---|
491 | + b*fHiGainSignal[khi]
|
---|
492 | + (a*a*a-a)*fHiGainSecondDeriv[klo]
|
---|
493 | + (b*b*b-b)*fHiGainSecondDeriv[khi];
|
---|
494 |
|
---|
495 | if (y > fAbMax)
|
---|
496 | {
|
---|
497 | fAbMax = y;
|
---|
498 | fAbMaxPos = x;
|
---|
499 | }
|
---|
500 |
|
---|
501 | }
|
---|
502 |
|
---|
503 | //
|
---|
504 | // Search for the absolute maximum from maxpos to maxpos+1 in steps of 0.2
|
---|
505 | //
|
---|
506 | if (fAbMaxPos > upper-0.1)
|
---|
507 | {
|
---|
508 |
|
---|
509 | upper = 1. + maxpos;
|
---|
510 | lower = (Float_t)maxpos;
|
---|
511 | x = lower;
|
---|
512 | a = 1.;
|
---|
513 | b = 0.;
|
---|
514 | khi = maxpos+1;
|
---|
515 | klo = maxpos;
|
---|
516 |
|
---|
517 | while (x<upper-0.3)
|
---|
518 | {
|
---|
519 |
|
---|
520 | x += step;
|
---|
521 | a -= step;
|
---|
522 | b += step;
|
---|
523 |
|
---|
524 | y = a*fHiGainSignal[klo]
|
---|
525 | + b*fHiGainSignal[khi]
|
---|
526 | + (a*a*a-a)*fHiGainSecondDeriv[klo]
|
---|
527 | + (b*b*b-b)*fHiGainSecondDeriv[khi];
|
---|
528 |
|
---|
529 | if (y > fAbMax)
|
---|
530 | {
|
---|
531 | fAbMax = y;
|
---|
532 | fAbMaxPos = x;
|
---|
533 | }
|
---|
534 | }
|
---|
535 | }
|
---|
536 | //
|
---|
537 | // Now, the time, abmax and khicont and klocont are set correctly within the previous precision.
|
---|
538 | // Try a better precision.
|
---|
539 | //
|
---|
540 | const Float_t up = fAbMaxPos+step - 3.0*fResolution;
|
---|
541 | const Float_t lo = fAbMaxPos-step + 3.0*fResolution;
|
---|
542 | const Float_t maxpossave = fAbMaxPos;
|
---|
543 |
|
---|
544 | x = fAbMaxPos;
|
---|
545 | a = upper - x;
|
---|
546 | b = x - lower;
|
---|
547 |
|
---|
548 | step = 2.*fResolution; // step size of 0.1 FADC slices
|
---|
549 |
|
---|
550 | while (x<up)
|
---|
551 | {
|
---|
552 |
|
---|
553 | x += step;
|
---|
554 | a -= step;
|
---|
555 | b += step;
|
---|
556 |
|
---|
557 | y = a*fHiGainSignal[klo]
|
---|
558 | + b*fHiGainSignal[khi]
|
---|
559 | + (a*a*a-a)*fHiGainSecondDeriv[klo]
|
---|
560 | + (b*b*b-b)*fHiGainSecondDeriv[khi];
|
---|
561 |
|
---|
562 | if (y > fAbMax)
|
---|
563 | {
|
---|
564 | fAbMax = y;
|
---|
565 | fAbMaxPos = x;
|
---|
566 | }
|
---|
567 | }
|
---|
568 |
|
---|
569 | //
|
---|
570 | // Second, try from time down to time-0.2 in steps of fResolution.
|
---|
571 | //
|
---|
572 | x = maxpossave;
|
---|
573 |
|
---|
574 | //
|
---|
575 | // Test the possibility that the absolute maximum has not been found between
|
---|
576 | // maxpos and maxpos+0.05, then we have to look between maxpos-0.05 and maxpos
|
---|
577 | // which requires new setting of klocont and khicont
|
---|
578 | //
|
---|
579 | if (x < lower + fResolution)
|
---|
580 | {
|
---|
581 | klo--;
|
---|
582 | khi--;
|
---|
583 | upper -= 1.;
|
---|
584 | lower -= 1.;
|
---|
585 | }
|
---|
586 |
|
---|
587 | a = upper - x;
|
---|
588 | b = x - lower;
|
---|
589 |
|
---|
590 | while (x>lo)
|
---|
591 | {
|
---|
592 |
|
---|
593 | x -= step;
|
---|
594 | a += step;
|
---|
595 | b -= step;
|
---|
596 |
|
---|
597 | y = a*fHiGainSignal[klo]
|
---|
598 | + b*fHiGainSignal[khi]
|
---|
599 | + (a*a*a-a)*fHiGainSecondDeriv[klo]
|
---|
600 | + (b*b*b-b)*fHiGainSecondDeriv[khi];
|
---|
601 |
|
---|
602 | if (y > fAbMax)
|
---|
603 | {
|
---|
604 | fAbMax = y;
|
---|
605 | fAbMaxPos = x;
|
---|
606 | }
|
---|
607 | }
|
---|
608 |
|
---|
609 | if (fExtractionType == kAmplitude)
|
---|
610 | {
|
---|
611 | time = fAbMaxPos + (Int_t)fHiGainFirst;
|
---|
612 | sum = fAbMax;
|
---|
613 | return;
|
---|
614 | }
|
---|
615 |
|
---|
616 | fHalfMax = fAbMax/2.;
|
---|
617 |
|
---|
618 | //
|
---|
619 | // Now, loop from the maximum bin leftward down in order to find the position of the half maximum.
|
---|
620 | // First, find the right FADC slice:
|
---|
621 | //
|
---|
622 | klo = maxpos;
|
---|
623 | while (klo > 0)
|
---|
624 | {
|
---|
625 | klo--;
|
---|
626 | if (fHiGainSignal[klo] < fHalfMax)
|
---|
627 | break;
|
---|
628 | }
|
---|
629 |
|
---|
630 | khi = klo+1;
|
---|
631 | //
|
---|
632 | // Loop from the beginning of the slice upwards to reach the fHalfMax:
|
---|
633 | // With means of bisection:
|
---|
634 | //
|
---|
635 | x = (Float_t)klo;
|
---|
636 | a = 1.;
|
---|
637 | b = 0.;
|
---|
638 |
|
---|
639 | step = 0.5;
|
---|
640 | Bool_t back = kFALSE;
|
---|
641 |
|
---|
642 | Int_t maxcnt = 20;
|
---|
643 | Int_t cnt = 0;
|
---|
644 |
|
---|
645 | while (TMath::Abs(y-fHalfMax) > fResolution)
|
---|
646 | {
|
---|
647 |
|
---|
648 | if (back)
|
---|
649 | {
|
---|
650 | x -= step;
|
---|
651 | a += step;
|
---|
652 | b -= step;
|
---|
653 | }
|
---|
654 | else
|
---|
655 | {
|
---|
656 | x += step;
|
---|
657 | a -= step;
|
---|
658 | b += step;
|
---|
659 | }
|
---|
660 |
|
---|
661 | y = a*fHiGainSignal[klo]
|
---|
662 | + b*fHiGainSignal[khi]
|
---|
663 | + (a*a*a-a)*fHiGainSecondDeriv[klo]
|
---|
664 | + (b*b*b-b)*fHiGainSecondDeriv[khi];
|
---|
665 |
|
---|
666 | if (y > fHalfMax)
|
---|
667 | back = kTRUE;
|
---|
668 | else
|
---|
669 | back = kFALSE;
|
---|
670 |
|
---|
671 | if (++cnt > maxcnt)
|
---|
672 | break;
|
---|
673 |
|
---|
674 | step /= 2.;
|
---|
675 | }
|
---|
676 |
|
---|
677 | time = (Float_t)fHiGainFirst + x;
|
---|
678 | //
|
---|
679 | // Now integrate the whole thing!
|
---|
680 | //
|
---|
681 |
|
---|
682 | Float_t start = fAbMaxPos - fRiseTimeHiGain;
|
---|
683 | Float_t last = fAbMaxPos + fFallTimeHiGain;
|
---|
684 |
|
---|
685 | const Int_t diff = int(last) - range;
|
---|
686 |
|
---|
687 | if (diff > 0)
|
---|
688 | {
|
---|
689 | last -= diff;
|
---|
690 | start -= diff;
|
---|
691 | }
|
---|
692 |
|
---|
693 | CalcIntegralHiGain(sum, start, last);
|
---|
694 | }
|
---|
695 |
|
---|
696 |
|
---|
697 | // --------------------------------------------------------------------------
|
---|
698 | //
|
---|
699 | // Calculates the arrival time and charge for each pixel
|
---|
700 | //
|
---|
701 | void MExtractTimeAndChargeSpline::FindTimeAndChargeLoGain(Byte_t *first, Float_t &sum, Float_t &dsum,
|
---|
702 | Float_t &time, Float_t &dtime,
|
---|
703 | Byte_t &sat, const MPedestalPix &ped, const Bool_t abflag)
|
---|
704 | {
|
---|
705 |
|
---|
706 | Int_t range = fLoGainLast - fLoGainFirst + 1;
|
---|
707 | const Byte_t *end = first + range;
|
---|
708 | Byte_t *p = first;
|
---|
709 |
|
---|
710 | const Float_t pedes = ped.GetPedestal();
|
---|
711 | const Float_t ABoffs = ped.GetPedestalABoffset();
|
---|
712 |
|
---|
713 | const Float_t pedmean[2] = { pedes + ABoffs, pedes - ABoffs };
|
---|
714 |
|
---|
715 | fAbMax = 0.;
|
---|
716 | fAbMaxPos = 0.;
|
---|
717 | Int_t maxpos = 0;
|
---|
718 | Int_t max = 0;
|
---|
719 |
|
---|
720 | //
|
---|
721 | // Check for saturation in all other slices
|
---|
722 | //
|
---|
723 | Int_t ids = fLoGainFirst;
|
---|
724 | Float_t *sample = fLoGainSignal.GetArray();
|
---|
725 | while (p<end)
|
---|
726 | {
|
---|
727 |
|
---|
728 | *sample++ = (Float_t)*p - pedmean[(ids++ + abflag) & 0x1];
|
---|
729 |
|
---|
730 | if (*p > max)
|
---|
731 | {
|
---|
732 | maxpos = ids-fLoGainFirst-1;
|
---|
733 | max = *p;
|
---|
734 | }
|
---|
735 |
|
---|
736 | if (*p++ >= fSaturationLimit)
|
---|
737 | sat++;
|
---|
738 | }
|
---|
739 |
|
---|
740 | fAbMax = fLoGainSignal[maxpos];
|
---|
741 |
|
---|
742 | Float_t pp;
|
---|
743 |
|
---|
744 | fLoGainSecondDeriv[0] = 0.;
|
---|
745 | fLoGainFirstDeriv[0] = 0.;
|
---|
746 |
|
---|
747 | for (Int_t i=1;i<range-1;i++)
|
---|
748 | {
|
---|
749 | pp = fLoGainSecondDeriv[i-1] + 4.;
|
---|
750 | fLoGainSecondDeriv[i] = -1.0/pp;
|
---|
751 | fLoGainFirstDeriv [i] = fLoGainSignal[i+1] - fLoGainSignal[i] - fLoGainSignal[i] + fLoGainSignal[i-1];
|
---|
752 | fLoGainFirstDeriv [i] = (6.0*fLoGainFirstDeriv[i]-fLoGainFirstDeriv[i-1])/pp;
|
---|
753 | }
|
---|
754 |
|
---|
755 | fLoGainSecondDeriv[range-1] = 0.;
|
---|
756 |
|
---|
757 | for (Int_t k=range-2;k>=0;k--)
|
---|
758 | fLoGainSecondDeriv[k] = fLoGainSecondDeriv[k]*fLoGainSecondDeriv[k+1] + fLoGainFirstDeriv[k];
|
---|
759 | for (Int_t k=range-2;k>=0;k--)
|
---|
760 | fLoGainSecondDeriv[k] /= 6.;
|
---|
761 |
|
---|
762 | if (IsNoiseCalculation())
|
---|
763 | {
|
---|
764 | if (fRandomIter == int(1./fResolution))
|
---|
765 | fRandomIter = 0;
|
---|
766 |
|
---|
767 | const Float_t nsx = fRandomIter * fResolution;
|
---|
768 |
|
---|
769 | if (fExtractionType == kAmplitude)
|
---|
770 | {
|
---|
771 | const Float_t b = nsx;
|
---|
772 | const Float_t a = 1. - nsx;
|
---|
773 |
|
---|
774 | sum = a*fLoGainSignal[1]
|
---|
775 | + b*fLoGainSignal[2]
|
---|
776 | + (a*a*a-a)*fLoGainSecondDeriv[1]
|
---|
777 | + (b*b*b-b)*fLoGainSecondDeriv[2];
|
---|
778 | }
|
---|
779 | else
|
---|
780 | {
|
---|
781 | Float_t start = 2. + nsx;
|
---|
782 | Float_t last = start + fRiseTimeLoGain + fFallTimeLoGain;
|
---|
783 |
|
---|
784 | if (int(last) > range)
|
---|
785 | {
|
---|
786 | const Int_t diff = range - int(last);
|
---|
787 | last -= diff;
|
---|
788 | start -= diff;
|
---|
789 | }
|
---|
790 |
|
---|
791 | CalcIntegralLoGain(sum, start, last);
|
---|
792 | }
|
---|
793 | fRandomIter++;
|
---|
794 | return;
|
---|
795 | }
|
---|
796 | //
|
---|
797 | // Allow no saturated slice
|
---|
798 | // and
|
---|
799 | // Don't start if the maxpos is too close to the limits.
|
---|
800 | //
|
---|
801 | if (sat || maxpos < TMath::Ceil(fRiseTimeLoGain) || maxpos > range-2)
|
---|
802 | {
|
---|
803 |
|
---|
804 | dtime = 1.0;
|
---|
805 | if (fExtractionType == kAmplitude)
|
---|
806 | {
|
---|
807 | time = (Float_t)(fLoGainFirst + maxpos);
|
---|
808 | sum = fAbMax;
|
---|
809 | return;
|
---|
810 | }
|
---|
811 |
|
---|
812 | if (maxpos > range-2)
|
---|
813 | CalcIntegralLoGain(sum, (Float_t)range - fRiseTimeLoGain - fFallTimeLoGain -1., (Float_t)range - 0.001);
|
---|
814 | else
|
---|
815 | CalcIntegralLoGain(sum, 0.001, fRiseTimeLoGain + fFallTimeLoGain + 1.);
|
---|
816 |
|
---|
817 | time = (Float_t)(fLoGainFirst + maxpos - 1);
|
---|
818 | return;
|
---|
819 | }
|
---|
820 |
|
---|
821 | dtime = fResolution;
|
---|
822 |
|
---|
823 | //
|
---|
824 | // Now find the maximum
|
---|
825 | //
|
---|
826 | Float_t step = 0.2; // start with step size of 1ns and loop again with the smaller one
|
---|
827 | Float_t lower = -1. + maxpos;
|
---|
828 | Float_t upper = (Float_t)maxpos;
|
---|
829 | fAbMaxPos = upper;
|
---|
830 | Float_t x = lower;
|
---|
831 | Float_t y = 0.;
|
---|
832 | Float_t a = 1.;
|
---|
833 | Float_t b = 0.;
|
---|
834 | Int_t klo = maxpos-1;
|
---|
835 | Int_t khi = maxpos;
|
---|
836 |
|
---|
837 | //
|
---|
838 | // Search for the maximum, starting in interval maxpos-1 in steps of 0.2 till maxpos-0.2.
|
---|
839 | // If no maximum is found, go to interval maxpos+1.
|
---|
840 | //
|
---|
841 | while ( x < upper - 0.3 )
|
---|
842 | {
|
---|
843 |
|
---|
844 | x += step;
|
---|
845 | a -= step;
|
---|
846 | b += step;
|
---|
847 |
|
---|
848 | y = a*fLoGainSignal[klo]
|
---|
849 | + b*fLoGainSignal[khi]
|
---|
850 | + (a*a*a-a)*fLoGainSecondDeriv[klo]
|
---|
851 | + (b*b*b-b)*fLoGainSecondDeriv[khi];
|
---|
852 |
|
---|
853 | if (y > fAbMax)
|
---|
854 | {
|
---|
855 | fAbMax = y;
|
---|
856 | fAbMaxPos = x;
|
---|
857 | }
|
---|
858 |
|
---|
859 | }
|
---|
860 |
|
---|
861 | //
|
---|
862 | // Test the possibility that the absolute maximum has not been found before the
|
---|
863 | // maxpos and search from maxpos to maxpos+1 in steps of 0.2
|
---|
864 | //
|
---|
865 | if (fAbMaxPos > upper-0.1)
|
---|
866 | {
|
---|
867 |
|
---|
868 | upper = 1. + maxpos;
|
---|
869 | lower = (Float_t)maxpos;
|
---|
870 | x = lower;
|
---|
871 | a = 1.;
|
---|
872 | b = 0.;
|
---|
873 | khi = maxpos+1;
|
---|
874 | klo = maxpos;
|
---|
875 |
|
---|
876 | while (x<upper-0.3)
|
---|
877 | {
|
---|
878 |
|
---|
879 | x += step;
|
---|
880 | a -= step;
|
---|
881 | b += step;
|
---|
882 |
|
---|
883 | y = a*fLoGainSignal[klo]
|
---|
884 | + b*fLoGainSignal[khi]
|
---|
885 | + (a*a*a-a)*fLoGainSecondDeriv[klo]
|
---|
886 | + (b*b*b-b)*fLoGainSecondDeriv[khi];
|
---|
887 |
|
---|
888 | if (y > fAbMax)
|
---|
889 | {
|
---|
890 | fAbMax = y;
|
---|
891 | fAbMaxPos = x;
|
---|
892 | }
|
---|
893 | }
|
---|
894 | }
|
---|
895 |
|
---|
896 |
|
---|
897 | //
|
---|
898 | // Now, the time, abmax and khicont and klocont are set correctly within the previous precision.
|
---|
899 | // Try a better precision.
|
---|
900 | //
|
---|
901 | const Float_t up = fAbMaxPos+step - 3.0*fResolution;
|
---|
902 | const Float_t lo = fAbMaxPos-step + 3.0*fResolution;
|
---|
903 | const Float_t maxpossave = fAbMaxPos;
|
---|
904 |
|
---|
905 | x = fAbMaxPos;
|
---|
906 | a = upper - x;
|
---|
907 | b = x - lower;
|
---|
908 |
|
---|
909 | step = 2.*fResolution; // step size of 0.1 FADC slice
|
---|
910 |
|
---|
911 | while (x<up)
|
---|
912 | {
|
---|
913 |
|
---|
914 | x += step;
|
---|
915 | a -= step;
|
---|
916 | b += step;
|
---|
917 |
|
---|
918 | y = a*fLoGainSignal[klo]
|
---|
919 | + b*fLoGainSignal[khi]
|
---|
920 | + (a*a*a-a)*fLoGainSecondDeriv[klo]
|
---|
921 | + (b*b*b-b)*fLoGainSecondDeriv[khi];
|
---|
922 |
|
---|
923 | if (y > fAbMax)
|
---|
924 | {
|
---|
925 | fAbMax = y;
|
---|
926 | fAbMaxPos = x;
|
---|
927 | }
|
---|
928 | }
|
---|
929 |
|
---|
930 | //
|
---|
931 | // Second, try from time down to time-0.2 in steps of 0.025.
|
---|
932 | //
|
---|
933 | x = maxpossave;
|
---|
934 |
|
---|
935 | //
|
---|
936 | // Test the possibility that the absolute maximum has not been found between
|
---|
937 | // maxpos and maxpos+0.05, then we have to look between maxpos-0.05 and maxpos
|
---|
938 | // which requires new setting of klocont and khicont
|
---|
939 | //
|
---|
940 | if (x < lower + fResolution)
|
---|
941 | {
|
---|
942 | klo--;
|
---|
943 | khi--;
|
---|
944 | upper -= 1.;
|
---|
945 | lower -= 1.;
|
---|
946 | }
|
---|
947 |
|
---|
948 | a = upper - x;
|
---|
949 | b = x - lower;
|
---|
950 |
|
---|
951 | while (x>lo)
|
---|
952 | {
|
---|
953 |
|
---|
954 | x -= step;
|
---|
955 | a += step;
|
---|
956 | b -= step;
|
---|
957 |
|
---|
958 | y = a*fLoGainSignal[klo]
|
---|
959 | + b*fLoGainSignal[khi]
|
---|
960 | + (a*a*a-a)*fLoGainSecondDeriv[klo]
|
---|
961 | + (b*b*b-b)*fLoGainSecondDeriv[khi];
|
---|
962 |
|
---|
963 | if (y > fAbMax)
|
---|
964 | {
|
---|
965 | fAbMax = y;
|
---|
966 | fAbMaxPos = x;
|
---|
967 | }
|
---|
968 | }
|
---|
969 |
|
---|
970 | if (fExtractionType == kAmplitude)
|
---|
971 | {
|
---|
972 | time = fAbMaxPos + (Int_t)fLoGainFirst;
|
---|
973 | sum = fAbMax;
|
---|
974 | return;
|
---|
975 | }
|
---|
976 |
|
---|
977 | fHalfMax = fAbMax/2.;
|
---|
978 |
|
---|
979 | //
|
---|
980 | // Now, loop from the maximum bin leftward down in order to find the position of the half maximum.
|
---|
981 | // First, find the right FADC slice:
|
---|
982 | //
|
---|
983 | klo = maxpos;
|
---|
984 | while (klo > 0)
|
---|
985 | {
|
---|
986 | klo--;
|
---|
987 | if (fLoGainSignal[klo] < fHalfMax)
|
---|
988 | break;
|
---|
989 | }
|
---|
990 |
|
---|
991 | khi = klo+1;
|
---|
992 | //
|
---|
993 | // Loop from the beginning of the slice upwards to reach the fHalfMax:
|
---|
994 | // With means of bisection:
|
---|
995 | //
|
---|
996 | x = (Float_t)klo;
|
---|
997 | a = 1.;
|
---|
998 | b = 0.;
|
---|
999 |
|
---|
1000 | step = 0.5;
|
---|
1001 | Bool_t back = kFALSE;
|
---|
1002 |
|
---|
1003 | Int_t maxcnt = 20;
|
---|
1004 | Int_t cnt = 0;
|
---|
1005 |
|
---|
1006 | while (TMath::Abs(y-fHalfMax) > fResolution)
|
---|
1007 | {
|
---|
1008 |
|
---|
1009 | if (back)
|
---|
1010 | {
|
---|
1011 | x -= step;
|
---|
1012 | a += step;
|
---|
1013 | b -= step;
|
---|
1014 | }
|
---|
1015 | else
|
---|
1016 | {
|
---|
1017 | x += step;
|
---|
1018 | a -= step;
|
---|
1019 | b += step;
|
---|
1020 | }
|
---|
1021 |
|
---|
1022 | y = a*fLoGainSignal[klo]
|
---|
1023 | + b*fLoGainSignal[khi]
|
---|
1024 | + (a*a*a-a)*fLoGainSecondDeriv[klo]
|
---|
1025 | + (b*b*b-b)*fLoGainSecondDeriv[khi];
|
---|
1026 |
|
---|
1027 | if (y > fHalfMax)
|
---|
1028 | back = kTRUE;
|
---|
1029 | else
|
---|
1030 | back = kFALSE;
|
---|
1031 |
|
---|
1032 | if (++cnt > maxcnt)
|
---|
1033 | break;
|
---|
1034 |
|
---|
1035 | step /= 2.;
|
---|
1036 | }
|
---|
1037 |
|
---|
1038 | time = x + (Int_t)fLoGainFirst;
|
---|
1039 |
|
---|
1040 | //
|
---|
1041 | // Now integrate the whole thing!
|
---|
1042 | //
|
---|
1043 | Float_t start = fAbMaxPos - fRiseTimeLoGain;
|
---|
1044 | Float_t last = fAbMaxPos + fFallTimeLoGain;
|
---|
1045 |
|
---|
1046 | const Int_t diff = int(last) - range;
|
---|
1047 |
|
---|
1048 | if (diff > 0)
|
---|
1049 | {
|
---|
1050 | last -= diff;
|
---|
1051 | start -= diff;
|
---|
1052 | }
|
---|
1053 | CalcIntegralLoGain(sum, start, last);
|
---|
1054 | }
|
---|
1055 |
|
---|
1056 | void MExtractTimeAndChargeSpline::CalcIntegralHiGain(Float_t &sum, Float_t start, Float_t last)
|
---|
1057 | {
|
---|
1058 |
|
---|
1059 | const Float_t step = 0.2;
|
---|
1060 |
|
---|
1061 | if (start < 0)
|
---|
1062 | {
|
---|
1063 | last -= start;
|
---|
1064 | start = 0.;
|
---|
1065 | }
|
---|
1066 |
|
---|
1067 | Int_t klo = int(start);
|
---|
1068 | Int_t khi = klo+1;
|
---|
1069 |
|
---|
1070 | Float_t lo = TMath::Floor(start);
|
---|
1071 | Float_t up = lo + 1.;
|
---|
1072 |
|
---|
1073 | const Int_t m = int((start-klo)/step);
|
---|
1074 | start = step*m + klo; // Correct start for the digitization due to resolution
|
---|
1075 |
|
---|
1076 | Float_t x = start;
|
---|
1077 | Float_t a = up-start;
|
---|
1078 | Float_t b = start-lo;
|
---|
1079 |
|
---|
1080 | while (1)
|
---|
1081 | {
|
---|
1082 |
|
---|
1083 | while (x<up)
|
---|
1084 | {
|
---|
1085 | x += step;
|
---|
1086 |
|
---|
1087 | if (x > last)
|
---|
1088 | {
|
---|
1089 | sum *= step;
|
---|
1090 | return;
|
---|
1091 | }
|
---|
1092 |
|
---|
1093 | a -= step;
|
---|
1094 | b += step;
|
---|
1095 |
|
---|
1096 | sum += a*fHiGainSignal[klo]
|
---|
1097 | + b*fHiGainSignal[khi]
|
---|
1098 | + (a*a*a-a)*fHiGainSecondDeriv[klo]
|
---|
1099 | + (b*b*b-b)*fHiGainSecondDeriv[khi];
|
---|
1100 | }
|
---|
1101 |
|
---|
1102 | up += 1.;
|
---|
1103 | lo += 1.;
|
---|
1104 | klo++;
|
---|
1105 | khi++;
|
---|
1106 | start += 1.;
|
---|
1107 | a = 1.;
|
---|
1108 | b = 0.;
|
---|
1109 | }
|
---|
1110 |
|
---|
1111 | }
|
---|
1112 | void MExtractTimeAndChargeSpline::CalcIntegralLoGain(Float_t &sum, Float_t start, Float_t last)
|
---|
1113 | {
|
---|
1114 |
|
---|
1115 | const Float_t step = 0.1;
|
---|
1116 |
|
---|
1117 | if (start < 0)
|
---|
1118 | {
|
---|
1119 | last -= start;
|
---|
1120 | start = 0.;
|
---|
1121 | }
|
---|
1122 |
|
---|
1123 | Int_t klo = int(start);
|
---|
1124 | Int_t khi = klo+1;
|
---|
1125 |
|
---|
1126 | Float_t lo = TMath::Floor(start);
|
---|
1127 | Float_t up = lo + 1.;
|
---|
1128 |
|
---|
1129 | const Int_t m = int((start-klo)/step);
|
---|
1130 | start = step*m + klo; // Correct start for the digitization due to resolution
|
---|
1131 |
|
---|
1132 | Float_t x = start;
|
---|
1133 | Float_t a = up-start;
|
---|
1134 | Float_t b = start-lo;
|
---|
1135 |
|
---|
1136 | while (1)
|
---|
1137 | {
|
---|
1138 |
|
---|
1139 | while (x<up)
|
---|
1140 | {
|
---|
1141 | x += step;
|
---|
1142 |
|
---|
1143 | if (x > last)
|
---|
1144 | {
|
---|
1145 | sum *= step;
|
---|
1146 | return;
|
---|
1147 | }
|
---|
1148 |
|
---|
1149 | a -= step;
|
---|
1150 | b += step;
|
---|
1151 |
|
---|
1152 | sum += a*fLoGainSignal[klo]
|
---|
1153 | + b*fLoGainSignal[khi]
|
---|
1154 | + (a*a*a-a)*fLoGainSecondDeriv[klo]
|
---|
1155 | + (b*b*b-b)*fLoGainSecondDeriv[khi];
|
---|
1156 |
|
---|
1157 | }
|
---|
1158 |
|
---|
1159 | up += 1.;
|
---|
1160 | lo += 1.;
|
---|
1161 | klo++;
|
---|
1162 | khi++;
|
---|
1163 | start += 1.;
|
---|
1164 | a = 1.;
|
---|
1165 | b = 0.;
|
---|
1166 | }
|
---|
1167 |
|
---|
1168 | }
|
---|
1169 |
|
---|
1170 |
|
---|
1171 |
|
---|
1172 |
|
---|
1173 | // --------------------------------------------------------------------------
|
---|
1174 | //
|
---|
1175 | // In addition to the resources of the base-class MExtractor:
|
---|
1176 | // Resolution
|
---|
1177 | // RiseTimeHiGain
|
---|
1178 | // FallTimeHiGain
|
---|
1179 | // LoGainStretch
|
---|
1180 | // ExtractionType: amplitude, integral
|
---|
1181 | //
|
---|
1182 | Int_t MExtractTimeAndChargeSpline::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
|
---|
1183 | {
|
---|
1184 |
|
---|
1185 | Bool_t rc = kFALSE;
|
---|
1186 |
|
---|
1187 | if (IsEnvDefined(env, prefix, "Resolution", print))
|
---|
1188 | {
|
---|
1189 | SetResolution(GetEnvValue(env, prefix, "Resolution",fResolution));
|
---|
1190 | rc = kTRUE;
|
---|
1191 | }
|
---|
1192 | if (IsEnvDefined(env, prefix, "RiseTimeHiGain", print))
|
---|
1193 | {
|
---|
1194 | SetRiseTimeHiGain(GetEnvValue(env, prefix, "RiseTimeHiGain", fRiseTimeHiGain));
|
---|
1195 | rc = kTRUE;
|
---|
1196 | }
|
---|
1197 | if (IsEnvDefined(env, prefix, "FallTimeHiGain", print))
|
---|
1198 | {
|
---|
1199 | SetFallTimeHiGain(GetEnvValue(env, prefix, "FallTimeHiGain", fFallTimeHiGain));
|
---|
1200 | rc = kTRUE;
|
---|
1201 | }
|
---|
1202 | if (IsEnvDefined(env, prefix, "LoGainStretch", print))
|
---|
1203 | {
|
---|
1204 | SetLoGainStretch(GetEnvValue(env, prefix, "LoGainStretch", fLoGainStretch));
|
---|
1205 | rc = kTRUE;
|
---|
1206 | }
|
---|
1207 |
|
---|
1208 | if (IsEnvDefined(env, prefix, "ExtractionType", print))
|
---|
1209 | {
|
---|
1210 | TString type = GetEnvValue(env, prefix, "ExtractionType", "");
|
---|
1211 | type.ToLower();
|
---|
1212 | type = type.Strip(TString::kBoth);
|
---|
1213 | if (type==(TString)"amplitude")
|
---|
1214 | SetChargeType(kAmplitude);
|
---|
1215 | if (type==(TString)"integral")
|
---|
1216 | SetChargeType(kIntegral);
|
---|
1217 | rc=kTRUE;
|
---|
1218 | }
|
---|
1219 |
|
---|
1220 | return MExtractTimeAndCharge::ReadEnv(env, prefix, print) ? kTRUE : rc;
|
---|
1221 |
|
---|
1222 | }
|
---|
1223 |
|
---|
1224 |
|
---|