Logo ROOT  
Reference Guide
TMultiDimFit.cxx
Go to the documentation of this file.
1// @(#)root/hist:$Id$
2// Author: Christian Holm Christensen 07/11/2000
3
4/** \class TMultiDimFit
5 \ingroup Hist
6
7 Multidimensional Fits in ROOT.
8 ## Overview
9 A common problem encountered in different fields of applied science is
10 to find an expression for one physical quantity in terms of several
11 others, which are directly measurable.
12
13 An example in high energy physics is the evaluation of the momentum of
14 a charged particle from the observation of its trajectory in a magnetic
15 field. The problem is to relate the momentum of the particle to the
16 observations, which may consists of positional measurements at
17 intervals along the particle trajectory.
18
19 The exact functional relationship between the measured quantities
20 (e.g., the space-points) and the dependent quantity (e.g., the
21 momentum) is in general not known, but one possible way of solving the
22 problem, is to find an expression which reliably approximates the
23 dependence of the momentum on the observations.
24
25 This explicit function of the observations can be obtained by a
26 <I>least squares</I> fitting procedure applied to a representative
27 sample of the data, for which the dependent quantity (e.g., momentum)
28 and the independent observations are known. The function can then be
29 used to compute the quantity of interest for new observations of the
30 independent variables.
31
32 This class <TT>TMultiDimFit</TT> implements such a procedure in
33 ROOT. It is largely based on the CERNLIB MUDIFI package [2].
34 Though the basic concepts are still sound, and
35 therefore kept, a few implementation details have changed, and this
36 class can take advantage of MINUIT [4] to improve the errors
37 of the fitting, thanks to the class TMinuit.
38
39 In [5] and [6] H. Wind demonstrates the utility
40 of this procedure in the context of tracking, magnetic field
41 parameterisation, and so on. The outline of the method used in this
42 class is based on Winds discussion, and I refer these two excellents
43 text for more information.
44
45 And example of usage is given in multidimfit.C.
46
47 ## The Method
48 Let \f$ D \f$ by the dependent quantity of interest, which depends smoothly
49 on the observable quantities \f$ x_1, \ldots, x_N \f$ which we'll denote by
50 \f$\mathbf{x}\f$. Given a training sample of \f$ M\f$ tuples of the form, (TMultiDimFit::AddRow)
51
52 \f[
53 \left(\mathbf{x}_j, D_j, E_j\right)\quad,
54 \f]
55 where \f$\mathbf{x}_j = (x_{1,j},\ldots,x_{N,j})\f$ are \f$ N\f$ independent
56 variables, \f$ D_j\f$ is the known, quantity dependent at \f$\mathbf{x}_j\f$ and \f$ E_j\f$ is
57 the square error in \f$ D_j\f$, the class will try to find the parameterization
58 \f[
59 D_p(\mathbf{x}) = \sum_{l=1}^{L} c_l \prod_{i=1}^{N} p_{li}\left(x_i\right)
60 = \sum_{l=1}^{L} c_l F_l(\mathbf{x})
61 \f]
62 such that
63
64 \f[
65 S \equiv \sum_{j=1}^{M} \left(D_j - D_p\left(\mathbf{x}_j\right)\right)^2
66 \f]
67 is minimal. Here \f$p_{li}(x_i)\f$ are monomials, or Chebyshev or Legendre
68 polynomials, labelled \f$l = 1, \ldots, L\f$, in each variable \f$ x_i\f$,\f$ i=1, \ldots, N\f$.
69
70 So what TMultiDimFit does, is to determine the number of terms \f$ L\f$, and then \f$ L\f$ terms
71 (or functions) \f$ F_l\f$, and the \f$ L\f$ coefficients \f$ c_l\f$, so that \f$ S\f$ is minimal
72 (TMultiDimFit::FindParameterization).
73
74 Of course it's more than a little unlikely that \f$ S\f$ will ever become
75 exact zero as a result of the procedure outlined below. Therefore, the
76 user is asked to provide a minimum relative error \f$ \epsilon\f$ (TMultiDimFit::SetMinRelativeError),
77 and \f$ S\f$ will be considered minimized when
78
79 \f[
80 R = \frac{S}{\sum_{j=1}^M D_j^2} < \epsilon
81 \f]
82 Optionally, the user may impose a functional expression by specifying
83 the powers of each variable in \f$ L\f$ specified functions \f$ F_1, \ldots,F_L\f$ (TMultiDimFit::SetPowers).
84 In that case, only the coefficients \f$ c_l\f$ is calculated by the class.
85
86 ## Limiting the Number of Terms
87 As always when dealing with fits, there's a real chance of *over fitting*. As is well-known, it's
88 always possible to fit an \f$ N-1\f$ polynomial in \f$ x\f$ to \f$ N\f$ points \f$ (x,y)\f$ with
89 \f$\chi^2 = 0\f$, but the polynomial is not likely to fit new data at all [1].
90 Therefore, the user is asked to provide an upper limit, \f$ L_{max}\f$ to the number of terms in
91 \f$ D_p\f$ (TMultiDimFit::SetMaxTerms).
92
93 However, since there's an infinite number of \f$ F_l\f$ to choose from, the
94 user is asked to give the maximum power. \f$ P_{max,i}\f$, of each variable
95 \f$ x_i\f$ to be considered in the minimization of \f$ S\f$ (TMultiDimFit::SetMaxPowers).
96
97 One way of obtaining values for the maximum power in variable \f$ i\f$, is
98 to perform a regular fit to the dependent quantity \f$ D\f$, using a
99 polynomial only in \f$ x_i\f$. The maximum power is \f$ P_{max,i}\f$ is then the
100 power that does not significantly improve the one-dimensional
101 least-square fit over \f$ x_i\f$ to \f$ D\f$ [5].
102
103 There are still a huge amount of possible choices for \f$ F_l\f$; in fact
104 there are \f$\prod_{i=1}^{N} (P_{max,i} + 1)\f$ possible
105 choices. Obviously we need to limit this. To this end, the user is
106 asked to set a *power control limit*, \f$ Q\f$ (TMultiDimFit::SetPowerLimit), and a function
107 \f$ F_l\f$ is only accepted if
108 \f[
109 Q_l = \sum_{i=1}^{N} \frac{P_{li}}{P_{max,i}} < Q
110 \f]
111 where \f$ P_{li}\f$ is the leading power of variable \f$ x_i\f$ in function \f$ F_l\f$ (TMultiDimFit::MakeCandidates).
112 So the number of functions increase with \f$ Q\f$ (1, 2 is fine, 5 is way out).
113
114 ## Gram-Schmidt Orthogonalisation</A>
115 To further reduce the number of functions in the final expression,
116 only those functions that significantly reduce \f$ S\f$ is chosen. What
117 `significant' means, is chosen by the user, and will be
118 discussed below (see [2.3](TMultiFimFit.html#sec:selectiondetail)).
119
120 The functions \f$ F_l\f$ are generally not orthogonal, which means one will
121 have to evaluate all possible \f$ F_l\f$'s over all data-points before
122 finding the most significant [1]. We can, however, do
123 better then that. By applying the *modified Gram-Schmidt
124 orthogonalisation* algorithm [5] [3] to the
125 functions \f$ F_l\f$, we can evaluate the contribution to the reduction of
126 \f$ S\f$ from each function in turn, and we may delay the actual inversion
127 of the curvature-matrix (TMultiDimFit::MakeGramSchmidt).
128
129 So we are let to consider an \f$ M\times L\f$ matrix \f$\mathsf{F}\f$, an
130 element of which is given by
131 \f[
132 f_{jl} = F_j\left(x_{1j} , x_{2j}, \ldots, x_{Nj}\right)
133 = F_l(\mathbf{x}_j)\, \quad\mbox{with}~j=1,2,\ldots,M,
134 \f]
135 where \f$ j\f$ labels the \f$ M\f$ rows in the training sample and \f$ l\f$ labels
136 \f$ L\f$ functions of \f$ N\f$ variables, and \f$ L \leq M\f$. That is, \f$ f_{jl}\f$ is
137 the term (or function) numbered \f$ l\f$ evaluated at the data point
138 \f$ j\f$. We have to normalise \f$\mathbf{x}_j\f$ to \f$ [-1,1]\f$ for this to
139 succeed [5] (TMultiDimFit::MakeNormalized). We then define a
140 matrix \f$\mathsf{W}\f$ of which the columns \f$\mathbf{w}_j\f$ are given by
141 \f{eqnarray*}{
142 \mathbf{w}_1 &=& \mathbf{f}_1 = F_1\left(\mathbf x_1\right)\\
143 \mathbf{w}_l &=& \mathbf{f}_l - \sum^{l-1}_{k=1} \frac{\mathbf{f}_l \bullet
144 \mathbf{w}_k}{\mathbf{w}_k^2}\mathbf{w}_k\,.
145 \f}
146 and \f$\mathbf{w}_{l}\f$ is the component of \f$\mathbf{f}_{l} \f$ orthogonal
147 to \f$\mathbf{w}_{1}, \ldots, \mathbf{w}_{l-1}\f$. Hence we obtain [3],
148 \f[
149 \mathbf{w}_k\bullet\mathbf{w}_l = 0\quad\mbox{if}~k \neq l\quad.
150 \f]
151 We now take as a new model \f$\mathsf{W}\mathbf{a}\f$. We thus want to
152 minimize
153 \f[
154 S\equiv \left(\mathbf{D} - \mathsf{W}\mathbf{a}\right)^2\quad,
155 \f]
156 where \f$\mathbf{D} = \left(D_1,\ldots,D_M\right)\f$ is a vector of the
157 dependent quantity in the sample. Differentiation with respect to
158 \f$ a_j\f$ gives, using [6], <a name="eq:dS2"></a>
159 \f[
160 \mathbf{D}\bullet\mathbf{w}_l - a_l\mathbf{w}_l^2 = 0
161 \f]
162 or
163 \f[
164 a_l = \frac{\mathbf{D}_l\bullet\mathbf{w}_l}{\mathbf{w}_l^2}
165 \f]
166 Let \f$ S_j\f$ be the sum of squares of residuals when taking \f$ j\f$ functions
167 into account. Then
168 \f[
169 S_l = \left[\mathbf{D} - \sum^l_{k=1} a_k\mathbf{w}_k\right]^2
170 = \mathbf{D}^2 - 2\mathbf{D} \sum^l_{k=1} a_k\mathbf{w}_k
171 + \sum^l_{k=1} a_k^2\mathbf{w}_k^2
172 \f]
173 Using [9], we see that
174 \f{eqnarray*}{
175 S_l &=& \mathbf{D}^2 - 2 \sum^l_{k=1} a_k^2\mathbf{w}_k^2 +
176 \sum^j_{k=1} a_k^2\mathbf{w}_k^2\nonumber\\
177 &=& \mathbf{D}^2 - \sum^l_{k=1} a_k^2\mathbf{w}_k^2\nonumber\\
178 &=& \mathbf{D}^2 - \sum^l_{k=1} \frac{\left(\mathbf D\bullet \mathbf
179 w_k\right)}{\mathbf w_k^2}
180 \f}
181 So for each new function \f$ F_l\f$ included in the model, we get a
182 reduction of the sum of squares of residuals of \f$a_l^2\mathbf{w}_l^2\f$,
183 where \f$\mathbf{w}_l\f$ is given by [4] and \f$ a_l\f$ by [9]. Thus, using
184 the Gram-Schmidt orthogonalisation, we
185 can decide if we want to include this function in the final model,
186 *before* the matrix inversion.
187
188 ## Function Selection Based on Residual
189 Supposing that \f$ L-1\f$ steps of the procedure have been performed, the
190 problem now is to consider the \f$L^{\mbox{th}}\f$ function.
191
192 The sum of squares of residuals can be written as
193 \f[
194 S_L = \textbf{D}^T\bullet\textbf{D} -
195 \sum^L_{l=1}a^2_l\left(\textbf{w}_l^T\bullet\textbf{w}_l\right)
196 \f]
197 where the relation [9] have been taken into account. The
198 contribution of the \f$L^{\mbox{th}}\f$ function to the reduction of S, is
199 given by
200 \f[
201 \Delta S_L = a^2_L\left(\textbf{w}_L^T\bullet\textbf{w}_L\right)
202 \f]
203 Two test are now applied to decide whether this \f$L^{\mbox{th}}\f$
204 function is to be included in the final expression, or not.
205
206 ## Test 1
207 Denoting by \f$ H_{L-1}\f$ the subspace spanned by \f$\textbf{w}_1,\ldots,\textbf{w}_{L-1}\f$
208 the function \d$\textbf{w}_L\d$ is by construction (see 4) the projection of the function
209 \f$ F_L\f$ onto the direction perpendicular to \f$ H_{L-1}\f$. Now, if the
210 length of \f$\textbf{w}_L\f$ (given by \f$\textbf{w}_L\bullet\textbf{w}_L\f$)
211 is very small compared to the length of \f$\textbf{f}_L\f$ this new
212 function can not contribute much to the reduction of the sum of
213 squares of residuals. The test consists then in calculating the angle
214 \f$ \theta \f$ between the two vectors \f$\textbf{w}_L\f$ \f$ \textbf {f}_L\f$
215 (see also figure 1) and requiring that it's
216 *greater* then a threshold value which the user must set (TMultiDimFit::SetMinAngle).
217
218 \image html multidimfit_img86.gif "Figure 1: (a) angle \f$\theta\f$ between \f$\textbf{w}_l\f$ and \f$\textbf{f}_L\f$, (b) angle \f$ \phi \f$ between \f$\textbf{w}_L\f$ and \f$\textbf{D}\f$"
219
220 ## Test 2
221 Let \f$\textbf{D}\f$ be the data vector to be fitted. As illustrated in
222 figure 1, the \f$L^{\mbox{th}}\f$ function \f$\textbf{w}_L\f$
223 will contribute significantly to the reduction of \f$ S\f$, if the angle
224 \f$\phi^\prime\f$ between \f$\textbf{w}_L\f$ and \f$\textbf{D}\f$ is smaller than
225 an upper limit \f$ \phi \f$, defined by the user (MultiDimFit::SetMaxAngle)
226
227 However, the method automatically readjusts the value of this angle
228 while fitting is in progress, in order to make the selection criteria
229 less and less difficult to be fulfilled. The result is that the
230 functions contributing most to the reduction of \f$ S\f$ are chosen first
231 (TMultiDimFit::TestFunction).
232
233 In case \f$ \phi \f$ isn't defined, an alternative method of
234 performing this second test is used: The \f$L^{\mbox{th}}\f$
235 function \f$\textbf{f}_L\f$ is accepted if (refer also to equation (13))
236 \f[
237 \Delta S_L > \frac{S_{L-1}}{L_{max}-L}
238 \f]
239 where \f$ S_{L-1}\f$ is the sum of the \f$ L-1\f$ first residuals from the
240 \f$ L-1\f$ functions previously accepted; and \f$ L_{max}\f$ is the total number
241 of functions allowed in the final expression of the fit (defined by
242 user).
243
244 From this we see, that by restricting \f$ L_{max}\f$ -- the number of
245 terms in the final model -- the fit is more difficult to perform,
246 since the above selection criteria is more limiting.
247
248 The more coefficients we evaluate, the more the sum of squares of
249 residuals \f$ S\f$ will be reduced. We can evaluate \f$ S\f$ before inverting
250 \f$\mathsf{B}\f$ as shown below.
251
252 ## Coefficients and Coefficient Errors
253 Having found a parameterization, that is the \f$ F_l\f$'s and \f$ L\f$, that
254 minimizes \f$ S\f$, we still need to determine the coefficients
255 \f$ c_l\f$. However, it's a feature of how we choose the significant
256 functions, that the evaluation of the \f$ c_l\f$'s becomes trivial [5]. To derive
257 \f$\mathbf{c}\f$, we first note that
258 equation (4) can be written as
259 \f[
260 \mathsf{F} = \mathsf{W}\mathsf{B}
261 \f]
262 where
263 \f{eqnarray*}{
264 b_{ij} = \frac{\mathbf{f}_j \bullet \mathbf{w}_i}{\mathbf{w}_i^2}
265 & \mbox{if} & i < j\\
266 1 & \mbox{if} & i = j\\
267 0 & \mbox{if} & i > j
268 \f}
269 Consequently, \f$\mathsf{B}\f$ is an upper triangle matrix, which can be
270 readily inverted. So we now evaluate
271 \f[
272 \mathsf{F}\mathsf{B}^{-1} = \mathsf{W}
273 \f]
274 The model \f$\mathsf{W}\mathbf{a}\f$ can therefore be written as
275 \f$(\mathsf{F}\mathsf{B}^{-1})\mathbf{a} = \mathsf{F}(\mathsf{B}^{-1}\mathbf{a})\,.\f$
276
277 The original model \f$\mathsf{F}\mathbf{c}\f$ is therefore identical with
278 this if
279 \f[
280 \mathbf{c} = \left(\mathsf{B}^{-1}\mathbf{a}\right) =
281 \left[\mathbf{a}^T\left(\mathsf{B}^{-1}\right)^T\right]^T\,.
282 \f]
283 The reason we use \f$\left(\mathsf{B}^{-1}\right)^T\f$ rather then
284 \f$\mathsf{B}^{-1}\f$ is to save storage, since \f$\left(\mathsf{B}^{-1}\right)^T\f$
285 can be stored in the same matrix as \f$\mathsf{B}\f$ (TMultiDimFit::MakeCoefficients).
286 The errors in the coefficients is calculated by inverting the curvature matrix
287 of the non-orthogonal functions \f$ f_{lj}\f$ [1] (TMultiDimFit::MakeCoefficientErrors).
288
289 ## Considerations
290 It's important to realize that the training sample should be
291 representative of the problem at hand, in particular along the borders
292 of the region of interest. This is because the algorithm presented
293 here, is a *interpolation*, rather then a *extrapolation* [5].
294
295 Also, the independent variables \f$ x_{i}\f$ need to be linear
296 independent, since the procedure will perform poorly if they are not
297 [5]. One can find an linear transformation from ones
298 original variables \f$ \xi_{i}\f$ to a set of linear independent variables
299 \f$ x_{i}\f$, using a *Principal Components Analysis* (see TPrincipal), and
300 then use the transformed variable as input to this class [5] [6].
301
302 H. Wind also outlines a method for parameterising a multidimensional
303 dependence over a multidimensional set of variables. An example
304 of the method from [5], is a follows (please refer to
305 [5] for a full discussion):
306
307 1. Define \f$\mathbf{P} = (P_1, \ldots, P_5)\f$ are the 5 dependent
308 quantities that define a track.
309 2. Compute, for \f$ M\f$ different values of \f$\mathbf{P}\f$, the tracks
310 through the magnetic field, and determine the corresponding
311 \f$\mathbf{x} = (x_1, \ldots, x_N)\f$.
312 3. Use the simulated observations to determine, with a simple
313 approximation, the values of \f$\mathbf{P}_j\f$. We call these values
314 \f$\mathbf{P}^\prime_j, j = 1, \ldots, M\f$.
315 4. Determine from \f$\mathbf{x}\f$ a set of at least five relevant
316 coordinates \f$\mathbf{x}^\prime\f$, using contrains, *or
317 alternative:*
318 5. Perform a Principal Component Analysis (using TPrincipal), and use
319 to get a linear transformation \f$\mathbf{x} \rightarrow \mathbf{x}^\prime\f$, so that
320 \f$\mathbf{x}^\prime\f$ are constrained and linear independent.
321 6. Perform a Principal Component Analysis on
322 \f$Q_i = P_i / P^\prime_i\, i = 1, \ldots, 5\f$, to get linear
323 indenpendent (among themselves, but not independent of \f$\mathbf{x}\f$) quantities
324 \f$\mathbf{Q}^\prime\f$
325 7. For each component \f$Q^\prime_i\f$ make a multidimensional fit,
326 using \f$\mathbf{x}^\prime\f$ as the variables, thus determining a set of
327 coefficients \f$\mathbf{c}_i\f$.
328
329 To process data, using this parameterisation, do
330 1. Test wether the observation \f$\mathbf{x}\f$ within the domain of
331 the parameterization, using the result from the Principal Component
332 Analysis.
333 2. Determine \f$\mathbf{P}^\prime\f$ as before.
334 3. Determine \f$\mathbf{x}^\prime\f$ as before.
335 4. Use the result of the fit to determine \f$\mathbf{Q}^\prime\f$.
336 5. Transform back to \f$\mathbf{P}\f$ from \f$\mathbf{Q}^\prime\f$, using
337 the result from the Principal Component Analysis.
338
339 ## Testing the parameterization
340 The class also provides functionality for testing the, over the
341 training sample, found parameterization (TMultiDimFit::Fit). This is done by passing
342 the class a test sample of \f$ M_t\f$ tuples of the form
343 \f$(\mathbf{x}_{t,j},D_{t,j}, E_{t,j})\f$, where \f$\mathbf{x}_{t,j}\f$ are the independent
344 variables, \f$ D_{t,j}\f$ the known, dependent quantity, and \f$ E_{t,j}\f$ is
345 the square error in \f$ D_{t,j}\f$ (TMultiDimFit::AddTestRow).
346
347 The parameterization is then evaluated at every \f$\mathbf{x}_t\f$ in the
348 test sample, and
349 \f[
350 S_t \equiv \sum_{j=1}^{M_t} \left(D_{t,j} -
351 D_p\left(\mathbf{x}_{t,j}\right)\right)^2
352 \f]
353 is evaluated. The relative error over the test sample
354 \f[
355 R_t = \frac{S_t}{\sum_{j=1}^{M_t} D_{t,j}^2}
356 \f]
357 should not be to low or high compared to \f$ R\f$ from the training
358 sample. Also, multiple correlation coefficient from both samples should
359 be fairly close, otherwise one of the samples is not representative of
360 the problem. A large difference in the reduced \f$ \chi^2\f$ over the two
361 samples indicate an over fit, and the maximum number of terms in the
362 parameterisation should be reduced.
363
364 It's possible to use [4] to further improve the fit, using the test sample.
365
366 Christian Holm
367
368 ## Bibliography
369 - <a name="bevington"></a> Philip R. Bevington and D. Keith Robinson. *Data Reduction and Error Analysis for
370 the Physical Sciences*. McGraw-Hill, 2 edition, 1992.
371 - <a name="mudifi"></a> R. Brun et al. *Long writeup DD/75-23*, CERN, 1980.
372 - Gene H. Golub and Charles F. van Loan. *Matrix Computations*.
373 John Hopkins University Press, Baltimore, 3 edition, 1996.
374 - <a name="minuit"></a>F. James. *Minuit*. Long writeup D506, CERN, 1998.
375 - <a name="wind72"></a>H. Wind. *Function parameterization*. Proceedings of the 1972 CERN Computing and Data Processing
376 School, volume 72-21 of Yellow report. CERN, 1972.
377 - <a name="wind81"></a>H. Wind. 1. principal component analysis, 2. pattern recognition for track
378 finding, 3. interpolation and functional representation. Yellow report EP/81-12, CERN, 1981.
379
380[1]: classTMultiDimFit.html#bevington
381[2]: classTMultiDimFit.html#mudifi
382[4]: classTMultiDimFit.html#minuit
383[5]: classTMultiDimFit.html#wind72
384[6]: classTMultiDimFit.html#wind81
385[9]: classTMultiDimFit.html#eq:dS2
386*/
387
388
389#include "Riostream.h"
390#include "TMultiDimFit.h"
391#include "TMath.h"
392#include "TH1.h"
393#include "TH2.h"
394#include "TROOT.h"
395#include "TBrowser.h"
396#include "TDecompChol.h"
397#include "TDatime.h"
398
399
400#define RADDEG (180. / TMath::Pi())
401#define DEGRAD (TMath::Pi() / 180.)
402#define HIST_XORIG 0
403#define HIST_DORIG 1
404#define HIST_XNORM 2
405#define HIST_DSHIF 3
406#define HIST_RX 4
407#define HIST_RD 5
408#define HIST_RTRAI 6
409#define HIST_RTEST 7
410#define PARAM_MAXSTUDY 1
411#define PARAM_SEVERAL 2
412#define PARAM_RELERR 3
413#define PARAM_MAXTERMS 4
414
415
416////////////////////////////////////////////////////////////////////////////////
417
418static void mdfHelper(int&, double*, double&, double*, int);
419
420////////////////////////////////////////////////////////////////////////////////
421
423
424//____________________________________________________________________
425// Static instance. Used with mdfHelper and TMinuit
427
428
429////////////////////////////////////////////////////////////////////////////////
430/// Empty CTOR. Do not use
431
433{
434 fMeanQuantity = 0;
435 fMaxQuantity = 0;
436 fMinQuantity = 0;
437 fSumSqQuantity = 0;
439
440 fNVariables = 0;
441 fSampleSize = 0;
442 fTestSampleSize = 0;
443
444 fMinAngle = 1;
445 fMaxAngle = 0;
446 fMaxTerms = 0;
448 fMaxPowers = 0;
449 fPowerLimit = 0;
450
451 fMaxFunctions = 0;
452 fFunctionCodes = 0;
453 fMaxStudy = 0;
454 fMaxFuncNV = 0;
455
456 fMaxPowersFinal = 0;
457 fPowers = 0;
458 fPowerIndex = 0;
459
460 fMaxResidual = 0;
461 fMinResidual = 0;
462 fMaxResidualRow = 0;
463 fMinResidualRow = 0;
464 fSumSqResidual = 0;
465
466 fNCoefficients = 0;
467 fRMS = 0;
468 fChi2 = 0;
470
471 fError = 0;
472 fTestError = 0;
473 fPrecision = 0;
474 fTestPrecision = 0;
477
478 fHistograms = 0;
479 fHistogramMask = 0;
480 fBinVarX = 100;
481 fBinVarY = 100;
482
483 fFitter = 0;
488
489}
490
491
492////////////////////////////////////////////////////////////////////////////////
493/// Constructor
494/// Second argument is the type of polynomials to use in
495/// parameterisation, one of:
496/// TMultiDimFit::kMonomials
497/// TMultiDimFit::kChebyshev
498/// TMultiDimFit::kLegendre
499///
500/// Options:
501/// K Compute (k)correlation matrix
502/// V Be verbose
503///
504/// Default is no options.
505///
506
509 Option_t *option)
510: TNamed("multidimfit","Multi-dimensional fit object"),
511fQuantity(dimension),
512fSqError(dimension),
513fVariables(dimension*100),
514fMeanVariables(dimension),
515fMaxVariables(dimension),
516fMinVariables(dimension)
517{
518 fgInstance = this;
519
520 fMeanQuantity = 0;
521 fMaxQuantity = 0;
522 fMinQuantity = 0;
523 fSumSqQuantity = 0;
525
526 fNVariables = dimension;
527 fSampleSize = 0;
528 fTestSampleSize = 0;
529
530 fMinAngle = 1;
531 fMaxAngle = 0;
532 fMaxTerms = 0;
533 fMinRelativeError = 0.01;
534 fMaxPowers = new Int_t[dimension];
535 fPowerLimit = 1;
536
537 fMaxFunctions = 0;
538 fFunctionCodes = 0;
539 fMaxStudy = 0;
540 fMaxFuncNV = 0;
541
542 fMaxPowersFinal = new Int_t[dimension];
543 fPowers = 0;
544 fPowerIndex = 0;
545
546 fMaxResidual = 0;
547 fMinResidual = 0;
548 fMaxResidualRow = 0;
549 fMinResidualRow = 0;
550 fSumSqResidual = 0;
551
552 fNCoefficients = 0;
553 fRMS = 0;
554 fChi2 = 0;
556
557 fError = 0;
558 fTestError = 0;
559 fPrecision = 0;
560 fTestPrecision = 0;
563
564 fHistograms = 0;
565 fHistogramMask = 0;
566 fBinVarX = 100;
567 fBinVarY = 100;
568
569 fFitter = 0;
570 fPolyType = type;
574 TString opt = option;
575 opt.ToLower();
576
577 if (opt.Contains("k")) fShowCorrelation = kTRUE;
578 if (opt.Contains("v")) fIsVerbose = kTRUE;
579}
580
581
582////////////////////////////////////////////////////////////////////////////////
583/// Destructor
584
586{
587 delete [] fPowers;
588 delete [] fMaxPowers;
589 delete [] fMaxPowersFinal;
590 delete [] fPowerIndex;
591 delete [] fFunctionCodes;
592 if (fHistograms) fHistograms->Clear("nodelete");
593 delete fHistograms;
594}
595
596
597////////////////////////////////////////////////////////////////////////////////
598/// Add a row consisting of fNVariables independent variables, the
599/// known, dependent quantity, and optionally, the square error in
600/// the dependent quantity, to the training sample to be used for the
601/// parameterization.
602/// The mean of the variables and quantity is calculated on the fly,
603/// as outlined in TPrincipal::AddRow.
604/// This sample should be representative of the problem at hand.
605/// Please note, that if no error is given Poisson statistics is
606/// assumed and the square error is set to the value of dependent
607/// quantity. See also the
608/// <a href="#TMultiDimFit:description">class description</a>
609
611{
612 if (!x)
613 return;
614
615 if (++fSampleSize == 1) {
616 fMeanQuantity = D;
617 fMaxQuantity = D;
618 fMinQuantity = D;
619 fSumSqQuantity = D * D;// G.Q. erratum on August 15th, 2008
620 }
621 else {
624 fSumSqQuantity += D * D;
625
626 if (D >= fMaxQuantity) fMaxQuantity = D;
627 if (D <= fMinQuantity) fMinQuantity = D;
628 }
629
630
631 // If the vector isn't big enough to hold the new data, then
632 // expand the vector by half it's size.
633 Int_t size = fQuantity.GetNrows();
634 if (fSampleSize > size) {
635 fQuantity.ResizeTo(size + size/2);
636 fSqError.ResizeTo(size + size/2);
637 }
638
639 // Store the value
640 fQuantity(fSampleSize-1) = D;
641 fSqError(fSampleSize-1) = (E == 0 ? D : E);
642
643 // Store data point in internal vector
644 // If the vector isn't big enough to hold the new data, then
645 // expand the vector by half it's size
646 size = fVariables.GetNrows();
647 if (fSampleSize * fNVariables > size)
648 fVariables.ResizeTo(size + size/2);
649
650
651 // Increment the data point counter
652 Int_t i,j;
653 for (i = 0; i < fNVariables; i++) {
654 if (fSampleSize == 1) {
655 fMeanVariables(i) = x[i];
656 fMaxVariables(i) = x[i];
657 fMinVariables(i) = x[i];
658 }
659 else {
660 fMeanVariables(i) *= 1 - 1./Double_t(fSampleSize);
662
663 // Update the maximum value for this component
664 if (x[i] >= fMaxVariables(i)) fMaxVariables(i) = x[i];
665
666 // Update the minimum value for this component
667 if (x[i] <= fMinVariables(i)) fMinVariables(i) = x[i];
668
669 }
670
671 // Store the data.
672 j = (fSampleSize-1) * fNVariables + i;
673 fVariables(j) = x[i];
674 }
675}
676
677
678////////////////////////////////////////////////////////////////////////////////
679/// Add a row consisting of fNVariables independent variables, the
680/// known, dependent quantity, and optionally, the square error in
681/// the dependent quantity, to the test sample to be used for the
682/// test of the parameterization.
683/// This sample needn't be representative of the problem at hand.
684/// Please note, that if no error is given Poisson statistics is
685/// assumed and the square error is set to the value of dependent
686/// quantity. See also the
687/// <a href="#TMultiDimFit:description">class description</a>
688
690{
691 if (fTestSampleSize++ == 0) {
695 }
696
697 // If the vector isn't big enough to hold the new data, then
698 // expand the vector by half it's size.
700 if (fTestSampleSize > size) {
701 fTestQuantity.ResizeTo(size + size/2);
702 fTestSqError.ResizeTo(size + size/2);
703 }
704
705 // Store the value
707 fTestSqError(fTestSampleSize-1) = (E == 0 ? D : E);
708
709 // Store data point in internal vector
710 // If the vector isn't big enough to hold the new data, then
711 // expand the vector by half it's size
712 size = fTestVariables.GetNrows();
713 if (fTestSampleSize * fNVariables > size)
714 fTestVariables.ResizeTo(size + size/2);
715
716
717 // Increment the data point counter
718 Int_t i,j;
719 for (i = 0; i < fNVariables; i++) {
720 j = fNVariables * (fTestSampleSize - 1) + i;
721 fTestVariables(j) = x[i];
722
723 if (x[i] > fMaxVariables(i))
724 Warning("AddTestRow", "variable %d (row: %d) too large: %f > %f",
725 i, fTestSampleSize, x[i], fMaxVariables(i));
726 if (x[i] < fMinVariables(i))
727 Warning("AddTestRow", "variable %d (row: %d) too small: %f < %f",
728 i, fTestSampleSize, x[i], fMinVariables(i));
729 }
730}
731
732
733////////////////////////////////////////////////////////////////////////////////
734/// Browse the TMultiDimFit object in the TBrowser.
735
737{
738 if (fHistograms) {
739 TIter next(fHistograms);
740 TH1* h = 0;
741 while ((h = (TH1*)next()))
742 b->Add(h,h->GetName());
743 }
744 if (fVariables.IsValid())
745 b->Add(&fVariables, "Variables (Training)");
746 if (fQuantity.IsValid())
747 b->Add(&fQuantity, "Quantity (Training)");
748 if (fSqError.IsValid())
749 b->Add(&fSqError, "Error (Training)");
751 b->Add(&fMeanVariables, "Mean of Variables (Training)");
753 b->Add(&fMaxVariables, "Mean of Variables (Training)");
755 b->Add(&fMinVariables, "Min of Variables (Training)");
757 b->Add(&fTestVariables, "Variables (Test)");
759 b->Add(&fTestQuantity, "Quantity (Test)");
760 if (fTestSqError.IsValid())
761 b->Add(&fTestSqError, "Error (Test)");
762 if (fFunctions.IsValid())
763 b->Add(&fFunctions, "Functions");
765 b->Add(&fCoefficients,"Coefficients");
767 b->Add(&fCoefficientsRMS,"Coefficients Errors");
769 b->Add(&fOrthFunctions, "Orthogonal Functions");
771 b->Add(&fOrthFunctionNorms, "Orthogonal Functions Norms");
772 if (fResiduals.IsValid())
773 b->Add(&fResiduals, "Residuals");
775 b->Add(&fOrthCoefficients,"Orthogonal Coefficients");
777 b->Add(&fOrthCurvatureMatrix,"Orthogonal curvature matrix");
779 b->Add(&fCorrelationMatrix,"Correlation Matrix");
780 if (fFitter)
781 b->Add(fFitter, fFitter->GetName());
782}
783
784
785////////////////////////////////////////////////////////////////////////////////
786/// Clear internal structures and variables
787
789{
790 Int_t i, j, n = fNVariables, m = fMaxFunctions;
791
792 // Training sample, dependent quantity
793 fQuantity.Zero();
794 fSqError.Zero();
795 fMeanQuantity = 0;
796 fMaxQuantity = 0;
797 fMinQuantity = 0;
798 fSumSqQuantity = 0;
800
801 // Training sample, independent variables
803 fNVariables = 0;
804 fSampleSize = 0;
808
809 // Test sample
813 fTestSampleSize = 0;
814
815 // Functions
817 //for (i = 0; i < fMaxTerms; i++) fPowerIndex[i] = 0;
818 //for (i = 0; i < fMaxTerms; i++) fFunctionCodes[i] = 0;
819 fMaxFunctions = 0;
820 fMaxStudy = 0;
823
824 // Control parameters
826 fMinAngle = 0;
827 fMaxAngle = 0;
828 fMaxTerms = 0;
829
830 // Powers
831 for (i = 0; i < n; i++) {
832 fMaxPowers[i] = 0;
833 fMaxPowersFinal[i] = 0;
834 for (j = 0; j < m; j++)
835 fPowers[i * n + j] = 0;
836 }
837 fPowerLimit = 0;
838
839 // Residuals
840 fMaxResidual = 0;
841 fMinResidual = 0;
842 fMaxResidualRow = 0;
843 fMinResidualRow = 0;
844 fSumSqResidual = 0;
845
846 // Fit
847 fNCoefficients = 0;
850 fRMS = 0;
852 fError = 0;
853 fTestError = 0;
854 fPrecision = 0;
855 fTestPrecision = 0;
856
857 // Coefficients
861 fHistograms->Clear(option);
862
863 // Options
867}
868
869
870////////////////////////////////////////////////////////////////////////////////
871/// Evaluate parameterization at point x. Optional argument coeff is
872/// a vector of coefficients for the parameterisation, fNCoefficients
873/// elements long.
874
875Double_t TMultiDimFit::Eval(const Double_t *x, const Double_t* coeff) const
876{
877 Double_t returnValue = fMeanQuantity;
878 Double_t term = 0;
879 Int_t i, j;
880
881 for (i = 0; i < fNCoefficients; i++) {
882 // Evaluate the ith term in the expansion
883 term = (coeff ? coeff[i] : fCoefficients(i));
884 for (j = 0; j < fNVariables; j++) {
885 // Evaluate the factor (polynomial) in the j-th variable.
886 Int_t p = fPowers[fPowerIndex[i] * fNVariables + j];
887 Double_t y = 1 + 2. / (fMaxVariables(j) - fMinVariables(j))
888 * (x[j] - fMaxVariables(j));
889 term *= EvalFactor(p,y);
890 }
891 // Add this term to the final result
892 returnValue += term;
893 }
894 return returnValue;
895}
896
897
898////////////////////////////////////////////////////////////////////////////////
899/// Evaluate parameterization error at point x. Optional argument coeff is
900/// a vector of coefficients for the parameterisation, fNCoefficients
901/// elements long.
902
904{
905 Double_t returnValue = 0;
906 Double_t term = 0;
907 Int_t i, j;
908
909 for (i = 0; i < fNCoefficients; i++) {
910 // std::cout << "Error coef " << i << " -> " << fCoefficientsRMS(i) << std::endl;
911 }
912 for (i = 0; i < fNCoefficients; i++) {
913 // Evaluate the ith term in the expansion
914 term = (coeff ? coeff[i] : fCoefficientsRMS(i));
915 for (j = 0; j < fNVariables; j++) {
916 // Evaluate the factor (polynomial) in the j-th variable.
917 Int_t p = fPowers[fPowerIndex[i] * fNVariables + j];
918 Double_t y = 1 + 2. / (fMaxVariables(j) - fMinVariables(j))
919 * (x[j] - fMaxVariables(j));
920 term *= EvalFactor(p,y);
921 // std::cout << "i,j " << i << ", " << j << " " << p << " " << y << " " << EvalFactor(p,y) << " " << term << std::endl;
922 }
923 // Add this term to the final result
924 returnValue += term*term;
925 // std::cout << " i = " << i << " value = " << returnValue << std::endl;
926 }
927 returnValue = sqrt(returnValue);
928 return returnValue;
929}
930
931
932////////////////////////////////////////////////////////////////////////////////
933/// PRIVATE METHOD:
934/// Calculate the control parameter from the passed powers
935
937{
938 Double_t s = 0;
939 Double_t epsilon = 1e-6; // a small number
940 for (Int_t i = 0; i < fNVariables; i++) {
941 if (fMaxPowers[i] != 1)
942 s += (epsilon + iv[i] - 1) / (epsilon + fMaxPowers[i] - 1);
943 }
944 return s;
945}
946
947////////////////////////////////////////////////////////////////////////////////
948/// PRIVATE METHOD:
949/// Evaluate function with power p at variable value x
950
952{
953 Int_t i = 0;
954 Double_t p1 = 1;
955 Double_t p2 = 0;
956 Double_t p3 = 0;
957 Double_t r = 0;
958
959 switch(p) {
960 case 1:
961 r = 1;
962 break;
963 case 2:
964 r = x;
965 break;
966 default:
967 p2 = x;
968 for (i = 3; i <= p; i++) {
969 p3 = p2 * x;
970 if (fPolyType == kLegendre)
971 p3 = ((2 * i - 3) * p2 * x - (i - 2) * p1) / (i - 1);
972 else if (fPolyType == kChebyshev)
973 p3 = 2 * x * p2 - p1;
974 p1 = p2;
975 p2 = p3;
976 }
977 r = p3;
978 }
979
980 return r;
981}
982
983
984////////////////////////////////////////////////////////////////////////////////
985/// Find the parameterization
986///
987/// Options:
988/// None so far
989///
990/// For detailed description of what this entails, please refer to the
991/// <a href="#TMultiDimFit:description">class description</a>
992
994{
1001}
1002
1003////////////////////////////////////////////////////////////////////////////////
1004/// Try to fit the found parameterisation to the test sample.
1005///
1006/// Options
1007/// M use Minuit to improve coefficients
1008///
1009/// Also, refer to
1010/// <a href="#TMultiDimFit:description">class description</a>
1011
1013{
1014 Int_t i, j;
1016 Double_t sumSqD = 0;
1017 Double_t sumD = 0;
1018 Double_t sumSqR = 0;
1019 Double_t sumR = 0;
1020
1021 // Calculate the residuals over the test sample
1022 for (i = 0; i < fTestSampleSize; i++) {
1023 for (j = 0; j < fNVariables; j++)
1024 x[j] = fTestVariables(i * fNVariables + j);
1025 Double_t res = fTestQuantity(i) - Eval(x);
1026 sumD += fTestQuantity(i);
1027 sumSqD += fTestQuantity(i) * fTestQuantity(i);
1028 sumR += res;
1029 sumSqR += res * res;
1031 ((TH1D*)fHistograms->FindObject("res_test"))->Fill(res);
1032 }
1033 Double_t dAvg = sumSqD - (sumD * sumD) / fTestSampleSize;
1034 Double_t rAvg = sumSqR - (sumR * sumR) / fTestSampleSize;
1035 fTestCorrelationCoeff = (dAvg - rAvg) / dAvg;
1036 fTestError = sumSqR;
1037 fTestPrecision = sumSqR / sumSqD;
1038
1039 TString opt(option);
1040 opt.ToLower();
1041
1042 if (!opt.Contains("m"))
1043 MakeChi2();
1044
1046 Warning("Fit", "test sample is very small");
1047
1048 if (!opt.Contains("m")) {
1049 Error("Fit", "invalid option");
1050 delete [] x;
1051 return;
1052 }
1053
1055 if (!fFitter) {
1056 Error("Fit", "Cannot create Fitter");
1057 delete [] x;
1058 return;
1059 }
1061
1062 const Int_t maxArgs = 16;
1063 Int_t args = 1;
1064 Double_t* arglist = new Double_t[maxArgs];
1065 arglist[0] = -1;
1066 fFitter->ExecuteCommand("SET PRINT",arglist,args);
1067
1068 for (i = 0; i < fNCoefficients; i++) {
1069 Double_t startVal = fCoefficients(i);
1070 Double_t startErr = fCoefficientsRMS(i);
1071 fFitter->SetParameter(i, Form("coeff%02d",i),
1072 startVal, startErr, 0, 0);
1073 }
1074
1075 // arglist[0] = 0;
1076 args = 1;
1077 // fFitter->ExecuteCommand("SET PRINT",arglist,args);
1078 fFitter->ExecuteCommand("MIGRAD",arglist,args);
1079
1080 for (i = 0; i < fNCoefficients; i++) {
1081 Double_t val = 0, err = 0, low = 0, high = 0;
1082 fFitter->GetParameter(i, Form("coeff%02d",i),
1083 val, err, low, high);
1084 fCoefficients(i) = val;
1085 fCoefficientsRMS(i) = err;
1086 }
1087 delete [] x;
1088}
1089
1090////////////////////////////////////////////////////////////////////////////////
1091/// Return the static instance.
1092
1094{
1095 return fgInstance;
1096}
1097
1098////////////////////////////////////////////////////////////////////////////////
1099/// PRIVATE METHOD:
1100/// Create list of candidate functions for the parameterisation. See
1101/// also
1102/// <a href="#TMultiDimFit:description">class description</a>
1103
1105{
1106 Int_t i = 0;
1107 Int_t j = 0;
1108 Int_t k = 0;
1109
1110 // The temporary array to store the powers in. We don't need to
1111 // initialize this array however.
1113 Int_t *powers = new Int_t[fMaxFuncNV];
1114
1115 // store of `control variables'
1116 Double_t* control = new Double_t[fMaxFunctions];
1117
1118 // We've better initialize the variables
1119 Int_t *iv = new Int_t[fNVariables];
1120 for (i = 0; i < fNVariables; i++)
1121 iv[i] = 1;
1122
1123 if (!fIsUserFunction) {
1124
1125 // Number of funcs selected
1126 Int_t numberFunctions = 0;
1127
1128 // Absolute max number of functions
1129 Int_t maxNumberFunctions = 1;
1130 for (i = 0; i < fNVariables; i++)
1131 maxNumberFunctions *= fMaxPowers[i];
1132
1133 while (kTRUE) {
1134 // Get the control value for this function
1135 Double_t s = EvalControl(iv);
1136
1137 if (s <= fPowerLimit) {
1138
1139 // Call over-loadable method Select, as to allow the user to
1140 // interfere with the selection of functions.
1141 if (Select(iv)) {
1142 numberFunctions++;
1143
1144 // If we've reached the user defined limit of how many
1145 // functions we can consider, break out of the loop
1146 if (numberFunctions > fMaxFunctions)
1147 break;
1148
1149 // Store the control value, so we can sort array of powers
1150 // later on
1151 control[numberFunctions-1] = Int_t(1.0e+6*s);
1152
1153 // Store the powers in powers array.
1154 for (i = 0; i < fNVariables; i++) {
1155 j = (numberFunctions - 1) * fNVariables + i;
1156 powers[j] = iv[i];
1157 }
1158 } // if (Select())
1159 } // if (s <= fPowerLimit)
1160
1161 for (i = 0; i < fNVariables; i++)
1162 if (iv[i] < fMaxPowers[i])
1163 break;
1164
1165 // If all variables have reached their maximum power, then we
1166 // break out of the loop
1167 if (i == fNVariables) {
1168 fMaxFunctions = numberFunctions;
1169 break;
1170 }
1171
1172 // Next power in variable i
1173 if (i < fNVariables) iv[i]++;
1174
1175 for (j = 0; j < i; j++)
1176 iv[j] = 1;
1177 } // while (kTRUE)
1178 }
1179 else {
1180 // In case the user gave an explicit function
1181 for (i = 0; i < fMaxFunctions; i++) {
1182 // Copy the powers to working arrays
1183 for (j = 0; j < fNVariables; j++) {
1184 powers[i * fNVariables + j] = fPowers[i * fNVariables + j];
1185 iv[j] = fPowers[i * fNVariables + j];
1186 }
1187
1188 control[i] = Int_t(1.0e+6*EvalControl(iv));
1189 }
1190 }
1191
1192 // Now we need to sort the powers according to least `control
1193 // variable'
1194 Int_t *order = new Int_t[fMaxFunctions];
1195 for (i = 0; i < fMaxFunctions; i++)
1196 order[i] = i;
1198 fPowers = new Int_t[fMaxFuncNV];
1199
1200 for (i = 0; i < fMaxFunctions; i++) {
1201 Double_t x = control[i];
1202 Int_t l = order[i];
1203 k = i;
1204
1205 for (j = i; j < fMaxFunctions; j++) {
1206 if (control[j] <= x) {
1207 x = control[j];
1208 l = order[j];
1209 k = j;
1210 }
1211 }
1212
1213 if (k != i) {
1214 control[k] = control[i];
1215 control[i] = x;
1216 order[k] = order[i];
1217 order[i] = l;
1218 }
1219 }
1220
1221 for (i = 0; i < fMaxFunctions; i++)
1222 for (j = 0; j < fNVariables; j++)
1223 fPowers[i * fNVariables + j] = powers[order[i] * fNVariables + j];
1224
1225 delete [] control;
1226 delete [] powers;
1227 delete [] order;
1228 delete [] iv;
1229}
1230
1231
1232////////////////////////////////////////////////////////////////////////////////
1233/// Calculate Chi square over either the test sample. The optional
1234/// argument coeff is a vector of coefficients to use in the
1235/// evaluation of the parameterisation. If coeff == 0, then the found
1236/// coefficients is used.
1237/// Used my MINUIT for fit (see TMultDimFit::Fit)
1238
1240{
1241 fChi2 = 0;
1242 Int_t i, j;
1244 for (i = 0; i < fTestSampleSize; i++) {
1245 // Get the stored point
1246 for (j = 0; j < fNVariables; j++)
1247 x[j] = fTestVariables(i * fNVariables + j);
1248
1249 // Evaluate function. Scale to shifted values
1250 Double_t f = Eval(x,coeff);
1251
1252 // Calculate contribution to Chic square
1253 fChi2 += 1. / TMath::Max(fTestSqError(i),1e-20)
1254 * (fTestQuantity(i) - f) * (fTestQuantity(i) - f);
1255 }
1256
1257 // Clean up
1258 delete [] x;
1259
1260 return fChi2;
1261}
1262
1263
1264////////////////////////////////////////////////////////////////////////////////
1265/// Generate the file <filename> with .C appended if argument doesn't
1266/// end in .cxx or .C. The contains the implementation of the
1267/// function:
1268///
1269/// Double_t <funcname>(Double_t *x)
1270///
1271/// which does the same as TMultiDimFit::Eval. Please refer to this
1272/// method.
1273///
1274/// Further, the static variables:
1275///
1276/// Int_t gNVariables
1277/// Int_t gNCoefficients
1278/// Double_t gDMean
1279/// Double_t gXMean[]
1280/// Double_t gXMin[]
1281/// Double_t gXMax[]
1282/// Double_t gCoefficient[]
1283/// Int_t gPower[]
1284///
1285/// are initialized. The only ROOT header file needed is Rtypes.h
1286///
1287/// See TMultiDimFit::MakeRealCode for a list of options
1288
1289void TMultiDimFit::MakeCode(const char* filename, Option_t *option)
1290{
1291
1292 TString outName(filename);
1293 if (!outName.EndsWith(".C") && !outName.EndsWith(".cxx"))
1294 outName += ".C";
1295
1296 MakeRealCode(outName.Data(),"",option);
1297}
1298
1299
1300
1301////////////////////////////////////////////////////////////////////////////////
1302/// PRIVATE METHOD:
1303/// Compute the errors on the coefficients. For this to be done, the
1304/// curvature matrix of the non-orthogonal functions, is computed.
1305
1307{
1308 Int_t i = 0;
1309 Int_t j = 0;
1310 Int_t k = 0;
1314
1315 TMatrixDSym curvatureMatrix(fNCoefficients);
1316
1317 // Build the curvature matrix
1318 for (i = 0; i < fNCoefficients; i++) {
1319 iF = TMatrixDRow(fFunctions,i);
1320 for (j = 0; j <= i; j++) {
1321 jF = TMatrixDRow(fFunctions,j);
1322 for (k = 0; k < fSampleSize; k++)
1323 curvatureMatrix(i,j) +=
1324 1 / TMath::Max(fSqError(k), 1e-20) * iF(k) * jF(k);
1325 curvatureMatrix(j,i) = curvatureMatrix(i,j);
1326 }
1327 }
1328
1329 // Calculate Chi Square
1330 fChi2 = 0;
1331 for (i = 0; i < fSampleSize; i++) {
1332 Double_t f = 0;
1333 for (j = 0; j < fNCoefficients; j++)
1334 f += fCoefficients(j) * fFunctions(j,i);
1335 fChi2 += 1. / TMath::Max(fSqError(i),1e-20) * (fQuantity(i) - f)
1336 * (fQuantity(i) - f);
1337 }
1338
1339 // Invert the curvature matrix
1340 const TVectorD diag = TMatrixDDiag_const(curvatureMatrix);
1341 curvatureMatrix.NormByDiag(diag);
1342
1343 TDecompChol chol(curvatureMatrix);
1344 if (!chol.Decompose())
1345 Error("MakeCoefficientErrors", "curvature matrix is singular");
1346 chol.Invert(curvatureMatrix);
1347
1348 curvatureMatrix.NormByDiag(diag);
1349
1350 for (i = 0; i < fNCoefficients; i++)
1351 fCoefficientsRMS(i) = TMath::Sqrt(curvatureMatrix(i,i));
1352}
1353
1354
1355////////////////////////////////////////////////////////////////////////////////
1356/// PRIVATE METHOD:
1357/// Invert the model matrix B, and compute final coefficients. For a
1358/// more thorough discussion of what this means, please refer to the
1359/// <a href="#TMultiDimFit:description">class description</a>
1360///
1361/// First we invert the lower triangle matrix fOrthCurvatureMatrix
1362/// and store the inverted matrix in the upper triangle.
1363
1365{
1366 Int_t i = 0, j = 0;
1367 Int_t col = 0, row = 0;
1368
1369 // Invert the B matrix
1370 for (col = 1; col < fNCoefficients; col++) {
1371 for (row = col - 1; row > -1; row--) {
1372 fOrthCurvatureMatrix(row,col) = 0;
1373 for (i = row; i <= col ; i++)
1374 fOrthCurvatureMatrix(row,col) -=
1376 * fOrthCurvatureMatrix(i,col);
1377 }
1378 }
1379
1380 // Compute the final coefficients
1382
1383 for (i = 0; i < fNCoefficients; i++) {
1384 Double_t sum = 0;
1385 for (j = i; j < fNCoefficients; j++)
1387 fCoefficients(i) = sum;
1388 }
1389
1390 // Compute the final residuals
1392 for (i = 0; i < fSampleSize; i++)
1393 fResiduals(i) = fQuantity(i);
1394
1395 for (i = 0; i < fNCoefficients; i++)
1396 for (j = 0; j < fSampleSize; j++)
1397 fResiduals(j) -= fCoefficients(i) * fFunctions(i,j);
1398
1399 // Compute the max and minimum, and squared sum of the evaluated
1400 // residuals
1401 fMinResidual = 10e10;
1402 fMaxResidual = -10e10;
1403 Double_t sqRes = 0;
1404 for (i = 0; i < fSampleSize; i++){
1405 sqRes += fResiduals(i) * fResiduals(i);
1406 if (fResiduals(i) <= fMinResidual) {
1408 fMinResidualRow = i;
1409 }
1410 if (fResiduals(i) >= fMaxResidual) {
1412 fMaxResidualRow = i;
1413 }
1414 }
1415
1418
1419 // If we use histograms, fill some more
1423 for (i = 0; i < fSampleSize; i++) {
1425 ((TH2D*)fHistograms->FindObject("res_d"))->Fill(fQuantity(i),
1426 fResiduals(i));
1428 ((TH1D*)fHistograms->FindObject("res_train"))->Fill(fResiduals(i));
1429
1431 for (j = 0; j < fNVariables; j++)
1432 ((TH2D*)fHistograms->FindObject(Form("res_x_%d",j)))
1433 ->Fill(fVariables(i * fNVariables + j),fResiduals(i));
1434 }
1435 } // If histograms
1436
1437}
1438
1439
1440////////////////////////////////////////////////////////////////////////////////
1441/// PRIVATE METHOD:
1442/// Compute the correlation matrix
1443
1445{
1446 if (!fShowCorrelation)
1447 return;
1448
1450
1451 Double_t d2 = 0;
1452 Double_t ddotXi = 0; // G.Q. needs to be reinitialized in the loop over i fNVariables
1453 Double_t xiNorm = 0; // G.Q. needs to be reinitialized in the loop over i fNVariables
1454 Double_t xidotXj = 0; // G.Q. needs to be reinitialized in the loop over j fNVariables
1455 Double_t xjNorm = 0; // G.Q. needs to be reinitialized in the loop over j fNVariables
1456
1457 Int_t i, j, k, l, m; // G.Q. added m variable
1458 for (i = 0; i < fSampleSize; i++)
1459 d2 += fQuantity(i) * fQuantity(i);
1460
1461 for (i = 0; i < fNVariables; i++) {
1462 ddotXi = 0.; // G.Q. reinitialisation
1463 xiNorm = 0.; // G.Q. reinitialisation
1464 for (j = 0; j< fSampleSize; j++) {
1465 // Index of sample j of variable i
1466 k = j * fNVariables + i;
1467 ddotXi += fQuantity(j) * (fVariables(k) - fMeanVariables(i));
1468 xiNorm += (fVariables(k) - fMeanVariables(i))
1469 * (fVariables(k) - fMeanVariables(i));
1470 }
1471 fCorrelationMatrix(i,0) = ddotXi / TMath::Sqrt(d2 * xiNorm);
1472
1473 for (j = 0; j < i; j++) {
1474 xidotXj = 0.; // G.Q. reinitialisation
1475 xjNorm = 0.; // G.Q. reinitialisation
1476 for (k = 0; k < fSampleSize; k++) {
1477 // Index of sample j of variable i
1478 // l = j * fNVariables + k; // G.Q.
1479 l = k * fNVariables + j; // G.Q.
1480 m = k * fNVariables + i; // G.Q.
1481 // G.Q. xidotXj += (fVariables(i) - fMeanVariables(i))
1482 // G.Q. * (fVariables(l) - fMeanVariables(j));
1483 xidotXj += (fVariables(m) - fMeanVariables(i))
1484 * (fVariables(l) - fMeanVariables(j)); // G.Q. modified index for Xi
1485 xjNorm += (fVariables(l) - fMeanVariables(j))
1486 * (fVariables(l) - fMeanVariables(j));
1487 }
1488 //fCorrelationMatrix(i+1,j) = xidotXj / TMath::Sqrt(xiNorm * xjNorm);
1489 fCorrelationMatrix(i,j+1) = xidotXj / TMath::Sqrt(xiNorm * xjNorm);
1490 }
1491 }
1492}
1493
1494
1495
1496////////////////////////////////////////////////////////////////////////////////
1497/// PRIVATE METHOD:
1498/// Make Gram-Schmidt orthogonalisation. The class description gives
1499/// a thorough account of this algorithm, as well as
1500/// references. Please refer to the
1501/// <a href="#TMultiDimFit:description">class description</a>
1502
1504{
1505
1506 // calculate w_i, that is, evaluate the current function at data
1507 // point i
1508 Double_t f2 = 0;
1511 Int_t j = 0;
1512 Int_t k = 0;
1513
1514 for (j = 0; j < fSampleSize; j++) {
1515 fFunctions(fNCoefficients, j) = 1;
1517 // First, however, we need to calculate f_fNCoefficients
1518 for (k = 0; k < fNVariables; k++) {
1519 Int_t p = fPowers[function * fNVariables + k];
1520 Double_t x = fVariables(j * fNVariables + k);
1522 }
1523
1524 // Calculate f dot f in f2
1526 // Assign to w_fNCoefficients f_fNCoefficients
1528 }
1529
1530 // the first column of w is equal to f
1531 for (j = 0; j < fNCoefficients; j++) {
1532 Double_t fdw = 0;
1533 // Calculate (f_fNCoefficients dot w_j) / w_j^2
1534 for (k = 0; k < fSampleSize; k++) {
1535 fdw += fFunctions(fNCoefficients, k) * fOrthFunctions(j,k)
1536 / fOrthFunctionNorms(j);
1537 }
1538
1540 // and subtract it from the current value of w_ij
1541 for (k = 0; k < fSampleSize; k++)
1543 }
1544
1545 for (j = 0; j < fSampleSize; j++) {
1546 // calculate squared length of w_fNCoefficients
1550
1551 // calculate D dot w_fNCoefficients in A
1554 }
1555
1556 // First test, but only if didn't user specify
1557 if (!fIsUserFunction)
1560 return 0;
1561
1562 // The result found by this code for the first residual is always
1563 // much less then the one found be MUDIFI. That's because it's
1564 // supposed to be. The cause is the improved precision of Double_t
1565 // over DOUBLE PRECISION!
1569
1570 // Calculate the residual from including this fNCoefficients.
1572
1573 return dResidur;
1574}
1575
1576
1577////////////////////////////////////////////////////////////////////////////////
1578/// Make histograms of the result of the analysis. This message
1579/// should be sent after having read all data points, but before
1580/// finding the parameterization
1581///
1582/// Options:
1583/// A All the below
1584/// X Original independent variables
1585/// D Original dependent variables
1586/// N Normalised independent variables
1587/// S Shifted dependent variables
1588/// R1 Residuals versus normalised independent variables
1589/// R2 Residuals versus dependent variable
1590/// R3 Residuals computed on training sample
1591/// R4 Residuals computed on test sample
1592///
1593/// For a description of these quantities, refer to
1594/// <a href="#TMultiDimFit:description">class description</a>
1595
1597{
1598 TString opt(option);
1599 opt.ToLower();
1600
1601 if (opt.Length() < 1)
1602 return;
1603
1604 if (!fHistograms)
1605 fHistograms = new TList;
1606
1607 // Counter variable
1608 Int_t i = 0;
1609
1610 // Histogram of original variables
1611 if (opt.Contains("x") || opt.Contains("a")) {
1613 for (i = 0; i < fNVariables; i++)
1614 if (!fHistograms->FindObject(Form("x_%d_orig",i)))
1615 fHistograms->Add(new TH1D(Form("x_%d_orig",i),
1616 Form("Original variable # %d",i),
1618 fMaxVariables(i)));
1619 }
1620
1621 // Histogram of original dependent variable
1622 if (opt.Contains("d") || opt.Contains("a")) {
1624 if (!fHistograms->FindObject("d_orig"))
1625 fHistograms->Add(new TH1D("d_orig", "Original Quantity",
1627 }
1628
1629 // Histograms of normalized variables
1630 if (opt.Contains("n") || opt.Contains("a")) {
1632 for (i = 0; i < fNVariables; i++)
1633 if (!fHistograms->FindObject(Form("x_%d_norm",i)))
1634 fHistograms->Add(new TH1D(Form("x_%d_norm",i),
1635 Form("Normalized variable # %d",i),
1636 fBinVarX, -1,1));
1637 }
1638
1639 // Histogram of shifted dependent variable
1640 if (opt.Contains("s") || opt.Contains("a")) {
1642 if (!fHistograms->FindObject("d_shifted"))
1643 fHistograms->Add(new TH1D("d_shifted", "Shifted Quantity",
1646 }
1647
1648 // Residual from training sample versus independent variables
1649 if (opt.Contains("r1") || opt.Contains("a")) {
1651 for (i = 0; i < fNVariables; i++)
1652 if (!fHistograms->FindObject(Form("res_x_%d",i)))
1653 fHistograms->Add(new TH2D(Form("res_x_%d",i),
1654 Form("Computed residual versus x_%d", i),
1655 fBinVarX, -1, 1,
1656 fBinVarY,
1659 }
1660
1661 // Residual from training sample versus. dependent variable
1662 if (opt.Contains("r2") || opt.Contains("a")) {
1664 if (!fHistograms->FindObject("res_d"))
1665 fHistograms->Add(new TH2D("res_d",
1666 "Computed residuals vs Quantity",
1667 fBinVarX,
1670 fBinVarY,
1673 }
1674
1675 // Residual from training sample
1676 if (opt.Contains("r3") || opt.Contains("a")) {
1678 if (!fHistograms->FindObject("res_train"))
1679 fHistograms->Add(new TH1D("res_train",
1680 "Computed residuals over training sample",
1683
1684 }
1685 if (opt.Contains("r4") || opt.Contains("a")) {
1687 if (!fHistograms->FindObject("res_test"))
1688 fHistograms->Add(new TH1D("res_test",
1689 "Distribution of residuals from test",
1692 }
1693}
1694
1695
1696////////////////////////////////////////////////////////////////////////////////
1697/// Generate the file <classname>MDF.cxx which contains the
1698/// implementation of the method:
1699///
1700/// Double_t <classname>::MDF(Double_t *x)
1701///
1702/// which does the same as TMultiDimFit::Eval. Please refer to this
1703/// method.
1704///
1705/// Further, the public static members:
1706///
1707/// Int_t <classname>::fgNVariables
1708/// Int_t <classname>::fgNCoefficients
1709/// Double_t <classname>::fgDMean
1710/// Double_t <classname>::fgXMean[] //[fgNVariables]
1711/// Double_t <classname>::fgXMin[] //[fgNVariables]
1712/// Double_t <classname>::fgXMax[] //[fgNVariables]
1713/// Double_t <classname>::fgCoefficient[] //[fgNCoeffficents]
1714/// Int_t <classname>::fgPower[] //[fgNCoeffficents*fgNVariables]
1715///
1716/// are initialized, and assumed to exist. The class declaration is
1717/// assumed to be in <classname>.h and assumed to be provided by the
1718/// user.
1719///
1720/// See TMultiDimFit::MakeRealCode for a list of options
1721///
1722/// The minimal class definition is:
1723///
1724/// class <classname> {
1725/// public:
1726/// Int_t <classname>::fgNVariables; // Number of variables
1727/// Int_t <classname>::fgNCoefficients; // Number of terms
1728/// Double_t <classname>::fgDMean; // Mean from training sample
1729/// Double_t <classname>::fgXMean[]; // Mean from training sample
1730/// Double_t <classname>::fgXMin[]; // Min from training sample
1731/// Double_t <classname>::fgXMax[]; // Max from training sample
1732/// Double_t <classname>::fgCoefficient[]; // Coefficients
1733/// Int_t <classname>::fgPower[]; // Function powers
1734///
1735/// Double_t Eval(Double_t *x);
1736/// };
1737///
1738/// Whether the method <classname>::Eval should be static or not, is
1739/// up to the user.
1740
1741void TMultiDimFit::MakeMethod(const Char_t* classname, Option_t* option)
1742{
1743 MakeRealCode(Form("%sMDF.cxx", classname), classname, option);
1744}
1745
1746
1747
1748////////////////////////////////////////////////////////////////////////////////
1749/// PRIVATE METHOD:
1750/// Normalize data to the interval [-1;1]. This is needed for the
1751/// classes method to work.
1752
1754{
1755 Int_t i = 0;
1756 Int_t j = 0;
1757 Int_t k = 0;
1758
1759 for (i = 0; i < fSampleSize; i++) {
1761 ((TH1D*)fHistograms->FindObject("d_orig"))->Fill(fQuantity(i));
1762
1765
1767 ((TH1D*)fHistograms->FindObject("d_shifted"))->Fill(fQuantity(i));
1768
1769 for (j = 0; j < fNVariables; j++) {
1770 Double_t range = 1. / (fMaxVariables(j) - fMinVariables(j));
1771 k = i * fNVariables + j;
1772
1773 // Fill histograms of original independent variables
1775 ((TH1D*)fHistograms->FindObject(Form("x_%d_orig",j)))
1776 ->Fill(fVariables(k));
1777
1778 // Normalise independent variables
1779 fVariables(k) = 1 + 2 * range * (fVariables(k) - fMaxVariables(j));
1780
1781 // Fill histograms of normalised independent variables
1783 ((TH1D*)fHistograms->FindObject(Form("x_%d_norm",j)))
1784 ->Fill(fVariables(k));
1785
1786 }
1787 }
1788 // Shift min and max of dependent variable
1791
1792 // Shift mean of independent variables
1793 for (i = 0; i < fNVariables; i++) {
1794 Double_t range = 1. / (fMaxVariables(i) - fMinVariables(i));
1795 fMeanVariables(i) = 1 + 2 * range * (fMeanVariables(i)
1796 - fMaxVariables(i));
1797 }
1798}
1799
1800
1801////////////////////////////////////////////////////////////////////////////////
1802/// PRIVATE METHOD:
1803/// Find the parameterization over the training sample. A full account
1804/// of the algorithm is given in the
1805/// <a href="#TMultiDimFit:description">class description</a>
1806
1808{
1809 Int_t i = -1;
1810 Int_t j = 0;
1811 Int_t k = 0;
1812 Int_t maxPass = 3;
1813 Int_t studied = 0;
1814 Double_t squareResidual = fSumSqAvgQuantity;
1815 fNCoefficients = 0;
1822 fFunctions = 1;
1823
1826 Int_t l;
1827 for (l=0;l<fMaxFunctions;l++) fFunctionCodes[l] = 0;
1828 for (l=0;l<fMaxTerms;l++) fPowerIndex[l] = 0;
1829
1830 if (fMaxAngle != 0) maxPass = 100;
1831 if (fIsUserFunction) maxPass = 1;
1832
1833 // Loop over the number of functions we want to study.
1834 // increment inspection counter
1835 while(kTRUE) {
1836
1837 // Reach user defined limit of studies
1838 if (studied++ >= fMaxStudy) {
1840 break;
1841 }
1842
1843 // Considered all functions several times
1844 if (k >= maxPass) {
1846 break;
1847 }
1848
1849 // increment function counter
1850 i++;
1851
1852 // If we've reached the end of the functions, restart pass
1853 if (i == fMaxFunctions) {
1854 if (fMaxAngle != 0)
1855 fMaxAngle += (90 - fMaxAngle) / 2;
1856 i = 0;
1857 studied--;
1858 k++;
1859 continue;
1860 }
1861 if (studied == 1)
1862 fFunctionCodes[i] = 0;
1863 else if (fFunctionCodes[i] >= 2)
1864 continue;
1865
1866 // Print a happy message
1867 if (fIsVerbose && studied == 1)
1868 std::cout << "Coeff SumSqRes Contrib Angle QM Func"
1869 << " Value W^2 Powers" << std::endl;
1870
1871 // Make the Gram-Schmidt
1872 Double_t dResidur = MakeGramSchmidt(i);
1873
1874 if (dResidur == 0) {
1875 // This function is no good!
1876 // First test is in MakeGramSchmidt
1877 fFunctionCodes[i] = 1;
1878 continue;
1879 }
1880
1881 // If user specified function, assume they know what they are doing
1882 if (!fIsUserFunction) {
1883 // Flag this function as considered
1884 fFunctionCodes[i] = 2;
1885
1886 // Test if this function contributes to the fit
1887 if (!TestFunction(squareResidual, dResidur)) {
1888 fFunctionCodes[i] = 1;
1889 continue;
1890 }
1891 }
1892
1893 // If we get to here, the function currently considered is
1894 // fNCoefficients, so we increment the counter
1895 // Flag this function as OK, and store and the number in the
1896 // index.
1897 fFunctionCodes[i] = 3;
1900
1901 // We add the current contribution to the sum of square of
1902 // residuals;
1903 squareResidual -= dResidur;
1904
1905
1906 // Calculate control parameter from this function
1907 for (j = 0; j < fNVariables; j++) {
1908 if (fNCoefficients == 1
1909 || fMaxPowersFinal[j] <= fPowers[i * fNVariables + j] - 1)
1910 fMaxPowersFinal[j] = fPowers[i * fNVariables + j] - 1;
1911 }
1913
1914 // Print the statistics about this function
1915 if (fIsVerbose) {
1916 std::cout << std::setw(5) << fNCoefficients << " "
1917 << std::setw(10) << std::setprecision(4) << squareResidual << " "
1918 << std::setw(10) << std::setprecision(4) << dResidur << " "
1919 << std::setw(7) << std::setprecision(3) << fMaxAngle << " "
1920 << std::setw(7) << std::setprecision(3) << s << " "
1921 << std::setw(5) << i << " "
1922 << std::setw(10) << std::setprecision(4)
1924 << std::setw(10) << std::setprecision(4)
1926 << std::flush;
1927 for (j = 0; j < fNVariables; j++)
1928 std::cout << " " << fPowers[i * fNVariables + j] - 1 << std::flush;
1929 std::cout << std::endl;
1930 }
1931
1932 if (fNCoefficients >= fMaxTerms /* && fIsVerbose */) {
1934 break;
1935 }
1936
1937 Double_t err = TMath::Sqrt(TMath::Max(1e-20,squareResidual) /
1939 if (err < fMinRelativeError) {
1941 break;
1942 }
1943
1944 }
1945
1946 fError = TMath::Max(1e-20,squareResidual);
1949}
1950
1951
1952////////////////////////////////////////////////////////////////////////////////
1953/// PRIVATE METHOD:
1954/// This is the method that actually generates the code for the
1955/// evaluation the parameterization on some point.
1956/// It's called by TMultiDimFit::MakeCode and TMultiDimFit::MakeMethod.
1957///
1958/// The options are: NONE so far
1959
1960void TMultiDimFit::MakeRealCode(const char *filename,
1961 const char *classname,
1962 Option_t *)
1963{
1964 Int_t i, j;
1965
1966 Bool_t isMethod = (classname[0] == '\0' ? kFALSE : kTRUE);
1967 const char *prefix = (isMethod ? Form("%s::", classname) : "");
1968 const char *cv_qual = (isMethod ? "" : "static ");
1969
1970 std::ofstream outFile(filename,std::ios::out|std::ios::trunc);
1971 if (!outFile) {
1972 Error("MakeRealCode","couldn't open output file '%s'",filename);
1973 return;
1974 }
1975
1976 if (fIsVerbose)
1977 std::cout << "Writing on file \"" << filename << "\" ... " << std::flush;
1978 //
1979 // Write header of file
1980 //
1981 // Emacs mode line ;-)
1982 outFile << "// -*- mode: c++ -*-" << std::endl;
1983 // Info about creator
1984 outFile << "// " << std::endl
1985 << "// File " << filename
1986 << " generated by TMultiDimFit::MakeRealCode" << std::endl;
1987 // Time stamp
1988 TDatime date;
1989 outFile << "// on " << date.AsString() << std::endl;
1990 // ROOT version info
1991 outFile << "// ROOT version " << gROOT->GetVersion()
1992 << std::endl << "//" << std::endl;
1993 // General information on the code
1994 outFile << "// This file contains the function " << std::endl
1995 << "//" << std::endl
1996 << "// double " << prefix << "MDF(double *x); " << std::endl
1997 << "//" << std::endl
1998 << "// For evaluating the parameterization obtained" << std::endl
1999 << "// from TMultiDimFit and the point x" << std::endl
2000 << "// " << std::endl
2001 << "// See TMultiDimFit class documentation for more "
2002 << "information " << std::endl << "// " << std::endl;
2003 // Header files
2004 if (isMethod)
2005 // If these are methods, we need the class header
2006 outFile << "#include \"" << classname << ".h\"" << std::endl;
2007
2008 //
2009 // Now for the data
2010 //
2011 outFile << "//" << std::endl
2012 << "// Static data variables" << std::endl
2013 << "//" << std::endl;
2014 outFile << cv_qual << "int " << prefix << "gNVariables = "
2015 << fNVariables << ";" << std::endl;
2016 outFile << cv_qual << "int " << prefix << "gNCoefficients = "
2017 << fNCoefficients << ";" << std::endl;
2018 outFile << cv_qual << "double " << prefix << "gDMean = "
2019 << fMeanQuantity << ";" << std::endl;
2020
2021 // Assignment to mean vector.
2022 outFile << "// Assignment to mean vector." << std::endl;
2023 outFile << cv_qual << "double " << prefix
2024 << "gXMean[] = {" << std::endl;
2025 for (i = 0; i < fNVariables; i++)
2026 outFile << (i != 0 ? ", " : " ") << fMeanVariables(i) << std::flush;
2027 outFile << " };" << std::endl << std::endl;
2028
2029 // Assignment to minimum vector.
2030 outFile << "// Assignment to minimum vector." << std::endl;
2031 outFile << cv_qual << "double " << prefix
2032 << "gXMin[] = {" << std::endl;
2033 for (i = 0; i < fNVariables; i++)
2034 outFile << (i != 0 ? ", " : " ") << fMinVariables(i) << std::flush;
2035 outFile << " };" << std::endl << std::endl;
2036
2037 // Assignment to maximum vector.
2038 outFile << "// Assignment to maximum vector." << std::endl;
2039 outFile << cv_qual << "double " << prefix
2040 << "gXMax[] = {" << std::endl;
2041 for (i = 0; i < fNVariables; i++)
2042 outFile << (i != 0 ? ", " : " ") << fMaxVariables(i) << std::flush;
2043 outFile << " };" << std::endl << std::endl;
2044
2045 // Assignment to coefficients vector.
2046 outFile << "// Assignment to coefficients vector." << std::endl;
2047 outFile << cv_qual << "double " << prefix
2048 << "gCoefficient[] = {" << std::flush;
2049 for (i = 0; i < fNCoefficients; i++)
2050 outFile << (i != 0 ? "," : "") << std::endl
2051 << " " << fCoefficients(i) << std::flush;
2052 outFile << std::endl << " };" << std::endl << std::endl;
2053
2054 // Assignment to error coefficients vector.
2055 outFile << "// Assignment to error coefficients vector." << std::endl;
2056 outFile << cv_qual << "double " << prefix
2057 << "gCoefficientRMS[] = {" << std::flush;
2058 for (i = 0; i < fNCoefficients; i++)
2059 outFile << (i != 0 ? "," : "") << std::endl
2060 << " " << fCoefficientsRMS(i) << std::flush;
2061 outFile << std::endl << " };" << std::endl << std::endl;
2062
2063 // Assignment to powers vector.
2064 outFile << "// Assignment to powers vector." << std::endl
2065 << "// The powers are stored row-wise, that is" << std::endl
2066 << "// p_ij = " << prefix
2067 << "gPower[i * NVariables + j];" << std::endl;
2068 outFile << cv_qual << "int " << prefix
2069 << "gPower[] = {" << std::flush;
2070 for (i = 0; i < fNCoefficients; i++) {
2071 for (j = 0; j < fNVariables; j++) {
2072 if (j != 0) outFile << std::flush << " ";
2073 else outFile << std::endl << " ";
2074 outFile << fPowers[fPowerIndex[i] * fNVariables + j]
2075 << (i == fNCoefficients - 1 && j == fNVariables - 1 ? "" : ",")
2076 << std::flush;
2077 }
2078 }
2079 outFile << std::endl << "};" << std::endl << std::endl;
2080
2081
2082 //
2083 // Finally we reach the function itself
2084 //
2085 outFile << "// " << std::endl
2086 << "// The "
2087 << (isMethod ? "method " : "function ")
2088 << " double " << prefix
2089 << "MDF(double *x)"
2090 << std::endl << "// " << std::endl;
2091 outFile << "double " << prefix
2092 << "MDF(double *x) {" << std::endl
2093 << " double returnValue = " << prefix << "gDMean;" << std::endl
2094 << " int i = 0, j = 0, k = 0;" << std::endl
2095 << " for (i = 0; i < " << prefix << "gNCoefficients ; i++) {"
2096 << std::endl
2097 << " // Evaluate the ith term in the expansion" << std::endl
2098 << " double term = " << prefix << "gCoefficient[i];"
2099 << std::endl
2100 << " for (j = 0; j < " << prefix << "gNVariables; j++) {"
2101 << std::endl
2102 << " // Evaluate the polynomial in the jth variable." << std::endl
2103 << " int power = "<< prefix << "gPower["
2104 << prefix << "gNVariables * i + j]; " << std::endl
2105 << " double p1 = 1, p2 = 0, p3 = 0, r = 0;" << std::endl
2106 << " double v = 1 + 2. / ("
2107 << prefix << "gXMax[j] - " << prefix
2108 << "gXMin[j]) * (x[j] - " << prefix << "gXMax[j]);" << std::endl
2109 << " // what is the power to use!" << std::endl
2110 << " switch(power) {" << std::endl
2111 << " case 1: r = 1; break; " << std::endl
2112 << " case 2: r = v; break; " << std::endl
2113 << " default: " << std::endl
2114 << " p2 = v; " << std::endl
2115 << " for (k = 3; k <= power; k++) { " << std::endl
2116 << " p3 = p2 * v;" << std::endl;
2117 if (fPolyType == kLegendre)
2118 outFile << " p3 = ((2 * i - 3) * p2 * v - (i - 2) * p1)"
2119 << " / (i - 1);" << std::endl;
2120 if (fPolyType == kChebyshev)
2121 outFile << " p3 = 2 * v * p2 - p1; " << std::endl;
2122 outFile << " p1 = p2; p2 = p3; " << std::endl << " }" << std::endl
2123 << " r = p3;" << std::endl << " }" << std::endl
2124 << " // multiply this term by the poly in the jth var" << std::endl
2125 << " term *= r; " << std::endl << " }" << std::endl
2126 << " // Add this term to the final result" << std::endl
2127 << " returnValue += term;" << std::endl << " }" << std::endl
2128 << " return returnValue;" << std::endl << "}" << std::endl << std::endl;
2129
2130 // EOF
2131 outFile << "// EOF for " << filename << std::endl;
2132
2133 // Close the file
2134 outFile.close();
2135
2136 if (fIsVerbose)
2137 std::cout << "done" << std::endl;
2138}
2139
2140
2141////////////////////////////////////////////////////////////////////////////////
2142/// Print statistics etc.
2143/// Options are
2144/// P Parameters
2145/// S Statistics
2146/// C Coefficients
2147/// R Result of parameterisation
2148/// F Result of fit
2149/// K Correlation Matrix
2150/// M Pretty print formula
2151///
2152
2154{
2155 Int_t i = 0;
2156 Int_t j = 0;
2157
2158 TString opt(option);
2159 opt.ToLower();
2160
2161 if (opt.Contains("p")) {
2162 // Print basic parameters for this object
2163 std::cout << "User parameters:" << std::endl
2164 << "----------------" << std::endl
2165 << " Variables: " << fNVariables << std::endl
2166 << " Data points: " << fSampleSize << std::endl
2167 << " Max Terms: " << fMaxTerms << std::endl
2168 << " Power Limit Parameter: " << fPowerLimit << std::endl
2169 << " Max functions: " << fMaxFunctions << std::endl
2170 << " Max functions to study: " << fMaxStudy << std::endl
2171 << " Max angle (optional): " << fMaxAngle << std::endl
2172 << " Min angle: " << fMinAngle << std::endl
2173 << " Relative Error accepted: " << fMinRelativeError << std::endl
2174 << " Maximum Powers: " << std::flush;
2175 for (i = 0; i < fNVariables; i++)
2176 std::cout << " " << fMaxPowers[i] - 1 << std::flush;
2177 std::cout << std::endl << std::endl
2178 << " Parameterisation will be done using " << std::flush;
2179 if (fPolyType == kChebyshev)
2180 std::cout << "Chebyshev polynomials" << std::endl;
2181 else if (fPolyType == kLegendre)
2182 std::cout << "Legendre polynomials" << std::endl;
2183 else
2184 std::cout << "Monomials" << std::endl;
2185 std::cout << std::endl;
2186 }
2187
2188 if (opt.Contains("s")) {
2189 // Print statistics for read data
2190 std::cout << "Sample statistics:" << std::endl
2191 << "------------------" << std::endl
2192 << " D" << std::flush;
2193 for (i = 0; i < fNVariables; i++)
2194 std::cout << " " << std::setw(10) << i+1 << std::flush;
2195 std::cout << std::endl << " Max: " << std::setw(10) << std::setprecision(7)
2196 << fMaxQuantity << std::flush;
2197 for (i = 0; i < fNVariables; i++)
2198 std::cout << " " << std::setw(10) << std::setprecision(4)
2199 << fMaxVariables(i) << std::flush;
2200 std::cout << std::endl << " Min: " << std::setw(10) << std::setprecision(7)
2201 << fMinQuantity << std::flush;
2202 for (i = 0; i < fNVariables; i++)
2203 std::cout << " " << std::setw(10) << std::setprecision(4)
2204 << fMinVariables(i) << std::flush;
2205 std::cout << std::endl << " Mean: " << std::setw(10) << std::setprecision(7)
2206 << fMeanQuantity << std::flush;
2207 for (i = 0; i < fNVariables; i++)
2208 std::cout << " " << std::setw(10) << std::setprecision(4)
2209 << fMeanVariables(i) << std::flush;
2210 std::cout << std::endl << " Function Sum Squares: " << fSumSqQuantity
2211 << std::endl << std::endl;
2212 }
2213
2214 if (opt.Contains("r")) {
2215 std::cout << "Results of Parameterisation:" << std::endl
2216 << "----------------------------" << std::endl
2217 << " Total reduction of square residuals "
2218 << fSumSqResidual << std::endl
2219 << " Relative precision obtained: "
2220 << fPrecision << std::endl
2221 << " Error obtained: "
2222 << fError << std::endl
2223 << " Multiple correlation coefficient: "
2224 << fCorrelationCoeff << std::endl
2225 << " Reduced Chi square over sample: "
2226 << fChi2 / (fSampleSize - fNCoefficients) << std::endl
2227 << " Maximum residual value: "
2228 << fMaxResidual << std::endl
2229 << " Minimum residual value: "
2230 << fMinResidual << std::endl
2231 << " Estimated root mean square: "
2232 << fRMS << std::endl
2233 << " Maximum powers used: " << std::flush;
2234 for (j = 0; j < fNVariables; j++)
2235 std::cout << fMaxPowersFinal[j] << " " << std::flush;
2236 std::cout << std::endl
2237 << " Function codes of candidate functions." << std::endl
2238 << " 1: considered,"
2239 << " 2: too little contribution,"
2240 << " 3: accepted." << std::flush;
2241 for (i = 0; i < fMaxFunctions; i++) {
2242 if (i % 60 == 0)
2243 std::cout << std::endl << " " << std::flush;
2244 else if (i % 10 == 0)
2245 std::cout << " " << std::flush;
2246 std::cout << fFunctionCodes[i];
2247 }
2248 std::cout << std::endl << " Loop over candidates stopped because " << std::flush;
2249 switch(fParameterisationCode){
2250 case PARAM_MAXSTUDY:
2251 std::cout << "max allowed studies reached" << std::endl; break;
2252 case PARAM_SEVERAL:
2253 std::cout << "all candidates considered several times" << std::endl; break;
2254 case PARAM_RELERR:
2255 std::cout << "wanted relative error obtained" << std::endl; break;
2256 case PARAM_MAXTERMS:
2257 std::cout << "max number of terms reached" << std::endl; break;
2258 default:
2259 std::cout << "some unknown reason" << std::endl;
2260 break;
2261 }
2262 std::cout << std::endl;
2263 }
2264
2265 if (opt.Contains("f")) {
2266 std::cout << "Results of Fit:" << std::endl
2267 << "---------------" << std::endl
2268 << " Test sample size: "
2269 << fTestSampleSize << std::endl
2270 << " Multiple correlation coefficient: "
2271 << fTestCorrelationCoeff << std::endl
2272 << " Relative precision obtained: "
2273 << fTestPrecision << std::endl
2274 << " Error obtained: "
2275 << fTestError << std::endl
2276 << " Reduced Chi square over sample: "
2277 << fChi2 / (fSampleSize - fNCoefficients) << std::endl
2278 << std::endl;
2279 if (fFitter) {
2280 fFitter->PrintResults(1,1);
2281 std::cout << std::endl;
2282 }
2283 }
2284
2285 if (opt.Contains("c")){
2286 std::cout << "Coefficients:" << std::endl
2287 << "-------------" << std::endl
2288 << " # Value Error Powers" << std::endl
2289 << " ---------------------------------------" << std::endl;
2290 for (i = 0; i < fNCoefficients; i++) {
2291 std::cout << " " << std::setw(3) << i << " "
2292 << std::setw(12) << fCoefficients(i) << " "
2293 << std::setw(12) << fCoefficientsRMS(i) << " " << std::flush;
2294 for (j = 0; j < fNVariables; j++)
2295 std::cout << " " << std::setw(3)
2296 << fPowers[fPowerIndex[i] * fNVariables + j] - 1 << std::flush;
2297 std::cout << std::endl;
2298 }
2299 std::cout << std::endl;
2300 }
2301 if (opt.Contains("k") && fCorrelationMatrix.IsValid()) {
2302 std::cout << "Correlation Matrix:" << std::endl
2303 << "-------------------";
2305 }
2306
2307 if (opt.Contains("m")) {
2308 std::cout << "Parameterization:" << std::endl
2309 << "-----------------" << std::endl
2310 << " Normalised variables: " << std::endl;
2311 for (i = 0; i < fNVariables; i++)
2312 std::cout << "\ty_" << i << "\t= 1 + 2 * (x_" << i << " - "
2313 << fMaxVariables(i) << ") / ("
2314 << fMaxVariables(i) << " - " << fMinVariables(i) << ")"
2315 << std::endl;
2316 std::cout << std::endl
2317 << " f(";
2318 for (i = 0; i < fNVariables; i++) {
2319 std::cout << "y_" << i;
2320 if (i != fNVariables-1) std::cout << ", ";
2321 }
2322 std::cout << ") = ";
2323 for (i = 0; i < fNCoefficients; i++) {
2324 if (i != 0)
2325 std::cout << std::endl << "\t" << (fCoefficients(i) < 0 ? "- " : "+ ")
2327 else
2328 std::cout << fCoefficients(i);
2329 for (j = 0; j < fNVariables; j++) {
2330 Int_t p = fPowers[fPowerIndex[i] * fNVariables + j];
2331 switch (p) {
2332 case 1: break;
2333 case 2: std::cout << " * y_" << j; break;
2334 default:
2335 switch(fPolyType) {
2336 case kLegendre: std::cout << " * L_" << p-1 << "(y_" << j << ")"; break;
2337 case kChebyshev: std::cout << " * C_" << p-1 << "(y_" << j << ")"; break;
2338 default: std::cout << " * y_" << j << "^" << p-1; break;
2339 }
2340 }
2341
2342 }
2343 }
2344 std::cout << std::endl;
2345 }
2346}
2347
2348
2349////////////////////////////////////////////////////////////////////////////////
2350/// Selection method. User can override this method for specialized
2351/// selection of acceptable functions in fit. Default is to select
2352/// all. This message is sent during the build-up of the function
2353/// candidates table once for each set of powers in
2354/// variables. Notice, that the argument array contains the powers
2355/// PLUS ONE. For example, to De select the function
2356/// f = x1^2 * x2^4 * x3^5,
2357/// this method should return kFALSE if given the argument
2358/// { 3, 4, 6 }
2359
2361{
2362 return kTRUE;
2363}
2364
2365////////////////////////////////////////////////////////////////////////////////
2366/// Set the max angle (in degrees) between the initial data vector to
2367/// be fitted, and the new candidate function to be included in the
2368/// fit. By default it is 0, which automatically chooses another
2369/// selection criteria. See also
2370/// <a href="#TMultiDimFit:description">class description</a>
2371
2373{
2374 if (ang >= 90 || ang < 0) {
2375 Warning("SetMaxAngle", "angle must be in [0,90)");
2376 return;
2377 }
2378
2379 fMaxAngle = ang;
2380}
2381
2382////////////////////////////////////////////////////////////////////////////////
2383/// Set the min angle (in degrees) between a new candidate function
2384/// and the subspace spanned by the previously accepted
2385/// functions. See also
2386/// <a href="#TMultiDimFit:description">class description</a>
2387
2389{
2390 if (ang > 90 || ang <= 0) {
2391 Warning("SetMinAngle", "angle must be in [0,90)");
2392 return;
2393 }
2394
2395 fMinAngle = ang;
2396
2397}
2398
2399
2400////////////////////////////////////////////////////////////////////////////////
2401/// Define a user function. The input array must be of the form
2402/// (p11, ..., p1N, ... ,pL1, ..., pLN)
2403/// Where N is the dimension of the data sample, L is the number of
2404/// terms (given in terms) and the first number, labels the term, the
2405/// second the variable. More information is given in the
2406/// <a href="#TMultiDimFit:description">class description</a>
2407
2408void TMultiDimFit::SetPowers(const Int_t* powers, Int_t terms)
2409{
2411 fMaxFunctions = terms;
2412 fMaxTerms = terms;
2413 fMaxStudy = terms;
2415 fPowers = new Int_t[fMaxFuncNV];
2416 Int_t i, j;
2417 for (i = 0; i < fMaxFunctions; i++)
2418 for(j = 0; j < fNVariables; j++)
2419 fPowers[i * fNVariables + j] = powers[i * fNVariables + j] + 1;
2420}
2421
2422////////////////////////////////////////////////////////////////////////////////
2423/// Set the user parameter for the function selection. The bigger the
2424/// limit, the more functions are used. The meaning of this variable
2425/// is defined in the
2426/// <a href="#TMultiDimFit:description">class description</a>
2427
2429{
2430 fPowerLimit = limit;
2431}
2432
2433////////////////////////////////////////////////////////////////////////////////
2434/// Set the maximum power to be considered in the fit for each
2435/// variable. See also
2436/// <a href="#TMultiDimFit:description">class description</a>
2437
2439{
2440 if (!powers)
2441 return;
2442
2443 for (Int_t i = 0; i < fNVariables; i++)
2444 fMaxPowers[i] = powers[i]+1;
2445}
2446
2447////////////////////////////////////////////////////////////////////////////////
2448/// Set the acceptable relative error for when sum of square
2449/// residuals is considered minimized. For a full account, refer to
2450/// the
2451/// <a href="#TMultiDimFit:description">class description</a>
2452
2454{
2455 fMinRelativeError = error;
2456}
2457
2458
2459////////////////////////////////////////////////////////////////////////////////
2460/// PRIVATE METHOD:
2461/// Test whether the currently considered function contributes to the
2462/// fit. See also
2463/// <a href="#TMultiDimFit:description">class description</a>
2464
2466 Double_t dResidur)
2467{
2468 if (fNCoefficients != 0) {
2469 // Now for the second test:
2470 if (fMaxAngle == 0) {
2471 // If the user hasn't supplied a max angle do the test as,
2472 if (dResidur <
2473 squareResidual / (fMaxTerms - fNCoefficients + 1 + 1E-10)) {
2474 return kFALSE;
2475 }
2476 }
2477 else {
2478 // If the user has provided a max angle, test if the calculated
2479 // angle is less then the max angle.
2480 if (TMath::Sqrt(dResidur/fSumSqAvgQuantity) <
2482 return kFALSE;
2483 }
2484 }
2485 }
2486 // If we get here, the function is OK
2487 return kTRUE;
2488}
2489
2490
2491////////////////////////////////////////////////////////////////////////////////
2492/// Helper function for doing the minimisation of Chi2 using Minuit
2493
2494void mdfHelper(int& /*npar*/, double* /*divs*/, double& chi2,
2495 double* coeffs, int /*flag*/)
2496{
2497 // Get pointer to current TMultiDimFit object.
2499 chi2 = mdf->MakeChi2(coeffs);
2500}
ROOT::R::TRInterface & r
Definition: Object.C:4
#define b(i)
Definition: RSha256.hxx:100
#define f(i)
Definition: RSha256.hxx:104
#define h(i)
Definition: RSha256.hxx:106
#define e(i)
Definition: RSha256.hxx:103
int Int_t
Definition: RtypesCore.h:43
char Char_t
Definition: RtypesCore.h:31
const Bool_t kFALSE
Definition: RtypesCore.h:90
double Double_t
Definition: RtypesCore.h:57
const Bool_t kTRUE
Definition: RtypesCore.h:89
const char Option_t
Definition: RtypesCore.h:64
#define ClassImp(name)
Definition: Rtypes.h:361
#define SETBIT(n, i)
Definition: Rtypes.h:84
#define TESTBIT(n, i)
Definition: Rtypes.h:86
int type
Definition: TGX11.cxx:120
double sqrt(double)
TMatrixTRow< Double_t > TMatrixDRow
TMatrixTDiag_const< Double_t > TMatrixDDiag_const
#define HIST_DSHIF
#define PARAM_RELERR
#define HIST_DORIG
#define HIST_XNORM
#define HIST_RX
#define HIST_XORIG
#define PARAM_SEVERAL
#define PARAM_MAXSTUDY
static void mdfHelper(int &, double *, double &, double *, int)
Helper function for doing the minimisation of Chi2 using Minuit.
#define HIST_RTEST
#define DEGRAD
#define HIST_RTRAI
#define HIST_RD
#define PARAM_MAXTERMS
#define gROOT
Definition: TROOT.h:406
char * Form(const char *fmt,...)
Using a TBrowser one can browse all ROOT objects.
Definition: TBrowser.h:37
This class stores the date and time with a precision of one second in an unsigned 32 bit word (950130...
Definition: TDatime.h:37
const char * AsString() const
Return the date & time as a string (ctime() format).
Definition: TDatime.cxx:101
Cholesky Decomposition class.
Definition: TDecompChol.h:25
Bool_t Invert(TMatrixDSym &inv)
For a symmetric matrix A(m,m), its inverse A_inv(m,m) is returned .
virtual Bool_t Decompose()
Matrix A is decomposed in component U so that A = U^T * U If the decomposition succeeds,...
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:614
The TH1 histogram class.
Definition: TH1.h:56
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:292
A doubly linked list.
Definition: TList.h:44
virtual void Add(TObject *obj)
Definition: TList.h:87
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:577
virtual void Clear(Option_t *option="")
Remove all objects from the list.
Definition: TList.cxx:401
void Print(Option_t *name="") const
Print the matrix as a table of elements.
virtual TMatrixTBase< Element > & Zero()
Set matrix elements to zero.
Bool_t IsValid() const
Definition: TMatrixTBase.h:147
virtual TMatrixTBase< Element > & NormByDiag(const TVectorT< Element > &v, Option_t *option="D")
option:
virtual TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t=-1)
Set size of the matrix to nrows x ncols New dynamic elements are created, the overlapping part of the...
Definition: TMatrixT.cxx:1213
Multidimensional Fits in ROOT.
Definition: TMultiDimFit.h:15
virtual void MakeCorrelation()
PRIVATE METHOD: Compute the correlation matrix.
TVectorD fCoefficients
Definition: TMultiDimFit.h:82
Double_t fPrecision
Definition: TMultiDimFit.h:90
Byte_t fHistogramMask
Definition: TMultiDimFit.h:97
virtual void MakeCoefficientErrors()
PRIVATE METHOD: Compute the errors on the coefficients.
TMatrixD fOrthCurvatureMatrix
Definition: TMultiDimFit.h:81
Bool_t fShowCorrelation
Definition: TMultiDimFit.h:104
Bool_t fIsUserFunction
Definition: TMultiDimFit.h:105
Int_t fNVariables
Definition: TMultiDimFit.h:37
Int_t fMaxResidualRow
Definition: TMultiDimFit.h:75
Double_t fMinQuantity
Definition: TMultiDimFit.h:32
virtual Double_t Eval(const Double_t *x, const Double_t *coeff=0) const
Evaluate parameterization at point x.
Double_t fSumSqQuantity
Definition: TMultiDimFit.h:33
TVectorD fMaxVariables
Definition: TMultiDimFit.h:39
Int_t fSampleSize
Definition: TMultiDimFit.h:42
Double_t fSumSqAvgQuantity
Definition: TMultiDimFit.h:34
Int_t fMaxTerms
Definition: TMultiDimFit.h:52
Int_t * fPowers
Definition: TMultiDimFit.h:69
virtual void MakeNormalized()
PRIVATE METHOD: Normalize data to the interval [-1;1].
Int_t * fMaxPowers
Definition: TMultiDimFit.h:54
Int_t fMaxStudy
Definition: TMultiDimFit.h:61
TVectorD fTestQuantity
Definition: TMultiDimFit.h:44
TVectorD fTestSqError
Definition: TMultiDimFit.h:45
Int_t fMaxFunctions
Definition: TMultiDimFit.h:59
Int_t fBinVarY
Definition: TMultiDimFit.h:99
static TMultiDimFit * Instance()
Return the static instance.
virtual ~TMultiDimFit()
Destructor.
virtual void MakeCode(const char *functionName="MDF", Option_t *option="")
Generate the file <filename> with .C appended if argument doesn't end in .cxx or ....
virtual void MakeCandidates()
PRIVATE METHOD: Create list of candidate functions for the parameterisation.
Double_t fTestError
Definition: TMultiDimFit.h:89
Double_t fMinResidual
Definition: TMultiDimFit.h:74
void SetPowerLimit(Double_t limit=1e-3)
Set the user parameter for the function selection.
virtual void FindParameterization(Option_t *option="")
Find the parameterization.
virtual Double_t MakeGramSchmidt(Int_t function)
PRIVATE METHOD: Make Gram-Schmidt orthogonalisation.
TVectorD fSqError
Definition: TMultiDimFit.h:29
virtual void AddRow(const Double_t *x, Double_t D, Double_t E=0)
Add a row consisting of fNVariables independent variables, the known, dependent quantity,...
TVectorD fQuantity
Definition: TMultiDimFit.h:28
TVectorD fTestVariables
Definition: TMultiDimFit.h:46
TVectorD fOrthCoefficients
Definition: TMultiDimFit.h:80
Int_t fTestSampleSize
Definition: TMultiDimFit.h:48
TVectorD fResiduals
Definition: TMultiDimFit.h:72
virtual Bool_t TestFunction(Double_t squareResidual, Double_t dResidur)
PRIVATE METHOD: Test whether the currently considered function contributes to the fit.
TVectorD fVariables
Definition: TMultiDimFit.h:36
Double_t fCorrelationCoeff
Definition: TMultiDimFit.h:92
void SetMaxPowers(const Int_t *powers)
Set the maximum power to be considered in the fit for each variable.
Int_t * fMaxPowersFinal
Definition: TMultiDimFit.h:68
Int_t fMinResidualRow
Definition: TMultiDimFit.h:76
virtual void SetPowers(const Int_t *powers, Int_t terms)
Define a user function.
Double_t fMaxQuantity
Definition: TMultiDimFit.h:31
virtual void MakeRealCode(const char *filename, const char *classname, Option_t *option="")
PRIVATE METHOD: This is the method that actually generates the code for the evaluation the parameteri...
Double_t fRMS
Definition: TMultiDimFit.h:84
TVirtualFitter * fFitter
Definition: TMultiDimFit.h:101
TMatrixD fFunctions
Definition: TMultiDimFit.h:58
virtual void AddTestRow(const Double_t *x, Double_t D, Double_t E=0)
Add a row consisting of fNVariables independent variables, the known, dependent quantity,...
Int_t fParameterisationCode
Definition: TMultiDimFit.h:86
Double_t fMaxAngle
Definition: TMultiDimFit.h:51
Int_t fNCoefficients
Definition: TMultiDimFit.h:79
virtual void Fit(Option_t *option="")
Try to fit the found parameterisation to the test sample.
TMultiDimFit()
Empty CTOR. Do not use.
Double_t fTestCorrelationCoeff
Definition: TMultiDimFit.h:94
Double_t fMinAngle
Definition: TMultiDimFit.h:50
void SetMinRelativeError(Double_t error)
Set the acceptable relative error for when sum of square residuals is considered minimized.
Double_t fSumSqResidual
Definition: TMultiDimFit.h:77
TVectorD fMeanVariables
Definition: TMultiDimFit.h:38
virtual void MakeMethod(const Char_t *className="MDF", Option_t *option="")
Generate the file <classname>MDF.cxx which contains the implementation of the method:
void SetMaxAngle(Double_t angle=0)
Set the max angle (in degrees) between the initial data vector to be fitted, and the new candidate fu...
Double_t fTestPrecision
Definition: TMultiDimFit.h:91
virtual void MakeCoefficients()
PRIVATE METHOD: Invert the model matrix B, and compute final coefficients.
TMatrixD fCorrelationMatrix
Definition: TMultiDimFit.h:93
Int_t * fPowerIndex
Definition: TMultiDimFit.h:70
TMatrixD fOrthFunctions
Definition: TMultiDimFit.h:64
virtual void Clear(Option_t *option="")
Clear internal structures and variables.
virtual Double_t EvalControl(const Int_t *powers) const
PRIVATE METHOD: Calculate the control parameter from the passed powers.
virtual void Browse(TBrowser *b)
Browse the TMultiDimFit object in the TBrowser.
virtual void MakeParameterization()
PRIVATE METHOD: Find the parameterization over the training sample.
EMDFPolyType fPolyType
Fit object (MINUIT)
Definition: TMultiDimFit.h:103
virtual Bool_t Select(const Int_t *iv)
Selection method.
virtual Double_t EvalError(const Double_t *x, const Double_t *coeff=0) const
Evaluate parameterization error at point x.
TList * fHistograms
Definition: TMultiDimFit.h:96
Double_t fChi2
Definition: TMultiDimFit.h:85
Double_t fMaxResidual
Definition: TMultiDimFit.h:73
Double_t fError
Definition: TMultiDimFit.h:88
Int_t fBinVarX
Definition: TMultiDimFit.h:98
virtual Double_t EvalFactor(Int_t p, Double_t x) const
PRIVATE METHOD: Evaluate function with power p at variable value x.
void SetMinAngle(Double_t angle=1)
Set the min angle (in degrees) between a new candidate function and the subspace spanned by the previ...
virtual void Print(Option_t *option="ps") const
Print statistics etc.
TVectorD fMinVariables
Definition: TMultiDimFit.h:40
TVectorD fOrthFunctionNorms
Definition: TMultiDimFit.h:65
Double_t fMeanQuantity
Definition: TMultiDimFit.h:30
TVectorD fCoefficientsRMS
Definition: TMultiDimFit.h:83
Double_t fPowerLimit
Definition: TMultiDimFit.h:55
static TMultiDimFit * fgInstance
Definition: TMultiDimFit.h:25
virtual void MakeHistograms(Option_t *option="A")
Make histograms of the result of the analysis.
Bool_t fIsVerbose
Definition: TMultiDimFit.h:106
Int_t fMaxFuncNV
Definition: TMultiDimFit.h:62
Double_t fMinRelativeError
Definition: TMultiDimFit.h:53
virtual Double_t MakeChi2(const Double_t *coeff=0)
Calculate Chi square over either the test sample.
Int_t * fFunctionCodes
Definition: TMultiDimFit.h:60
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:877
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:891
Basic string class.
Definition: TString.h:131
Ssiz_t Length() const
Definition: TString.h:405
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1125
Bool_t EndsWith(const char *pat, ECaseCompare cmp=kExact) const
Return true if string ends with the specified string.
Definition: TString.cxx:2177
const char * Data() const
Definition: TString.h:364
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:619
TVectorT< Element > & Zero()
Set vector elements to zero.
Definition: TVectorT.cxx:454
TVectorT< Element > & ResizeTo(Int_t lwb, Int_t upb)
Resize the vector to [lwb:upb] .
Definition: TVectorT.cxx:295
Int_t GetNrows() const
Definition: TVectorT.h:75
Bool_t IsValid() const
Definition: TVectorT.h:83
virtual void PrintResults(Int_t level, Double_t amin) const =0
virtual void SetFCN(void(*fcn)(Int_t &, Double_t *, Double_t &f, Double_t *, Int_t))
To set the address of the minimization objective function called by the native compiler (see function...
virtual Int_t SetParameter(Int_t ipar, const char *parname, Double_t value, Double_t verr, Double_t vlow, Double_t vhigh)=0
virtual Int_t ExecuteCommand(const char *command, Double_t *args, Int_t nargs)=0
virtual Double_t GetParameter(Int_t ipar) const =0
static TVirtualFitter * Fitter(TObject *obj, Int_t maxpar=25)
Static function returning a pointer to the current fitter.
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
const Int_t n
Definition: legend1.C:16
void function(const Char_t *name_, T fun, const Char_t *docstring=0)
Definition: RExports.h:151
static constexpr double s
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:212
constexpr Double_t E()
Base of natural log:
Definition: TMath.h:97
Double_t Sqrt(Double_t x)
Definition: TMath.h:681
Double_t Cos(Double_t)
Definition: TMath.h:631
Double_t Sin(Double_t)
Definition: TMath.h:627
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
auto * m
Definition: textangle.C:8
auto * l
Definition: textangle.C:4
static long int sum(long int i)
Definition: Factory.cxx:2275
REAL epsilon
Definition: triangle.c:617