laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
Lau1DCubicSpline.cc
Go to the documentation of this file.
1 
2 /*
3 Copyright 2015 University of Warwick
4 
5 Licensed under the Apache License, Version 2.0 (the "License");
6 you may not use this file except in compliance with the License.
7 You may obtain a copy of the License at
8 
9  http://www.apache.org/licenses/LICENSE-2.0
10 
11 Unless required by applicable law or agreed to in writing, software
12 distributed under the License is distributed on an "AS IS" BASIS,
13 WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 See the License for the specific language governing permissions and
15 limitations under the License.
16 */
17 
18 /*
19 Laura++ package authors:
20 John Back
21 Paul Harrison
22 Thomas Latham
23 */
24 
29 #include "Lau1DCubicSpline.hh"
30 
31 #include <cmath>
32 #include <cstdlib>
33 #include <iostream>
34 
35 #include <TH2.h>
36 #include <TSystem.h>
37 
38 Lau1DCubicSpline::Lau1DCubicSpline( const std::vector<Double_t>& xs,
39  const std::vector<Double_t>& ys,
40  LauSplineType type,
41  LauSplineBoundaryType leftBound,
42  LauSplineBoundaryType rightBound,
43  Double_t dydx0,
44  Double_t dydxn ) :
45  nKnots_( xs.size() ),
46  x_( xs ),
47  y_( ys ),
48  type_( type ),
49  leftBound_( leftBound ),
50  rightBound_( rightBound ),
51  dydx0_( dydx0 ),
52  dydxn_( dydxn )
53 {
54  init();
55 }
56 
58 {
59 }
60 
61 Double_t Lau1DCubicSpline::evaluate( Double_t x ) const
62 {
63  // do not attempt to extrapolate the spline
64  if ( x < x_[0] || x > x_[nKnots_ - 1] ) {
65  std::cout << "WARNING in Lau1DCubicSpline::evaluate : function is only defined between "
66  << x_[0] << " and " << x_[nKnots_ - 1] << std::endl;
67  std::cout << " value at " << x << " returned as 0"
68  << std::endl;
69  return 0.;
70  }
71 
72  // first determine which 'cell' of the spline x is in
73  // cell i runs from knot i to knot i+1
74  Int_t cell( 0 );
75 
76  while ( x > x_[cell + 1] ) {
77  ++cell;
78  }
79 
80  // obtain x- and y-values of the neighbouring knots
81  Double_t xLow = x_[cell];
82  Double_t xHigh = x_[cell + 1];
83  Double_t yLow = y_[cell];
84  Double_t yHigh = y_[cell + 1];
85 
87  return yHigh * ( x - xLow ) / ( xHigh - xLow ) + yLow * ( xHigh - x ) / ( xHigh - xLow );
88  }
89 
90  // obtain t, the normalised x-coordinate within the cell,
91  // and the coefficients a and b, which are defined in cell i as:
92  //
93  // a_i = k_i *(x_i+1 - x_i) - (y_i+1 - y_i),
94  // b_i = -k_i+1*(x_i+1 - x_i) + (y_i+1 - y_i)
95  //
96  // where k_i is (by construction) the first derivative at knot i
97  Double_t t = ( x - xLow ) / ( xHigh - xLow );
98  Double_t a = dydx_[cell] * ( xHigh - xLow ) - ( yHigh - yLow );
99  Double_t b = -1. * dydx_[cell + 1] * ( xHigh - xLow ) + ( yHigh - yLow );
100 
101  Double_t retVal = ( 1 - t ) * yLow + t * yHigh + t * ( 1 - t ) * ( a * ( 1 - t ) + b * t );
102 
103  return retVal;
104 }
105 
106 void Lau1DCubicSpline::updateYValues( const std::vector<Double_t>& ys )
107 {
108  y_ = ys;
109  this->calcDerivatives();
110 }
111 
113 {
114  if ( type_ != type ) {
115  type_ = type;
116  this->calcDerivatives();
117  }
118 }
119 
121  LauSplineBoundaryType rightBound,
122  Double_t dydx0,
123  Double_t dydxn )
124 {
125  Bool_t updateDerivatives( kFALSE );
126 
127  if ( leftBound_ != leftBound || rightBound_ != rightBound ) {
128  leftBound_ = leftBound;
129  rightBound_ = rightBound;
130  updateDerivatives = kTRUE;
131  }
132 
133  if ( dydx0_ != dydx0 ) {
134  dydx0_ = dydx0;
136  updateDerivatives = kTRUE;
137  }
138 
139  if ( dydxn_ != dydxn ) {
140  dydxn_ = dydxn;
142  updateDerivatives = kTRUE;
143  }
144 
145  if ( updateDerivatives ) {
146  this->calcDerivatives();
147  }
148 }
149 
151 {
152  if ( y_.size() != x_.size() ) {
153  std::cout << "ERROR in Lau1DCubicSpline::init : The number of y-values given does not match the number of x-values"
154  << std::endl;
155  std::cout << " Found " << y_.size() << ", expected "
156  << x_.size() << std::endl;
157  gSystem->Exit( EXIT_FAILURE );
158  }
159 
160  dydx_.insert( dydx_.begin(), nKnots_, 0. );
161 
162  a_.insert( a_.begin(), nKnots_, 0. );
163  b_.insert( b_.begin(), nKnots_, 0. );
164  c_.insert( c_.begin(), nKnots_, 0. );
165  d_.insert( d_.begin(), nKnots_, 0. );
166 
167  this->calcDerivatives();
168 }
169 
171 {
172  switch ( type_ ) {
174  this->calcDerivativesStandard();
175  break;
177  this->calcDerivativesAkima();
178  break;
180  //derivatives not needed for linear interpolation
181  break;
182  }
183 }
184 
186 {
187  // derivatives are determined such that the second derivative is continuous at internal knots
188 
189  // derivatives, k_i, are the solutions to a set of linear equations of the form:
190  // a_i * k_i-1 + b_i * k+i + c_i * k_i+1 = d_i with a_0 = 0, c_n-1 = 0
191 
192  // this is solved using the tridiagonal matrix algorithm as on en.wikipedia.org/wiki/Tridiagonal_matrix_algorithm
193 
194  // first and last equations give boundary conditions
195  // - for natural boundary, require f''(x) = 0 at end knot
196  // - for 'not a knot' boundary, require f'''(x) continuous at second knot
197  // - for clamped boundary, require predefined value of f'(x) at end knot
198 
199  // non-zero values of a_0 and c_n-1 would give cyclic boundary conditions
200  a_[0] = 0.;
201  c_[nKnots_ - 1] = 0.;
202 
203  // set left boundary condition
205  b_[0] = 2. / ( x_[1] - x_[0] );
206  c_[0] = 1. / ( x_[1] - x_[0] );
207  d_[0] = 3. * ( y_[1] - y_[0] ) / ( ( x_[1] - x_[0] ) * ( x_[1] - x_[0] ) );
208 
209  } else if ( leftBound_ == Lau1DCubicSpline::NotAKnot ) {
210  // define the width, h, and the 'slope', delta, of the first cell
211  Double_t h1( x_[1] - x_[0] ), h2( x_[2] - x_[0] );
212  Double_t delta1( ( y_[1] - y_[0] ) / h1 ), delta2( ( y_[2] - y_[1] ) / h2 );
213 
214  // these coefficients can be determined by requiring f'''_0(x_1) = f'''_1(x_1)
215  // the requirement f''_0(x_1) = f''_1(x_1) has been used to remove the dependence on k_2
216  b_[0] = h2;
217  c_[0] = h1 + h2;
218  d_[0] = delta1 * ( 2. * h2 * h2 + 3. * h1 * h2 ) / ( h1 + h2 ) +
219  delta2 * 5. * h1 * h1 / ( h1 + h2 );
220 
221  } else { //Clamped
222  b_[0] = 1.;
223  c_[0] = 0.;
224  d_[0] = dydx0_;
225  }
226 
227  // set right boundary condition
229  a_[nKnots_ - 1] = 1. / ( x_[nKnots_ - 1] - x_[nKnots_ - 2] );
230  b_[nKnots_ - 1] = 2. / ( x_[nKnots_ - 1] - x_[nKnots_ - 2] );
231  d_[nKnots_ - 1] = 3. * ( y_[nKnots_ - 1] - y_[nKnots_ - 2] ) /
232  ( ( x_[nKnots_ - 1] - x_[nKnots_ - 2] ) *
233  ( x_[nKnots_ - 1] - x_[nKnots_ - 2] ) );
234 
235  } else if ( rightBound_ == Lau1DCubicSpline::NotAKnot ) {
236  // define the width, h, and the 'slope', delta, of the last cell
237  Double_t hnm1( x_[nKnots_ - 1] - x_[nKnots_ - 2] ), hnm2( x_[nKnots_ - 2] - x_[nKnots_ - 3] );
238  Double_t deltanm1( ( y_[nKnots_ - 1] - y_[nKnots_ - 2] ) / hnm1 ),
239  deltanm2( ( y_[nKnots_ - 2] - y_[nKnots_ - 3] ) / hnm2 );
240 
241  // these coefficients can be determined by requiring f'''_n-3(x_n-2) = f'''_n-2(x_n-2)
242  // the requirement f''_n-3(x_n-2) = f''_n-2(x_n-2) has been used to remove
243  // the dependence on k_n-3
244  a_[nKnots_ - 1] = hnm2 + hnm1;
245  b_[nKnots_ - 1] = hnm1;
246  d_[nKnots_ - 1] = deltanm2 * hnm1 * hnm1 / ( hnm2 + hnm1 ) +
247  deltanm1 * ( 2. * hnm2 * hnm2 + 3. * hnm2 * hnm1 ) / ( hnm2 + hnm1 );
248 
249  } else { //Clamped
250  a_[nKnots_ - 1] = 0.;
251  b_[nKnots_ - 1] = 1.;
252  d_[nKnots_ - 1] = dydxn_;
253  }
254 
255  // the remaining equations ensure that f_i-1''(x_i) = f''_i(x_i) for all internal knots
256  for ( UInt_t i = 1; i < nKnots_ - 1; ++i ) {
257  a_[i] = 1. / ( x_[i] - x_[i - 1] );
258  b_[i] = 2. / ( x_[i] - x_[i - 1] ) + 2. / ( x_[i + 1] - x_[i] );
259  c_[i] = 1. / ( x_[i + 1] - x_[i] );
260  d_[i] = 3. * ( y_[i] - y_[i - 1] ) / ( ( x_[i] - x_[i - 1] ) * ( x_[i] - x_[i - 1] ) ) +
261  3. * ( y_[i + 1] - y_[i] ) / ( ( x_[i + 1] - x_[i] ) * ( x_[i + 1] - x_[i] ) );
262  }
263 
264  // forward sweep - replace c_i and d_i with
265  //
266  // c'_i = c_i / b_i for i = 0
267  // c_i / (b_i - a_i * c_i-1) otherwise
268  //
269  // and
270  //
271  // d'_i = d_i / b_i for i = 0
272  // (d_i - a_i * d_i-1) / (b_i - a_i * c_i-1) otherwise
273 
274  c_[0] /= b_[0];
275  d_[0] /= b_[0];
276 
277  for ( UInt_t i = 1; i < nKnots_ - 1; ++i ) {
278  c_[i] = c_[i] / ( b_[i] - a_[i] * c_[i - 1] );
279  d_[i] = ( d_[i] - a_[i] * d_[i - 1] ) / ( b_[i] - a_[i] * c_[i - 1] );
280  }
281 
282  d_[nKnots_ - 1] = ( d_[nKnots_ - 1] - a_[nKnots_ - 1] * d_[nKnots_ - 2] ) /
283  ( b_[nKnots_ - 1] - a_[nKnots_ - 1] * c_[nKnots_ - 2] );
284 
285  // back substitution - calculate k_i = dy/dx at each knot
286  //
287  // k_i = d'_i for i = n-1
288  // d'_i - c'_i * k_i+1 otherwise
289 
290  dydx_[nKnots_ - 1] = d_[nKnots_ - 1];
291 
292  for ( Int_t i = nKnots_ - 2; i >= 0; --i ) {
293  dydx_[i] = d_[i] - c_[i] * dydx_[i + 1];
294  }
295 }
296 
298 {
299  //derivatives are calculated according to the Akima method
300  // J.ACM vol. 17 no. 4 pp 589-602
301 
302  Double_t am1( 0. ), an( 0. ), anp1( 0. );
303 
304  // a[i] is the slope of the segment from point i-1 to point i
305  //
306  // n.b. segment 0 is before point 0 and segment n is after point n-1
307  // internal segments are numbered 1 - n-1
308  for ( UInt_t i = 1; i < nKnots_; ++i ) {
309  a_[i] = ( y_[i] - y_[i - 1] ) / ( x_[i] - x_[i - 1] );
310  }
311 
312  // calculate slopes between additional points on each end according to the Akima method
313  // method assumes that the additional points follow a quadratic defined by the last three points
314  // this leads to the relations a[2] - a[1] = a[1] - a[0] = a[0] - a[-1] and a[n-1] - a[n-2] = a[n] - a[n-1] = a[n+1] - a[n]
315  a_[0] = 2 * a_[1] - a_[2];
316  am1 = 2 * a_[0] - a_[1];
317 
318  an = 2 * a_[nKnots_ - 1] - a_[nKnots_ - 2];
319  anp1 = 2 * an - a_[nKnots_ - 1];
320 
321  // b[i] is the weight of a[i] towards dydx[i]
322  // c[i] is the weight of a[i+1] towards dydx[i]
323  // See Appendix A of J.ACM vol. 17 no. 4 pp 589-602 for a justification of these weights
324  b_[0] = TMath::Abs( a_[1] - a_[0] );
325  b_[1] = TMath::Abs( a_[2] - a_[1] );
326 
327  c_[0] = TMath::Abs( a_[0] - am1 );
328  c_[1] = TMath::Abs( a_[1] - a_[0] );
329 
330  for ( UInt_t i = 2; i < nKnots_ - 2; ++i ) {
331  b_[i] = TMath::Abs( a_[i + 2] - a_[i + 1] );
332  c_[i] = TMath::Abs( a_[i] - a_[i - 1] );
333  }
334 
335  b_[nKnots_ - 2] = TMath::Abs( an - a_[nKnots_ - 1] );
336  b_[nKnots_ - 1] = TMath::Abs( anp1 - an );
337 
338  c_[nKnots_ - 2] = TMath::Abs( a_[nKnots_ - 2] - a_[nKnots_ - 3] );
339  c_[nKnots_ - 1] = TMath::Abs( a_[nKnots_ - 1] - a_[nKnots_ - 2] );
340 
341  // dy/dx calculated as the weighted average of a[i] and a[i+1]:
342  // dy/dx_i = ( | a_i+2 - a_i+1 | a_i + | a_i - a_i-1 | a_i+1 ) / ( | a_i+2 - a_i+1 | + | a_i - a_i-1 | )
343  // in the special case a_i-1 == a_i != a_i+1 == a_i+2 this function is undefined so dy/dx is then defined as (a_i + a_i+1) / 2
344  for ( UInt_t i = 0; i < nKnots_ - 2; ++i ) {
345  if ( b_[i] == 0 && c_[i] == 0 ) {
346  dydx_[i] = ( a_[i] + a_[i + 1] ) / 2.;
347  } else {
348  dydx_[i] = ( b_[i] * a_[i] + c_[i] * a_[i + 1] ) / ( b_[i] + c_[i] );
349  }
350  }
351 
352  if ( b_[nKnots_ - 1] == 0 && c_[nKnots_ - 1] == 0 ) {
353  dydx_[nKnots_ - 1] = ( a_[nKnots_ - 1] + an ) / 2.;
354  } else {
355  dydx_[nKnots_ - 1] = ( b_[nKnots_ - 1] * a_[nKnots_ - 1] + c_[nKnots_ - 1] * an ) /
356  ( b_[nKnots_ - 1] + c_[nKnots_ - 1] );
357  }
358 }
void calcDerivativesStandard()
Calculate the first derivatives according to the standard method.
Double_t dydxn_
The gradient at the right boundary for a clamped spline.
void init()
Initialise the class.
Double_t dydx0_
The gradient at the left boundary for a clamped spline.
std::vector< Double_t > y_
The y-value at each knot.
std::vector< Double_t > c_
The 'c' coefficients used to determine the derivatives.
void calcDerivatives()
Calculate the first derivative at each knot.
void updateBoundaryConditions(LauSplineBoundaryType leftBound, LauSplineBoundaryType rightBound, Double_t dydx0=0.0, Double_t dydxn=0.0)
Update the boundary conditions for the spline.
void calcDerivativesAkima()
Calculate the first derivatives according to the Akima method.
void updateType(LauSplineType type)
Update the type of interpolation to perform.
const UInt_t nKnots_
The number of knots in the spline.
std::vector< Double_t > a_
The 'a' coefficients used to determine the derivatives.
void updateYValues(const std::vector< Double_t > &ys)
Update the y-values of the knots.
LauSplineBoundaryType leftBound_
The left-hand boundary condition on the spline.
LauSplineType
Define the allowed interpolation types.
File containing declaration of Lau1DCubicSpline class.
std::vector< Double_t > x_
The x-value at each knot.
LauSplineType type_
The type of interpolation to be performed.
std::vector< Double_t > b_
The 'b' coefficients used to determine the derivatives.
Lau1DCubicSpline(const std::vector< Double_t > &xs, const std::vector< Double_t > &ys, LauSplineType type=Lau1DCubicSpline::StandardSpline, LauSplineBoundaryType leftBound=Lau1DCubicSpline::NotAKnot, LauSplineBoundaryType rightBound=Lau1DCubicSpline::NotAKnot, Double_t dydx0=0.0, Double_t dydxn=0.0)
Constructor.
LauSplineBoundaryType
Define the allowed boundary condition types.
LauSplineBoundaryType rightBound_
The right-hand boundary condition on the spline.
virtual ~Lau1DCubicSpline()
Destructor.
Double_t evaluate(Double_t x) const
Evaluate the function at given point.
std::vector< Double_t > dydx_
The first derivative at each knot.
std::vector< Double_t > d_
The 'd' coefficients used to determine the derivatives.