g4tools  5.4.0
spline
Go to the documentation of this file.
1 // Copyright (C) 2010, Guy Barrand. All rights reserved.
2 // See the file tools.license for terms.
3 
4 #ifndef tools_spline
5 #define tools_spline
6 
7 // From Federico Carminati code found in root-6.08.06/TSpline.h, TSpline.cxx.
8 
9 #include "mnmx"
10 #include <cstddef>
11 #include <vector>
12 #include <ostream>
13 #include <cmath>
14 
15 namespace tools {
16 namespace spline {
17 
18 class base_poly {
19 public:
20  base_poly():fX(0),fY(0) {}
21  base_poly(double x,double y):fX(x),fY(y) {}
22  virtual ~base_poly(){}
23 public:
24  base_poly(base_poly const &a_from):fX(a_from.fX),fY(a_from.fY) {}
25  base_poly& operator=(base_poly const &a_from) {
26  if(this==&a_from) return *this;
27  fX = a_from.fX;
28  fY = a_from.fY;
29  return *this;
30  }
31 public:
32  const double& X() const {return fX;}
33  const double& Y() const {return fY;}
34  double &X() {return fX;}
35  double &Y() {return fY;}
36 protected:
37  double fX; // abscissa
38  double fY; // constant term
39 };
40 
41 class cubic_poly : public base_poly {
42 public:
43  cubic_poly():fB(0), fC(0), fD(0) {}
44  cubic_poly(double x, double y, double b, double c, double d):base_poly(x,y), fB(b), fC(c), fD(d) {}
45 public:
46  cubic_poly(cubic_poly const &a_from)
47  :base_poly(a_from), fB(a_from.fB), fC(a_from.fC), fD(a_from.fD) {}
48  cubic_poly& operator=(cubic_poly const &a_from) {
49  if(this==&a_from) return *this;
50  base_poly::operator=(a_from);
51  fB = a_from.fB;
52  fC = a_from.fC;
53  fD = a_from.fD;
54  return *this;
55  }
56 public:
57  double &B() {return fB;}
58  double &C() {return fC;}
59  double &D() {return fD;}
60  double eval(double x) const {double dx=x-fX;return (fY+dx*(fB+dx*(fC+dx*fD)));}
61 protected:
62  double fB; // first order expansion coefficient : fB*1! is the first derivative at x
63  double fC; // second order expansion coefficient : fC*2! is the second derivative at x
64  double fD; // third order expansion coefficient : fD*3! is the third derivative at x
65 };
66 
67 class quintic_poly : public base_poly {
68 public:
69  quintic_poly():fB(0), fC(0), fD(0), fE(0), fF(0) {}
70  quintic_poly(double x, double y, double b, double c, double d, double e, double f)
71  :base_poly(x,y), fB(b), fC(c), fD(d), fE(e), fF(f) {}
72 public:
73  quintic_poly(quintic_poly const &a_from)
74  :base_poly(a_from)
75  ,fB(a_from.fB),fC(a_from.fC),fD(a_from.fD),fE(a_from.fE),fF(a_from.fF) {}
77  if(this==&a_from) return *this;
78  base_poly::operator=(a_from);
79  fB = a_from.fB;
80  fC = a_from.fC;
81  fD = a_from.fD;
82  fE = a_from.fE;
83  fF = a_from.fF;
84  return *this;
85  }
86 public:
87  double &B() {return fB;}
88  double &C() {return fC;}
89  double &D() {return fD;}
90  double &E() {return fE;}
91  double &F() {return fF;}
92  double eval(double x) const {double dx=x-fX;return (fY+dx*(fB+dx*(fC+dx*(fD+dx*(fE+dx*fF)))));}
93 protected:
94  double fB; // first order expansion coefficient : fB*1! is the first derivative at x
95  double fC; // second order expansion coefficient : fC*2! is the second derivative at x
96  double fD; // third order expansion coefficient : fD*3! is the third derivative at x
97  double fE; // fourth order expansion coefficient : fE*4! is the fourth derivative at x
98  double fF; // fifth order expansion coefficient : fF*5! is the fifth derivative at x
99 };
100 
104 class base_spline {
105 protected:
106  base_spline(std::ostream& a_out):m_out(a_out), fDelta(-1), fXmin(0), fXmax(0), fNp(0), fKstep(false) {}
107 public:
108  base_spline(std::ostream& a_out,double delta, double xmin, double xmax, size_t np, bool step)
109  :m_out(a_out),fDelta(delta), fXmin(xmin),fXmax(xmax), fNp(np), fKstep(step)
110  {}
111  virtual ~base_spline() {}
112 protected:
113  base_spline(const base_spline& a_from)
114  :m_out(a_from.m_out)
115  ,fDelta(a_from.fDelta),fXmin(a_from.fXmin),fXmax(a_from.fXmax),fNp(a_from.fNp),fKstep(a_from.fKstep) {}
117  if(this==&a_from) return *this;
118  fDelta=a_from.fDelta;
119  fXmin=a_from.fXmin;
120  fXmax=a_from.fXmax;
121  fNp=a_from.fNp;
122  fKstep=a_from.fKstep;
123  return *this;
124  }
125 protected:
126  std::ostream& m_out;
127  double fDelta; // Distance between equidistant knots
128  double fXmin; // Minimum value of abscissa
129  double fXmax; // Maximum value of abscissa
130  size_t fNp; // Number of knots
131  bool fKstep; // True of equidistant knots
132 };
133 
134 
136 // //
137 // cubic //
138 // //
139 // Class to create third splines to interpolate knots //
140 // Arbitrary conditions can be introduced for first and second //
141 // derivatives at beginning and ending points //
142 // //
144 
145 class cubic : public base_spline {
146 protected:
147  cubic(std::ostream& a_out) : base_spline(a_out) , fPoly(0), fValBeg(0), fValEnd(0), fBegCond(-1), fEndCond(-1) {}
148 public:
149  cubic(std::ostream& a_out,size_t a_n,double a_x[], double a_y[], double a_valbeg = 0, double a_valend = 0)
150  :base_spline(a_out,-1,0,0,a_n,false)
151  ,fValBeg(a_valbeg), fValEnd(a_valend), fBegCond(0), fEndCond(0)
152  {
153  if(!a_n) {
154  m_out << "tools::spline::cubic : a_np is null." << std::endl;
155  return;
156  }
157  fXmin = a_x[0];
158  fXmax = a_x[a_n-1];
159  fPoly.resize(a_n);
160  for (size_t i=0; i<a_n; ++i) {
161  fPoly[i].X() = a_x[i];
162  fPoly[i].Y() = a_y[i];
163  }
164  build_coeff(); // Build the spline coefficients
165  }
166 public:
167  cubic(const cubic& a_from)
168  :base_spline(a_from)
169  ,fPoly(a_from.fPoly),fValBeg(a_from.fValBeg),fValEnd(a_from.fValEnd),fBegCond(a_from.fBegCond),fEndCond(a_from.fEndCond)
170  {}
171  cubic& operator=(const cubic& a_from) {
172  if(this==&a_from) return *this;
173  base_spline::operator=(a_from);
174  fPoly = a_from.fPoly;
175  fValBeg=a_from.fValBeg;
176  fValEnd=a_from.fValEnd;
177  fBegCond=a_from.fBegCond;
178  fEndCond=a_from.fEndCond;
179  return *this;
180  }
181 public:
182  double eval(double x) const {
183  if(!fNp) return 0;
184  // Eval this spline at x
185  size_t klow = find_x(x);
186  if ( (fNp > 1) && (klow >= (fNp-1))) klow = fNp-2; //see: https://savannah.cern.ch/bugs/?71651
187  return fPoly[klow].eval(x);
188  }
189 protected:
190  template<typename T>
191  static int TMath_Nint(T x) {
192  // Round to nearest integer. Rounds half integers to the nearest even integer.
193  int i;
194  if (x >= 0) {
195  i = int(x + 0.5);
196  if ( i & 1 && x + 0.5 == T(i) ) i--;
197  } else {
198  i = int(x - 0.5);
199  if ( i & 1 && x - 0.5 == T(i) ) i++;
200  }
201  return i;
202  }
203  static int TMath_FloorNint(double x) { return TMath_Nint(::floor(x)); }
204 
205  size_t find_x(double x) const {
206  int klow=0, khig=fNp-1;
207  //
208  // If out of boundaries, extrapolate
209  // It may be badly wrong
210  if(x<=fXmin) klow=0;
211  else if(x>=fXmax) klow=khig;
212  else {
213  if(fKstep) { // Equidistant knots, use histogramming :
214  klow = TMath_FloorNint((x-fXmin)/fDelta);
215  // Correction for rounding errors
216  if (x < fPoly[klow].X())
217  klow = mx<int>(klow-1,0);
218  else if (klow < khig) {
219  if (x > fPoly[klow+1].X()) ++klow;
220  }
221  } else {
222  int khalf;
223  //
224  // Non equidistant knots, binary search
225  while((khig-klow)>1) {
226  khalf = (klow+khig)/2;
227  if(x>fPoly[khalf].X()) klow=khalf;
228  else khig=khalf;
229  }
230  //
231  // This could be removed, sanity check
232  if( (x<fPoly[klow].X()) || (fPoly[klow+1].X()<x) ) {
233  m_out << "tools::spline::cubic::find_x : Binary search failed"
234  << " x(" << klow << ") = " << fPoly[klow].X() << " < x= " << x
235  << " < x(" << klow+1 << ") = " << fPoly[klow+1].X() << "."
236  << "." << std::endl;
237  }
238  }
239  }
240  return klow;
241  }
242 
243  void build_coeff() {
273 
274  int j, l;
275  double divdf1,divdf3,dtau,g=0;
276  // ***** a tridiagonal linear system for the unknown slopes s(i) of
277  // f at tau(i), i=1,...,n, is generated and then solved by gauss elim-
278  // ination, with s(i) ending up in c(2,i), all i.
279  // c(3,.) and c(4,.) are used initially for temporary storage.
280  l = fNp-1;
281  // compute first differences of x sequence and store in C also,
282  // compute first divided difference of data and store in D.
283  {for (size_t m=1; m<fNp ; ++m) {
284  fPoly[m].C() = fPoly[m].X() - fPoly[m-1].X();
285  fPoly[m].D() = (fPoly[m].Y() - fPoly[m-1].Y())/fPoly[m].C();
286  }}
287  // construct first equation from the boundary condition, of the form
288  // D[0]*s[0] + C[0]*s[1] = B[0]
289  if(fBegCond==0) {
290  if(fNp == 2) {
291  // no condition at left end and n = 2.
292  fPoly[0].D() = 1.;
293  fPoly[0].C() = 1.;
294  fPoly[0].B() = 2.*fPoly[1].D();
295  } else {
296  // not-a-knot condition at left end and n .gt. 2.
297  fPoly[0].D() = fPoly[2].C();
298  fPoly[0].C() = fPoly[1].C() + fPoly[2].C();
299  fPoly[0].B() = ((fPoly[1].C()+2.*fPoly[0].C())*fPoly[1].D()*fPoly[2].C()+
300  fPoly[1].C()*fPoly[1].C()*fPoly[2].D())/fPoly[0].C();
301  }
302  } else if (fBegCond==1) {
303  // slope prescribed at left end.
304  fPoly[0].B() = fValBeg;
305  fPoly[0].D() = 1.;
306  fPoly[0].C() = 0.;
307  } else if (fBegCond==2) {
308  // second derivative prescribed at left end.
309  fPoly[0].D() = 2.;
310  fPoly[0].C() = 1.;
311  fPoly[0].B() = 3.*fPoly[1].D() - fPoly[1].C()/2.*fValBeg;
312  }
313  bool forward_gauss_elimination = true;
314  if(fNp > 2) {
315  // if there are interior knots, generate the corresp. equations and car-
316  // ry out the forward pass of gauss elimination, after which the m-th
317  // equation reads D[m]*s[m] + C[m]*s[m+1] = B[m].
318  {for (int m=1; m<l; ++m) {
319  g = -fPoly[m+1].C()/fPoly[m-1].D();
320  fPoly[m].B() = g*fPoly[m-1].B() + 3.*(fPoly[m].C()*fPoly[m+1].D()+fPoly[m+1].C()*fPoly[m].D());
321  fPoly[m].D() = g*fPoly[m-1].C() + 2.*(fPoly[m].C() + fPoly[m+1].C());
322  }}
323  // construct last equation from the second boundary condition, of the form
324  // (-g*D[n-2])*s[n-2] + D[n-1]*s[n-1] = B[n-1]
325  // if slope is prescribed at right end, one can go directly to back-
326  // substitution, since c array happens to be set up just right for it
327  // at this point.
328  if(fEndCond == 0) {
329  if (fNp > 3 || fBegCond != 0) {
330  // not-a-knot and n .ge. 3, and either n.gt.3 or also not-a-knot at
331  // left end point.
332  g = fPoly[fNp-2].C() + fPoly[fNp-1].C();
333  fPoly[fNp-1].B() = ((fPoly[fNp-1].C()+2.*g)*fPoly[fNp-1].D()*fPoly[fNp-2].C()
334  + fPoly[fNp-1].C()*fPoly[fNp-1].C()*(fPoly[fNp-2].Y()-fPoly[fNp-3].Y())/fPoly[fNp-2].C())/g;
335  g = -g/fPoly[fNp-2].D();
336  fPoly[fNp-1].D() = fPoly[fNp-2].C();
337  } else {
338  // either (n=3 and not-a-knot also at left) or (n=2 and not not-a-
339  // knot at left end point).
340  fPoly[fNp-1].B() = 2.*fPoly[fNp-1].D();
341  fPoly[fNp-1].D() = 1.;
342  g = -1./fPoly[fNp-2].D();
343  }
344  } else if (fEndCond == 1) {
345  fPoly[fNp-1].B() = fValEnd;
346  forward_gauss_elimination = false;
347  } else if (fEndCond == 2) {
348  // second derivative prescribed at right endpoint.
349  fPoly[fNp-1].B() = 3.*fPoly[fNp-1].D() + fPoly[fNp-1].C()/2.*fValEnd;
350  fPoly[fNp-1].D() = 2.;
351  g = -1./fPoly[fNp-2].D();
352  }
353  } else {
354  if(fEndCond == 0) {
355  if (fBegCond > 0) {
356  // either (n=3 and not-a-knot also at left) or (n=2 and not not-a-
357  // knot at left end point).
358  fPoly[fNp-1].B() = 2.*fPoly[fNp-1].D();
359  fPoly[fNp-1].D() = 1.;
360  g = -1./fPoly[fNp-2].D();
361  } else {
362  // not-a-knot at right endpoint and at left endpoint and n = 2.
363  fPoly[fNp-1].B() = fPoly[fNp-1].D();
364  forward_gauss_elimination = false;
365  }
366  } else if(fEndCond == 1) {
367  fPoly[fNp-1].B() = fValEnd;
368  forward_gauss_elimination = false;
369  } else if(fEndCond == 2) {
370  // second derivative prescribed at right endpoint.
371  fPoly[fNp-1].B() = 3.*fPoly[fNp-1].D() + fPoly[fNp-1].C()/2.*fValEnd;
372  fPoly[fNp-1].D() = 2.;
373  g = -1./fPoly[fNp-2].D();
374  }
375  }
376  // complete forward pass of gauss elimination.
377  if(forward_gauss_elimination) {
378  fPoly[fNp-1].D() = g*fPoly[fNp-2].C() + fPoly[fNp-1].D();
379  fPoly[fNp-1].B() = (g*fPoly[fNp-2].B() + fPoly[fNp-1].B())/fPoly[fNp-1].D();
380  }
381  // carry out back substitution
382  j = l-1;
383  do {
384  fPoly[j].B() = (fPoly[j].B() - fPoly[j].C()*fPoly[j+1].B())/fPoly[j].D();
385  --j;
386  } while (j>=0);
387  // ****** generate cubic coefficients in each interval, i.e., the deriv.s
388  // at its left endpoint, from value and slope at its endpoints.
389  for (size_t i=1; i<fNp; ++i) {
390  dtau = fPoly[i].C();
391  divdf1 = (fPoly[i].Y() - fPoly[i-1].Y())/dtau;
392  divdf3 = fPoly[i-1].B() + fPoly[i].B() - 2.*divdf1;
393  fPoly[i-1].C() = (divdf1 - fPoly[i-1].B() - divdf3)/dtau;
394  fPoly[i-1].D() = (divdf3/dtau)/dtau;
395  }
396  }
397 
398 protected:
399  std::vector<cubic_poly> fPoly; //[fNp] Array of polynomial terms
400  double fValBeg; // Initial value of first or second derivative
401  double fValEnd; // End value of first or second derivative
402  int fBegCond; // 0=no beg cond, 1=first derivative, 2=second derivative
403  int fEndCond; // 0=no end cond, 1=first derivative, 2=second derivative
404 };
405 
407 // //
408 // quintic //
409 // //
410 // Class to create quintic natural splines to interpolate knots //
411 // Arbitrary conditions can be introduced for first and second //
412 // derivatives using double knots (see build_coeff) for more on this. //
413 // Double knots are automatically introduced at ending points //
414 // //
416 class quintic : public base_spline {
417 protected:
418  quintic(std::ostream& a_out):base_spline(a_out),fPoly() {}
419 public:
420  quintic(std::ostream& a_out,size_t a_n ,double a_x[], double a_y[])
421  :base_spline(a_out,-1,0,0,a_n,false) {
422  if(!a_n) {
423  m_out << "tools::spline::quintic : a_np is null." << std::endl;
424  return;
425  }
426  fXmin = a_x[0];
427  fXmax = a_x[a_n-1];
428  fPoly.resize(fNp);
429  for (size_t i=0; i<a_n; ++i) {
430  fPoly[i].X() = a_x[i];
431  fPoly[i].Y() = a_y[i];
432  }
433  build_coeff(); // Build the spline coefficients.
434  }
435 public:
436  quintic(const quintic& a_from):base_spline(a_from),fPoly(a_from.fPoly) {}
437  quintic& operator=(const quintic& a_from) {
438  if(this==&a_from) return *this;
439  base_spline::operator=(a_from);
440  fPoly = a_from.fPoly;
441  return *this;
442  }
443 public:
444  double eval(double x) const {if(!fNp) return 0;size_t klow=find_x(x);return fPoly[klow].eval(x);}
445 
446 protected:
447  size_t find_x(double x) const {
448  int klow=0;
449 
450  // If out of boundaries, extrapolate
451  // It may be badly wrong
452  if(x<=fXmin) klow=0;
453  else if(x>=fXmax) klow=fNp-1;
454  else {
455  if(fKstep) { // Equidistant knots, use histogramming :
456  klow = mn<int>(int((x-fXmin)/fDelta),fNp-1);
457  } else {
458  int khig=fNp-1;
459  int khalf;
460  // Non equidistant knots, binary search
461  while((khig-klow)>1) {
462  khalf = (klow+khig)/2;
463  if(x>fPoly[khalf].X()) klow=khalf;
464  else khig=khalf;
465  }
466  }
467 
468  // This could be removed, sanity check
469  if( (x<fPoly[klow].X()) || (fPoly[klow+1].X()<x) ) {
470  m_out << "tools::spline::quintic::find_x : Binary search failed"
471  << " x(" << klow << ") = " << fPoly[klow].X() << " < x= " << x
472  << " < x(" << klow+1<< ") = " << fPoly[klow+1].X() << "."
473  << std::endl;
474  }
475  }
476  return klow;
477  }
478 
479  void build_coeff() {
545 
546  size_t i, m;
547  double pqqr, p, q, r, s, t, u, v,
548  b1, p2, p3, q2, q3, r2, pq, pr, qr;
549 
550  if (fNp <= 2) return;
551 
552  // coefficients of a positive definite, pentadiagonal matrix,
553  // stored in D, E, F from 1 to n-3.
554  m = fNp-2;
555  q = fPoly[1].X()-fPoly[0].X();
556  r = fPoly[2].X()-fPoly[1].X();
557  q2 = q*q;
558  r2 = r*r;
559  qr = q+r;
560  fPoly[0].D() = fPoly[0].E() = 0;
561  if (q) fPoly[1].D() = q*6.*q2/(qr*qr);
562  else fPoly[1].D() = 0;
563 
564  if (m > 1) {
565  for (i = 1; i < m; ++i) {
566  p = q;
567  q = r;
568  r = fPoly[i+2].X()-fPoly[i+1].X();
569  p2 = q2;
570  q2 = r2;
571  r2 = r*r;
572  pq = qr;
573  qr = q+r;
574  if (q) {
575  q3 = q2*q;
576  pr = p*r;
577  pqqr = pq*qr;
578  fPoly[i+1].D() = q3*6./(qr*qr);
579  fPoly[i].D() += (q+q)*(pr*15.*pr+(p+r)*q
580  *(pr* 20.+q2*7.)+q2*
581  ((p2+r2)*8.+pr*21.+q2+q2))/(pqqr*pqqr);
582  fPoly[i-1].D() += q3*6./(pq*pq);
583  fPoly[i].E() = q2*(p*qr+pq*3.*(qr+r+r))/(pqqr*qr);
584  fPoly[i-1].E() += q2*(r*pq+qr*3.*(pq+p+p))/(pqqr*pq);
585  fPoly[i-1].F() = q3/pqqr;
586  } else
587  fPoly[i+1].D() = fPoly[i].E() = fPoly[i-1].F() = 0;
588  }
589  }
590  if (r) fPoly[m-1].D() += r*6.*r2/(qr*qr);
591 
592  // First and second order divided differences of the given function
593  // values, stored in b from 2 to n and in c from 3 to n
594  // respectively. care is taken of double and triple knots.
595  for (i = 1; i < fNp; ++i) {
596  if (fPoly[i].X() != fPoly[i-1].X()) {
597  fPoly[i].B() =
598  (fPoly[i].Y()-fPoly[i-1].Y())/(fPoly[i].X()-fPoly[i-1].X());
599  } else {
600  fPoly[i].B() = fPoly[i].Y();
601  fPoly[i].Y() = fPoly[i-1].Y();
602  }
603  }
604  for (i = 2; i < fNp; ++i) {
605  if (fPoly[i].X() != fPoly[i-2].X()) {
606  fPoly[i].C() =
607  (fPoly[i].B()-fPoly[i-1].B())/(fPoly[i].X()-fPoly[i-2].X());
608  } else {
609  fPoly[i].C() = fPoly[i].B()*.5;
610  fPoly[i].B() = fPoly[i-1].B();
611  }
612  }
613 
614  // Solve the linear system with c(i+2) - c(i+1) as right-hand side.
615  if (m > 1) {
616  p=fPoly[0].C()=fPoly[m-1].E()=fPoly[0].F()
617  =fPoly[m-2].F()=fPoly[m-1].F()=0;
618  fPoly[1].C() = fPoly[3].C()-fPoly[2].C();
619  fPoly[1].D() = 1./fPoly[1].D();
620 
621  if (m > 2) {
622  for (i = 2; i < m; ++i) {
623  q = fPoly[i-1].D()*fPoly[i-1].E();
624  fPoly[i].D() = 1./(fPoly[i].D()-p*fPoly[i-2].F()-q*fPoly[i-1].E());
625  fPoly[i].E() -= q*fPoly[i-1].F();
626  fPoly[i].C() = fPoly[i+2].C()-fPoly[i+1].C()-p*fPoly[i-2].C()
627  -q*fPoly[i-1].C();
628  p = fPoly[i-1].D()*fPoly[i-1].F();
629  }
630  }
631  }
632 
633  fPoly[fNp-2].C() = fPoly[fNp-1].C() = 0;
634  if (fNp > 3)
635  for (i=fNp-3; i > 0; --i)
636  fPoly[i].C() = (fPoly[i].C()-fPoly[i].E()*fPoly[i+1].C()
637  -fPoly[i].F()*fPoly[i+2].C())*fPoly[i].D();
638 
639  // Integrate the third derivative of s(x)
640  m = fNp-1;
641  q = fPoly[1].X()-fPoly[0].X();
642  r = fPoly[2].X()-fPoly[1].X();
643  b1 = fPoly[1].B();
644  q3 = q*q*q;
645  qr = q+r;
646  if (qr) {
647  v = fPoly[1].C()/qr;
648  t = v;
649  } else
650  v = t = 0;
651  if (q) fPoly[0].F() = v/q;
652  else fPoly[0].F() = 0;
653  for (i = 1; i < m; ++i) {
654  p = q;
655  q = r;
656  if (i != m-1) r = fPoly[i+2].X()-fPoly[i+1].X();
657  else r = 0;
658  p3 = q3;
659  q3 = q*q*q;
660  pq = qr;
661  qr = q+r;
662  s = t;
663  if (qr) t = (fPoly[i+1].C()-fPoly[i].C())/qr;
664  else t = 0;
665  u = v;
666  v = t-s;
667  if (pq) {
668  fPoly[i].F() = fPoly[i-1].F();
669  if (q) fPoly[i].F() = v/q;
670  fPoly[i].E() = s*5.;
671  fPoly[i].D() = (fPoly[i].C()-q*s)*10;
672  fPoly[i].C() =
673  fPoly[i].D()*(p-q)+(fPoly[i+1].B()-fPoly[i].B()+(u-fPoly[i].E())*
674  p3-(v+fPoly[i].E())*q3)/pq;
675  fPoly[i].B() = (p*(fPoly[i+1].B()-v*q3)+q*(fPoly[i].B()-u*p3))/pq-p
676  *q*(fPoly[i].D()+fPoly[i].E()*(q-p));
677  } else {
678  fPoly[i].C() = fPoly[i-1].C();
679  fPoly[i].D() = fPoly[i].E() = fPoly[i].F() = 0;
680  }
681  }
682 
683  // End points x(1) and x(n)
684  p = fPoly[1].X()-fPoly[0].X();
685  s = fPoly[0].F()*p*p*p;
686  fPoly[0].E() = fPoly[0].D() = 0;
687  fPoly[0].C() = fPoly[1].C()-s*10;
688  fPoly[0].B() = b1-(fPoly[0].C()+s)*p;
689 
690  q = fPoly[fNp-1].X()-fPoly[fNp-2].X();
691  t = fPoly[fNp-2].F()*q*q*q;
692  fPoly[fNp-1].E() = fPoly[fNp-1].D() = 0;
693  fPoly[fNp-1].C() = fPoly[fNp-2].C()+t*10;
694  fPoly[fNp-1].B() += (fPoly[fNp-1].C()-t)*q;
695  }
696 protected:
697  std::vector<quintic_poly> fPoly; //[fNp] Array of polynomial terms
698 };
699 
700 }}
701 
702 #endif
tools::spline::cubic::find_x
size_t find_x(double x) const
Definition: spline:205
tools::spline::cubic_poly::operator=
cubic_poly & operator=(cubic_poly const &a_from)
Definition: spline:48
tools::spline::quintic_poly::quintic_poly
quintic_poly()
Definition: spline:69
tools::spline::cubic_poly
Definition: spline:41
tools::spline::base_spline::fDelta
double fDelta
Definition: spline:127
tools::spline::cubic::eval
double eval(double x) const
Definition: spline:182
tools::spline::cubic::cubic
cubic(std::ostream &a_out)
Definition: spline:147
tools::spline::cubic_poly::fB
double fB
Definition: spline:62
tools::spline::cubic_poly::C
double & C()
Definition: spline:58
tools::spline::cubic::fValEnd
double fValEnd
Definition: spline:401
tools::spline::base_spline::m_out
std::ostream & m_out
Definition: spline:126
tools::spline::base_poly::operator=
base_poly & operator=(base_poly const &a_from)
Definition: spline:25
tools::spline::cubic::TMath_FloorNint
static int TMath_FloorNint(double x)
Definition: spline:203
tools::spline::quintic_poly::E
double & E()
Definition: spline:90
tools::spline::base_poly::base_poly
base_poly()
Definition: spline:20
tools::spline::quintic::quintic
quintic(std::ostream &a_out)
Definition: spline:418
tools::spline::base_poly::X
double & X()
Definition: spline:34
tools::spline::quintic_poly::D
double & D()
Definition: spline:89
tools::spline::base_spline::operator=
base_spline & operator=(const base_spline &a_from)
Definition: spline:116
tools::spline::quintic::eval
double eval(double x) const
Definition: spline:444
tools::spline::base_spline::~base_spline
virtual ~base_spline()
Definition: spline:111
tools::spline::quintic_poly::quintic_poly
quintic_poly(quintic_poly const &a_from)
Definition: spline:73
tools::spline::quintic_poly::fF
double fF
Definition: spline:98
tools::spline::cubic_poly::cubic_poly
cubic_poly(cubic_poly const &a_from)
Definition: spline:46
tools::spline::base_spline::fNp
size_t fNp
Definition: spline:130
tools::spline::quintic_poly::fD
double fD
Definition: spline:96
tools::spline::cubic_poly::B
double & B()
Definition: spline:57
tools::spline::cubic::fPoly
std::vector< cubic_poly > fPoly
Definition: spline:399
tools::spline::quintic::quintic
quintic(const quintic &a_from)
Definition: spline:436
tools::spline::base_poly::X
const double & X() const
Definition: spline:32
tools::spline::cubic::TMath_Nint
static int TMath_Nint(T x)
Definition: spline:191
tools::spline::base_spline::base_spline
base_spline(std::ostream &a_out, double delta, double xmin, double xmax, size_t np, bool step)
Definition: spline:108
tools::spline::base_poly
Definition: spline:18
tools::spline::base_spline::base_spline
base_spline(std::ostream &a_out)
Definition: spline:106
tools::spline::base_poly::fX
double fX
Definition: spline:37
tools::spline::cubic_poly::fD
double fD
Definition: spline:64
tools::spline::cubic::operator=
cubic & operator=(const cubic &a_from)
Definition: spline:171
tools::spline::cubic::cubic
cubic(const cubic &a_from)
Definition: spline:167
tools::spline::quintic_poly::B
double & B()
Definition: spline:87
tools::spline::quintic::quintic
quintic(std::ostream &a_out, size_t a_n, double a_x[], double a_y[])
Definition: spline:420
tools::spline::base_spline::fKstep
bool fKstep
Definition: spline:131
tools::spline::base_poly::Y
double & Y()
Definition: spline:35
tools::spline::base_poly::Y
const double & Y() const
Definition: spline:33
tools::spline::quintic::operator=
quintic & operator=(const quintic &a_from)
Definition: spline:437
tools::spline::cubic::fValBeg
double fValBeg
Definition: spline:400
tools::spline::cubic::build_coeff
void build_coeff()
Definition: spline:243
tools::spline::cubic_poly::cubic_poly
cubic_poly()
Definition: spline:43
tools::spline::cubic_poly::cubic_poly
cubic_poly(double x, double y, double b, double c, double d)
Definition: spline:44
tools::spline::base_spline
Definition: spline:104
tools::spline::base_spline::fXmin
double fXmin
Definition: spline:128
tools::spline::quintic_poly::fB
double fB
Definition: spline:94
tools::spline::quintic_poly::operator=
quintic_poly & operator=(quintic_poly const &a_from)
Definition: spline:76
tools::spline::base_poly::base_poly
base_poly(base_poly const &a_from)
Definition: spline:24
mnmx
tools::spline::cubic_poly::eval
double eval(double x) const
Definition: spline:60
tools::spline::base_spline::fXmax
double fXmax
Definition: spline:129
tools::spline::base_spline::base_spline
base_spline(const base_spline &a_from)
Definition: spline:113
tools::spline::base_poly::base_poly
base_poly(double x, double y)
Definition: spline:21
tools
inlined C code : ///////////////////////////////////
Definition: aida_ntuple:26
tools::spline::quintic::fPoly
std::vector< quintic_poly > fPoly
Definition: spline:697
tools::spline::cubic::fEndCond
int fEndCond
Definition: spline:403
tools::spline::quintic_poly::fC
double fC
Definition: spline:95
tools::spline::quintic_poly::fE
double fE
Definition: spline:97
tools::spline::quintic::build_coeff
void build_coeff()
Definition: spline:479
tools::spline::cubic_poly::D
double & D()
Definition: spline:59
tools::spline::quintic_poly::eval
double eval(double x) const
Definition: spline:92
tools::spline::quintic_poly::quintic_poly
quintic_poly(double x, double y, double b, double c, double d, double e, double f)
Definition: spline:70
tools::spline::cubic_poly::fC
double fC
Definition: spline:63
tools::spline::cubic::cubic
cubic(std::ostream &a_out, size_t a_n, double a_x[], double a_y[], double a_valbeg=0, double a_valend=0)
Definition: spline:149
tools::spline::cubic::fBegCond
int fBegCond
Definition: spline:402
tools::spline::quintic_poly
Definition: spline:67
tools::spline::base_poly::~base_poly
virtual ~base_poly()
Definition: spline:22
tools::spline::base_poly::fY
double fY
Definition: spline:38
tools::spline::cubic
Definition: spline:145
tools::spline::quintic::find_x
size_t find_x(double x) const
Definition: spline:447
tools::spline::quintic_poly::F
double & F()
Definition: spline:91
tools::spline::quintic_poly::C
double & C()
Definition: spline:88
tools::spline::quintic
Definition: spline:416