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