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