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