GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TPulseAnalyzer.h
Go to the documentation of this file.
1#ifndef TPULSE_ANALYZER_H
2#define TPULSE_ANALYZER_H
3
4/** \addtogroup Fitting Fitting & Analysis
5 * @{
6 */
7
8#include "TFragment.h"
9#include "TGRSIFunctions.h"
10
11#include <vector>
12#include <cstdlib>
13#include <cstdio>
14#include <cstring>
15#include <cmath>
16
17#include "TNamed.h"
18#include "Rtypes.h"
19#include "TH1.h"
20#include "TF1.h"
21#include "TMath.h"
22#include "TGraph.h"
23
24////////////////////////////////////////////////////////////////////////////////
25/// \class TPulseAnalyzer
26///
27/// Mostly a direct port of SFU code
28/// I have stripped out some surplus and encapsulated it, but I havent changed much
29/// I'm sure there is more that can be stripped and reformated but it is working currently.
30///
31////////////////////////////////////////////////////////////////////////////////
32
34private:
35 struct WaveFormPar {
36 // parameters for baseline
38 double baseline; // baseline
39 double baselineStDev; // baseline variance
40 double baselinefin; // baseline
41 double baselineStDevfin; // baseline variance
42 bool bflag; // flag for baseline determiniation
43 // paramaters for oscillation
44 double baseamp; // amplitude of baseline oscillations
45 double basefreq; // tick frequency of baseline oscillations
46 double basephase; // flag for baseline determination
47 int osciflag; // flag for baseline determination
48 // paremeters for exclusion zone
49 double max; // max value of the filtered waveform
50 double tmax; // x position of the max value of the filtered waveform
51 double baselineMin; // max crossing point for exclusion zone
52 double baselineMax; // max crossing point for exclusion zone
53 int temax; // x position, exclusion zone upper limit
54 int temin; // x position, exclusion zone lower limit
55 double afit, bfit; // parameters for the line which fits risetime above temax
56 int mflag; // flag for tmax found
57 int teflag; // flag for exclusion zone determined
58 double t0; // required for compilation of map.c - check that it works
59 double t0_local;
60 // new stuff necessary for compiliation of Kris' waveform analyzer changes
61 double b0;
62 double b1;
63 long double s0;
64 long double s1;
65 long double s2;
66 double t90;
67 double t50;
68 double t30;
69 double t10;
70 double t10t90;
76 int thigh;
77 double sig2noise;
78 double amplitude; // amplitude from sili fits
79 double tauDecay;
80 double tauRise;
81 };
82
83 struct LinePar {
84 double slope;
85 double intercept;
86 double chisq;
87 double ndf;
88 };
89
90 struct ParPar {
91 double constant;
92 double linear;
93 double quadratic;
94 double chisq;
95 double ndf;
96 };
97
98 struct SinPar {
99 double A;
100 double t0;
101 double C;
102 };
103
104 struct ShapePar {
105 double chisq;
106 int ndf;
107 int type;
108 // decay constants for the fits
109 long double t[5]; // NOLINT(*-avoid-c-arrays)
110 // associated aplitudes for the decay constants
111 long double am[5]; // NOLINT(*-avoid-c-arrays)
112 double rf[5]; // NOLINT(*-avoid-c-arrays)
113
114 // new stuff necessary for compiliation of Kris' waveform analyzer changes
115 long double chisq_ex;
116 long double chisq_f;
118 int ndf_f;
119 };
120
121public:
123 explicit TPulseAnalyzer(const TFragment& fragment, double = 0);
124 explicit TPulseAnalyzer(const std::vector<Short_t>& wave, double = 0, std::string name = "");
126 TPulseAnalyzer(TPulseAnalyzer&&) noexcept = default;
127 TPulseAnalyzer& operator=(const TPulseAnalyzer&) = default;
128 TPulseAnalyzer& operator=(TPulseAnalyzer&&) noexcept = default;
129 virtual ~TPulseAnalyzer();
130
131 void SetData(const TFragment& fragment, double = 0);
132 void SetData(const std::vector<Short_t>& wave, double = 0);
133 void Clear(Option_t* opt = "");
134 bool IsSet() const { return set; }
135
136 inline double Get_wpar_T0() { return cWpar->t0; }
137 inline double Get_wpar_baselinefin() { return cWpar->baselinefin; }
138 inline double Get_wpar_amplitude() { return cWpar->amplitude; }
139 inline double Get_wpar_decay() { return cWpar->tauDecay; }
140 inline double Get_wpar_rise() { return cWpar->tauRise; }
141
142 bool SiliShapePrepare(double tauDecay, double tauRise);
143 bool GetSiliShape(double tauDecay, double tauRise); // Added for Spice, parameters to be found : t0 and Amplitude
144 bool GetSiliShapeTF1(double tauDecay, double tauRise, double baseline, double basefreq = 0); // Added for Spice, full slow fit to establish parameters, needs initial estimates
145 inline bool Get_bflag() { return cWpar->bflag; }
146 TF1 Getsilifit();
147 double GetsiliSmirnov();
148 void Drawsilifit();
149
150 static double SiLiFitFunction(double* i, double* p);
151
152 double fit_newT0();
153 double fit_rf(double = 2 * 8.48409);
154 double get_sig2noise();
155 int16_t good_baseline();
156 void print_WavePar();
157 TH1I* GetWaveHist();
158 TGraph* GetWaveGraph();
159 static int fNameIter;
160 void DrawWave();
161 void DrawT0fit();
162 void DrawRFFit();
163
164 // CsI functions:
165 double CsIPID();
166 double CsIt0();
167 void DrawCsIExclusion();
168 void DrawCsIFit();
169 int GetCsIChiSq();
170 int GetCsIFitType();
171
172private:
173 bool set{false};
175 int cN{0};
176 // TFragment* frag;
177 std::vector<Short_t> cWavebuffer;
178 SinPar* spar{nullptr};
179 ShapePar* shpar{nullptr};
180 WaveFormPar* csiTestWpar[4]{nullptr}; // NOLINT(*-avoid-c-arrays)
181 ShapePar* csiTestShpar[4]{nullptr}; // NOLINT(*-avoid-c-arrays)
182
183 std::string fName;
184
185 // pulse fitting parameters
186 int FILTER; // integration region for noise reduction (in samples)
187 int T0RANGE; // tick range over which baseline is calulated
189
190 // linear equation dataholders
192 long double lineq_matrix[20][20]; // NOLINT(*-avoid-c-arrays)
193 long double lineq_vector[20]; // NOLINT(*-avoid-c-arrays)
194 long double lineq_solution[20]; // NOLINT(*-avoid-c-arrays)
195 long double copy_matrix[20][20]; // NOLINT(*-avoid-c-arrays)
196
197 // CsI functions
198 void GetCsIExclusionZone();
199 double GetCsITau(int, ShapePar*);
200 double GetCsIt0(ShapePar*, WaveFormPar*);
201 int GetCsIShape();
202 int FitCsIShape(int, ShapePar*, WaveFormPar*);
203
204 bool CsISet;
205 double EPS;
206
207 void SetCsI(bool option = true) { CsISet = option; }
208 bool CsIIsSet() const { return CsISet; }
209
210 // internal methods
211 int solve_lin_eq();
212 long double determinant(int);
213
214 int fit_parabola(int, int, ParPar*);
215 int fit_smooth_parabola(int, int, double, ParPar*);
216 int fit_line(int, int, LinePar*);
217 double get_linear_T0();
218 double get_parabolic_T0();
219 double get_smooth_T0();
220
221 void get_baseline();
222 void get_baseline_fin();
223 void get_tmax();
224
225 double get_tfrac(double, double, double);
226 void get_t10();
227 void get_t30();
228 void get_t50();
229 void get_t90();
230
231 double get_sin_par(double);
232
233 void GetQuickPara();
234
235 // bad chi squares for debugging
236 const static int BADCHISQ_SMOOTH_T0 = -1024 - 2; // smooth_t0 gives bad result
237 const static int BADCHISQ_PAR_T0 = -1024 - 3; // parabolic_t0 gives bad result
238 const static int BADCHISQ_LIN_T0 = -1024 - 4; // linear_t0 gives bad result
239 const static int BADCHISQ_MAT = -1024 - 5; // matrix for fit is not invertable
240 const static int BADCHISQ_FAIL_DIRECT = -1024 - 9; // generic fit failure
241 const static int BAD_EXCLUSION_ZONE = -1024 - 10; // fails finding exclusion zone
242 const static int BAD_MAX = -1024 - 11;
243 const static int BAD_BASELINE_RANGE = -1024 - 12;
244 const static int BAD_PID = -1024; // PID value returned when the fit fails
245
246 // new definitions for Kris' changes to the waveform analyzer
247 const static int PIN_BASELINE_RANGE = 16; // minimum ticks before max for a valid signal
248 const static int MAX_SAMPLES = 4096;
249
250 const static int CSI_BASELINE_RANGE = 50; // baseline range in channels
251 const static int NOISE_LEVEL_CSI = 100; // noise level for CsI
252 const static int NSHAPE = 5; // number of trial functions for waveform fit
253
254 const static int BADCHISQ_T0 = -1024 - 7;
255 const static int BADCHISQ_NEG = -1024 - 1;
256 const static int BADCHISQ_AMPL = -1024 - 6;
257
258 /// \cond CLASSIMP
259 ClassDef(TPulseAnalyzer, 4) // NOLINT(readability-else-after-return)
260 /// \endcond
261};
262/*! @} */
263#endif
long double lineq_solution[20]
std::vector< Short_t > cWavebuffer
double fit_rf(double=2 *8.48409)
long double determinant(int)
bool CsIIsSet() const
static const int BADCHISQ_PAR_T0
static const int BADCHISQ_LIN_T0
bool GetSiliShapeTF1(double tauDecay, double tauRise, double baseline, double basefreq=0)
static const int CSI_BASELINE_RANGE
static const int BAD_EXCLUSION_ZONE
double GetCsIt0(ShapePar *, WaveFormPar *)
double get_tfrac(double, double, double)
TPulseAnalyzer(const TPulseAnalyzer &)=default
static const int NSHAPE
void SetCsI(bool option=true)
static const int BADCHISQ_FAIL_DIRECT
void Clear(Option_t *opt="")
bool IsSet() const
double GetCsITau(int, ShapePar *)
static const int BAD_BASELINE_RANGE
static double SiLiFitFunction(double *i, double *p)
double Get_wpar_decay()
long double copy_matrix[20][20]
ShapePar * csiTestShpar[4]
static const int BADCHISQ_SMOOTH_T0
long double lineq_vector[20]
bool SiliShapePrepare(double tauDecay, double tauRise)
int fit_parabola(int, int, ParPar *)
int fit_line(int, int, LinePar *)
WaveFormPar * cWpar
TPulseAnalyzer(TPulseAnalyzer &&) noexcept=default
static int fNameIter
std::string fName
static const int BADCHISQ_AMPL
int fit_smooth_parabola(int, int, double, ParPar *)
double Get_wpar_rise()
static const int MAX_SAMPLES
static const int BADCHISQ_MAT
static const int PIN_BASELINE_RANGE
static const int BAD_MAX
WaveFormPar * csiTestWpar[4]
int16_t good_baseline()
long double lineq_matrix[20][20]
int FitCsIShape(int, ShapePar *, WaveFormPar *)
double get_parabolic_T0()
static const int NOISE_LEVEL_CSI
static const int BADCHISQ_NEG
static const int BADCHISQ_T0
void SetData(const TFragment &fragment, double=0)
TGraph * GetWaveGraph()
bool GetSiliShape(double tauDecay, double tauRise)
double Get_wpar_amplitude()
double Get_wpar_baselinefin()
double get_sin_par(double)
static const int BAD_PID