GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TLMFitter.h
Go to the documentation of this file.
1// Author: Ryan Dunlop 11/15
2
3#ifndef TLMFITTER_H
4#define TLMFITTER_H
5
6/** \addtogroup Fitting Fitting & Analysis
7 * @{
8 */
9
10#include "Globals.h"
11
12#include <vector>
13#include <complex>
14#include <iostream>
15#include <fstream>
16
17#include "Rtypes.h"
18#include "TObject.h"
19#include "TH1.h"
20#include "TF1.h"
21
22/////////////////////////////////////////////////////////////////
23///
24/// \class TLMFitter
25///
26/// This Class can be used to fit weighted-poisson distributed
27/// data. It is originally from Numerical recipes, and adapted
28/// by G.F Grinyer. It is based on the non-linear
29/// Levenberg-Marquardt minimization algorithm.
30///
31/////////////////////////////////////////////////////////////////
32
33// Vector Types
34
35// Overloaded complex operations to handle mixed float and double
36// This takes care of e.g. 1.0/z, z complex<float>
37
38template <class T>
39class NRVec {
40private:
41 int nn; // size of array. upper index is nn-1
42 T* v;
43
44public:
45 NRVec();
46 explicit NRVec(int n); // Zero-based array
47 NRVec(const T& a, int n); // initialize to constant value
48 NRVec(const T* a, int n); // Initialize to array
49 NRVec(const NRVec& rhs); // Copy constructor
50 NRVec(NRVec&& rhs) noexcept = default; // Move constructor
51 NRVec& operator=(const NRVec& rhs); // copy assignment
52 NRVec& operator=(NRVec&& rhs) noexcept = default; // move assignment
53 NRVec& operator=(const T& a); // assign a to every element
54 inline T& operator[](int i); // i'th element
55 inline const T& operator[](int i) const;
56 inline int size() const;
57 ~NRVec();
58};
59
60template <class T>
61NRVec<T>::NRVec() : nn(0), v(0)
62{
63}
64
65template <class T>
66NRVec<T>::NRVec(int n) : nn(n), v(new T[n]())
67{
68 for(int i = 0; i < n; ++i) {
69 v[i] = 0.0;
70 }
71}
72
73template <class T>
74NRVec<T>::NRVec(const T& a, int n) : nn(n), v(new T[n]())
75{
76 for(int i = 0; i < n; i++) {
77 v[i] = a;
78 }
79}
80
81template <class T>
82NRVec<T>::NRVec(const T* a, int n) : nn(n), v(new T[n]())
83{
84 for(int i = 0; i < n; i++) {
85 v[i] = *a++;
86 }
87}
88
89template <class T>
90NRVec<T>::NRVec(const NRVec<T>& rhs) : nn(rhs.nn), v(new T[nn]())
91{
92 for(int i = 0; i < nn; i++) {
93 v[i] = rhs[i];
94 }
95}
96
97template <class T>
99// postcondition: normal assignment via copying has been performed;
100// if vector and rhs were different sizes, vector
101// has been resized to match the size of rhs
102{
103 if(this != &rhs) {
104 if(nn != rhs.nn) {
105 delete[] v;
106 nn = rhs.nn;
107 v = new T[nn];
108 }
109 for(int i = 0; i < nn; i++) {
110 v[i] = rhs[i];
111 }
112 }
113 return *this;
114}
115
116template <class T>
117NRVec<T>& NRVec<T>::operator=(const T& a) // assign a to every element
118{
119 for(int i = 0; i < nn; i++) {
120 v[i] = a;
121 }
122 return *this;
123}
124
125template <class T>
126inline T& NRVec<T>::operator[](const int i) // subscripting
127{
128 return v[i];
129}
130
131template <class T>
132inline const T& NRVec<T>::operator[](const int i) const // subscripting
133{
134 return v[i];
135}
136
137template <class T>
138inline int NRVec<T>::size() const
139{
140 return nn;
141}
142
143template <class T>
145{
146 if(v != nullptr) {
147 delete[] (v);
148 }
149}
150
151template <class T>
152class NRMat {
153private:
154 int nn;
155 int mm;
156 T** v;
157
158public:
159 NRMat();
160 NRMat(int n, int m); // Zero-based array
161 NRMat(const T& a, int n, int m); // Initialize to constant
162 NRMat(const T* a, int n, int m); // Initialize to array
163 NRMat(const NRMat& rhs); // Copy constructor
164 NRMat(NRMat&& rhs) noexcept = default; // Move constructor
165 NRMat& operator=(const NRMat& rhs); // copy assignment
166 NRMat& operator=(NRMat&& rhs) noexcept = default; // move assignment
167 NRMat& operator=(const T& a); // assign a to every element
168 inline T* operator[](int i); // subscripting: pointer to row i
169 inline const T* operator[](int i) const;
170 inline int nrows() const;
171 inline int ncols() const;
172 ~NRMat();
173};
174
175template <class T>
176NRMat<T>::NRMat() : nn(0), mm(0), v(0)
177{
178}
179
180template <class T>
181NRMat<T>::NRMat(int n, int m) : nn(n), mm(m), v(new T*[n]())
182{
183 v[0] = new T[m * n]();
184 for(int i = 1; i < n; i++) {
185 v[i] = v[i - 1] + m;
186 }
187
188 for(int i = 0; i < n; ++i) {
189 for(int j = 0; j < m; ++j) {
190 v[i][j] = 0.0;
191 }
192 }
193}
194
195template <class T>
196NRMat<T>::NRMat(const T& a, int n, int m) : nn(n), mm(m), v(new T*[n]())
197{
198 v[0] = new T[m * n]();
199 for(int i = 1; i < n; i++) {
200 v[i] = v[i - 1] + m;
201 }
202 for(int i = 0; i < n; i++) {
203 for(int j = 0; j < m; j++) {
204 v[i][j] = a;
205 }
206 }
207}
208
209template <class T>
210NRMat<T>::NRMat(const T* a, int n, int m) : nn(n), mm(m), v(new T*[n]())
211{
212 v[0] = new T[m * n]();
213 for(int i = 1; i < n; i++) {
214 v[i] = v[i - 1] + m;
215 }
216 for(int i = 0; i < n; i++) {
217 for(int j = 0; j < m; j++) {
218 v[i][j] = *a++;
219 }
220 }
221}
222
223template <class T>
224NRMat<T>::NRMat(const NRMat& rhs) : nn(rhs.nn), mm(rhs.mm), v(new T*[nn]())
225{
226 v[0] = new T[mm * nn]();
227 for(int i = 1; i < nn; i++) {
228 v[i] = v[i - 1] + mm;
229 }
230 for(int i = 0; i < nn; i++) {
231 for(int j = 0; j < mm; j++) {
232 v[i][j] = rhs[i][j];
233 }
234 }
235}
236
237template <class T>
239// postcondition: normal assignment via copying has been performed;
240// if matrix and rhs were different sizes, matrix
241// has been resized to match the size of rhs
242{
243 if(this != &rhs) {
244 if(nn != rhs.nn || mm != rhs.mm) {
245 delete[] v[0];
246 delete[] v;
247 nn = rhs.nn;
248 mm = rhs.mm;
249 v = new T*[nn];
250 v[0] = new T[mm * nn];
251 }
252 for(int i = 1; i < nn; i++) {
253 v[i] = v[i - 1] + mm;
254 }
255 for(int i = 0; i < nn; i++) {
256 for(int j = 0; j < mm; j++) {
257 v[i][j] = rhs[i][j];
258 }
259 }
260 }
261 return *this;
262}
263
264template <class T>
265NRMat<T>& NRMat<T>::operator=(const T& a) // assign a to every element
266{
267 for(int i = 0; i < nn; i++) {
268 for(int j = 0; j < mm; j++) {
269 v[i][j] = a;
270 }
271 }
272 return *this;
273}
274
275template <class T>
276inline T* NRMat<T>::operator[](const int i) // subscripting: pointer to row i
277{
278 return v[i];
279}
280
281template <class T>
282inline const T* NRMat<T>::operator[](const int i) const
283{
284 return v[i];
285}
286
287template <class T>
288inline int NRMat<T>::nrows() const
289{
290 return nn;
291}
292
293template <class T>
294inline int NRMat<T>::ncols() const
295{
296 return mm;
297}
298
299template <class T>
301{
302 if(v != nullptr) {
303 delete[] (v[0]);
304 delete[] (v);
305 }
306}
307
308template <class T>
309class NRMat3d {
310private:
311 int nn;
312 int mm;
313 int kk;
314 T*** v;
315
316public:
317 NRMat3d();
318 NRMat3d(int n, int m, int k);
319 inline T** operator[](int i); // subscripting: pointer to row i
320 inline const T* const* operator[](int i) const;
321 inline int dim1() const;
322 inline int dim2() const;
323 inline int dim3() const;
324 NRMat3d(const NRMat3d&) = default;
325 NRMat3d(NRMat3d&&) noexcept = default;
326 NRMat3d& operator=(const NRMat3d&) = default;
327 NRMat3d& operator=(NRMat3d&&) noexcept = default;
328 ~NRMat3d();
329};
330
331template <class T>
332NRMat3d<T>::NRMat3d() : nn(0), mm(0), kk(0), v(0)
333{
334}
335
336template <class T>
337NRMat3d<T>::NRMat3d(int n, int m, int k) : nn(n), mm(m), kk(k), v(new T**[n])
338{
339 v[0] = new T*[n * m];
340 v[0][0] = new T[n * m * k];
341 for(int j = 1; j < m; j++) {
342 v[0][j] = v[0][j - 1] + k;
343 }
344 for(int i = 1; i < n; i++) {
345 v[i] = v[i - 1] + m;
346 v[i][0] = v[i - 1][0] + m * k;
347 for(int j = 1; j < m; j++) {
348 v[i][j] = v[i][j - 1] + k;
349 }
350 }
351}
352
353template <class T>
354inline T** NRMat3d<T>::operator[](const int i) // subscripting: pointer to row i
355{
356 return v[i];
357}
358
359template <class T>
360inline const T* const* NRMat3d<T>::operator[](const int i) const
361{
362 return v[i];
363}
364
365template <class T>
366inline int NRMat3d<T>::dim1() const
367{
368 return nn;
369}
370
371template <class T>
372inline int NRMat3d<T>::dim2() const
373{
374 return mm;
375}
376
377template <class T>
378inline int NRMat3d<T>::dim3() const
379{
380 return kk;
381}
382
383template <class T>
385{
386 if(v != 0) {
387 delete[] (v[0][0]);
388 delete[] (v[0]);
389 delete[] (v);
390 }
391}
392
397
399using Vec_I_CHR = const NRVec<char>;
402
407
409using Vec_I_INT = const NRVec<int>;
412
417
422
427
429using Vec_I_SP = const NRVec<float>;
432
437
442
447
448// Matrix Types
449
454
456using Mat_I_CHR = const NRMat<char>;
459
464
466using Mat_I_INT = const NRMat<int>;
469
474
479
484
486using Mat_I_SP = const NRMat<float>;
489
494
499
504
505// 3D Matrix Types
506
511
512// Miscellaneous Types
513
517
518template <class T>
519inline T SQR(const T a)
520{
521 return a * a;
522}
523
524template <class T>
525inline T MAX(const T& a, const T& b)
526{
527 return b > a ? (b) : (a);
528}
529
530inline float MAX(const double& a, const float& b)
531{
532 return b > a ? (b) : static_cast<float>(a);
533}
534
535inline float MAX(const float& a, const double& b)
536{
537 return b > a ? static_cast<float>(b) : (a);
538}
539
540template <class T>
541inline T SIGN(const T& a, const T& b)
542{
543 return b >= 0 ? (a >= 0 ? a : -a) : (a >= 0 ? -a : a);
544}
545
546inline float SIGN(const float& a, const double& b)
547{
548 return b >= 0 ? (a >= 0 ? a : -a) : (a >= 0 ? -a : a);
549}
550
551inline float SIGN(const double& a, const float& b)
552{
553 return b >= 0 ? static_cast<float>(a >= 0 ? a : -a) : static_cast<float>(a >= 0 ? -a : a);
554}
555
556template <class T>
557inline void SWAP(T& a, T& b)
558{
559 T dum = a;
560 a = b;
561 b = dum;
562}
563
564inline std::complex<float> operator+(const double& a, const std::complex<float>& b)
565{
566 return static_cast<float>(a) + b;
567}
568
569inline std::complex<float> operator+(const std::complex<float>& a, const double& b)
570{
571 return a + static_cast<float>(b);
572}
573
574inline std::complex<float> operator-(const double& a, const std::complex<float>& b)
575{
576 return static_cast<float>(a) - b;
577}
578
579inline std::complex<float> operator-(const std::complex<float>& a, const double& b)
580{
581 return a - static_cast<float>(b);
582}
583
584inline std::complex<float> operator*(const double& a, const std::complex<float>& b)
585{
586 return static_cast<float>(a) * b;
587}
588
589inline std::complex<float> operator*(const std::complex<float>& a, const double& b)
590{
591 return a * static_cast<float>(b);
592}
593
594inline std::complex<float> operator/(const double& a, const std::complex<float>& b)
595{
596 return static_cast<float>(a) / b;
597}
598
599inline std::complex<float> operator/(const std::complex<float>& a, const double& b)
600{
601 return a / static_cast<float>(b);
602}
603
604class TLMFitter : public TObject {
605public:
606 TLMFitter() = default;
607 TLMFitter(const TLMFitter&) = default;
608 TLMFitter(TLMFitter&&) noexcept = default;
609 TLMFitter& operator=(const TLMFitter&) = default;
610 TLMFitter& operator=(TLMFitter&&) noexcept = default;
611 ~TLMFitter() = default;
612
613private:
615 TH1* fHist{nullptr};
616 TF1* fFunction{nullptr};
618 int fRangeMin{0};
619 int fRangeMax{0};
620
621public:
622 template <class T>
623 class NRVec;
624 template <class T>
625 class NRMat;
626 template <class T>
627 class NRMat3d;
628
629 void Fit(TH1* hist, TF1* func);
630
631protected:
632 void SetFitterRange(int min, int max)
633 {
634 fRangeMin = min;
635 fRangeMax = max;
636 }
637
638 inline void nrerror(const std::string& error_text)
639 // Numerical Recipes standard error handler
640 {
641 std::cerr << "Numerical Recipes run-time error..." << std::endl;
642 std::cerr << error_text << std::endl;
643 std::cerr << "...now exiting to system..." << std::endl;
644 exit(1);
645 }
646
647 void funcs(const double& x, Vec_IO_double& a, double& y, Vec_O_double& dyda);
649 Mat_O_double& alpha, double& chisq, Vec_I_double& W, double& alamda);
650
652 Vec_O_double& beta, double& chisq, Vec_I_double& W, double& chisqexp);
654 void covsrt(Mat_IO_double& covar, Vec_I_BOOL& ia, int mfit);
656 Vec_double& dyda, int chisqnumber, const double& bin_width, Vec_double& yfit, const int& bin);
657
658public:
659 /// \cond CLASSIMP
660 ClassDefOverride(TLMFitter, 1) // NOLINT(readability-else-after-return)
661 /// \endcond
662};
663/*! @} */
664#endif // TLMFitter_H
TH1D * hist
Definition UserFillObj.h:3
T *** v
Definition TLMFitter.h:314
NRMat3d(NRMat3d &&) noexcept=default
NRMat3d(const NRMat3d &)=default
int nn
Definition TLMFitter.h:154
T ** v
Definition TLMFitter.h:156
NRMat(NRMat &&rhs) noexcept=default
NRMat & operator=(NRMat &&rhs) noexcept=default
int mm
Definition TLMFitter.h:155
int nn
Definition TLMFitter.h:41
NRVec & operator=(NRVec &&rhs) noexcept=default
T * v
Definition TLMFitter.h:42
NRVec(NRVec &&rhs) noexcept=default
TH1 * fHist
Definition TLMFitter.h:615
void nrerror(const std::string &error_text)
Definition TLMFitter.h:638
void covsrt(Mat_IO_double &covar, Vec_I_BOOL &ia, int mfit)
TLMFitter(const TLMFitter &)=default
TLMFitter(TLMFitter &&) noexcept=default
int integrator(Vec_I_double &x, Vec_I_double &y, Vec_double &sig, Vec_I_double &W, Vec_IO_double &a, Vec_double &dyda, int chisqnumber, const double &bin_width, Vec_double &yfit, const int &bin)
void SetFitterRange(int min, int max)
Definition TLMFitter.h:632
int fIntegrationSteps
Definition TLMFitter.h:614
int fInitChi2Number
Definition TLMFitter.h:617
int fRangeMin
Definition TLMFitter.h:618
TLMFitter()=default
int fRangeMax
Definition TLMFitter.h:619
void mrqcof(Vec_I_double &x, Vec_I_double &y, Vec_double &sig, Vec_IO_double &a, Vec_I_BOOL &ia, Mat_O_double &alpha, Vec_O_double &beta, double &chisq, Vec_I_double &W, double &chisqexp)
TF1 * fFunction
Definition TLMFitter.h:616
void mrqmin(Vec_I_double &x, Vec_I_double &y, Vec_double &sig, Vec_IO_double &a, Vec_I_BOOL &ia, Mat_O_double &covar, Mat_O_double &alpha, double &chisq, Vec_I_double &W, double &alamda)
*******************************************************************‍/
void gaussj(Mat_IO_double &a, Mat_IO_double &b)
void funcs(const double &x, Vec_IO_double &a, double &y, Vec_O_double &dyda)
Definition TLMFitter.cxx:23
void Fit(TH1 *hist, TF1 *func)
Definition TLMFitter.cxx:35
T SQR(const T a)
Definition TLMFitter.h:519
T & operator[](int i)
Definition TLMFitter.h:126
int nrows() const
Definition TLMFitter.h:288
~NRMat()
Definition TLMFitter.h:300
int ncols() const
Definition TLMFitter.h:294
T ** operator[](int i)
Definition TLMFitter.h:354
std::complex< float > operator+(const double &a, const std::complex< float > &b)
Definition TLMFitter.h:564
int size() const
Definition TLMFitter.h:138
T * operator[](int i)
Definition TLMFitter.h:276
int dim1() const
Definition TLMFitter.h:366
NRVec()
Definition TLMFitter.h:61
void SWAP(T &a, T &b)
Definition TLMFitter.h:557
std::complex< float > operator-(const double &a, const std::complex< float > &b)
Definition TLMFitter.h:574
NRMat & operator=(const NRMat &rhs)
Definition TLMFitter.h:238
T MAX(const T &a, const T &b)
Definition TLMFitter.h:525
int dim3() const
Definition TLMFitter.h:378
std::complex< float > operator/(const double &a, const std::complex< float > &b)
Definition TLMFitter.h:594
NRVec & operator=(const NRVec &rhs)
Definition TLMFitter.h:98
T SIGN(const T &a, const T &b)
Definition TLMFitter.h:541
int dim2() const
Definition TLMFitter.h:372
~NRVec()
Definition TLMFitter.h:144
std::complex< float > operator*(const double &a, const std::complex< float > &b)
Definition TLMFitter.h:584
NRMat()
Definition TLMFitter.h:176