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