ROOT logo
// @(#)root/hist:$Id$
// Author: Federico Carminati   28/02/2000

/*************************************************************************
 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/

//////////////////////////////////////////////////////////////////////////
//                                                                      //
// TSpline                                                              //
//                                                                      //
// Base class for spline implementation containing the Draw/Paint       //
// methods                                                              //
//                                                                      //
//////////////////////////////////////////////////////////////////////////

#include "TROOT.h"
#include "TSpline.h"
#include "TVirtualPad.h"
#include "TH1.h"
#include "TF1.h"
#include "TSystem.h"
#include "Riostream.h"
#include "TClass.h"
#include "TMath.h"

ClassImp(TSplinePoly)
ClassImp(TSplinePoly3)
ClassImp(TSplinePoly5)
ClassImp(TSpline3)
ClassImp(TSpline5)
ClassImp(TSpline)


//______________________________________________________________________________
TSpline::TSpline(const TSpline &sp) :
  TNamed(sp),
  TAttLine(sp),
  TAttFill(sp),
  TAttMarker(sp),
  fDelta(sp.fDelta),
  fXmin(sp.fXmin),
  fXmax(sp.fXmax),
  fNp(sp.fNp),
  fKstep(sp.fKstep),
  fHistogram(0),
  fGraph(0),
  fNpx(sp.fNpx)
{
   //copy constructor
}


//______________________________________________________________________________
TSpline::~TSpline()
{
   //destructor
   if(fHistogram) delete fHistogram;
   if(fGraph) delete fGraph;
}


//______________________________________________________________________________
TSpline& TSpline::operator=(const TSpline &sp)
{
   //assignment operator
   if(this!=&sp) {
      TNamed::operator=(sp);
      TAttLine::operator=(sp);
      TAttFill::operator=(sp);
      TAttMarker::operator=(sp);
      fDelta=sp.fDelta;
      fXmin=sp.fXmin;
      fXmax=sp.fXmax;
      fNp=sp.fNp;
      fKstep=sp.fKstep;
      fHistogram=0;
      fGraph=0;
      fNpx=sp.fNpx;
   }
   return *this;
}


//______________________________________________________________________________
void TSpline::Draw(Option_t *option)
{
   // Draw this function with its current attributes.
   //
   // Possible option values are:
   //   "SAME"  superimpose on top of existing picture
   //   "L"     connect all computed points with a straight line
   //   "C"     connect all computed points with a smooth curve.
   //   "P"     add a polymarker at each knot
   //
   // Note that the default value is "L". Therefore to draw on top
   // of an existing picture, specify option "LSAME"

   TString opt = option;
   opt.ToLower();
   if (gPad && !opt.Contains("same")) gPad->Clear();

   AppendPad(option);
}


//______________________________________________________________________________
Int_t TSpline::DistancetoPrimitive(Int_t px, Int_t py)
{
   // Compute distance from point px,py to a spline.

   if (!fHistogram) return 999;
   return fHistogram->DistancetoPrimitive(px, py);
}


//______________________________________________________________________________
void TSpline::ExecuteEvent(Int_t event, Int_t px, Int_t py)
{
   // Execute action corresponding to one event.

   if (!fHistogram) return;
   fHistogram->ExecuteEvent(event, px, py);
}


//______________________________________________________________________________
void TSpline::Paint(Option_t *option)
{
   // Paint this function with its current attributes.

   Int_t i;
   Double_t xv;

   TString opt = option;
   opt.ToLower();
   Double_t xmin, xmax, pmin, pmax;
   pmin = gPad->PadtoX(gPad->GetUxmin());
   pmax = gPad->PadtoX(gPad->GetUxmax());
   xmin = fXmin;
   xmax = fXmax;
   if (opt.Contains("same")) {
      if (xmax < pmin) return;  // Otto: completely outside
      if (xmin > pmax) return;
      if (xmin < pmin) xmin = pmin;
      if (xmax > pmax) xmax = pmax;
   } else {
      gPad->Clear();
   }

   //  Create a temporary histogram and fill each channel with the function value
   if (fHistogram)
      if ((!gPad->GetLogx() && fHistogram->TestBit(TH1::kLogX)) ||
          (gPad->GetLogx() && !fHistogram->TestBit(TH1::kLogX)))
         { delete fHistogram; fHistogram = 0;}

   if (fHistogram) {
      //if (xmin != fXmin || xmax != fXmax)
      fHistogram->GetXaxis()->SetLimits(xmin,xmax);
   } else {
   //      if logx, we must bin in logx and not in x !!!
   //      otherwise if several decades, one gets crazy results
      if (xmin > 0 && gPad->GetLogx()) {
         Double_t *xbins  = new Double_t[fNpx+1];
         Double_t xlogmin = TMath::Log10(xmin);
         Double_t xlogmax = TMath::Log10(xmax);
         Double_t dlogx   = (xlogmax-xlogmin)/((Double_t)fNpx);
         for (i=0;i<=fNpx;i++) {
            xbins[i] = gPad->PadtoX(xlogmin+ i*dlogx);
         }
         fHistogram = new TH1F("Spline",GetTitle(),fNpx,xbins);
         fHistogram->SetBit(TH1::kLogX);
         delete [] xbins;
      } else {
         fHistogram = new TH1F("Spline",GetTitle(),fNpx,xmin,xmax);
      }
      if (!fHistogram) return;
      fHistogram->SetDirectory(0);
   }
   for (i=1;i<=fNpx;i++) {
      xv = fHistogram->GetBinCenter(i);
      fHistogram->SetBinContent(i,this->Eval(xv));
   }

   // Copy Function attributes to histogram attributes
   fHistogram->SetBit(TH1::kNoStats);
   fHistogram->SetLineColor(GetLineColor());
   fHistogram->SetLineStyle(GetLineStyle());
   fHistogram->SetLineWidth(GetLineWidth());
   fHistogram->SetFillColor(GetFillColor());
   fHistogram->SetFillStyle(GetFillStyle());
   fHistogram->SetMarkerColor(GetMarkerColor());
   fHistogram->SetMarkerStyle(GetMarkerStyle());
   fHistogram->SetMarkerSize(GetMarkerSize());

   //  Draw the histogram
   //  but first strip off the 'p' option if any
   char *o = (char *) opt.Data();
   Int_t j=0;
   i=0;
   Bool_t graph=kFALSE;
   do
      if(o[i]=='p') graph=kTRUE ; else o[j++]=o[i];
   while(o[i++]);
   if (opt.Length() == 0 ) fHistogram->Paint("lf");
   else if (opt == "same") fHistogram->Paint("lfsame");
   else                    fHistogram->Paint(opt.Data());

   // Think about the graph, if demanded
   if(graph) {
      if(!fGraph) {
         Double_t *xx = new Double_t[fNp];
         Double_t *yy = new Double_t[fNp];
         for(i=0; i<fNp; ++i)
            GetKnot(i,xx[i],yy[i]);
         fGraph=new TGraph(fNp,xx,yy);
         delete [] xx;
         delete [] yy;
      }
      fGraph->SetMarkerColor(GetMarkerColor());
      fGraph->SetMarkerStyle(GetMarkerStyle());
      fGraph->SetMarkerSize(GetMarkerSize());
      fGraph->Paint("p");
   }
}


//______________________________________________________________________________
void TSpline::Streamer(TBuffer &R__b)
{
   // Stream an object of class TSpline.

   if (R__b.IsReading()) {
      UInt_t R__s, R__c;
      Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
      if (R__v > 1) {
         R__b.ReadClassBuffer(TSpline::Class(), this, R__v, R__s, R__c);
         return;
      }
      //====process old versions before automatic schema evolution
      TNamed::Streamer(R__b);
      TAttLine::Streamer(R__b);
      TAttFill::Streamer(R__b);
      TAttMarker::Streamer(R__b);

      fNp = 0;
      /*
      R__b >> fDelta;
      R__b >> fXmin;
      R__b >> fXmax;
      R__b >> fNp;
      R__b >> fKstep;
      R__b >> fHistogram;
      R__b >> fGraph;
      R__b >> fNpx;
      */
      R__b.CheckByteCount(R__s, R__c, TSpline::IsA());
      //====end of old versions

   } else {
      R__b.WriteClassBuffer(TSpline::Class(),this);
   }
}


//////////////////////////////////////////////////////////////////////////
//                                                                      //
// TSplinePoly                                                          //
//                                                                      //
// Base class for TSpline knot                                          //
//                                                                      //
//////////////////////////////////////////////////////////////////////////


//______________________________________________________________________________
TSplinePoly &TSplinePoly::operator=(TSplinePoly const &other)
{
   //assignment operator
   if(this != &other) {
      TObject::operator=(other);
      CopyPoly(other);
   }
   return *this;
}

//______________________________________________________________________________
void TSplinePoly::CopyPoly(TSplinePoly const &other)
{
   //utility called by the copy constructors and = operator
   fX = other.fX;
   fY = other.fY;
}

//////////////////////////////////////////////////////////////////////////
//                                                                      //
// TSplinePoly3                                                         //
//                                                                      //
// Class for TSpline3 knot                                              //
//                                                                      //
//////////////////////////////////////////////////////////////////////////


//______________________________________________________________________________
TSplinePoly3 &TSplinePoly3::operator=(TSplinePoly3 const &other)
{
   //assignment operator
   if(this != &other) {
      TSplinePoly::operator=(other);
      CopyPoly(other);
   }
   return *this;
}

//______________________________________________________________________________
void TSplinePoly3::CopyPoly(TSplinePoly3 const &other)
{
   //utility called by the copy constructors and = operator
   fB = other.fB;
   fC = other.fC;
   fD = other.fD;
}

//////////////////////////////////////////////////////////////////////////
//                                                                      //
// TSplinePoly5                                                         //
//                                                                      //
// Class for TSpline5 knot                                              //
//                                                                      //
//////////////////////////////////////////////////////////////////////////


//______________________________________________________________________________
TSplinePoly5 &TSplinePoly5::operator=(TSplinePoly5 const &other)
{
   //assignment operator
   if(this != &other) {
      TSplinePoly::operator=(other);
      CopyPoly(other);
   }
   return *this;
}

//______________________________________________________________________________
void TSplinePoly5::CopyPoly(TSplinePoly5 const &other)
{
   //utility called by the copy constructors and = operator
   fB = other.fB;
   fC = other.fC;
   fD = other.fD;
   fE = other.fE;
   fF = other.fF;
}



//////////////////////////////////////////////////////////////////////////
//                                                                      //
// TSpline3                                                             //
//                                                                      //
// Class to create third splines to interpolate knots                   //
// Arbitrary conditions can be introduced for first and second          //
// derivatives at beginning and ending points                           //
//                                                                      //
//////////////////////////////////////////////////////////////////////////


//______________________________________________________________________________
TSpline3::TSpline3(const char *title,
                   Double_t x[], Double_t y[], Int_t n, const char *opt,
                   Double_t valbeg, Double_t valend) :
  TSpline(title,-1,x[0],x[n-1],n,kFALSE),
  fValBeg(valbeg), fValEnd(valend), fBegCond(0), fEndCond(0)
{
   // Third spline creator given an array of
   // arbitrary knots in increasing abscissa order and
   // possibly end point conditions

   fName="Spline3";

   // Set endpoint conditions
   if(opt) SetCond(opt);

   // Create the plynomial terms and fill
   // them with node information
   fPoly = new TSplinePoly3[n];
   for (Int_t i=0; i<n; ++i) {
      fPoly[i].X() = x[i];
      fPoly[i].Y() = y[i];
   }

   // Build the spline coefficients
   BuildCoeff();
}

//______________________________________________________________________________
TSpline3::TSpline3(const char *title,
                   Double_t xmin, Double_t xmax,
                   Double_t y[], Int_t n, const char *opt,
                   Double_t valbeg, Double_t valend) :
  TSpline(title,(xmax-xmin)/(n-1), xmin, xmax, n, kTRUE),
  fValBeg(valbeg), fValEnd(valend),
  fBegCond(0), fEndCond(0)
{
   // Third spline creator given an array of
   // arbitrary function values on equidistant n abscissa
   // values from xmin to xmax and possibly end point conditions

   fName="Spline3";

   // Set endpoint conditions
   if(opt) SetCond(opt);

   // Create the plynomial terms and fill
   // them with node information
   fPoly = new TSplinePoly3[n];
   for (Int_t i=0; i<n; ++i) {
      fPoly[i].X() = fXmin+i*fDelta;
      fPoly[i].Y() = y[i];
   }

   // Build the spline coefficients
   BuildCoeff();
}


//______________________________________________________________________________
TSpline3::TSpline3(const char *title,
                   Double_t x[], const TF1 *func, Int_t n, const char *opt,
                   Double_t valbeg, Double_t valend) :
  TSpline(title,-1, x[0], x[n-1], n, kFALSE),
  fValBeg(valbeg), fValEnd(valend),
  fBegCond(0), fEndCond(0)
{
   // Third spline creator given an array of
   // arbitrary abscissas in increasing order and a function
   // to interpolate and possibly end point conditions

   fName="Spline3";

   // Set endpoint conditions
   if(opt) SetCond(opt);

   // Create the plynomial terms and fill
   // them with node information
   fPoly = new TSplinePoly3[n];
   for (Int_t i=0; i<n; ++i) {
      fPoly[i].X() = x[i];
      fPoly[i].Y() = ((TF1*)func)->Eval(x[i]);
   }

   // Build the spline coefficients
   BuildCoeff();
}


//______________________________________________________________________________
TSpline3::TSpline3(const char *title,
                   Double_t xmin, Double_t xmax,
                   const TF1 *func, Int_t n, const char *opt,
                   Double_t valbeg, Double_t valend) :
  TSpline(title,(xmax-xmin)/(n-1), xmin, xmax, n, kTRUE),
  fValBeg(valbeg), fValEnd(valend),
  fBegCond(0), fEndCond(0)
{
   // Third spline creator given a function to be
   // evaluated on n equidistand abscissa points between xmin
   // and xmax and possibly end point conditions

   fName="Spline3";

   // Set endpoint conditions
   if(opt) SetCond(opt);

   // Create the plynomial terms and fill
   // them with node information
   fPoly = new TSplinePoly3[n];
   //when func is null we return. In this case it is assumed that the spline
   //points will be given later via SetPoint and SetPointCoeff
   if (!func) {fKstep = kFALSE; fDelta = -1; return;}
   for (Int_t i=0; i<n; ++i) {
      Double_t x=fXmin+i*fDelta;
      fPoly[i].X() = x;
      fPoly[i].Y() = ((TF1*)func)->Eval(x);
   }

   // Build the spline coefficients
   BuildCoeff();
}


//______________________________________________________________________________
TSpline3::TSpline3(const char *title,
                   const TGraph *g, const char *opt,
                   Double_t valbeg, Double_t valend) :
  TSpline(title,-1,0,0,g->GetN(),kFALSE),
  fValBeg(valbeg), fValEnd(valend),
  fBegCond(0), fEndCond(0)
{
   // Third spline creator given a TGraph with
   // abscissa in increasing order and possibly end
   // point conditions

   fName="Spline3";

   // Set endpoint conditions
   if(opt) SetCond(opt);

   // Create the plynomial terms and fill
   // them with node information
   fPoly = new TSplinePoly3[fNp];
   for (Int_t i=0; i<fNp; ++i) {
      Double_t xx, yy;
      g->GetPoint(i,xx,yy);
      fPoly[i].X()=xx;
      fPoly[i].Y()=yy;
   }
   fXmin = fPoly[0].X();
   fXmax = fPoly[fNp-1].X();

   // Build the spline coefficients
   BuildCoeff();
}


//______________________________________________________________________________
TSpline3::TSpline3(const TH1 *h, const char *opt,
                   Double_t valbeg, Double_t valend) :
  TSpline(h->GetTitle(),-1,0,0,h->GetNbinsX(),kFALSE),
  fValBeg(valbeg), fValEnd(valend),
  fBegCond(0), fEndCond(0)
{
   // Third spline creator given a TH1 

   fName=h->GetName();

   // Set endpoint conditions
   if(opt) SetCond(opt);

   // Create the plynomial terms and fill
   // them with node information
   fPoly = new TSplinePoly3[fNp];
   for (Int_t i=0; i<fNp; ++i) {
      fPoly[i].X()=h->GetXaxis()->GetBinCenter(i+1);
      fPoly[i].Y()=h->GetBinContent(i+1);
   }
   fXmin = fPoly[0].X();
   fXmax = fPoly[fNp-1].X();

   // Build the spline coefficients
   BuildCoeff();
}


//______________________________________________________________________________
TSpline3::TSpline3(const TSpline3& sp3) :
  TSpline(sp3),
  fPoly(0),
  fValBeg(sp3.fValBeg),
  fValEnd(sp3.fValEnd),
  fBegCond(sp3.fBegCond),
  fEndCond(sp3.fEndCond)
{
   //copy constructor
   if (fNp > 0) fPoly = new TSplinePoly3[fNp];
   for (Int_t i=0; i<fNp; ++i) 
      fPoly[i] = sp3.fPoly[i];
}


//______________________________________________________________________________
TSpline3& TSpline3::operator=(const TSpline3& sp3)
{
   //assignment operator
   if(this!=&sp3) {
      TSpline::operator=(sp3);
      fPoly= 0;
      if (fNp > 0) fPoly = new TSplinePoly3[fNp];
      for (Int_t i=0; i<fNp; ++i) 
         fPoly[i] = sp3.fPoly[i];
      
      fValBeg=sp3.fValBeg;
      fValEnd=sp3.fValEnd;
      fBegCond=sp3.fBegCond;
      fEndCond=sp3.fEndCond;
   }
   return *this;
}


//______________________________________________________________________________
void TSpline3::SetCond(const char *opt)
{
   // Check the boundary conditions

   const char *b1 = strstr(opt,"b1");
   const char *e1 = strstr(opt,"e1");
   const char *b2 = strstr(opt,"b2");
   const char *e2 = strstr(opt,"e2");
   if (b1 && b2)
      Error("SetCond","Cannot specify first and second derivative at first point");
   if (e1 && e2)
      Error("SetCond","Cannot specify first and second derivative at last point");
   if (b1) fBegCond=1;
   else if (b2) fBegCond=2;
   if (e1) fEndCond=1;
   else if (e2) fEndCond=2;
}


//______________________________________________________________________________
void TSpline3::Test()
{
   // Test method for TSpline5
   //
   //   n          number of data points.
   //   m          2*m-1 is order of spline.
   //                 m = 2 always for third spline.
   //   nn,nm1,mm,
   //   mm1,i,k,
   //   j,jj       temporary integer variables.
   //   z,p        temporary double precision variables.
   //   x[n]       the sequence of knots.
   //   y[n]       the prescribed function values at the knots.
   //   a[200][4]  two dimensional array whose columns are
   //                 the computed spline coefficients
   //   diff[3]    maximum values of differences of values and
   //                 derivatives to right and left of knots.
   //   com[3]     maximum values of coefficients.
   //
   //
   //   test of TSpline3 with nonequidistant knots and
   //      equidistant knots follows.

   Double_t hx;
   Double_t diff[3];
   Double_t a[800], c[4];
   Int_t i, j, k, m, n;
   Double_t x[200], y[200], z;
   Int_t jj, mm;
   Int_t mm1, nm1;
   Double_t com[3];
   printf("1         TEST OF TSpline3 WITH NONEQUIDISTANT KNOTS\n");
   n = 5;
   x[0] = -3;
   x[1] = -1;
   x[2] = 0;
   x[3] = 3;
   x[4] = 4;
   y[0] = 7;
   y[1] = 11;
   y[2] = 26;
   y[3] = 56;
   y[4] = 29;
   m = 2;
   mm = m << 1;
   mm1 = mm-1;
   printf("\n-N = %3d    M =%2d\n",n,m);
   TSpline3 *spline = new TSpline3("Test",x,y,n);
   for (i = 0; i < n; ++i)
      spline->GetCoeff(i,hx, a[i],a[i+200],a[i+400],a[i+600]);
   delete spline;
   for (i = 0; i < mm1; ++i) diff[i] = com[i] = 0;
   for (k = 0; k < n; ++k) {
      for (i = 0; i < mm; ++i) c[i] = a[k+i*200];
      printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
      printf("%12.8f\n",x[k]);
      if (k == n-1) {
         printf("%16.8f\n",c[0]);
      } else {
         for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
         printf("\n");
         for (i = 0; i < mm1; ++i)
            if ((z=TMath::Abs(a[k+i*200])) > com[i]) com[i] = z;
         z = x[k+1]-x[k];
         for (i = 1; i < mm; ++i)
            for (jj = i; jj < mm; ++jj) {
               j = mm+i-jj;
               c[j-2] = c[j-1]*z+c[j-2];
            }
         for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
         printf("\n");
         for (i = 0; i < mm1; ++i)
            if (!(k >= n-2 && i != 0))
               if((z = TMath::Abs(c[i]-a[k+1+i*200]))
                  > diff[i]) diff[i] = z;
      }
   }
   printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
   for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
   printf("\n");
   printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
   if (TMath::Abs(c[0]) > com[0])
      com[0] = TMath::Abs(c[0]);
   for (i = 0; i < mm1; ++i) printf("%16.8f",com[i]);
   printf("\n");
   m = 2;
   for (n = 10; n <= 100; n += 10) {
      mm = m << 1;
      mm1 = mm-1;
      nm1 = n-1;
      for (i = 0; i < nm1; i += 2) {
         x[i] = i+1;
         x[i+1] = i+2;
         y[i] = 1;
         y[i+1] = 0;
      }
      if (n % 2 != 0) {
         x[n-1] = n;
         y[n-1] = 1;
      }
      printf("\n-N = %3d    M =%2d\n",n,m);
      spline = new TSpline3("Test",x,y,n);
      for (i = 0; i < n; ++i)
         spline->GetCoeff(i,hx,a[i],a[i+200],a[i+400],a[i+600]);
      delete spline;
      for (i = 0; i < mm1; ++i)
         diff[i] = com[i] = 0;
      for (k = 0; k < n; ++k) {
         for (i = 0; i < mm; ++i)
            c[i] = a[k+i*200];
         if (n < 11) {
            printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
            printf("%12.8f\n",x[k]);
            if (k == n-1) printf("%16.8f\n",c[0]);
         }
         if (k == n-1) break;
         if (n <= 10) {
            for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
            printf("\n");
         }
         for (i = 0; i < mm1; ++i)
         if ((z=TMath::Abs(a[k+i*200])) > com[i])
            com[i] = z;
         z = x[k+1]-x[k];
         for (i = 1; i < mm; ++i)
            for (jj = i; jj < mm; ++jj) {
               j = mm+i-jj;
               c[j-2] = c[j-1]*z+c[j-2];
            }
         if (n <= 10) {
            for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
            printf("\n");
         }
         for (i = 0; i < mm1; ++i)
         if (!(k >= n-2 && i != 0))
            if ((z = TMath::Abs(c[i]-a[k+1+i*200]))
               > diff[i]) diff[i] = z;
      }
      printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
      for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
      printf("\n");
      printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
      if (TMath::Abs(c[0]) > com[0])
         com[0] = TMath::Abs(c[0]);
      for (i = 0; i < mm1; ++i) printf("%16.8E",com[i]);
         printf("\n");
   }
}


//______________________________________________________________________________
Int_t TSpline3::FindX(Double_t x) const
{
   // Find X

   Int_t klow=0, khig=fNp-1;
   //
   // If out of boundaries, extrapolate
   // It may be badly wrong
   if(x<=fXmin) klow=0;
   else if(x>=fXmax) klow=khig;
   else {
      if(fKstep) {
         //
         // Equidistant knots, use histogramming
         klow = TMath::FloorNint((x-fXmin)/fDelta);
         // Correction for rounding errors
         if (x < fPoly[klow].X())
            klow = TMath::Max(klow-1,0);
         else if (klow < khig) {
            if (x > fPoly[klow+1].X()) ++klow;
         }
      } else {
         Int_t khalf;
         //
         // Non equidistant knots, binary search
         while(khig-klow>1)
            if(x>fPoly[khalf=(klow+khig)/2].X())
               klow=khalf;
            else
               khig=khalf;
         //
         // This could be removed, sanity check
         if(!(fPoly[klow].X()<=x && x<=fPoly[klow+1].X()))
            Error("Eval",
                  "Binary search failed x(%d) = %f < x= %f < x(%d) = %f\n",
                  klow,fPoly[klow].X(),x,klow+1,fPoly[klow+1].X());
      }
   }
   return klow;
}


//______________________________________________________________________________
Double_t TSpline3::Eval(Double_t x) const
{
   // Eval this spline at x

   Int_t klow=FindX(x);
   if (klow >= fNp-1 && fNp > 1) klow = fNp-2; //see: https://savannah.cern.ch/bugs/?71651
   return fPoly[klow].Eval(x);
}


//______________________________________________________________________________
Double_t TSpline3::Derivative(Double_t x) const
{
   // Derivative

   Int_t klow=FindX(x);
   if (klow >= fNp-1) klow = fNp-2; //see: https://savannah.cern.ch/bugs/?71651
   return fPoly[klow].Derivative(x);
}


//______________________________________________________________________________
void TSpline3::SaveAs(const char *filename, Option_t * /*option*/) const
{
   // write this spline as a C++ function that can be executed without ROOT
   // the name of the function is the name of the file up to the "." if any

   //open the file
   ofstream *f = new ofstream(filename,ios::out);
   if (f == 0 || gSystem->AccessPathName(filename,kWritePermission)) {
      Error("SaveAs","Cannot open file:%s\n",filename);
      return;
   }

   //write the function name and the spline constants
   char buffer[512];
   Int_t nch = strlen(filename);
   snprintf(buffer,512,"double %s",filename);
   char *dot = strstr(buffer,".");
   if (dot) *dot = 0;
   strlcat(buffer,"(double x) {\n",512);
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   const int fNp = %d, fKstep = %d;\n",fNp,fKstep);
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   const double fDelta = %g, fXmin = %g, fXmax = %g;\n",fDelta,fXmin,fXmax);
   nch = strlen(buffer); f->write(buffer,nch);

   //write the spline coefficients
   //array fX
   snprintf(buffer,512,"   const double fX[%d] = {",fNp);
   nch = strlen(buffer); f->write(buffer,nch);
   buffer[0] = 0;
   Int_t i;
   char numb[20];
   for (i=0;i<fNp;i++) {
      snprintf(numb,20," %g,",fPoly[i].X());
      nch = strlen(numb);
      if (i == fNp-1) numb[nch-1]=0;
      strlcat(buffer,numb,512);
      if (i%5 == 4 || i == fNp-1) {
         nch = strlen(buffer); f->write(buffer,nch);
         if (i != fNp-1) snprintf(buffer,512,"\n                       ");
      }
   }
   snprintf(buffer,512," };\n");
   nch = strlen(buffer); f->write(buffer,nch);
   //array fY
   snprintf(buffer,512,"   const double fY[%d] = {",fNp);
   nch = strlen(buffer); f->write(buffer,nch);
   buffer[0] = 0;
   for (i=0;i<fNp;i++) {
      snprintf(numb,20," %g,",fPoly[i].Y());
      nch = strlen(numb);
      if (i == fNp-1) numb[nch-1]=0;
      strlcat(buffer,numb,512);
      if (i%5 == 4 || i == fNp-1) {
         nch = strlen(buffer); f->write(buffer,nch);
         if (i != fNp-1) snprintf(buffer,512,"\n                       ");
      }
   }
   snprintf(buffer,512," };\n");
   nch = strlen(buffer); f->write(buffer,nch);
   //array fB
   snprintf(buffer,512,"   const double fB[%d] = {",fNp);
   nch = strlen(buffer); f->write(buffer,nch);
   buffer[0] = 0;
   for (i=0;i<fNp;i++) {
      snprintf(numb,20," %g,",fPoly[i].B());
      nch = strlen(numb);
      if (i == fNp-1) numb[nch-1]=0;
      strlcat(buffer,numb,512);
      if (i%5 == 4 || i == fNp-1) {
         nch = strlen(buffer); f->write(buffer,nch);
         if (i != fNp-1) snprintf(buffer,512,"\n                       ");
      }
   }
   snprintf(buffer,512," };\n");
   nch = strlen(buffer); f->write(buffer,nch);
   //array fC
   snprintf(buffer,512,"   const double fC[%d] = {",fNp);
   nch = strlen(buffer); f->write(buffer,nch);
   buffer[0] = 0;
   for (i=0;i<fNp;i++) {
      snprintf(numb,20," %g,",fPoly[i].C());
      nch = strlen(numb);
      if (i == fNp-1) numb[nch-1]=0;
      strlcat(buffer,numb,512);
      if (i%5 == 4 || i == fNp-1) {
         nch = strlen(buffer); f->write(buffer,nch);
         if (i != fNp-1) snprintf(buffer,512,"\n                       ");
      }
   }
   snprintf(buffer,512," };\n");
   nch = strlen(buffer); f->write(buffer,nch);
    //array fD
   snprintf(buffer,512,"   const double fD[%d] = {",fNp);
   nch = strlen(buffer); f->write(buffer,nch);
   buffer[0] = 0;
   for (i=0;i<fNp;i++) {
      snprintf(numb,20," %g,",fPoly[i].D());
      nch = strlen(numb);
      if (i == fNp-1) numb[nch-1]=0;
      strlcat(buffer,numb,512);
      if (i%5 == 4 || i == fNp-1) {
         nch = strlen(buffer); f->write(buffer,nch);
         if (i != fNp-1) snprintf(buffer,512,"\n                       ");
      }
   }
   snprintf(buffer,512," };\n");
   nch = strlen(buffer); f->write(buffer,nch);

   //generate code for the spline evaluation
   snprintf(buffer,512,"   int klow=0;\n");
   nch = strlen(buffer); f->write(buffer,nch);

   snprintf(buffer,512,"   // If out of boundaries, extrapolate. It may be badly wrong\n");
   snprintf(buffer,512,"   if(x<=fXmin) klow=0;\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   else if(x>=fXmax) klow=fNp-1;\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   else {\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"     if(fKstep) {\n");
   nch = strlen(buffer); f->write(buffer,nch);

   snprintf(buffer,512,"       // Equidistant knots, use histogramming\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"       klow = int((x-fXmin)/fDelta);\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"       if (klow < fNp-1) klow = fNp-1;\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"     } else {\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"       int khig=fNp-1, khalf;\n");
   nch = strlen(buffer); f->write(buffer,nch);

   snprintf(buffer,512,"       // Non equidistant knots, binary search\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"       while(khig-klow>1)\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"         if(x>fX[khalf=(klow+khig)/2]) klow=khalf;\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"         else khig=khalf;\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"     }\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   }\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   // Evaluate now\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   double dx=x-fX[klow];\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   return (fY[klow]+dx*(fB[klow]+dx*(fC[klow]+dx*fD[klow])));\n");
   nch = strlen(buffer); f->write(buffer,nch);

   //close file
   f->write("}\n",2);

   if (f) { f->close(); delete f;}
}


//______________________________________________________________________________
void TSpline3::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
{
    // Save primitive as a C++ statement(s) on output stream out

   char quote = '"';
   out<<"   "<<endl;
   if (gROOT->ClassSaved(TSpline3::Class())) {
      out<<"   ";
   } else {
      out<<"   TSpline3 *";
   }
   out<<"spline3 = new TSpline3("<<quote<<GetTitle()<<quote<<","
      <<fXmin<<","<<fXmax<<",(TF1*)0,"<<fNp<<","<<quote<<quote<<","
      <<fValBeg<<","<<fValEnd<<");"<<endl;
   out<<"   spline3->SetName("<<quote<<GetName()<<quote<<");"<<endl;

   SaveFillAttributes(out,"spline3",0,1001);
   SaveLineAttributes(out,"spline3",1,1,1);
   SaveMarkerAttributes(out,"spline3",1,1,1);
   if (fNpx != 100) out<<"   spline3->SetNpx("<<fNpx<<");"<<endl;

   for (Int_t i=0;i<fNp;i++) {
      out<<"   spline3->SetPoint("<<i<<","<<fPoly[i].X()<<","<<fPoly[i].Y()<<");"<<endl;
      out<<"   spline3->SetPointCoeff("<<i<<","<<fPoly[i].B()<<","<<fPoly[i].C()<<","<<fPoly[i].D()<<");"<<endl;
   }
   out<<"   spline3->Draw("<<quote<<option<<quote<<");"<<endl;
}


//______________________________________________________________________________
void TSpline3::SetPoint(Int_t i, Double_t x, Double_t y)
{
   //set point number i.
   
   if (i < 0 || i >= fNp) return;
   fPoly[i].X()= x;
   fPoly[i].Y()= y;
}

//______________________________________________________________________________
void TSpline3::SetPointCoeff(Int_t i, Double_t b, Double_t c, Double_t d)
{
   // set point coefficient number i
 
   if (i < 0 || i >= fNp) return;
   fPoly[i].B()= b;
   fPoly[i].C()= c;
   fPoly[i].D()= d;
}

//______________________________________________________________________________
void TSpline3::BuildCoeff()
{
   //      subroutine cubspl ( tau, c, n, ibcbeg, ibcend )
   //  from  * a practical guide to splines *  by c. de boor
   //     ************************  input  ***************************
   //     n = number of data points. assumed to be .ge. 2.
   //     (tau(i), c(1,i), i=1,...,n) = abscissae and ordinates of the
   //        data points. tau is assumed to be strictly increasing.
   //     ibcbeg, ibcend = boundary condition indicators, and
   //     c(2,1), c(2,n) = boundary condition information. specifically,
   //        ibcbeg = 0  means no boundary condition at tau(1) is given.
   //           in this case, the not-a-knot condition is used, i.e. the
   //           jump in the third derivative across tau(2) is forced to
   //           zero, thus the first and the second cubic polynomial pieces
   //           are made to coincide.)
   //        ibcbeg = 1  means that the slope at tau(1) is made to equal
   //           c(2,1), supplied by input.
   //        ibcbeg = 2  means that the second derivative at tau(1) is
   //           made to equal c(2,1), supplied by input.
   //        ibcend = 0, 1, or 2 has analogous meaning concerning the
   //           boundary condition at tau(n), with the additional infor-
   //           mation taken from c(2,n).
   //     ***********************  output  **************************
   //     c(j,i), j=1,...,4; i=1,...,l (= n-1) = the polynomial coefficients
   //        of the cubic interpolating spline with interior knots (or
   //        joints) tau(2), ..., tau(n-1). precisely, in the interval
   //        (tau(i), tau(i+1)), the spline f is given by
   //           f(x) = c(1,i)+h*(c(2,i)+h*(c(3,i)+h*c(4,i)/3.)/2.)
   //        where h = x - tau(i). the function program *ppvalu* may be
   //        used to evaluate f or its derivatives from tau,c, l = n-1,
   //        and k=4.

   Int_t i, j, l, m;
   Double_t   divdf1,divdf3,dtau,g=0;
   //***** a tridiagonal linear system for the unknown slopes s(i) of
   //  f  at tau(i), i=1,...,n, is generated and then solved by gauss elim-
   //  ination, with s(i) ending up in c(2,i), all i.
   //     c(3,.) and c(4,.) are used initially for temporary storage.
   l = fNp-1;
   // compute first differences of x sequence and store in C also,
   // compute first divided difference of data and store in D.
   for (m=1; m<fNp ; ++m) {
      fPoly[m].C() = fPoly[m].X() - fPoly[m-1].X();
      fPoly[m].D() = (fPoly[m].Y() - fPoly[m-1].Y())/fPoly[m].C();
   }
   // construct first equation from the boundary condition, of the form
   //             D[0]*s[0] + C[0]*s[1] = B[0]
   if(fBegCond==0) {
      if(fNp == 2) {
         //     no condition at left end and n = 2.
         fPoly[0].D() = 1.;
         fPoly[0].C() = 1.;
         fPoly[0].B() = 2.*fPoly[1].D();
      } else {
         //     not-a-knot condition at left end and n .gt. 2.
         fPoly[0].D() = fPoly[2].C();
         fPoly[0].C() = fPoly[1].C() + fPoly[2].C();
         fPoly[0].B() =((fPoly[1].C()+2.*fPoly[0].C())*fPoly[1].D()*fPoly[2].C()+fPoly[1].C()*fPoly[1].C()*fPoly[2].D())/fPoly[0].C();
      }
   } else if (fBegCond==1) {
      //     slope prescribed at left end.
      fPoly[0].B() = fValBeg;
      fPoly[0].D() = 1.;
      fPoly[0].C() = 0.;
   } else if (fBegCond==2) {
      //     second derivative prescribed at left end.
      fPoly[0].D() = 2.;
      fPoly[0].C() = 1.;
      fPoly[0].B() = 3.*fPoly[1].D() - fPoly[1].C()/2.*fValBeg;
   }
   if(fNp > 2) {
      //  if there are interior knots, generate the corresp. equations and car-
      //  ry out the forward pass of gauss elimination, after which the m-th
      //  equation reads    D[m]*s[m] + C[m]*s[m+1] = B[m].
      for (m=1; m<l; ++m) {
         g = -fPoly[m+1].C()/fPoly[m-1].D();
         fPoly[m].B() = g*fPoly[m-1].B() + 3.*(fPoly[m].C()*fPoly[m+1].D()+fPoly[m+1].C()*fPoly[m].D());
         fPoly[m].D() = g*fPoly[m-1].C() + 2.*(fPoly[m].C() + fPoly[m+1].C());
      }
      // construct last equation from the second boundary condition, of the form
      //           (-g*D[n-2])*s[n-2] + D[n-1]*s[n-1] = B[n-1]
      //     if slope is prescribed at right end, one can go directly to back-
      //     substitution, since c array happens to be set up just right for it
      //     at this point.
      if(fEndCond == 0) {
         if (fNp > 3 || fBegCond != 0) {
            //     not-a-knot and n .ge. 3, and either n.gt.3 or  also not-a-knot at
            //     left end point.
            g = fPoly[fNp-2].C() + fPoly[fNp-1].C();
            fPoly[fNp-1].B() = ((fPoly[fNp-1].C()+2.*g)*fPoly[fNp-1].D()*fPoly[fNp-2].C()
                         + fPoly[fNp-1].C()*fPoly[fNp-1].C()*(fPoly[fNp-2].Y()-fPoly[fNp-3].Y())/fPoly[fNp-2].C())/g;
            g = -g/fPoly[fNp-2].D();
            fPoly[fNp-1].D() = fPoly[fNp-2].C();
         } else {
            //     either (n=3 and not-a-knot also at left) or (n=2 and not not-a-
            //     knot at left end point).
            fPoly[fNp-1].B() = 2.*fPoly[fNp-1].D();
            fPoly[fNp-1].D() = 1.;
            g = -1./fPoly[fNp-2].D();
         }
      } else if (fEndCond == 1) {
         fPoly[fNp-1].B() = fValEnd;
         goto L30;
      } else if (fEndCond == 2) {
         //     second derivative prescribed at right endpoint.
         fPoly[fNp-1].B() = 3.*fPoly[fNp-1].D() + fPoly[fNp-1].C()/2.*fValEnd;
         fPoly[fNp-1].D() = 2.;
         g = -1./fPoly[fNp-2].D();
      }
   } else {
      if(fEndCond == 0) {
         if (fBegCond > 0) {
            //     either (n=3 and not-a-knot also at left) or (n=2 and not not-a-
            //     knot at left end point).
            fPoly[fNp-1].B() = 2.*fPoly[fNp-1].D();
            fPoly[fNp-1].D() = 1.;
            g = -1./fPoly[fNp-2].D();
         } else {
            //     not-a-knot at right endpoint and at left endpoint and n = 2.
            fPoly[fNp-1].B() = fPoly[fNp-1].D();
            goto L30;
         }
      } else if(fEndCond == 1) {
         fPoly[fNp-1].B() = fValEnd;
         goto L30;
      } else if(fEndCond == 2) {
         //     second derivative prescribed at right endpoint.
         fPoly[fNp-1].B() = 3.*fPoly[fNp-1].D() + fPoly[fNp-1].C()/2.*fValEnd;
         fPoly[fNp-1].D() = 2.;
         g = -1./fPoly[fNp-2].D();
      }
   }
   // complete forward pass of gauss elimination.
   fPoly[fNp-1].D() = g*fPoly[fNp-2].C() + fPoly[fNp-1].D();
   fPoly[fNp-1].B() = (g*fPoly[fNp-2].B() + fPoly[fNp-1].B())/fPoly[fNp-1].D();
   // carry out back substitution
L30: j = l-1;
   do {
      fPoly[j].B() = (fPoly[j].B() - fPoly[j].C()*fPoly[j+1].B())/fPoly[j].D();
      --j;
   }  while (j>=0);
   //****** generate cubic coefficients in each interval, i.e., the deriv.s
   //  at its left endpoint, from value and slope at its endpoints.
   for (i=1; i<fNp; ++i) {
      dtau = fPoly[i].C();
      divdf1 = (fPoly[i].Y() - fPoly[i-1].Y())/dtau;
      divdf3 = fPoly[i-1].B() + fPoly[i].B() - 2.*divdf1;
      fPoly[i-1].C() = (divdf1 - fPoly[i-1].B() - divdf3)/dtau;
      fPoly[i-1].D() = (divdf3/dtau)/dtau;
   }
}


//______________________________________________________________________________
void TSpline3::Streamer(TBuffer &R__b)
{
   // Stream an object of class TSpline3.

   if (R__b.IsReading()) {
      UInt_t R__s, R__c;
      Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
      if (R__v > 1) {
         R__b.ReadClassBuffer(TSpline3::Class(), this, R__v, R__s, R__c);
         return;
      }
      //====process old versions before automatic schema evolution
      TSpline::Streamer(R__b);
      if (fNp > 0) {
         fPoly = new TSplinePoly3[fNp];
         for(Int_t i=0; i<fNp; ++i) {
            fPoly[i].Streamer(R__b);
         }
      }
      //      R__b >> fPoly;
      R__b >> fValBeg;
      R__b >> fValEnd;
      R__b >> fBegCond;
      R__b >> fEndCond;
   } else {
      R__b.WriteClassBuffer(TSpline3::Class(),this);
   }
}


//////////////////////////////////////////////////////////////////////////
//                                                                      //
// TSpline5                                                             //
//                                                                      //
// Class to create quintic natural splines to interpolate knots         //
// Arbitrary conditions can be introduced for first and second          //
// derivatives using double knots (see BuildCoeff) for more on this.    //
// Double knots are automatically introduced at ending points           //
//                                                                      //
//////////////////////////////////////////////////////////////////////////


//______________________________________________________________________________
TSpline5::TSpline5(const char *title,
                   Double_t x[], Double_t y[], Int_t n,
                   const char *opt, Double_t b1, Double_t e1,
                   Double_t b2, Double_t e2) :
  TSpline(title,-1, x[0], x[n-1], n, kFALSE)
{
   // Quintic natural spline creator given an array of
   // arbitrary knots in increasing abscissa order and
   // possibly end point conditions

   Int_t beg, end;
   const char *cb1, *ce1, *cb2, *ce2;
   fName="Spline5";

   // Check endpoint conditions
   BoundaryConditions(opt,beg,end,cb1,ce1,cb2,ce2);

   // Create the plynomial terms and fill
   // them with node information
   fPoly = new TSplinePoly5[fNp];
   for (Int_t i=0; i<n; ++i) {
      fPoly[i+beg].X() = x[i];
      fPoly[i+beg].Y() = y[i];
   }

   // Set the double knots at boundaries
   SetBoundaries(b1,e1,b2,e2,cb1,ce1,cb2,ce2);

   // Build the spline coefficients
   BuildCoeff();
}


//______________________________________________________________________________
TSpline5::TSpline5(const char *title,
                   Double_t xmin, Double_t xmax,
                   Double_t y[], Int_t n,
                   const char *opt, Double_t b1, Double_t e1,
                   Double_t b2, Double_t e2) :
  TSpline(title,(xmax-xmin)/(n-1), xmin, xmax, n, kTRUE)
{
   // Quintic natural spline creator given an array of
   // arbitrary function values on equidistant n abscissa
   // values from xmin to xmax and possibly end point conditions

   Int_t beg, end;
   const char *cb1, *ce1, *cb2, *ce2;
   fName="Spline5";

   // Check endpoint conditions
   BoundaryConditions(opt,beg,end,cb1,ce1,cb2,ce2);

   // Create the plynomial terms and fill
   // them with node information
   fPoly = new TSplinePoly5[fNp];
   for (Int_t i=0; i<n; ++i) {
      fPoly[i+beg].X() = fXmin+i*fDelta;
      fPoly[i+beg].Y() = y[i];
   }

   // Set the double knots at boundaries
   SetBoundaries(b1,e1,b2,e2,cb1,ce1,cb2,ce2);

   // Build the spline coefficients
   BuildCoeff();
}


//______________________________________________________________________________
TSpline5::TSpline5(const char *title,
                   Double_t x[], const TF1 *func, Int_t n,
                   const char *opt, Double_t b1, Double_t e1,
                   Double_t b2, Double_t e2) :
  TSpline(title,-1, x[0], x[n-1], n, kFALSE)
{
   // Quintic natural spline creator given an array of
   // arbitrary abscissas in increasing order and a function
   // to interpolate and possibly end point conditions

   Int_t beg, end;
   const char *cb1, *ce1, *cb2, *ce2;
   fName="Spline5";

   // Check endpoint conditions
   BoundaryConditions(opt,beg,end,cb1,ce1,cb2,ce2);

   // Create the plynomial terms and fill
   // them with node information
   fPoly = new TSplinePoly5[fNp];
   for (Int_t i=0; i<n; i++) {
      fPoly[i+beg].X() = x[i];
      fPoly[i+beg].Y() = ((TF1*)func)->Eval(x[i]);
   }

   // Set the double knots at boundaries
   SetBoundaries(b1,e1,b2,e2,cb1,ce1,cb2,ce2);

   // Build the spline coefficients
   BuildCoeff();
}


//______________________________________________________________________________
TSpline5::TSpline5(const char *title,
                   Double_t xmin, Double_t xmax,
                   const TF1 *func, Int_t n,
                   const char *opt, Double_t b1, Double_t e1,
                   Double_t b2, Double_t e2) :
  TSpline(title,(xmax-xmin)/(n-1), xmin, xmax, n, kTRUE)
{
   // Quintic natural spline creator given a function to be
   // evaluated on n equidistand abscissa points between xmin
   // and xmax and possibly end point conditions

   Int_t beg, end;
   const char *cb1, *ce1, *cb2, *ce2;
   fName="Spline5";

   // Check endpoint conditions
   BoundaryConditions(opt,beg,end,cb1,ce1,cb2,ce2);

   // Create the plynomial terms and fill
   // them with node information
   fPoly = new TSplinePoly5[fNp];
   for (Int_t i=0; i<n; ++i) {
      Double_t x=fXmin+i*fDelta;
      fPoly[i+beg].X() = x;
      if (func) fPoly[i+beg].Y() = ((TF1*)func)->Eval(x);
   }
   if (!func) {fDelta = -1; fKstep = kFALSE;}

   // Set the double knots at boundaries
   SetBoundaries(b1,e1,b2,e2,cb1,ce1,cb2,ce2);

   // Build the spline coefficients
   if (func) BuildCoeff();
}


//______________________________________________________________________________
TSpline5::TSpline5(const char *title,
                   const TGraph *g,
                   const char *opt, Double_t b1, Double_t e1,
                   Double_t b2, Double_t e2) :
  TSpline(title,-1,0,0,g->GetN(),kFALSE)
{
   // Quintic natural spline creator given a TGraph with
   // abscissa in increasing order and possibly end
   // point conditions

   Int_t beg, end;
   const char *cb1, *ce1, *cb2, *ce2;
   fName="Spline5";

   // Check endpoint conditions
   BoundaryConditions(opt,beg,end,cb1,ce1,cb2,ce2);

   // Create the plynomial terms and fill
   // them with node information
   fPoly = new TSplinePoly5[fNp];
   for (Int_t i=0; i<fNp-beg; ++i) {
      Double_t xx, yy;
      g->GetPoint(i,xx,yy);
      fPoly[i+beg].X()=xx;
      fPoly[i+beg].Y()=yy;
   }

   // Set the double knots at boundaries
   SetBoundaries(b1,e1,b2,e2,cb1,ce1,cb2,ce2);
   fXmin = fPoly[0].X();
   fXmax = fPoly[fNp-1].X();

   // Build the spline coefficients
   BuildCoeff();
}


//______________________________________________________________________________
TSpline5::TSpline5(const TH1 *h,
                   const char *opt, Double_t b1, Double_t e1,
                   Double_t b2, Double_t e2) :
  TSpline(h->GetTitle(),-1,0,0,h->GetNbinsX(),kFALSE)
{
   // Quintic natural spline creator given a TH1

   Int_t beg, end;
   const char *cb1, *ce1, *cb2, *ce2;
   fName=h->GetName();

   // Check endpoint conditions
   BoundaryConditions(opt,beg,end,cb1,ce1,cb2,ce2);

   // Create the plynomial terms and fill
   // them with node information
   fPoly = new TSplinePoly5[fNp];
   for (Int_t i=0; i<fNp-beg; ++i) {
      fPoly[i+beg].X()=h->GetXaxis()->GetBinCenter(i+1);
      fPoly[i+beg].Y()=h->GetBinContent(i+1);
   }

   // Set the double knots at boundaries
   SetBoundaries(b1,e1,b2,e2,cb1,ce1,cb2,ce2);
   fXmin = fPoly[0].X();
   fXmax = fPoly[fNp-1].X();

   // Build the spline coefficients
   BuildCoeff();
}


//______________________________________________________________________________
TSpline5::TSpline5(const TSpline5& sp5) :
  TSpline(sp5),
  fPoly(0)
{
   //copy constructor
   if (fNp > 0) fPoly = new TSplinePoly5[fNp];
   for (Int_t i=0; i<fNp; ++i) {
      fPoly[i] = sp5.fPoly[i];
   }
}


//______________________________________________________________________________
TSpline5& TSpline5::operator=(const TSpline5& sp5)
{
   //assignment operator
   if(this!=&sp5) {
      TSpline::operator=(sp5);
      fPoly=0;
      if (fNp > 0) fPoly = new TSplinePoly5[fNp];
      for (Int_t i=0; i<fNp; ++i) {
         fPoly[i] = sp5.fPoly[i];
      }
   }
   return *this;
}


//______________________________________________________________________________
void TSpline5::BoundaryConditions(const char *opt,Int_t &beg,Int_t &end,
                                  const char *&cb1,const char *&ce1,
                                  const char *&cb2,const char *&ce2)
{
   // Check the boundary conditions and the
   // amount of extra double knots needed

   cb1=ce1=cb2=ce2=0;
   beg=end=0;
   if(opt) {
      cb1 = strstr(opt,"b1");
      ce1 = strstr(opt,"e1");
      cb2 = strstr(opt,"b2");
      ce2 = strstr(opt,"e2");
      if(cb2) {
         fNp=fNp+2;
         beg=2;
      } else if(cb1) {
         fNp=fNp+1;
         beg=1;
      }
      if(ce2) {
         fNp=fNp+2;
         end=2;
      } else if(ce1) {
         fNp=fNp+1;
         end=1;
      }
   }
}


//______________________________________________________________________________
void TSpline5::SetBoundaries(Double_t b1, Double_t e1, Double_t b2, Double_t e2,
                             const char *cb1, const char *ce1, const char *cb2,
                             const char *ce2)
{
   // Set the boundary conditions at double/triple knots

   if(cb2) {

      // Second derivative at the beginning
      fPoly[0].X() = fPoly[1].X() = fPoly[2].X();
      fPoly[0].Y() = fPoly[2].Y();
      fPoly[2].Y()=b2;

      // If first derivative not given, we take the finite
      // difference from first and second point... not recommended
      if(cb1)
         fPoly[1].Y()=b1;
      else
         fPoly[1].Y()=(fPoly[3].Y()-fPoly[0].Y())/(fPoly[3].X()-fPoly[2].X());
   } else if(cb1) {

      // First derivative at the end
      fPoly[0].X() = fPoly[1].X();
      fPoly[0].Y() = fPoly[1].Y();
      fPoly[1].Y()=b1;
   }
   if(ce2) {

      // Second derivative at the end
      fPoly[fNp-1].X() = fPoly[fNp-2].X() = fPoly[fNp-3].X();
      fPoly[fNp-1].Y()=e2;

      // If first derivative not given, we take the finite
      // difference from first and second point... not recommended
      if(ce1)
         fPoly[fNp-2].Y()=e1;
      else
         fPoly[fNp-2].Y()=
         (fPoly[fNp-3].Y()-fPoly[fNp-4].Y())
         /(fPoly[fNp-3].X()-fPoly[fNp-4].X());
   } else if(ce1) {

      // First derivative at the end
      fPoly[fNp-1].X() = fPoly[fNp-2].X();
      fPoly[fNp-1].Y()=e1;
   }
}


//______________________________________________________________________________
Int_t TSpline5::FindX(Double_t x) const
{
   // Find X

   Int_t klow=0;

   // If out of boundaries, extrapolate
   // It may be badly wrong
   if(x<=fXmin) klow=0;
   else if(x>=fXmax) klow=fNp-1;
   else {
      if(fKstep) {

         // Equidistant knots, use histogramming
         klow = TMath::Min(Int_t((x-fXmin)/fDelta),fNp-1);
      } else {
         Int_t khig=fNp-1, khalf;

         // Non equidistant knots, binary search
         while(khig-klow>1)
            if(x>fPoly[khalf=(klow+khig)/2].X())
               klow=khalf;
            else
               khig=khalf;
      }

      // This could be removed, sanity check
      if(!(fPoly[klow].X()<=x && x<=fPoly[klow+1].X()))
         Error("Eval",
               "Binary search failed x(%d) = %f < x(%d) = %f\n",
                klow,fPoly[klow].X(),klow+1,fPoly[klow+1].X());
   }
   return klow;
}


//______________________________________________________________________________
Double_t TSpline5::Eval(Double_t x) const
{
   // Eval this spline at x

   Int_t klow=FindX(x);
   return fPoly[klow].Eval(x);
}


//______________________________________________________________________________
Double_t TSpline5::Derivative(Double_t x) const
{
   // Derivative

   Int_t klow=FindX(x);
   return fPoly[klow].Derivative(x);
}


//______________________________________________________________________________
void TSpline5::SaveAs(const char *filename, Option_t * /*option*/) const
{
   // write this spline as a C++ function that can be executed without ROOT
   // the name of the function is the name of the file up to the "." if any

   //open the file
   ofstream *f = new ofstream(filename,ios::out);
   if (f == 0 || gSystem->AccessPathName(filename,kWritePermission)) {
      Error("SaveAs","Cannot open file:%s\n",filename);
      return;
   }

   //write the function name and the spline constants
   char buffer[512];
   Int_t nch = strlen(filename);
   snprintf(buffer,512,"double %s",filename);
   char *dot = strstr(buffer,".");
   if (dot) *dot = 0;
   strlcat(buffer,"(double x) {\n",512);
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   const int fNp = %d, fKstep = %d;\n",fNp,fKstep);
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   const double fDelta = %g, fXmin = %g, fXmax = %g;\n",fDelta,fXmin,fXmax);
   nch = strlen(buffer); f->write(buffer,nch);

   //write the spline coefficients
   //array fX
   snprintf(buffer,512,"   const double fX[%d] = {",fNp);
   nch = strlen(buffer); f->write(buffer,nch);
   buffer[0] = 0;
   Int_t i;
   char numb[20];
   for (i=0;i<fNp;i++) {
      snprintf(numb,20," %g,",fPoly[i].X());
      nch = strlen(numb);
      if (i == fNp-1) numb[nch-1]=0;
      strlcat(buffer,numb,512);
      if (i%5 == 4 || i == fNp-1) {
         nch = strlen(buffer); f->write(buffer,nch);
         if (i != fNp-1) snprintf(buffer,512,"\n                       ");
      }
   }
   snprintf(buffer,512," };\n");
   nch = strlen(buffer); f->write(buffer,nch);
   //array fY
   snprintf(buffer,512,"   const double fY[%d] = {",fNp);
   nch = strlen(buffer); f->write(buffer,nch);
   buffer[0] = 0;
   for (i=0;i<fNp;i++) {
      snprintf(numb,20," %g,",fPoly[i].Y());
      nch = strlen(numb);
      if (i == fNp-1) numb[nch-1]=0;
      strlcat(buffer,numb,512);
      if (i%5 == 4 || i == fNp-1) {
         nch = strlen(buffer); f->write(buffer,nch);
         if (i != fNp-1) snprintf(buffer,512,"\n                       ");
      }
   }
   snprintf(buffer,512," };\n");
   nch = strlen(buffer); f->write(buffer,nch);
   //array fB
   snprintf(buffer,512,"   const double fB[%d] = {",fNp);
   nch = strlen(buffer); f->write(buffer,nch);
   buffer[0] = 0;
   for (i=0;i<fNp;i++) {
      snprintf(numb,20," %g,",fPoly[i].B());
      nch = strlen(numb);
      if (i == fNp-1) numb[nch-1]=0;
      strlcat(buffer,numb,512);
      if (i%5 == 4 || i == fNp-1) {
         nch = strlen(buffer); f->write(buffer,nch);
         if (i != fNp-1) snprintf(buffer,512,"\n                       ");
      }
   }
   snprintf(buffer,512," };\n");
   nch = strlen(buffer); f->write(buffer,nch);
   //array fC
   snprintf(buffer,512,"   const double fC[%d] = {",fNp);
   nch = strlen(buffer); f->write(buffer,nch);
   buffer[0] = 0;
   for (i=0;i<fNp;i++) {
      snprintf(numb,20," %g,",fPoly[i].C());
      nch = strlen(numb);
      if (i == fNp-1) numb[nch-1]=0;
      strlcat(buffer,numb,512);
      if (i%5 == 4 || i == fNp-1) {
         nch = strlen(buffer); f->write(buffer,nch);
         if (i != fNp-1) snprintf(buffer,512,"\n                       ");
      }
   }
   snprintf(buffer,512," };\n");
   nch = strlen(buffer); f->write(buffer,nch);
    //array fD
   snprintf(buffer,512,"   const double fD[%d] = {",fNp);
   nch = strlen(buffer); f->write(buffer,nch);
   buffer[0] = 0;
   for (i=0;i<fNp;i++) {
      snprintf(numb,20," %g,",fPoly[i].D());
      nch = strlen(numb);
      if (i == fNp-1) numb[nch-1]=0;
      strlcat(buffer,numb,512);
      if (i%5 == 4 || i == fNp-1) {
         nch = strlen(buffer); f->write(buffer,nch);
         if (i != fNp-1) snprintf(buffer,512,"\n                       ");
      }
   }
   snprintf(buffer,512," };\n");
   nch = strlen(buffer); f->write(buffer,nch);
    //array fE
   snprintf(buffer,512,"   const double fE[%d] = {",fNp);
   nch = strlen(buffer); f->write(buffer,nch);
   buffer[0] = 0;
   for (i=0;i<fNp;i++) {
      snprintf(numb,20," %g,",fPoly[i].E());
      nch = strlen(numb);
      if (i == fNp-1) numb[nch-1]=0;
      strlcat(buffer,numb,512);
      if (i%5 == 4 || i == fNp-1) {
         nch = strlen(buffer); f->write(buffer,nch);
         if (i != fNp-1) snprintf(buffer,512,"\n                       ");
      }
   }
   snprintf(buffer,512," };\n");
   nch = strlen(buffer); f->write(buffer,nch);
    //array fF
   snprintf(buffer,512,"   const double fF[%d] = {",fNp);
   nch = strlen(buffer); f->write(buffer,nch);
   buffer[0] = 0;
   for (i=0;i<fNp;i++) {
      snprintf(numb,20," %g,",fPoly[i].F());
      nch = strlen(numb);
      if (i == fNp-1) numb[nch-1]=0;
      strlcat(buffer,numb,512);
      if (i%5 == 4 || i == fNp-1) {
         nch = strlen(buffer); f->write(buffer,nch);
         if (i != fNp-1) snprintf(buffer,512,"\n                       ");
      }
   }
   snprintf(buffer,512," };\n");
   nch = strlen(buffer); f->write(buffer,nch);

   //generate code for the spline evaluation
   snprintf(buffer,512,"   int klow=0;\n");
   nch = strlen(buffer); f->write(buffer,nch);

   snprintf(buffer,512,"   // If out of boundaries, extrapolate. It may be badly wrong\n");
   snprintf(buffer,512,"   if(x<=fXmin) klow=0;\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   else if(x>=fXmax) klow=fNp-1;\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   else {\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"     if(fKstep) {\n");
   nch = strlen(buffer); f->write(buffer,nch);

   snprintf(buffer,512,"       // Equidistant knots, use histogramming\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"       klow = int((x-fXmin)/fDelta);\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"       if (klow < fNp-1) klow = fNp-1;\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"     } else {\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"       int khig=fNp-1, khalf;\n");
   nch = strlen(buffer); f->write(buffer,nch);

   snprintf(buffer,512,"       // Non equidistant knots, binary search\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"       while(khig-klow>1)\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"         if(x>fX[khalf=(klow+khig)/2]) klow=khalf;\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"         else khig=khalf;\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"     }\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   }\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   // Evaluate now\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   double dx=x-fX[klow];\n");
   nch = strlen(buffer); f->write(buffer,nch);
   snprintf(buffer,512,"   return (fY[klow]+dx*(fB[klow]+dx*(fC[klow]+dx*(fD[klow]+dx*(fE[klow]+dx*fF[klow])))));\n");
   nch = strlen(buffer); f->write(buffer,nch);

   //close file
   f->write("}\n",2);

   if (f) { f->close(); delete f;}
}


//______________________________________________________________________________
void TSpline5::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
{
    // Save primitive as a C++ statement(s) on output stream out

   char quote = '"';
   out<<"   "<<endl;
   if (gROOT->ClassSaved(TSpline5::Class())) {
      out<<"   ";
   } else {
      out<<"   TSpline5 *";
   }
   Double_t b1 = fPoly[1].Y();
   Double_t e1 = fPoly[fNp-1].Y();
   Double_t b2 = fPoly[2].Y();
   Double_t e2 = fPoly[fNp-1].Y();
   out<<"spline5 = new TSpline5("<<quote<<GetTitle()<<quote<<","
      <<fXmin<<","<<fXmax<<",(TF1*)0,"<<fNp<<","<<quote<<quote<<","
      <<b1<<","<<e1<<","<<b2<<","<<e2<<");"<<endl;
   out<<"   spline5->SetName("<<quote<<GetName()<<quote<<");"<<endl;

   SaveFillAttributes(out,"spline5",0,1001);
   SaveLineAttributes(out,"spline5",1,1,1);
   SaveMarkerAttributes(out,"spline5",1,1,1);
   if (fNpx != 100) out<<"   spline5->SetNpx("<<fNpx<<");"<<endl;

   for (Int_t i=0;i<fNp;i++) {
      out<<"   spline5->SetPoint("<<i<<","<<fPoly[i].X()<<","<<fPoly[i].Y()<<");"<<endl;
      out<<"   spline5->SetPointCoeff("<<i<<","<<fPoly[i].B()<<","<<fPoly[i].C()<<","<<fPoly[i].D()<<","<<fPoly[i].E()<<","<<fPoly[i].F()<<");"<<endl;
   }
   out<<"   spline5->Draw("<<quote<<option<<quote<<");"<<endl;
}


//______________________________________________________________________________
void TSpline5::SetPoint(Int_t i, Double_t x, Double_t y)
{
   //set point number i.
   
   
   if (i < 0 || i >= fNp) return;
   fPoly[i].X()= x;
   fPoly[i].Y()= y;
}

//______________________________________________________________________________
void TSpline5::SetPointCoeff(Int_t i, Double_t b, Double_t c, Double_t d,
                             Double_t e, Double_t f)
{
   // set point coefficient number i
   
   if (i < 0 || i >= fNp) return;
   fPoly[i].B()= b;
   fPoly[i].C()= c;
   fPoly[i].D()= d;
   fPoly[i].E()= e;
   fPoly[i].F()= f;
}



//______________________________________________________________________________
void TSpline5::BuildCoeff()
{
   //     algorithm 600, collected algorithms from acm.
   //     algorithm appeared in acm-trans. math. software, vol.9, no. 2,
   //     jun., 1983, p. 258-259.
   //
   //     TSpline5 computes the coefficients of a quintic natural quintic spli
   //     s(x) with knots x(i) interpolating there to given function values:
   //               s(x(i)) = y(i)  for i = 1,2, ..., n.
   //     in each interval (x(i),x(i+1)) the spline function s(xx) is a
   //     polynomial of fifth degree:
   //     s(xx) = ((((f(i)*p+e(i))*p+d(i))*p+c(i))*p+b(i))*p+y(i)    (*)
   //           = ((((-f(i)*q+e(i+1))*q-d(i+1))*q+c(i+1))*q-b(i+1))*q+y(i+1)
   //     where  p = xx - x(i)  and  q = x(i+1) - xx.
   //     (note the first subscript in the second expression.)
   //     the different polynomials are pieced together so that s(x) and
   //     its derivatives up to s"" are continuous.
   //
   //        input:
   //
   //     n          number of data points, (at least three, i.e. n > 2)
   //     x(1:n)     the strictly increasing or decreasing sequence of
   //                knots.  the spacing must be such that the fifth power
   //                of x(i+1) - x(i) can be formed without overflow or
   //                underflow of exponents.
   //     y(1:n)     the prescribed function values at the knots.
   //
   //        output:
   //
   //     b,c,d,e,f  the computed spline coefficients as in (*).
   //         (1:n)  specifically
   //                b(i) = s'(x(i)), c(i) = s"(x(i))/2, d(i) = s"'(x(i))/6,
   //                e(i) = s""(x(i))/24,  f(i) = s""'(x(i))/120.
   //                f(n) is neither used nor altered.  the five arrays
   //                b,c,d,e,f must always be distinct.
   //
   //        option:
   //
   //     it is possible to specify values for the first and second
   //     derivatives of the spline function at arbitrarily many knots.
   //     this is done by relaxing the requirement that the sequence of
   //     knots be strictly increasing or decreasing.  specifically:
   //
   //     if x(j) = x(j+1) then s(x(j)) = y(j) and s'(x(j)) = y(j+1),
   //     if x(j) = x(j+1) = x(j+2) then in addition s"(x(j)) = y(j+2).
   //
   //     note that s""(x) is discontinuous at a double knot and, in
   //     addition, s"'(x) is discontinuous at a triple knot.  the
   //     subroutine assigns y(i) to y(i+1) in these cases and also to
   //     y(i+2) at a triple knot.  the representation (*) remains
   //     valid in each open interval (x(i),x(i+1)).  at a double knot,
   //     x(j) = x(j+1), the output coefficients have the following values:
   //       y(j) = s(x(j))          = y(j+1)
   //       b(j) = s'(x(j))         = b(j+1)
   //       c(j) = s"(x(j))/2       = c(j+1)
   //       d(j) = s"'(x(j))/6      = d(j+1)
   //       e(j) = s""(x(j)-0)/24     e(j+1) = s""(x(j)+0)/24
   //       f(j) = s""'(x(j)-0)/120   f(j+1) = s""'(x(j)+0)/120
   //     at a triple knot, x(j) = x(j+1) = x(j+2), the output
   //     coefficients have the following values:
   //       y(j) = s(x(j))         = y(j+1)    = y(j+2)
   //       b(j) = s'(x(j))        = b(j+1)    = b(j+2)
   //       c(j) = s"(x(j))/2      = c(j+1)    = c(j+2)
   //       d(j) = s"'((x(j)-0)/6    d(j+1) = 0  d(j+2) = s"'(x(j)+0)/6
   //       e(j) = s""(x(j)-0)/24    e(j+1) = 0  e(j+2) = s""(x(j)+0)/24
   //       f(j) = s""'(x(j)-0)/120  f(j+1) = 0  f(j+2) = s""'(x(j)+0)/120

   Int_t i, m;
   Double_t pqqr, p, q, r, s, t, u, v,
      b1, p2, p3, q2, q3, r2, pq, pr, qr;

   if (fNp <= 2) {
      return;
   }

   //     coefficients of a positive definite, pentadiagonal matrix,
   //     stored in D, E, F from 1 to n-3.
   m = fNp-2;
   q = fPoly[1].X()-fPoly[0].X();
   r = fPoly[2].X()-fPoly[1].X();
   q2 = q*q;
   r2 = r*r;
   qr = q+r;
   fPoly[0].D() = fPoly[0].E() = 0;
   if (q) fPoly[1].D() = q*6.*q2/(qr*qr);
   else fPoly[1].D() = 0;

   if (m > 1) {
      for (i = 1; i < m; ++i) {
         p = q;
         q = r;
         r = fPoly[i+2].X()-fPoly[i+1].X();
         p2 = q2;
         q2 = r2;
         r2 = r*r;
         pq = qr;
         qr = q+r;
         if (q) {
            q3 = q2*q;
            pr = p*r;
            pqqr = pq*qr;
            fPoly[i+1].D() = q3*6./(qr*qr);
            fPoly[i].D() += (q+q)*(pr*15.*pr+(p+r)*q
                                 *(pr* 20.+q2*7.)+q2*
                                 ((p2+r2)*8.+pr*21.+q2+q2))/(pqqr*pqqr);
            fPoly[i-1].D() += q3*6./(pq*pq);
            fPoly[i].E() = q2*(p*qr+pq*3.*(qr+r+r))/(pqqr*qr);
            fPoly[i-1].E() += q2*(r*pq+qr*3.*(pq+p+p))/(pqqr*pq);
            fPoly[i-1].F() = q3/pqqr;
         } else
            fPoly[i+1].D() = fPoly[i].E() = fPoly[i-1].F() = 0;
      }
   }
   if (r) fPoly[m-1].D() += r*6.*r2/(qr*qr);

   //     First and second order divided differences of the given function
   //     values, stored in b from 2 to n and in c from 3 to n
   //     respectively. care is taken of double and triple knots.
   for (i = 1; i < fNp; ++i) {
      if (fPoly[i].X() != fPoly[i-1].X()) {
         fPoly[i].B() =
            (fPoly[i].Y()-fPoly[i-1].Y())/(fPoly[i].X()-fPoly[i-1].X());
      } else {
         fPoly[i].B() = fPoly[i].Y();
         fPoly[i].Y() = fPoly[i-1].Y();
      }
   }
   for (i = 2; i < fNp; ++i) {
      if (fPoly[i].X() != fPoly[i-2].X()) {
         fPoly[i].C() =
            (fPoly[i].B()-fPoly[i-1].B())/(fPoly[i].X()-fPoly[i-2].X());
      } else {
         fPoly[i].C() = fPoly[i].B()*.5;
         fPoly[i].B() = fPoly[i-1].B();
      }
   }

   //     Solve the linear system with c(i+2) - c(i+1) as right-hand side. */
   if (m > 1) {
      p=fPoly[0].C()=fPoly[m-1].E()=fPoly[0].F()
         =fPoly[m-2].F()=fPoly[m-1].F()=0;
      fPoly[1].C() = fPoly[3].C()-fPoly[2].C();
      fPoly[1].D() = 1./fPoly[1].D();

      if (m > 2) {
         for (i = 2; i < m; ++i) {
            q = fPoly[i-1].D()*fPoly[i-1].E();
            fPoly[i].D() = 1./(fPoly[i].D()-p*fPoly[i-2].F()-q*fPoly[i-1].E());
            fPoly[i].E() -= q*fPoly[i-1].F();
            fPoly[i].C() = fPoly[i+2].C()-fPoly[i+1].C()-p*fPoly[i-2].C()
                           -q*fPoly[i-1].C();
            p = fPoly[i-1].D()*fPoly[i-1].F();
         }
      }
   }

   fPoly[fNp-2].C() = fPoly[fNp-1].C() = 0;
   if (fNp > 3)
      for (i=fNp-3; i > 0; --i)
         fPoly[i].C() = (fPoly[i].C()-fPoly[i].E()*fPoly[i+1].C()
                        -fPoly[i].F()*fPoly[i+2].C())*fPoly[i].D();

   //     Integrate the third derivative of s(x)
   m = fNp-1;
   q = fPoly[1].X()-fPoly[0].X();
   r = fPoly[2].X()-fPoly[1].X();
   b1 = fPoly[1].B();
   q3 = q*q*q;
   qr = q+r;
   if (qr) {
      v = fPoly[1].C()/qr;
      t = v;
   } else
      v = t = 0;
   if (q) fPoly[0].F() = v/q;
   else fPoly[0].F() = 0;
   for (i = 1; i < m; ++i) {
      p = q;
      q = r;
      if (i != m-1) r = fPoly[i+2].X()-fPoly[i+1].X();
      else r = 0;
      p3 = q3;
      q3 = q*q*q;
      pq = qr;
      qr = q+r;
      s = t;
      if (qr) t = (fPoly[i+1].C()-fPoly[i].C())/qr;
      else t = 0;
      u = v;
      v = t-s;
      if (pq) {
         fPoly[i].F() = fPoly[i-1].F();
         if (q) fPoly[i].F() = v/q;
         fPoly[i].E() = s*5.;
         fPoly[i].D() = (fPoly[i].C()-q*s)*10;
         fPoly[i].C() =
         fPoly[i].D()*(p-q)+(fPoly[i+1].B()-fPoly[i].B()+(u-fPoly[i].E())*
                            p3-(v+fPoly[i].E())*q3)/pq;
         fPoly[i].B() = (p*(fPoly[i+1].B()-v*q3)+q*(fPoly[i].B()-u*p3))/pq-p
         *q*(fPoly[i].D()+fPoly[i].E()*(q-p));
      } else {
         fPoly[i].C() = fPoly[i-1].C();
         fPoly[i].D() = fPoly[i].E() = fPoly[i].F() = 0;
      }
   }

   //     End points x(1) and x(n)
   p = fPoly[1].X()-fPoly[0].X();
   s = fPoly[0].F()*p*p*p;
   fPoly[0].E() = fPoly[0].D() = 0;
   fPoly[0].C() = fPoly[1].C()-s*10;
   fPoly[0].B() = b1-(fPoly[0].C()+s)*p;

   q = fPoly[fNp-1].X()-fPoly[fNp-2].X();
   t = fPoly[fNp-2].F()*q*q*q;
   fPoly[fNp-1].E() = fPoly[fNp-1].D() = 0;
   fPoly[fNp-1].C() = fPoly[fNp-2].C()+t*10;
   fPoly[fNp-1].B() += (fPoly[fNp-1].C()-t)*q;
}


//______________________________________________________________________________
void TSpline5::Test()
{
   // Test method for TSpline5
   //
   //
   //   n          number of data points.
   //   m          2*m-1 is order of spline.
   //                 m = 3 always for quintic spline.
   //   nn,nm1,mm,
   //   mm1,i,k,
   //   j,jj       temporary integer variables.
   //   z,p        temporary double precision variables.
   //   x[n]       the sequence of knots.
   //   y[n]       the prescribed function values at the knots.
   //   a[200][6]  two dimensional array whose columns are
   //                 the computed spline coefficients
   //   diff[5]    maximum values of differences of values and
   //                 derivatives to right and left of knots.
   //   com[5]     maximum values of coefficients.
   //
   //
   //   test of TSpline5 with nonequidistant knots and
   //      equidistant knots follows.

   Double_t hx;
   Double_t diff[5];
   Double_t a[1200], c[6];
   Int_t i, j, k, m, n;
   Double_t p, x[200], y[200], z;
   Int_t jj, mm, nn;
   Int_t mm1, nm1;
   Double_t com[5];

   printf("1         TEST OF TSpline5 WITH NONEQUIDISTANT KNOTS\n");
   n = 5;
   x[0] = -3;
   x[1] = -1;
   x[2] = 0;
   x[3] = 3;
   x[4] = 4;
   y[0] = 7;
   y[1] = 11;
   y[2] = 26;
   y[3] = 56;
   y[4] = 29;
   m = 3;
   mm = m << 1;
   mm1 = mm-1;
   printf("\n-N = %3d    M =%2d\n",n,m);
   TSpline5 *spline = new TSpline5("Test",x,y,n);
   for (i = 0; i < n; ++i)
      spline->GetCoeff(i,hx, a[i],a[i+200],a[i+400],
                       a[i+600],a[i+800],a[i+1000]);
   delete spline;
   for (i = 0; i < mm1; ++i) diff[i] = com[i] = 0;
   for (k = 0; k < n; ++k) {
      for (i = 0; i < mm; ++i) c[i] = a[k+i*200];
      printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
      printf("%12.8f\n",x[k]);
      if (k == n-1) {
         printf("%16.8f\n",c[0]);
      } else {
         for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
         printf("\n");
         for (i = 0; i < mm1; ++i)
            if ((z=TMath::Abs(a[k+i*200])) > com[i]) com[i] = z;
         z = x[k+1]-x[k];
         for (i = 1; i < mm; ++i)
            for (jj = i; jj < mm; ++jj) {
               j = mm+i-jj;
               c[j-2] = c[j-1]*z+c[j-2];
            }
         for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
         printf("\n");
         for (i = 0; i < mm1; ++i)
            if (!(k >= n-2 && i != 0))
               if((z = TMath::Abs(c[i]-a[k+1+i*200]))
                  > diff[i]) diff[i] = z;
      }
   }
   printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
   for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
   printf("\n");
   printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
   if (TMath::Abs(c[0]) > com[0])
      com[0] = TMath::Abs(c[0]);
   for (i = 0; i < mm1; ++i) printf("%16.8f",com[i]);
   printf("\n");
   m = 3;
   for (n = 10; n <= 100; n += 10) {
      mm = m << 1;
      mm1 = mm-1;
      nm1 = n-1;
      for (i = 0; i < nm1; i += 2) {
         x[i] = i+1;
         x[i+1] = i+2;
         y[i] = 1;
         y[i+1] = 0;
      }
      if (n % 2 != 0) {
         x[n-1] = n;
         y[n-1] = 1;
      }
      printf("\n-N = %3d    M =%2d\n",n,m);
      spline = new TSpline5("Test",x,y,n);
      for (i = 0; i < n; ++i)
         spline->GetCoeff(i,hx,a[i],a[i+200],a[i+400],
                          a[i+600],a[i+800],a[i+1000]);
      delete spline;
      for (i = 0; i < mm1; ++i)
         diff[i] = com[i] = 0;
      for (k = 0; k < n; ++k) {
         for (i = 0; i < mm; ++i)
            c[i] = a[k+i*200];
         if (n < 11) {
            printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
            printf("%12.8f\n",x[k]);
            if (k == n-1) printf("%16.8f\n",c[0]);
         }
         if (k == n-1) break;
         if (n <= 10) {
            for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
            printf("\n");
         }
         for (i = 0; i < mm1; ++i)
            if ((z=TMath::Abs(a[k+i*200])) > com[i])
               com[i] = z;
         z = x[k+1]-x[k];
         for (i = 1; i < mm; ++i)
            for (jj = i; jj < mm; ++jj) {
               j = mm+i-jj;
               c[j-2] = c[j-1]*z+c[j-2];
            }
         if (n <= 10) {
            for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
            printf("\n");
         }
         for (i = 0; i < mm1; ++i)
            if (!(k >= n-2 && i != 0))
               if ((z = TMath::Abs(c[i]-a[k+1+i*200]))
                  > diff[i]) diff[i] = z;
      }
      printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
      for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
      printf("\n");
      printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
      if (TMath::Abs(c[0]) > com[0])
         com[0] = TMath::Abs(c[0]);
      for (i = 0; i < mm1; ++i) printf("%16.8E",com[i]);
      printf("\n");
   }

   //     Test of TSpline5 with nonequidistant double knots follows
   printf("1  TEST OF TSpline5 WITH NONEQUIDISTANT DOUBLE KNOTS\n");
   n = 5;
   x[0] = -3;
   x[1] = -3;
   x[2] = -1;
   x[3] = -1;
   x[4] = 0;
   x[5] = 0;
   x[6] = 3;
   x[7] = 3;
   x[8] = 4;
   x[9] = 4;
   y[0] = 7;
   y[1] = 2;
   y[2] = 11;
   y[3] = 15;
   y[4] = 26;
   y[5] = 10;
   y[6] = 56;
   y[7] = -27;
   y[8] = 29;
   y[9] = -30;
   m = 3;
   nn = n << 1;
   mm = m << 1;
   mm1 = mm-1;
   printf("-N = %3d    M =%2d\n",n,m);
   spline = new TSpline5("Test",x,y,nn);
   for (i = 0; i < nn; ++i)
      spline->GetCoeff(i,hx,a[i],a[i+200],a[i+400],
                       a[i+600],a[i+800],a[i+1000]);
   delete spline;
   for (i = 0; i < mm1; ++i)
      diff[i] = com[i] = 0;
   for (k = 0; k < nn; ++k) {
      for (i = 0; i < mm; ++i)
         c[i] = a[k+i*200];
      printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
      printf("%12.8f\n",x[k]);
      if (k == nn-1) {
         printf("%16.8f\n",c[0]);
         break;
      }
      for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
      printf("\n");
      for (i = 0; i < mm1; ++i)
         if ((z=TMath::Abs(a[k+i*200])) > com[i]) com[i] = z;
      z = x[k+1]-x[k];
      for (i = 1; i < mm; ++i)
         for (jj = i; jj < mm; ++jj) {
            j = mm+i-jj;
            c[j-2] = c[j-1]*z+c[j-2];
         }
      for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
      printf("\n");
      for (i = 0; i < mm1; ++i)
         if (!(k >= nn-2 && i != 0))
            if ((z = TMath::Abs(c[i]-a[k+1+i*200]))
               > diff[i]) diff[i] = z;
   }
   printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
   for (i = 1; i <= mm1; ++i) {
      printf("%18.9E",diff[i-1]);
   }
   printf("\n");
   if (TMath::Abs(c[0]) > com[0])
      com[0] = TMath::Abs(c[0]);
   printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
   for (i = 0; i < mm1; ++i) printf("%16.8f",com[i]);
   printf("\n");
   m = 3;
   for (n = 10; n <= 100; n += 10) {
      nn = n << 1;
      mm = m << 1;
      mm1 = mm-1;
      p = 0;
      for (i = 0; i < n; ++i) {
         p += TMath::Abs(TMath::Sin(i+1));
         x[(i << 1)] = p;
         x[(i << 1)+1] = p;
         y[(i << 1)] = TMath::Cos(i+1)-.5;
         y[(i << 1)+1] = TMath::Cos((i << 1)+2)-.5;
      }
      printf("-N = %3d    M =%2d\n",n,m);
      spline = new TSpline5("Test",x,y,nn);
      for (i = 0; i < nn; ++i)
         spline->GetCoeff(i,hx,a[i],a[i+200],a[i+400],
                          a[i+600],a[i+800],a[i+1000]);
      delete spline;
      for (i = 0; i < mm1; ++i)
         diff[i] = com[i] = 0;
      for (k = 0; k < nn; ++k) {
         for (i = 0; i < mm; ++i)
            c[i] = a[k+i*200];
         if (n < 11) {
            printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
            printf("%12.8f\n",x[k]);
            if (k == nn-1) printf("%16.8f\n",c[0]);
         }
         if (k == nn-1) break;
         if (n <= 10) {
            for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
            printf("\n");
         }
         for (i = 0; i < mm1; ++i)
            if ((z=TMath::Abs(a[k+i*200])) > com[i]) com[i] = z;
         z = x[k+1]-x[k];
         for (i = 1; i < mm; ++i) {
            for (jj = i; jj < mm; ++jj) {
               j = mm+i-jj;
               c[j-2] = c[j-1]*z+c[j-2];
            }
         }
         if (n <= 10) {
            for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
            printf("\n");
         }
         for (i = 0; i < mm1; ++i)
            if (!(k >= nn-2 && i != 0))
               if ((z = TMath::Abs(c[i]-a[k+1+i*200]))
                  > diff[i]) diff[i] = z;
      }
      printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
      for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
      printf("\n");
      printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
      if (TMath::Abs(c[0]) > com[0])
         com[0] = TMath::Abs(c[0]);
      for (i = 0; i < mm1; ++i) printf("%18.9E",com[i]);
      printf("\n");
   }

   //     test of TSpline5 with nonequidistant knots, one double knot,
   //        one triple knot, follows.
   printf("1         TEST OF TSpline5 WITH NONEQUIDISTANT KNOTS,\n");
   printf("             ONE DOUBLE, ONE TRIPLE KNOT\n");
   n = 8;
   x[0] = -3;
   x[1] = -1;
   x[2] = -1;
   x[3] = 0;
   x[4] = 3;
   x[5] = 3;
   x[6] = 3;
   x[7] = 4;
   y[0] = 7;
   y[1] = 11;
   y[2] = 15;
   y[3] = 26;
   y[4] = 56;
   y[5] = -30;
   y[6] = -7;
   y[7] = 29;
   m = 3;
   mm = m << 1;
   mm1 = mm-1;
   printf("-N = %3d    M =%2d\n",n,m);
   spline=new TSpline5("Test",x,y,n);
   for (i = 0; i < n; ++i)
      spline->GetCoeff(i,hx,a[i],a[i+200],a[i+400],
                       a[i+600],a[i+800],a[i+1000]);
   delete spline;
   for (i = 0; i < mm1; ++i)
      diff[i] = com[i] = 0;
   for (k = 0; k < n; ++k) {
      for (i = 0; i < mm; ++i)
         c[i] = a[k+i*200];
      printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
      printf("%12.8f\n",x[k]);
      if (k == n-1) {
         printf("%16.8f\n",c[0]);
         break;
      }
      for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
      printf("\n");
      for (i = 0; i < mm1; ++i)
         if ((z=TMath::Abs(a[k+i*200])) > com[i]) com[i] = z;
      z = x[k+1]-x[k];
      for (i = 1; i < mm; ++i)
         for (jj = i; jj < mm; ++jj) {
            j = mm+i-jj;
            c[j-2] = c[j-1]*z+c[j-2];
         }
      for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
      printf("\n");
      for (i = 0; i < mm1; ++i)
         if (!(k >= n-2 && i != 0))
            if ((z = TMath::Abs(c[i]-a[k+1+i*200]))
               > diff[i]) diff[i] = z;
   }
   printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
   for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
   printf("\n");
   printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
   if (TMath::Abs(c[0]) > com[0])
      com[0] = TMath::Abs(c[0]);
   for (i = 0; i < mm1; ++i) printf("%16.8f",com[i]);
   printf("\n");

   //     Test of TSpline5 with nonequidistant knots, two double knots,
   //        one triple knot,follows.
   printf("1         TEST OF TSpline5 WITH NONEQUIDISTANT KNOTS,\n");
   printf("             TWO DOUBLE, ONE TRIPLE KNOT\n");
   n = 10;
   x[0] = 0;
   x[1] = 2;
   x[2] = 2;
   x[3] = 3;
   x[4] = 3;
   x[5] = 3;
   x[6] = 5;
   x[7] = 8;
   x[8] = 9;
   x[9] = 9;
   y[0] = 163;
   y[1] = 237;
   y[2] = -127;
   y[3] = 119;
   y[4] = -65;
   y[5] = 192;
   y[6] = 293;
   y[7] = 326;
   y[8] = 0;
   y[9] = -414;
   m = 3;
   mm = m << 1;
   mm1 = mm-1;
   printf("-N = %3d    M =%2d\n",n,m);
   spline = new TSpline5("Test",x,y,n);
   for (i = 0; i < n; ++i)
      spline->GetCoeff(i,hx,a[i],a[i+200],a[i+400],
                       a[i+600],a[i+800],a[i+1000]);
   delete spline;
   for (i = 0; i < mm1; ++i)
      diff[i] = com[i] = 0;
   for (k = 0; k < n; ++k) {
      for (i = 0; i < mm; ++i)
         c[i] = a[k+i*200];
      printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
      printf("%12.8f\n",x[k]);
      if (k == n-1) {
         printf("%16.8f\n",c[0]);
         break;
      }
      for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
      printf("\n");
      for (i = 0; i < mm1; ++i)
         if ((z=TMath::Abs(a[k+i*200])) > com[i]) com[i] = z;
      z = x[k+1]-x[k];
      for (i = 1; i < mm; ++i)
         for (jj = i; jj < mm; ++jj) {
            j = mm+i-jj;
            c[j-2] = c[j-1]*z+c[j-2];
         }
      for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
      printf("\n");
      for (i = 0; i < mm1; ++i)
         if (!(k >= n-2 && i != 0))
            if((z = TMath::Abs(c[i]-a[k+1+i*200]))
               > diff[i]) diff[i] = z;
   }
   printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
   for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
   printf("\n");
   printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
   if (TMath::Abs(c[0]) > com[0])
      com[0] = TMath::Abs(c[0]);
   for (i = 0; i < mm1; ++i) printf("%16.8f",com[i]);
   printf("\n");
}


//______________________________________________________________________________
void TSpline5::Streamer(TBuffer &R__b)
{
   // Stream an object of class TSpline5.

   if (R__b.IsReading()) {
      UInt_t R__s, R__c;
      Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
      if (R__v > 1) {
         R__b.ReadClassBuffer(TSpline5::Class(), this, R__v, R__s, R__c);
         return;
      }
      //====process old versions before automatic schema evolution
      TSpline::Streamer(R__b);
      if (fNp > 0) {
         fPoly = new TSplinePoly5[fNp];
         for(Int_t i=0; i<fNp; ++i) {
            fPoly[i].Streamer(R__b);
         }
      }
      //      R__b >> fPoly;
   } else {
      R__b.WriteClassBuffer(TSpline5::Class(),this);
   }
}
 TSpline.cxx:1
 TSpline.cxx:2
 TSpline.cxx:3
 TSpline.cxx:4
 TSpline.cxx:5
 TSpline.cxx:6
 TSpline.cxx:7
 TSpline.cxx:8
 TSpline.cxx:9
 TSpline.cxx:10
 TSpline.cxx:11
 TSpline.cxx:12
 TSpline.cxx:13
 TSpline.cxx:14
 TSpline.cxx:15
 TSpline.cxx:16
 TSpline.cxx:17
 TSpline.cxx:18
 TSpline.cxx:19
 TSpline.cxx:20
 TSpline.cxx:21
 TSpline.cxx:22
 TSpline.cxx:23
 TSpline.cxx:24
 TSpline.cxx:25
 TSpline.cxx:26
 TSpline.cxx:27
 TSpline.cxx:28
 TSpline.cxx:29
 TSpline.cxx:30
 TSpline.cxx:31
 TSpline.cxx:32
 TSpline.cxx:33
 TSpline.cxx:34
 TSpline.cxx:35
 TSpline.cxx:36
 TSpline.cxx:37
 TSpline.cxx:38
 TSpline.cxx:39
 TSpline.cxx:40
 TSpline.cxx:41
 TSpline.cxx:42
 TSpline.cxx:43
 TSpline.cxx:44
 TSpline.cxx:45
 TSpline.cxx:46
 TSpline.cxx:47
 TSpline.cxx:48
 TSpline.cxx:49
 TSpline.cxx:50
 TSpline.cxx:51
 TSpline.cxx:52
 TSpline.cxx:53
 TSpline.cxx:54
 TSpline.cxx:55
 TSpline.cxx:56
 TSpline.cxx:57
 TSpline.cxx:58
 TSpline.cxx:59
 TSpline.cxx:60
 TSpline.cxx:61
 TSpline.cxx:62
 TSpline.cxx:63
 TSpline.cxx:64
 TSpline.cxx:65
 TSpline.cxx:66
 TSpline.cxx:67
 TSpline.cxx:68
 TSpline.cxx:69
 TSpline.cxx:70
 TSpline.cxx:71
 TSpline.cxx:72
 TSpline.cxx:73
 TSpline.cxx:74
 TSpline.cxx:75
 TSpline.cxx:76
 TSpline.cxx:77
 TSpline.cxx:78
 TSpline.cxx:79
 TSpline.cxx:80
 TSpline.cxx:81
 TSpline.cxx:82
 TSpline.cxx:83
 TSpline.cxx:84
 TSpline.cxx:85
 TSpline.cxx:86
 TSpline.cxx:87
 TSpline.cxx:88
 TSpline.cxx:89
 TSpline.cxx:90
 TSpline.cxx:91
 TSpline.cxx:92
 TSpline.cxx:93
 TSpline.cxx:94
 TSpline.cxx:95
 TSpline.cxx:96
 TSpline.cxx:97
 TSpline.cxx:98
 TSpline.cxx:99
 TSpline.cxx:100
 TSpline.cxx:101
 TSpline.cxx:102
 TSpline.cxx:103
 TSpline.cxx:104
 TSpline.cxx:105
 TSpline.cxx:106
 TSpline.cxx:107
 TSpline.cxx:108
 TSpline.cxx:109
 TSpline.cxx:110
 TSpline.cxx:111
 TSpline.cxx:112
 TSpline.cxx:113
 TSpline.cxx:114
 TSpline.cxx:115
 TSpline.cxx:116
 TSpline.cxx:117
 TSpline.cxx:118
 TSpline.cxx:119
 TSpline.cxx:120
 TSpline.cxx:121
 TSpline.cxx:122
 TSpline.cxx:123
 TSpline.cxx:124
 TSpline.cxx:125
 TSpline.cxx:126
 TSpline.cxx:127
 TSpline.cxx:128
 TSpline.cxx:129
 TSpline.cxx:130
 TSpline.cxx:131
 TSpline.cxx:132
 TSpline.cxx:133
 TSpline.cxx:134
 TSpline.cxx:135
 TSpline.cxx:136
 TSpline.cxx:137
 TSpline.cxx:138
 TSpline.cxx:139
 TSpline.cxx:140
 TSpline.cxx:141
 TSpline.cxx:142
 TSpline.cxx:143
 TSpline.cxx:144
 TSpline.cxx:145
 TSpline.cxx:146
 TSpline.cxx:147
 TSpline.cxx:148
 TSpline.cxx:149
 TSpline.cxx:150
 TSpline.cxx:151
 TSpline.cxx:152
 TSpline.cxx:153
 TSpline.cxx:154
 TSpline.cxx:155
 TSpline.cxx:156
 TSpline.cxx:157
 TSpline.cxx:158
 TSpline.cxx:159
 TSpline.cxx:160
 TSpline.cxx:161
 TSpline.cxx:162
 TSpline.cxx:163
 TSpline.cxx:164
 TSpline.cxx:165
 TSpline.cxx:166
 TSpline.cxx:167
 TSpline.cxx:168
 TSpline.cxx:169
 TSpline.cxx:170
 TSpline.cxx:171
 TSpline.cxx:172
 TSpline.cxx:173
 TSpline.cxx:174
 TSpline.cxx:175
 TSpline.cxx:176
 TSpline.cxx:177
 TSpline.cxx:178
 TSpline.cxx:179
 TSpline.cxx:180
 TSpline.cxx:181
 TSpline.cxx:182
 TSpline.cxx:183
 TSpline.cxx:184
 TSpline.cxx:185
 TSpline.cxx:186
 TSpline.cxx:187
 TSpline.cxx:188
 TSpline.cxx:189
 TSpline.cxx:190
 TSpline.cxx:191
 TSpline.cxx:192
 TSpline.cxx:193
 TSpline.cxx:194
 TSpline.cxx:195
 TSpline.cxx:196
 TSpline.cxx:197
 TSpline.cxx:198
 TSpline.cxx:199
 TSpline.cxx:200
 TSpline.cxx:201
 TSpline.cxx:202
 TSpline.cxx:203
 TSpline.cxx:204
 TSpline.cxx:205
 TSpline.cxx:206
 TSpline.cxx:207
 TSpline.cxx:208
 TSpline.cxx:209
 TSpline.cxx:210
 TSpline.cxx:211
 TSpline.cxx:212
 TSpline.cxx:213
 TSpline.cxx:214
 TSpline.cxx:215
 TSpline.cxx:216
 TSpline.cxx:217
 TSpline.cxx:218
 TSpline.cxx:219
 TSpline.cxx:220
 TSpline.cxx:221
 TSpline.cxx:222
 TSpline.cxx:223
 TSpline.cxx:224
 TSpline.cxx:225
 TSpline.cxx:226
 TSpline.cxx:227
 TSpline.cxx:228
 TSpline.cxx:229
 TSpline.cxx:230
 TSpline.cxx:231
 TSpline.cxx:232
 TSpline.cxx:233
 TSpline.cxx:234
 TSpline.cxx:235
 TSpline.cxx:236
 TSpline.cxx:237
 TSpline.cxx:238
 TSpline.cxx:239
 TSpline.cxx:240
 TSpline.cxx:241
 TSpline.cxx:242
 TSpline.cxx:243
 TSpline.cxx:244
 TSpline.cxx:245
 TSpline.cxx:246
 TSpline.cxx:247
 TSpline.cxx:248
 TSpline.cxx:249
 TSpline.cxx:250
 TSpline.cxx:251
 TSpline.cxx:252
 TSpline.cxx:253
 TSpline.cxx:254
 TSpline.cxx:255
 TSpline.cxx:256
 TSpline.cxx:257
 TSpline.cxx:258
 TSpline.cxx:259
 TSpline.cxx:260
 TSpline.cxx:261
 TSpline.cxx:262
 TSpline.cxx:263
 TSpline.cxx:264
 TSpline.cxx:265
 TSpline.cxx:266
 TSpline.cxx:267
 TSpline.cxx:268
 TSpline.cxx:269
 TSpline.cxx:270
 TSpline.cxx:271
 TSpline.cxx:272
 TSpline.cxx:273
 TSpline.cxx:274
 TSpline.cxx:275
 TSpline.cxx:276
 TSpline.cxx:277
 TSpline.cxx:278
 TSpline.cxx:279
 TSpline.cxx:280
 TSpline.cxx:281
 TSpline.cxx:282
 TSpline.cxx:283
 TSpline.cxx:284
 TSpline.cxx:285
 TSpline.cxx:286
 TSpline.cxx:287
 TSpline.cxx:288
 TSpline.cxx:289
 TSpline.cxx:290
 TSpline.cxx:291
 TSpline.cxx:292
 TSpline.cxx:293
 TSpline.cxx:294
 TSpline.cxx:295
 TSpline.cxx:296
 TSpline.cxx:297
 TSpline.cxx:298
 TSpline.cxx:299
 TSpline.cxx:300
 TSpline.cxx:301
 TSpline.cxx:302
 TSpline.cxx:303
 TSpline.cxx:304
 TSpline.cxx:305
 TSpline.cxx:306
 TSpline.cxx:307
 TSpline.cxx:308
 TSpline.cxx:309
 TSpline.cxx:310
 TSpline.cxx:311
 TSpline.cxx:312
 TSpline.cxx:313
 TSpline.cxx:314
 TSpline.cxx:315
 TSpline.cxx:316
 TSpline.cxx:317
 TSpline.cxx:318
 TSpline.cxx:319
 TSpline.cxx:320
 TSpline.cxx:321
 TSpline.cxx:322
 TSpline.cxx:323
 TSpline.cxx:324
 TSpline.cxx:325
 TSpline.cxx:326
 TSpline.cxx:327
 TSpline.cxx:328
 TSpline.cxx:329
 TSpline.cxx:330
 TSpline.cxx:331
 TSpline.cxx:332
 TSpline.cxx:333
 TSpline.cxx:334
 TSpline.cxx:335
 TSpline.cxx:336
 TSpline.cxx:337
 TSpline.cxx:338
 TSpline.cxx:339
 TSpline.cxx:340
 TSpline.cxx:341
 TSpline.cxx:342
 TSpline.cxx:343
 TSpline.cxx:344
 TSpline.cxx:345
 TSpline.cxx:346
 TSpline.cxx:347
 TSpline.cxx:348
 TSpline.cxx:349
 TSpline.cxx:350
 TSpline.cxx:351
 TSpline.cxx:352
 TSpline.cxx:353
 TSpline.cxx:354
 TSpline.cxx:355
 TSpline.cxx:356
 TSpline.cxx:357
 TSpline.cxx:358
 TSpline.cxx:359
 TSpline.cxx:360
 TSpline.cxx:361
 TSpline.cxx:362
 TSpline.cxx:363
 TSpline.cxx:364
 TSpline.cxx:365
 TSpline.cxx:366
 TSpline.cxx:367
 TSpline.cxx:368
 TSpline.cxx:369
 TSpline.cxx:370
 TSpline.cxx:371
 TSpline.cxx:372
 TSpline.cxx:373
 TSpline.cxx:374
 TSpline.cxx:375
 TSpline.cxx:376
 TSpline.cxx:377
 TSpline.cxx:378
 TSpline.cxx:379
 TSpline.cxx:380
 TSpline.cxx:381
 TSpline.cxx:382
 TSpline.cxx:383
 TSpline.cxx:384
 TSpline.cxx:385
 TSpline.cxx:386
 TSpline.cxx:387
 TSpline.cxx:388
 TSpline.cxx:389
 TSpline.cxx:390
 TSpline.cxx:391
 TSpline.cxx:392
 TSpline.cxx:393
 TSpline.cxx:394
 TSpline.cxx:395
 TSpline.cxx:396
 TSpline.cxx:397
 TSpline.cxx:398
 TSpline.cxx:399
 TSpline.cxx:400
 TSpline.cxx:401
 TSpline.cxx:402
 TSpline.cxx:403
 TSpline.cxx:404
 TSpline.cxx:405
 TSpline.cxx:406
 TSpline.cxx:407
 TSpline.cxx:408
 TSpline.cxx:409
 TSpline.cxx:410
 TSpline.cxx:411
 TSpline.cxx:412
 TSpline.cxx:413
 TSpline.cxx:414
 TSpline.cxx:415
 TSpline.cxx:416
 TSpline.cxx:417
 TSpline.cxx:418
 TSpline.cxx:419
 TSpline.cxx:420
 TSpline.cxx:421
 TSpline.cxx:422
 TSpline.cxx:423
 TSpline.cxx:424
 TSpline.cxx:425
 TSpline.cxx:426
 TSpline.cxx:427
 TSpline.cxx:428
 TSpline.cxx:429
 TSpline.cxx:430
 TSpline.cxx:431
 TSpline.cxx:432
 TSpline.cxx:433
 TSpline.cxx:434
 TSpline.cxx:435
 TSpline.cxx:436
 TSpline.cxx:437
 TSpline.cxx:438
 TSpline.cxx:439
 TSpline.cxx:440
 TSpline.cxx:441
 TSpline.cxx:442
 TSpline.cxx:443
 TSpline.cxx:444
 TSpline.cxx:445
 TSpline.cxx:446
 TSpline.cxx:447
 TSpline.cxx:448
 TSpline.cxx:449
 TSpline.cxx:450
 TSpline.cxx:451
 TSpline.cxx:452
 TSpline.cxx:453
 TSpline.cxx:454
 TSpline.cxx:455
 TSpline.cxx:456
 TSpline.cxx:457
 TSpline.cxx:458
 TSpline.cxx:459
 TSpline.cxx:460
 TSpline.cxx:461
 TSpline.cxx:462
 TSpline.cxx:463
 TSpline.cxx:464
 TSpline.cxx:465
 TSpline.cxx:466
 TSpline.cxx:467
 TSpline.cxx:468
 TSpline.cxx:469
 TSpline.cxx:470
 TSpline.cxx:471
 TSpline.cxx:472
 TSpline.cxx:473
 TSpline.cxx:474
 TSpline.cxx:475
 TSpline.cxx:476
 TSpline.cxx:477
 TSpline.cxx:478
 TSpline.cxx:479
 TSpline.cxx:480
 TSpline.cxx:481
 TSpline.cxx:482
 TSpline.cxx:483
 TSpline.cxx:484
 TSpline.cxx:485
 TSpline.cxx:486
 TSpline.cxx:487
 TSpline.cxx:488
 TSpline.cxx:489
 TSpline.cxx:490
 TSpline.cxx:491
 TSpline.cxx:492
 TSpline.cxx:493
 TSpline.cxx:494
 TSpline.cxx:495
 TSpline.cxx:496
 TSpline.cxx:497
 TSpline.cxx:498
 TSpline.cxx:499
 TSpline.cxx:500
 TSpline.cxx:501
 TSpline.cxx:502
 TSpline.cxx:503
 TSpline.cxx:504
 TSpline.cxx:505
 TSpline.cxx:506
 TSpline.cxx:507
 TSpline.cxx:508
 TSpline.cxx:509
 TSpline.cxx:510
 TSpline.cxx:511
 TSpline.cxx:512
 TSpline.cxx:513
 TSpline.cxx:514
 TSpline.cxx:515
 TSpline.cxx:516
 TSpline.cxx:517
 TSpline.cxx:518
 TSpline.cxx:519
 TSpline.cxx:520
 TSpline.cxx:521
 TSpline.cxx:522
 TSpline.cxx:523
 TSpline.cxx:524
 TSpline.cxx:525
 TSpline.cxx:526
 TSpline.cxx:527
 TSpline.cxx:528
 TSpline.cxx:529
 TSpline.cxx:530
 TSpline.cxx:531
 TSpline.cxx:532
 TSpline.cxx:533
 TSpline.cxx:534
 TSpline.cxx:535
 TSpline.cxx:536
 TSpline.cxx:537
 TSpline.cxx:538
 TSpline.cxx:539
 TSpline.cxx:540
 TSpline.cxx:541
 TSpline.cxx:542
 TSpline.cxx:543
 TSpline.cxx:544
 TSpline.cxx:545
 TSpline.cxx:546
 TSpline.cxx:547
 TSpline.cxx:548
 TSpline.cxx:549
 TSpline.cxx:550
 TSpline.cxx:551
 TSpline.cxx:552
 TSpline.cxx:553
 TSpline.cxx:554
 TSpline.cxx:555
 TSpline.cxx:556
 TSpline.cxx:557
 TSpline.cxx:558
 TSpline.cxx:559
 TSpline.cxx:560
 TSpline.cxx:561
 TSpline.cxx:562
 TSpline.cxx:563
 TSpline.cxx:564
 TSpline.cxx:565
 TSpline.cxx:566
 TSpline.cxx:567
 TSpline.cxx:568
 TSpline.cxx:569
 TSpline.cxx:570
 TSpline.cxx:571
 TSpline.cxx:572
 TSpline.cxx:573
 TSpline.cxx:574
 TSpline.cxx:575
 TSpline.cxx:576
 TSpline.cxx:577
 TSpline.cxx:578
 TSpline.cxx:579
 TSpline.cxx:580
 TSpline.cxx:581
 TSpline.cxx:582
 TSpline.cxx:583
 TSpline.cxx:584
 TSpline.cxx:585
 TSpline.cxx:586
 TSpline.cxx:587
 TSpline.cxx:588
 TSpline.cxx:589
 TSpline.cxx:590
 TSpline.cxx:591
 TSpline.cxx:592
 TSpline.cxx:593
 TSpline.cxx:594
 TSpline.cxx:595
 TSpline.cxx:596
 TSpline.cxx:597
 TSpline.cxx:598
 TSpline.cxx:599
 TSpline.cxx:600
 TSpline.cxx:601
 TSpline.cxx:602
 TSpline.cxx:603
 TSpline.cxx:604
 TSpline.cxx:605
 TSpline.cxx:606
 TSpline.cxx:607
 TSpline.cxx:608
 TSpline.cxx:609
 TSpline.cxx:610
 TSpline.cxx:611
 TSpline.cxx:612
 TSpline.cxx:613
 TSpline.cxx:614
 TSpline.cxx:615
 TSpline.cxx:616
 TSpline.cxx:617
 TSpline.cxx:618
 TSpline.cxx:619
 TSpline.cxx:620
 TSpline.cxx:621
 TSpline.cxx:622
 TSpline.cxx:623
 TSpline.cxx:624
 TSpline.cxx:625
 TSpline.cxx:626
 TSpline.cxx:627
 TSpline.cxx:628
 TSpline.cxx:629
 TSpline.cxx:630
 TSpline.cxx:631
 TSpline.cxx:632
 TSpline.cxx:633
 TSpline.cxx:634
 TSpline.cxx:635
 TSpline.cxx:636
 TSpline.cxx:637
 TSpline.cxx:638
 TSpline.cxx:639
 TSpline.cxx:640
 TSpline.cxx:641
 TSpline.cxx:642
 TSpline.cxx:643
 TSpline.cxx:644
 TSpline.cxx:645
 TSpline.cxx:646
 TSpline.cxx:647
 TSpline.cxx:648
 TSpline.cxx:649
 TSpline.cxx:650
 TSpline.cxx:651
 TSpline.cxx:652
 TSpline.cxx:653
 TSpline.cxx:654
 TSpline.cxx:655
 TSpline.cxx:656
 TSpline.cxx:657
 TSpline.cxx:658
 TSpline.cxx:659
 TSpline.cxx:660
 TSpline.cxx:661
 TSpline.cxx:662
 TSpline.cxx:663
 TSpline.cxx:664
 TSpline.cxx:665
 TSpline.cxx:666
 TSpline.cxx:667
 TSpline.cxx:668
 TSpline.cxx:669
 TSpline.cxx:670
 TSpline.cxx:671
 TSpline.cxx:672
 TSpline.cxx:673
 TSpline.cxx:674
 TSpline.cxx:675
 TSpline.cxx:676
 TSpline.cxx:677
 TSpline.cxx:678
 TSpline.cxx:679
 TSpline.cxx:680
 TSpline.cxx:681
 TSpline.cxx:682
 TSpline.cxx:683
 TSpline.cxx:684
 TSpline.cxx:685
 TSpline.cxx:686
 TSpline.cxx:687
 TSpline.cxx:688
 TSpline.cxx:689
 TSpline.cxx:690
 TSpline.cxx:691
 TSpline.cxx:692
 TSpline.cxx:693
 TSpline.cxx:694
 TSpline.cxx:695
 TSpline.cxx:696
 TSpline.cxx:697
 TSpline.cxx:698
 TSpline.cxx:699
 TSpline.cxx:700
 TSpline.cxx:701
 TSpline.cxx:702
 TSpline.cxx:703
 TSpline.cxx:704
 TSpline.cxx:705
 TSpline.cxx:706
 TSpline.cxx:707
 TSpline.cxx:708
 TSpline.cxx:709
 TSpline.cxx:710
 TSpline.cxx:711
 TSpline.cxx:712
 TSpline.cxx:713
 TSpline.cxx:714
 TSpline.cxx:715
 TSpline.cxx:716
 TSpline.cxx:717
 TSpline.cxx:718
 TSpline.cxx:719
 TSpline.cxx:720
 TSpline.cxx:721
 TSpline.cxx:722
 TSpline.cxx:723
 TSpline.cxx:724
 TSpline.cxx:725
 TSpline.cxx:726
 TSpline.cxx:727
 TSpline.cxx:728
 TSpline.cxx:729
 TSpline.cxx:730
 TSpline.cxx:731
 TSpline.cxx:732
 TSpline.cxx:733
 TSpline.cxx:734
 TSpline.cxx:735
 TSpline.cxx:736
 TSpline.cxx:737
 TSpline.cxx:738
 TSpline.cxx:739
 TSpline.cxx:740
 TSpline.cxx:741
 TSpline.cxx:742
 TSpline.cxx:743
 TSpline.cxx:744
 TSpline.cxx:745
 TSpline.cxx:746
 TSpline.cxx:747
 TSpline.cxx:748
 TSpline.cxx:749
 TSpline.cxx:750
 TSpline.cxx:751
 TSpline.cxx:752
 TSpline.cxx:753
 TSpline.cxx:754
 TSpline.cxx:755
 TSpline.cxx:756
 TSpline.cxx:757
 TSpline.cxx:758
 TSpline.cxx:759
 TSpline.cxx:760
 TSpline.cxx:761
 TSpline.cxx:762
 TSpline.cxx:763
 TSpline.cxx:764
 TSpline.cxx:765
 TSpline.cxx:766
 TSpline.cxx:767
 TSpline.cxx:768
 TSpline.cxx:769
 TSpline.cxx:770
 TSpline.cxx:771
 TSpline.cxx:772
 TSpline.cxx:773
 TSpline.cxx:774
 TSpline.cxx:775
 TSpline.cxx:776
 TSpline.cxx:777
 TSpline.cxx:778
 TSpline.cxx:779
 TSpline.cxx:780
 TSpline.cxx:781
 TSpline.cxx:782
 TSpline.cxx:783
 TSpline.cxx:784
 TSpline.cxx:785
 TSpline.cxx:786
 TSpline.cxx:787
 TSpline.cxx:788
 TSpline.cxx:789
 TSpline.cxx:790
 TSpline.cxx:791
 TSpline.cxx:792
 TSpline.cxx:793
 TSpline.cxx:794
 TSpline.cxx:795
 TSpline.cxx:796
 TSpline.cxx:797
 TSpline.cxx:798
 TSpline.cxx:799
 TSpline.cxx:800
 TSpline.cxx:801
 TSpline.cxx:802
 TSpline.cxx:803
 TSpline.cxx:804
 TSpline.cxx:805
 TSpline.cxx:806
 TSpline.cxx:807
 TSpline.cxx:808
 TSpline.cxx:809
 TSpline.cxx:810
 TSpline.cxx:811
 TSpline.cxx:812
 TSpline.cxx:813
 TSpline.cxx:814
 TSpline.cxx:815
 TSpline.cxx:816
 TSpline.cxx:817
 TSpline.cxx:818
 TSpline.cxx:819
 TSpline.cxx:820
 TSpline.cxx:821
 TSpline.cxx:822
 TSpline.cxx:823
 TSpline.cxx:824
 TSpline.cxx:825
 TSpline.cxx:826
 TSpline.cxx:827
 TSpline.cxx:828
 TSpline.cxx:829
 TSpline.cxx:830
 TSpline.cxx:831
 TSpline.cxx:832
 TSpline.cxx:833
 TSpline.cxx:834
 TSpline.cxx:835
 TSpline.cxx:836
 TSpline.cxx:837
 TSpline.cxx:838
 TSpline.cxx:839
 TSpline.cxx:840
 TSpline.cxx:841
 TSpline.cxx:842
 TSpline.cxx:843
 TSpline.cxx:844
 TSpline.cxx:845
 TSpline.cxx:846
 TSpline.cxx:847
 TSpline.cxx:848
 TSpline.cxx:849
 TSpline.cxx:850
 TSpline.cxx:851
 TSpline.cxx:852
 TSpline.cxx:853
 TSpline.cxx:854
 TSpline.cxx:855
 TSpline.cxx:856
 TSpline.cxx:857
 TSpline.cxx:858
 TSpline.cxx:859
 TSpline.cxx:860
 TSpline.cxx:861
 TSpline.cxx:862
 TSpline.cxx:863
 TSpline.cxx:864
 TSpline.cxx:865
 TSpline.cxx:866
 TSpline.cxx:867
 TSpline.cxx:868
 TSpline.cxx:869
 TSpline.cxx:870
 TSpline.cxx:871
 TSpline.cxx:872
 TSpline.cxx:873
 TSpline.cxx:874
 TSpline.cxx:875
 TSpline.cxx:876
 TSpline.cxx:877
 TSpline.cxx:878
 TSpline.cxx:879
 TSpline.cxx:880
 TSpline.cxx:881
 TSpline.cxx:882
 TSpline.cxx:883
 TSpline.cxx:884
 TSpline.cxx:885
 TSpline.cxx:886
 TSpline.cxx:887
 TSpline.cxx:888
 TSpline.cxx:889
 TSpline.cxx:890
 TSpline.cxx:891
 TSpline.cxx:892
 TSpline.cxx:893
 TSpline.cxx:894
 TSpline.cxx:895
 TSpline.cxx:896
 TSpline.cxx:897
 TSpline.cxx:898
 TSpline.cxx:899
 TSpline.cxx:900
 TSpline.cxx:901
 TSpline.cxx:902
 TSpline.cxx:903
 TSpline.cxx:904
 TSpline.cxx:905
 TSpline.cxx:906
 TSpline.cxx:907
 TSpline.cxx:908
 TSpline.cxx:909
 TSpline.cxx:910
 TSpline.cxx:911
 TSpline.cxx:912
 TSpline.cxx:913
 TSpline.cxx:914
 TSpline.cxx:915
 TSpline.cxx:916
 TSpline.cxx:917
 TSpline.cxx:918
 TSpline.cxx:919
 TSpline.cxx:920
 TSpline.cxx:921
 TSpline.cxx:922
 TSpline.cxx:923
 TSpline.cxx:924
 TSpline.cxx:925
 TSpline.cxx:926
 TSpline.cxx:927
 TSpline.cxx:928
 TSpline.cxx:929
 TSpline.cxx:930
 TSpline.cxx:931
 TSpline.cxx:932
 TSpline.cxx:933
 TSpline.cxx:934
 TSpline.cxx:935
 TSpline.cxx:936
 TSpline.cxx:937
 TSpline.cxx:938
 TSpline.cxx:939
 TSpline.cxx:940
 TSpline.cxx:941
 TSpline.cxx:942
 TSpline.cxx:943
 TSpline.cxx:944
 TSpline.cxx:945
 TSpline.cxx:946
 TSpline.cxx:947
 TSpline.cxx:948
 TSpline.cxx:949
 TSpline.cxx:950
 TSpline.cxx:951
 TSpline.cxx:952
 TSpline.cxx:953
 TSpline.cxx:954
 TSpline.cxx:955
 TSpline.cxx:956
 TSpline.cxx:957
 TSpline.cxx:958
 TSpline.cxx:959
 TSpline.cxx:960
 TSpline.cxx:961
 TSpline.cxx:962
 TSpline.cxx:963
 TSpline.cxx:964
 TSpline.cxx:965
 TSpline.cxx:966
 TSpline.cxx:967
 TSpline.cxx:968
 TSpline.cxx:969
 TSpline.cxx:970
 TSpline.cxx:971
 TSpline.cxx:972
 TSpline.cxx:973
 TSpline.cxx:974
 TSpline.cxx:975
 TSpline.cxx:976
 TSpline.cxx:977
 TSpline.cxx:978
 TSpline.cxx:979
 TSpline.cxx:980
 TSpline.cxx:981
 TSpline.cxx:982
 TSpline.cxx:983
 TSpline.cxx:984
 TSpline.cxx:985
 TSpline.cxx:986
 TSpline.cxx:987
 TSpline.cxx:988
 TSpline.cxx:989
 TSpline.cxx:990
 TSpline.cxx:991
 TSpline.cxx:992
 TSpline.cxx:993
 TSpline.cxx:994
 TSpline.cxx:995
 TSpline.cxx:996
 TSpline.cxx:997
 TSpline.cxx:998
 TSpline.cxx:999
 TSpline.cxx:1000
 TSpline.cxx:1001
 TSpline.cxx:1002
 TSpline.cxx:1003
 TSpline.cxx:1004
 TSpline.cxx:1005
 TSpline.cxx:1006
 TSpline.cxx:1007
 TSpline.cxx:1008
 TSpline.cxx:1009
 TSpline.cxx:1010
 TSpline.cxx:1011
 TSpline.cxx:1012
 TSpline.cxx:1013
 TSpline.cxx:1014
 TSpline.cxx:1015
 TSpline.cxx:1016
 TSpline.cxx:1017
 TSpline.cxx:1018
 TSpline.cxx:1019
 TSpline.cxx:1020
 TSpline.cxx:1021
 TSpline.cxx:1022
 TSpline.cxx:1023
 TSpline.cxx:1024
 TSpline.cxx:1025
 TSpline.cxx:1026
 TSpline.cxx:1027
 TSpline.cxx:1028
 TSpline.cxx:1029
 TSpline.cxx:1030
 TSpline.cxx:1031
 TSpline.cxx:1032
 TSpline.cxx:1033
 TSpline.cxx:1034
 TSpline.cxx:1035
 TSpline.cxx:1036
 TSpline.cxx:1037
 TSpline.cxx:1038
 TSpline.cxx:1039
 TSpline.cxx:1040
 TSpline.cxx:1041
 TSpline.cxx:1042
 TSpline.cxx:1043
 TSpline.cxx:1044
 TSpline.cxx:1045
 TSpline.cxx:1046
 TSpline.cxx:1047
 TSpline.cxx:1048
 TSpline.cxx:1049
 TSpline.cxx:1050
 TSpline.cxx:1051
 TSpline.cxx:1052
 TSpline.cxx:1053
 TSpline.cxx:1054
 TSpline.cxx:1055
 TSpline.cxx:1056
 TSpline.cxx:1057
 TSpline.cxx:1058
 TSpline.cxx:1059
 TSpline.cxx:1060
 TSpline.cxx:1061
 TSpline.cxx:1062
 TSpline.cxx:1063
 TSpline.cxx:1064
 TSpline.cxx:1065
 TSpline.cxx:1066
 TSpline.cxx:1067
 TSpline.cxx:1068
 TSpline.cxx:1069
 TSpline.cxx:1070
 TSpline.cxx:1071
 TSpline.cxx:1072
 TSpline.cxx:1073
 TSpline.cxx:1074
 TSpline.cxx:1075
 TSpline.cxx:1076
 TSpline.cxx:1077
 TSpline.cxx:1078
 TSpline.cxx:1079
 TSpline.cxx:1080
 TSpline.cxx:1081
 TSpline.cxx:1082
 TSpline.cxx:1083
 TSpline.cxx:1084
 TSpline.cxx:1085
 TSpline.cxx:1086
 TSpline.cxx:1087
 TSpline.cxx:1088
 TSpline.cxx:1089
 TSpline.cxx:1090
 TSpline.cxx:1091
 TSpline.cxx:1092
 TSpline.cxx:1093
 TSpline.cxx:1094
 TSpline.cxx:1095
 TSpline.cxx:1096
 TSpline.cxx:1097
 TSpline.cxx:1098
 TSpline.cxx:1099
 TSpline.cxx:1100
 TSpline.cxx:1101
 TSpline.cxx:1102
 TSpline.cxx:1103
 TSpline.cxx:1104
 TSpline.cxx:1105
 TSpline.cxx:1106
 TSpline.cxx:1107
 TSpline.cxx:1108
 TSpline.cxx:1109
 TSpline.cxx:1110
 TSpline.cxx:1111
 TSpline.cxx:1112
 TSpline.cxx:1113
 TSpline.cxx:1114
 TSpline.cxx:1115
 TSpline.cxx:1116
 TSpline.cxx:1117
 TSpline.cxx:1118
 TSpline.cxx:1119
 TSpline.cxx:1120
 TSpline.cxx:1121
 TSpline.cxx:1122
 TSpline.cxx:1123
 TSpline.cxx:1124
 TSpline.cxx:1125
 TSpline.cxx:1126
 TSpline.cxx:1127
 TSpline.cxx:1128
 TSpline.cxx:1129
 TSpline.cxx:1130
 TSpline.cxx:1131
 TSpline.cxx:1132
 TSpline.cxx:1133
 TSpline.cxx:1134
 TSpline.cxx:1135
 TSpline.cxx:1136
 TSpline.cxx:1137
 TSpline.cxx:1138
 TSpline.cxx:1139
 TSpline.cxx:1140
 TSpline.cxx:1141
 TSpline.cxx:1142
 TSpline.cxx:1143
 TSpline.cxx:1144
 TSpline.cxx:1145
 TSpline.cxx:1146
 TSpline.cxx:1147
 TSpline.cxx:1148
 TSpline.cxx:1149
 TSpline.cxx:1150
 TSpline.cxx:1151
 TSpline.cxx:1152
 TSpline.cxx:1153
 TSpline.cxx:1154
 TSpline.cxx:1155
 TSpline.cxx:1156
 TSpline.cxx:1157
 TSpline.cxx:1158
 TSpline.cxx:1159
 TSpline.cxx:1160
 TSpline.cxx:1161
 TSpline.cxx:1162
 TSpline.cxx:1163
 TSpline.cxx:1164
 TSpline.cxx:1165
 TSpline.cxx:1166
 TSpline.cxx:1167
 TSpline.cxx:1168
 TSpline.cxx:1169
 TSpline.cxx:1170
 TSpline.cxx:1171
 TSpline.cxx:1172
 TSpline.cxx:1173
 TSpline.cxx:1174
 TSpline.cxx:1175
 TSpline.cxx:1176
 TSpline.cxx:1177
 TSpline.cxx:1178
 TSpline.cxx:1179
 TSpline.cxx:1180
 TSpline.cxx:1181
 TSpline.cxx:1182
 TSpline.cxx:1183
 TSpline.cxx:1184
 TSpline.cxx:1185
 TSpline.cxx:1186
 TSpline.cxx:1187
 TSpline.cxx:1188
 TSpline.cxx:1189
 TSpline.cxx:1190
 TSpline.cxx:1191
 TSpline.cxx:1192
 TSpline.cxx:1193
 TSpline.cxx:1194
 TSpline.cxx:1195
 TSpline.cxx:1196
 TSpline.cxx:1197
 TSpline.cxx:1198
 TSpline.cxx:1199
 TSpline.cxx:1200
 TSpline.cxx:1201
 TSpline.cxx:1202
 TSpline.cxx:1203
 TSpline.cxx:1204
 TSpline.cxx:1205
 TSpline.cxx:1206
 TSpline.cxx:1207
 TSpline.cxx:1208
 TSpline.cxx:1209
 TSpline.cxx:1210
 TSpline.cxx:1211
 TSpline.cxx:1212
 TSpline.cxx:1213
 TSpline.cxx:1214
 TSpline.cxx:1215
 TSpline.cxx:1216
 TSpline.cxx:1217
 TSpline.cxx:1218
 TSpline.cxx:1219
 TSpline.cxx:1220
 TSpline.cxx:1221
 TSpline.cxx:1222
 TSpline.cxx:1223
 TSpline.cxx:1224
 TSpline.cxx:1225
 TSpline.cxx:1226
 TSpline.cxx:1227
 TSpline.cxx:1228
 TSpline.cxx:1229
 TSpline.cxx:1230
 TSpline.cxx:1231
 TSpline.cxx:1232
 TSpline.cxx:1233
 TSpline.cxx:1234
 TSpline.cxx:1235
 TSpline.cxx:1236
 TSpline.cxx:1237
 TSpline.cxx:1238
 TSpline.cxx:1239
 TSpline.cxx:1240
 TSpline.cxx:1241
 TSpline.cxx:1242
 TSpline.cxx:1243
 TSpline.cxx:1244
 TSpline.cxx:1245
 TSpline.cxx:1246
 TSpline.cxx:1247
 TSpline.cxx:1248
 TSpline.cxx:1249
 TSpline.cxx:1250
 TSpline.cxx:1251
 TSpline.cxx:1252
 TSpline.cxx:1253
 TSpline.cxx:1254
 TSpline.cxx:1255
 TSpline.cxx:1256
 TSpline.cxx:1257
 TSpline.cxx:1258
 TSpline.cxx:1259
 TSpline.cxx:1260
 TSpline.cxx:1261
 TSpline.cxx:1262
 TSpline.cxx:1263
 TSpline.cxx:1264
 TSpline.cxx:1265
 TSpline.cxx:1266
 TSpline.cxx:1267
 TSpline.cxx:1268
 TSpline.cxx:1269
 TSpline.cxx:1270
 TSpline.cxx:1271
 TSpline.cxx:1272
 TSpline.cxx:1273
 TSpline.cxx:1274
 TSpline.cxx:1275
 TSpline.cxx:1276
 TSpline.cxx:1277
 TSpline.cxx:1278
 TSpline.cxx:1279
 TSpline.cxx:1280
 TSpline.cxx:1281
 TSpline.cxx:1282
 TSpline.cxx:1283
 TSpline.cxx:1284
 TSpline.cxx:1285
 TSpline.cxx:1286
 TSpline.cxx:1287
 TSpline.cxx:1288
 TSpline.cxx:1289
 TSpline.cxx:1290
 TSpline.cxx:1291
 TSpline.cxx:1292
 TSpline.cxx:1293
 TSpline.cxx:1294
 TSpline.cxx:1295
 TSpline.cxx:1296
 TSpline.cxx:1297
 TSpline.cxx:1298
 TSpline.cxx:1299
 TSpline.cxx:1300
 TSpline.cxx:1301
 TSpline.cxx:1302
 TSpline.cxx:1303
 TSpline.cxx:1304
 TSpline.cxx:1305
 TSpline.cxx:1306
 TSpline.cxx:1307
 TSpline.cxx:1308
 TSpline.cxx:1309
 TSpline.cxx:1310
 TSpline.cxx:1311
 TSpline.cxx:1312
 TSpline.cxx:1313
 TSpline.cxx:1314
 TSpline.cxx:1315
 TSpline.cxx:1316
 TSpline.cxx:1317
 TSpline.cxx:1318
 TSpline.cxx:1319
 TSpline.cxx:1320
 TSpline.cxx:1321
 TSpline.cxx:1322
 TSpline.cxx:1323
 TSpline.cxx:1324
 TSpline.cxx:1325
 TSpline.cxx:1326
 TSpline.cxx:1327
 TSpline.cxx:1328
 TSpline.cxx:1329
 TSpline.cxx:1330
 TSpline.cxx:1331
 TSpline.cxx:1332
 TSpline.cxx:1333
 TSpline.cxx:1334
 TSpline.cxx:1335
 TSpline.cxx:1336
 TSpline.cxx:1337
 TSpline.cxx:1338
 TSpline.cxx:1339
 TSpline.cxx:1340
 TSpline.cxx:1341
 TSpline.cxx:1342
 TSpline.cxx:1343
 TSpline.cxx:1344
 TSpline.cxx:1345
 TSpline.cxx:1346
 TSpline.cxx:1347
 TSpline.cxx:1348
 TSpline.cxx:1349
 TSpline.cxx:1350
 TSpline.cxx:1351
 TSpline.cxx:1352
 TSpline.cxx:1353
 TSpline.cxx:1354
 TSpline.cxx:1355
 TSpline.cxx:1356
 TSpline.cxx:1357
 TSpline.cxx:1358
 TSpline.cxx:1359
 TSpline.cxx:1360
 TSpline.cxx:1361
 TSpline.cxx:1362
 TSpline.cxx:1363
 TSpline.cxx:1364
 TSpline.cxx:1365
 TSpline.cxx:1366
 TSpline.cxx:1367
 TSpline.cxx:1368
 TSpline.cxx:1369
 TSpline.cxx:1370
 TSpline.cxx:1371
 TSpline.cxx:1372
 TSpline.cxx:1373
 TSpline.cxx:1374
 TSpline.cxx:1375
 TSpline.cxx:1376
 TSpline.cxx:1377
 TSpline.cxx:1378
 TSpline.cxx:1379
 TSpline.cxx:1380
 TSpline.cxx:1381
 TSpline.cxx:1382
 TSpline.cxx:1383
 TSpline.cxx:1384
 TSpline.cxx:1385
 TSpline.cxx:1386
 TSpline.cxx:1387
 TSpline.cxx:1388
 TSpline.cxx:1389
 TSpline.cxx:1390
 TSpline.cxx:1391
 TSpline.cxx:1392
 TSpline.cxx:1393
 TSpline.cxx:1394
 TSpline.cxx:1395
 TSpline.cxx:1396
 TSpline.cxx:1397
 TSpline.cxx:1398
 TSpline.cxx:1399
 TSpline.cxx:1400
 TSpline.cxx:1401
 TSpline.cxx:1402
 TSpline.cxx:1403
 TSpline.cxx:1404
 TSpline.cxx:1405
 TSpline.cxx:1406
 TSpline.cxx:1407
 TSpline.cxx:1408
 TSpline.cxx:1409
 TSpline.cxx:1410
 TSpline.cxx:1411
 TSpline.cxx:1412
 TSpline.cxx:1413
 TSpline.cxx:1414
 TSpline.cxx:1415
 TSpline.cxx:1416
 TSpline.cxx:1417
 TSpline.cxx:1418
 TSpline.cxx:1419
 TSpline.cxx:1420
 TSpline.cxx:1421
 TSpline.cxx:1422
 TSpline.cxx:1423
 TSpline.cxx:1424
 TSpline.cxx:1425
 TSpline.cxx:1426
 TSpline.cxx:1427
 TSpline.cxx:1428
 TSpline.cxx:1429
 TSpline.cxx:1430
 TSpline.cxx:1431
 TSpline.cxx:1432
 TSpline.cxx:1433
 TSpline.cxx:1434
 TSpline.cxx:1435
 TSpline.cxx:1436
 TSpline.cxx:1437
 TSpline.cxx:1438
 TSpline.cxx:1439
 TSpline.cxx:1440
 TSpline.cxx:1441
 TSpline.cxx:1442
 TSpline.cxx:1443
 TSpline.cxx:1444
 TSpline.cxx:1445
 TSpline.cxx:1446
 TSpline.cxx:1447
 TSpline.cxx:1448
 TSpline.cxx:1449
 TSpline.cxx:1450
 TSpline.cxx:1451
 TSpline.cxx:1452
 TSpline.cxx:1453
 TSpline.cxx:1454
 TSpline.cxx:1455
 TSpline.cxx:1456
 TSpline.cxx:1457
 TSpline.cxx:1458
 TSpline.cxx:1459
 TSpline.cxx:1460
 TSpline.cxx:1461
 TSpline.cxx:1462
 TSpline.cxx:1463
 TSpline.cxx:1464
 TSpline.cxx:1465
 TSpline.cxx:1466
 TSpline.cxx:1467
 TSpline.cxx:1468
 TSpline.cxx:1469
 TSpline.cxx:1470
 TSpline.cxx:1471
 TSpline.cxx:1472
 TSpline.cxx:1473
 TSpline.cxx:1474
 TSpline.cxx:1475
 TSpline.cxx:1476
 TSpline.cxx:1477
 TSpline.cxx:1478
 TSpline.cxx:1479
 TSpline.cxx:1480
 TSpline.cxx:1481
 TSpline.cxx:1482
 TSpline.cxx:1483
 TSpline.cxx:1484
 TSpline.cxx:1485
 TSpline.cxx:1486
 TSpline.cxx:1487
 TSpline.cxx:1488
 TSpline.cxx:1489
 TSpline.cxx:1490
 TSpline.cxx:1491
 TSpline.cxx:1492
 TSpline.cxx:1493
 TSpline.cxx:1494
 TSpline.cxx:1495
 TSpline.cxx:1496
 TSpline.cxx:1497
 TSpline.cxx:1498
 TSpline.cxx:1499
 TSpline.cxx:1500
 TSpline.cxx:1501
 TSpline.cxx:1502
 TSpline.cxx:1503
 TSpline.cxx:1504
 TSpline.cxx:1505
 TSpline.cxx:1506
 TSpline.cxx:1507
 TSpline.cxx:1508
 TSpline.cxx:1509
 TSpline.cxx:1510
 TSpline.cxx:1511
 TSpline.cxx:1512
 TSpline.cxx:1513
 TSpline.cxx:1514
 TSpline.cxx:1515
 TSpline.cxx:1516
 TSpline.cxx:1517
 TSpline.cxx:1518
 TSpline.cxx:1519
 TSpline.cxx:1520
 TSpline.cxx:1521
 TSpline.cxx:1522
 TSpline.cxx:1523
 TSpline.cxx:1524
 TSpline.cxx:1525
 TSpline.cxx:1526
 TSpline.cxx:1527
 TSpline.cxx:1528
 TSpline.cxx:1529
 TSpline.cxx:1530
 TSpline.cxx:1531
 TSpline.cxx:1532
 TSpline.cxx:1533
 TSpline.cxx:1534
 TSpline.cxx:1535
 TSpline.cxx:1536
 TSpline.cxx:1537
 TSpline.cxx:1538
 TSpline.cxx:1539
 TSpline.cxx:1540
 TSpline.cxx:1541
 TSpline.cxx:1542
 TSpline.cxx:1543
 TSpline.cxx:1544
 TSpline.cxx:1545
 TSpline.cxx:1546
 TSpline.cxx:1547
 TSpline.cxx:1548
 TSpline.cxx:1549
 TSpline.cxx:1550
 TSpline.cxx:1551
 TSpline.cxx:1552
 TSpline.cxx:1553
 TSpline.cxx:1554
 TSpline.cxx:1555
 TSpline.cxx:1556
 TSpline.cxx:1557
 TSpline.cxx:1558
 TSpline.cxx:1559
 TSpline.cxx:1560
 TSpline.cxx:1561
 TSpline.cxx:1562
 TSpline.cxx:1563
 TSpline.cxx:1564
 TSpline.cxx:1565
 TSpline.cxx:1566
 TSpline.cxx:1567
 TSpline.cxx:1568
 TSpline.cxx:1569
 TSpline.cxx:1570
 TSpline.cxx:1571
 TSpline.cxx:1572
 TSpline.cxx:1573
 TSpline.cxx:1574
 TSpline.cxx:1575
 TSpline.cxx:1576
 TSpline.cxx:1577
 TSpline.cxx:1578
 TSpline.cxx:1579
 TSpline.cxx:1580
 TSpline.cxx:1581
 TSpline.cxx:1582
 TSpline.cxx:1583
 TSpline.cxx:1584
 TSpline.cxx:1585
 TSpline.cxx:1586
 TSpline.cxx:1587
 TSpline.cxx:1588
 TSpline.cxx:1589
 TSpline.cxx:1590
 TSpline.cxx:1591
 TSpline.cxx:1592
 TSpline.cxx:1593
 TSpline.cxx:1594
 TSpline.cxx:1595
 TSpline.cxx:1596
 TSpline.cxx:1597
 TSpline.cxx:1598
 TSpline.cxx:1599
 TSpline.cxx:1600
 TSpline.cxx:1601
 TSpline.cxx:1602
 TSpline.cxx:1603
 TSpline.cxx:1604
 TSpline.cxx:1605
 TSpline.cxx:1606
 TSpline.cxx:1607
 TSpline.cxx:1608
 TSpline.cxx:1609
 TSpline.cxx:1610
 TSpline.cxx:1611
 TSpline.cxx:1612
 TSpline.cxx:1613
 TSpline.cxx:1614
 TSpline.cxx:1615
 TSpline.cxx:1616
 TSpline.cxx:1617
 TSpline.cxx:1618
 TSpline.cxx:1619
 TSpline.cxx:1620
 TSpline.cxx:1621
 TSpline.cxx:1622
 TSpline.cxx:1623
 TSpline.cxx:1624
 TSpline.cxx:1625
 TSpline.cxx:1626
 TSpline.cxx:1627
 TSpline.cxx:1628
 TSpline.cxx:1629
 TSpline.cxx:1630
 TSpline.cxx:1631
 TSpline.cxx:1632
 TSpline.cxx:1633
 TSpline.cxx:1634
 TSpline.cxx:1635
 TSpline.cxx:1636
 TSpline.cxx:1637
 TSpline.cxx:1638
 TSpline.cxx:1639
 TSpline.cxx:1640
 TSpline.cxx:1641
 TSpline.cxx:1642
 TSpline.cxx:1643
 TSpline.cxx:1644
 TSpline.cxx:1645
 TSpline.cxx:1646
 TSpline.cxx:1647
 TSpline.cxx:1648
 TSpline.cxx:1649
 TSpline.cxx:1650
 TSpline.cxx:1651
 TSpline.cxx:1652
 TSpline.cxx:1653
 TSpline.cxx:1654
 TSpline.cxx:1655
 TSpline.cxx:1656
 TSpline.cxx:1657
 TSpline.cxx:1658
 TSpline.cxx:1659
 TSpline.cxx:1660
 TSpline.cxx:1661
 TSpline.cxx:1662
 TSpline.cxx:1663
 TSpline.cxx:1664
 TSpline.cxx:1665
 TSpline.cxx:1666
 TSpline.cxx:1667
 TSpline.cxx:1668
 TSpline.cxx:1669
 TSpline.cxx:1670
 TSpline.cxx:1671
 TSpline.cxx:1672
 TSpline.cxx:1673
 TSpline.cxx:1674
 TSpline.cxx:1675
 TSpline.cxx:1676
 TSpline.cxx:1677
 TSpline.cxx:1678
 TSpline.cxx:1679
 TSpline.cxx:1680
 TSpline.cxx:1681
 TSpline.cxx:1682
 TSpline.cxx:1683
 TSpline.cxx:1684
 TSpline.cxx:1685
 TSpline.cxx:1686
 TSpline.cxx:1687
 TSpline.cxx:1688
 TSpline.cxx:1689
 TSpline.cxx:1690
 TSpline.cxx:1691
 TSpline.cxx:1692
 TSpline.cxx:1693
 TSpline.cxx:1694
 TSpline.cxx:1695
 TSpline.cxx:1696
 TSpline.cxx:1697
 TSpline.cxx:1698
 TSpline.cxx:1699
 TSpline.cxx:1700
 TSpline.cxx:1701
 TSpline.cxx:1702
 TSpline.cxx:1703
 TSpline.cxx:1704
 TSpline.cxx:1705
 TSpline.cxx:1706
 TSpline.cxx:1707
 TSpline.cxx:1708
 TSpline.cxx:1709
 TSpline.cxx:1710
 TSpline.cxx:1711
 TSpline.cxx:1712
 TSpline.cxx:1713
 TSpline.cxx:1714
 TSpline.cxx:1715
 TSpline.cxx:1716
 TSpline.cxx:1717
 TSpline.cxx:1718
 TSpline.cxx:1719
 TSpline.cxx:1720
 TSpline.cxx:1721
 TSpline.cxx:1722
 TSpline.cxx:1723
 TSpline.cxx:1724
 TSpline.cxx:1725
 TSpline.cxx:1726
 TSpline.cxx:1727
 TSpline.cxx:1728
 TSpline.cxx:1729
 TSpline.cxx:1730
 TSpline.cxx:1731
 TSpline.cxx:1732
 TSpline.cxx:1733
 TSpline.cxx:1734
 TSpline.cxx:1735
 TSpline.cxx:1736
 TSpline.cxx:1737
 TSpline.cxx:1738
 TSpline.cxx:1739
 TSpline.cxx:1740
 TSpline.cxx:1741
 TSpline.cxx:1742
 TSpline.cxx:1743
 TSpline.cxx:1744
 TSpline.cxx:1745
 TSpline.cxx:1746
 TSpline.cxx:1747
 TSpline.cxx:1748
 TSpline.cxx:1749
 TSpline.cxx:1750
 TSpline.cxx:1751
 TSpline.cxx:1752
 TSpline.cxx:1753
 TSpline.cxx:1754
 TSpline.cxx:1755
 TSpline.cxx:1756
 TSpline.cxx:1757
 TSpline.cxx:1758
 TSpline.cxx:1759
 TSpline.cxx:1760
 TSpline.cxx:1761
 TSpline.cxx:1762
 TSpline.cxx:1763
 TSpline.cxx:1764
 TSpline.cxx:1765
 TSpline.cxx:1766
 TSpline.cxx:1767
 TSpline.cxx:1768
 TSpline.cxx:1769
 TSpline.cxx:1770
 TSpline.cxx:1771
 TSpline.cxx:1772
 TSpline.cxx:1773
 TSpline.cxx:1774
 TSpline.cxx:1775
 TSpline.cxx:1776
 TSpline.cxx:1777
 TSpline.cxx:1778
 TSpline.cxx:1779
 TSpline.cxx:1780
 TSpline.cxx:1781
 TSpline.cxx:1782
 TSpline.cxx:1783
 TSpline.cxx:1784
 TSpline.cxx:1785
 TSpline.cxx:1786
 TSpline.cxx:1787
 TSpline.cxx:1788
 TSpline.cxx:1789
 TSpline.cxx:1790
 TSpline.cxx:1791
 TSpline.cxx:1792
 TSpline.cxx:1793
 TSpline.cxx:1794
 TSpline.cxx:1795
 TSpline.cxx:1796
 TSpline.cxx:1797
 TSpline.cxx:1798
 TSpline.cxx:1799
 TSpline.cxx:1800
 TSpline.cxx:1801
 TSpline.cxx:1802
 TSpline.cxx:1803
 TSpline.cxx:1804
 TSpline.cxx:1805
 TSpline.cxx:1806
 TSpline.cxx:1807
 TSpline.cxx:1808
 TSpline.cxx:1809
 TSpline.cxx:1810
 TSpline.cxx:1811
 TSpline.cxx:1812
 TSpline.cxx:1813
 TSpline.cxx:1814
 TSpline.cxx:1815
 TSpline.cxx:1816
 TSpline.cxx:1817
 TSpline.cxx:1818
 TSpline.cxx:1819
 TSpline.cxx:1820
 TSpline.cxx:1821
 TSpline.cxx:1822
 TSpline.cxx:1823
 TSpline.cxx:1824
 TSpline.cxx:1825
 TSpline.cxx:1826
 TSpline.cxx:1827
 TSpline.cxx:1828
 TSpline.cxx:1829
 TSpline.cxx:1830
 TSpline.cxx:1831
 TSpline.cxx:1832
 TSpline.cxx:1833
 TSpline.cxx:1834
 TSpline.cxx:1835
 TSpline.cxx:1836
 TSpline.cxx:1837
 TSpline.cxx:1838
 TSpline.cxx:1839
 TSpline.cxx:1840
 TSpline.cxx:1841
 TSpline.cxx:1842
 TSpline.cxx:1843
 TSpline.cxx:1844
 TSpline.cxx:1845
 TSpline.cxx:1846
 TSpline.cxx:1847
 TSpline.cxx:1848
 TSpline.cxx:1849
 TSpline.cxx:1850
 TSpline.cxx:1851
 TSpline.cxx:1852
 TSpline.cxx:1853
 TSpline.cxx:1854
 TSpline.cxx:1855
 TSpline.cxx:1856
 TSpline.cxx:1857
 TSpline.cxx:1858
 TSpline.cxx:1859
 TSpline.cxx:1860
 TSpline.cxx:1861
 TSpline.cxx:1862
 TSpline.cxx:1863
 TSpline.cxx:1864
 TSpline.cxx:1865
 TSpline.cxx:1866
 TSpline.cxx:1867
 TSpline.cxx:1868
 TSpline.cxx:1869
 TSpline.cxx:1870
 TSpline.cxx:1871
 TSpline.cxx:1872
 TSpline.cxx:1873
 TSpline.cxx:1874
 TSpline.cxx:1875
 TSpline.cxx:1876
 TSpline.cxx:1877
 TSpline.cxx:1878
 TSpline.cxx:1879
 TSpline.cxx:1880
 TSpline.cxx:1881
 TSpline.cxx:1882
 TSpline.cxx:1883
 TSpline.cxx:1884
 TSpline.cxx:1885
 TSpline.cxx:1886
 TSpline.cxx:1887
 TSpline.cxx:1888
 TSpline.cxx:1889
 TSpline.cxx:1890
 TSpline.cxx:1891
 TSpline.cxx:1892
 TSpline.cxx:1893
 TSpline.cxx:1894
 TSpline.cxx:1895
 TSpline.cxx:1896
 TSpline.cxx:1897
 TSpline.cxx:1898
 TSpline.cxx:1899
 TSpline.cxx:1900
 TSpline.cxx:1901
 TSpline.cxx:1902
 TSpline.cxx:1903
 TSpline.cxx:1904
 TSpline.cxx:1905
 TSpline.cxx:1906
 TSpline.cxx:1907
 TSpline.cxx:1908
 TSpline.cxx:1909
 TSpline.cxx:1910
 TSpline.cxx:1911
 TSpline.cxx:1912
 TSpline.cxx:1913
 TSpline.cxx:1914
 TSpline.cxx:1915
 TSpline.cxx:1916
 TSpline.cxx:1917
 TSpline.cxx:1918
 TSpline.cxx:1919
 TSpline.cxx:1920
 TSpline.cxx:1921
 TSpline.cxx:1922
 TSpline.cxx:1923
 TSpline.cxx:1924
 TSpline.cxx:1925
 TSpline.cxx:1926
 TSpline.cxx:1927
 TSpline.cxx:1928
 TSpline.cxx:1929
 TSpline.cxx:1930
 TSpline.cxx:1931
 TSpline.cxx:1932
 TSpline.cxx:1933
 TSpline.cxx:1934
 TSpline.cxx:1935
 TSpline.cxx:1936
 TSpline.cxx:1937
 TSpline.cxx:1938
 TSpline.cxx:1939
 TSpline.cxx:1940
 TSpline.cxx:1941
 TSpline.cxx:1942
 TSpline.cxx:1943
 TSpline.cxx:1944
 TSpline.cxx:1945
 TSpline.cxx:1946
 TSpline.cxx:1947
 TSpline.cxx:1948
 TSpline.cxx:1949
 TSpline.cxx:1950
 TSpline.cxx:1951
 TSpline.cxx:1952
 TSpline.cxx:1953
 TSpline.cxx:1954
 TSpline.cxx:1955
 TSpline.cxx:1956
 TSpline.cxx:1957
 TSpline.cxx:1958
 TSpline.cxx:1959
 TSpline.cxx:1960
 TSpline.cxx:1961
 TSpline.cxx:1962
 TSpline.cxx:1963
 TSpline.cxx:1964
 TSpline.cxx:1965
 TSpline.cxx:1966
 TSpline.cxx:1967
 TSpline.cxx:1968
 TSpline.cxx:1969
 TSpline.cxx:1970
 TSpline.cxx:1971
 TSpline.cxx:1972
 TSpline.cxx:1973
 TSpline.cxx:1974
 TSpline.cxx:1975
 TSpline.cxx:1976
 TSpline.cxx:1977
 TSpline.cxx:1978
 TSpline.cxx:1979
 TSpline.cxx:1980
 TSpline.cxx:1981
 TSpline.cxx:1982
 TSpline.cxx:1983
 TSpline.cxx:1984
 TSpline.cxx:1985
 TSpline.cxx:1986
 TSpline.cxx:1987
 TSpline.cxx:1988
 TSpline.cxx:1989
 TSpline.cxx:1990
 TSpline.cxx:1991
 TSpline.cxx:1992
 TSpline.cxx:1993
 TSpline.cxx:1994
 TSpline.cxx:1995
 TSpline.cxx:1996
 TSpline.cxx:1997
 TSpline.cxx:1998
 TSpline.cxx:1999
 TSpline.cxx:2000
 TSpline.cxx:2001
 TSpline.cxx:2002
 TSpline.cxx:2003
 TSpline.cxx:2004
 TSpline.cxx:2005
 TSpline.cxx:2006
 TSpline.cxx:2007
 TSpline.cxx:2008
 TSpline.cxx:2009
 TSpline.cxx:2010
 TSpline.cxx:2011
 TSpline.cxx:2012
 TSpline.cxx:2013
 TSpline.cxx:2014
 TSpline.cxx:2015
 TSpline.cxx:2016
 TSpline.cxx:2017
 TSpline.cxx:2018
 TSpline.cxx:2019
 TSpline.cxx:2020
 TSpline.cxx:2021
 TSpline.cxx:2022
 TSpline.cxx:2023
 TSpline.cxx:2024
 TSpline.cxx:2025
 TSpline.cxx:2026
 TSpline.cxx:2027
 TSpline.cxx:2028
 TSpline.cxx:2029
 TSpline.cxx:2030
 TSpline.cxx:2031
 TSpline.cxx:2032
 TSpline.cxx:2033
 TSpline.cxx:2034
 TSpline.cxx:2035
 TSpline.cxx:2036
 TSpline.cxx:2037
 TSpline.cxx:2038
 TSpline.cxx:2039
 TSpline.cxx:2040
 TSpline.cxx:2041
 TSpline.cxx:2042
 TSpline.cxx:2043
 TSpline.cxx:2044
 TSpline.cxx:2045
 TSpline.cxx:2046
 TSpline.cxx:2047
 TSpline.cxx:2048
 TSpline.cxx:2049
 TSpline.cxx:2050
 TSpline.cxx:2051
 TSpline.cxx:2052
 TSpline.cxx:2053
 TSpline.cxx:2054
 TSpline.cxx:2055
 TSpline.cxx:2056
 TSpline.cxx:2057
 TSpline.cxx:2058
 TSpline.cxx:2059
 TSpline.cxx:2060
 TSpline.cxx:2061
 TSpline.cxx:2062
 TSpline.cxx:2063
 TSpline.cxx:2064
 TSpline.cxx:2065
 TSpline.cxx:2066
 TSpline.cxx:2067
 TSpline.cxx:2068
 TSpline.cxx:2069
 TSpline.cxx:2070
 TSpline.cxx:2071
 TSpline.cxx:2072
 TSpline.cxx:2073
 TSpline.cxx:2074
 TSpline.cxx:2075
 TSpline.cxx:2076
 TSpline.cxx:2077
 TSpline.cxx:2078
 TSpline.cxx:2079
 TSpline.cxx:2080
 TSpline.cxx:2081
 TSpline.cxx:2082
 TSpline.cxx:2083
 TSpline.cxx:2084
 TSpline.cxx:2085
 TSpline.cxx:2086
 TSpline.cxx:2087
 TSpline.cxx:2088
 TSpline.cxx:2089
 TSpline.cxx:2090
 TSpline.cxx:2091
 TSpline.cxx:2092
 TSpline.cxx:2093
 TSpline.cxx:2094
 TSpline.cxx:2095
 TSpline.cxx:2096
 TSpline.cxx:2097
 TSpline.cxx:2098
 TSpline.cxx:2099
 TSpline.cxx:2100
 TSpline.cxx:2101
 TSpline.cxx:2102
 TSpline.cxx:2103
 TSpline.cxx:2104
 TSpline.cxx:2105
 TSpline.cxx:2106
 TSpline.cxx:2107
 TSpline.cxx:2108
 TSpline.cxx:2109
 TSpline.cxx:2110
 TSpline.cxx:2111
 TSpline.cxx:2112
 TSpline.cxx:2113
 TSpline.cxx:2114
 TSpline.cxx:2115
 TSpline.cxx:2116
 TSpline.cxx:2117
 TSpline.cxx:2118
 TSpline.cxx:2119
 TSpline.cxx:2120
 TSpline.cxx:2121
 TSpline.cxx:2122
 TSpline.cxx:2123
 TSpline.cxx:2124
 TSpline.cxx:2125
 TSpline.cxx:2126
 TSpline.cxx:2127
 TSpline.cxx:2128
 TSpline.cxx:2129
 TSpline.cxx:2130
 TSpline.cxx:2131
 TSpline.cxx:2132
 TSpline.cxx:2133
 TSpline.cxx:2134
 TSpline.cxx:2135
 TSpline.cxx:2136
 TSpline.cxx:2137
 TSpline.cxx:2138
 TSpline.cxx:2139
 TSpline.cxx:2140
 TSpline.cxx:2141
 TSpline.cxx:2142
 TSpline.cxx:2143
 TSpline.cxx:2144
 TSpline.cxx:2145
 TSpline.cxx:2146
 TSpline.cxx:2147
 TSpline.cxx:2148
 TSpline.cxx:2149
 TSpline.cxx:2150
 TSpline.cxx:2151
 TSpline.cxx:2152
 TSpline.cxx:2153
 TSpline.cxx:2154
 TSpline.cxx:2155
 TSpline.cxx:2156
 TSpline.cxx:2157
 TSpline.cxx:2158
 TSpline.cxx:2159
 TSpline.cxx:2160
 TSpline.cxx:2161
 TSpline.cxx:2162
 TSpline.cxx:2163
 TSpline.cxx:2164
 TSpline.cxx:2165
 TSpline.cxx:2166
 TSpline.cxx:2167
 TSpline.cxx:2168
 TSpline.cxx:2169
 TSpline.cxx:2170
 TSpline.cxx:2171
 TSpline.cxx:2172
 TSpline.cxx:2173
 TSpline.cxx:2174
 TSpline.cxx:2175
 TSpline.cxx:2176
 TSpline.cxx:2177
 TSpline.cxx:2178
 TSpline.cxx:2179
 TSpline.cxx:2180
 TSpline.cxx:2181
 TSpline.cxx:2182
 TSpline.cxx:2183
 TSpline.cxx:2184
 TSpline.cxx:2185
 TSpline.cxx:2186
 TSpline.cxx:2187
 TSpline.cxx:2188
 TSpline.cxx:2189
 TSpline.cxx:2190
 TSpline.cxx:2191
 TSpline.cxx:2192
 TSpline.cxx:2193
 TSpline.cxx:2194
 TSpline.cxx:2195
 TSpline.cxx:2196
 TSpline.cxx:2197
 TSpline.cxx:2198
 TSpline.cxx:2199
 TSpline.cxx:2200
 TSpline.cxx:2201
 TSpline.cxx:2202
 TSpline.cxx:2203
 TSpline.cxx:2204
 TSpline.cxx:2205
 TSpline.cxx:2206
 TSpline.cxx:2207
 TSpline.cxx:2208
 TSpline.cxx:2209
 TSpline.cxx:2210
 TSpline.cxx:2211
 TSpline.cxx:2212
 TSpline.cxx:2213
 TSpline.cxx:2214
 TSpline.cxx:2215
 TSpline.cxx:2216
 TSpline.cxx:2217
 TSpline.cxx:2218
 TSpline.cxx:2219
 TSpline.cxx:2220
 TSpline.cxx:2221
 TSpline.cxx:2222
 TSpline.cxx:2223
 TSpline.cxx:2224
 TSpline.cxx:2225
 TSpline.cxx:2226
 TSpline.cxx:2227
 TSpline.cxx:2228
 TSpline.cxx:2229
 TSpline.cxx:2230
 TSpline.cxx:2231
 TSpline.cxx:2232
 TSpline.cxx:2233
 TSpline.cxx:2234
 TSpline.cxx:2235
 TSpline.cxx:2236
 TSpline.cxx:2237
 TSpline.cxx:2238
 TSpline.cxx:2239
 TSpline.cxx:2240
 TSpline.cxx:2241
 TSpline.cxx:2242
 TSpline.cxx:2243
 TSpline.cxx:2244
 TSpline.cxx:2245
 TSpline.cxx:2246
 TSpline.cxx:2247
 TSpline.cxx:2248
 TSpline.cxx:2249
 TSpline.cxx:2250
 TSpline.cxx:2251
 TSpline.cxx:2252
 TSpline.cxx:2253
 TSpline.cxx:2254
 TSpline.cxx:2255
 TSpline.cxx:2256
 TSpline.cxx:2257
 TSpline.cxx:2258
 TSpline.cxx:2259
 TSpline.cxx:2260
 TSpline.cxx:2261
 TSpline.cxx:2262
 TSpline.cxx:2263
 TSpline.cxx:2264
 TSpline.cxx:2265
 TSpline.cxx:2266
 TSpline.cxx:2267
 TSpline.cxx:2268
 TSpline.cxx:2269
 TSpline.cxx:2270
 TSpline.cxx:2271
 TSpline.cxx:2272
 TSpline.cxx:2273
 TSpline.cxx:2274
 TSpline.cxx:2275
 TSpline.cxx:2276
 TSpline.cxx:2277
 TSpline.cxx:2278
 TSpline.cxx:2279
 TSpline.cxx:2280
 TSpline.cxx:2281
 TSpline.cxx:2282
 TSpline.cxx:2283
 TSpline.cxx:2284
 TSpline.cxx:2285
 TSpline.cxx:2286
 TSpline.cxx:2287
 TSpline.cxx:2288
 TSpline.cxx:2289
 TSpline.cxx:2290
 TSpline.cxx:2291
 TSpline.cxx:2292
 TSpline.cxx:2293
 TSpline.cxx:2294
 TSpline.cxx:2295
 TSpline.cxx:2296
 TSpline.cxx:2297
 TSpline.cxx:2298
 TSpline.cxx:2299
 TSpline.cxx:2300
 TSpline.cxx:2301
 TSpline.cxx:2302
 TSpline.cxx:2303
 TSpline.cxx:2304
 TSpline.cxx:2305
 TSpline.cxx:2306
 TSpline.cxx:2307
 TSpline.cxx:2308
 TSpline.cxx:2309
 TSpline.cxx:2310
 TSpline.cxx:2311
 TSpline.cxx:2312
 TSpline.cxx:2313
 TSpline.cxx:2314
 TSpline.cxx:2315
 TSpline.cxx:2316
 TSpline.cxx:2317
 TSpline.cxx:2318
 TSpline.cxx:2319
 TSpline.cxx:2320
 TSpline.cxx:2321
 TSpline.cxx:2322
 TSpline.cxx:2323
 TSpline.cxx:2324
 TSpline.cxx:2325
 TSpline.cxx:2326
 TSpline.cxx:2327
 TSpline.cxx:2328
 TSpline.cxx:2329
 TSpline.cxx:2330
 TSpline.cxx:2331
 TSpline.cxx:2332
 TSpline.cxx:2333
 TSpline.cxx:2334
 TSpline.cxx:2335
 TSpline.cxx:2336
 TSpline.cxx:2337
 TSpline.cxx:2338
 TSpline.cxx:2339
 TSpline.cxx:2340
 TSpline.cxx:2341
 TSpline.cxx:2342
 TSpline.cxx:2343
 TSpline.cxx:2344
 TSpline.cxx:2345
 TSpline.cxx:2346
 TSpline.cxx:2347
 TSpline.cxx:2348
 TSpline.cxx:2349
 TSpline.cxx:2350
 TSpline.cxx:2351
 TSpline.cxx:2352
 TSpline.cxx:2353
 TSpline.cxx:2354
 TSpline.cxx:2355
 TSpline.cxx:2356
 TSpline.cxx:2357
 TSpline.cxx:2358
 TSpline.cxx:2359
 TSpline.cxx:2360
 TSpline.cxx:2361
 TSpline.cxx:2362
 TSpline.cxx:2363
 TSpline.cxx:2364
 TSpline.cxx:2365
 TSpline.cxx:2366
 TSpline.cxx:2367
 TSpline.cxx:2368
 TSpline.cxx:2369
 TSpline.cxx:2370
 TSpline.cxx:2371
 TSpline.cxx:2372
 TSpline.cxx:2373
 TSpline.cxx:2374
 TSpline.cxx:2375
 TSpline.cxx:2376
 TSpline.cxx:2377
 TSpline.cxx:2378
 TSpline.cxx:2379
 TSpline.cxx:2380
 TSpline.cxx:2381
 TSpline.cxx:2382
 TSpline.cxx:2383
 TSpline.cxx:2384
 TSpline.cxx:2385
 TSpline.cxx:2386
 TSpline.cxx:2387
 TSpline.cxx:2388
 TSpline.cxx:2389
 TSpline.cxx:2390
 TSpline.cxx:2391
 TSpline.cxx:2392
 TSpline.cxx:2393
 TSpline.cxx:2394
 TSpline.cxx:2395
 TSpline.cxx:2396
 TSpline.cxx:2397
 TSpline.cxx:2398
 TSpline.cxx:2399
 TSpline.cxx:2400
 TSpline.cxx:2401
 TSpline.cxx:2402
 TSpline.cxx:2403
 TSpline.cxx:2404
 TSpline.cxx:2405
 TSpline.cxx:2406
 TSpline.cxx:2407
 TSpline.cxx:2408
 TSpline.cxx:2409
 TSpline.cxx:2410
 TSpline.cxx:2411
 TSpline.cxx:2412
 TSpline.cxx:2413
 TSpline.cxx:2414
 TSpline.cxx:2415
 TSpline.cxx:2416
 TSpline.cxx:2417
 TSpline.cxx:2418
 TSpline.cxx:2419
 TSpline.cxx:2420
 TSpline.cxx:2421
 TSpline.cxx:2422
 TSpline.cxx:2423
 TSpline.cxx:2424
 TSpline.cxx:2425
 TSpline.cxx:2426
 TSpline.cxx:2427
 TSpline.cxx:2428
 TSpline.cxx:2429
 TSpline.cxx:2430
 TSpline.cxx:2431
 TSpline.cxx:2432
 TSpline.cxx:2433
 TSpline.cxx:2434
 TSpline.cxx:2435
 TSpline.cxx:2436
 TSpline.cxx:2437
 TSpline.cxx:2438
 TSpline.cxx:2439
 TSpline.cxx:2440
 TSpline.cxx:2441
 TSpline.cxx:2442
 TSpline.cxx:2443
 TSpline.cxx:2444
 TSpline.cxx:2445
 TSpline.cxx:2446
 TSpline.cxx:2447
 TSpline.cxx:2448
 TSpline.cxx:2449
 TSpline.cxx:2450
 TSpline.cxx:2451
 TSpline.cxx:2452
 TSpline.cxx:2453
 TSpline.cxx:2454
 TSpline.cxx:2455
 TSpline.cxx:2456
 TSpline.cxx:2457
 TSpline.cxx:2458
 TSpline.cxx:2459
 TSpline.cxx:2460
 TSpline.cxx:2461
 TSpline.cxx:2462
 TSpline.cxx:2463
 TSpline.cxx:2464
 TSpline.cxx:2465
 TSpline.cxx:2466
 TSpline.cxx:2467
 TSpline.cxx:2468
 TSpline.cxx:2469
 TSpline.cxx:2470
 TSpline.cxx:2471
 TSpline.cxx:2472
 TSpline.cxx:2473
 TSpline.cxx:2474
 TSpline.cxx:2475
 TSpline.cxx:2476
 TSpline.cxx:2477
 TSpline.cxx:2478
 TSpline.cxx:2479
 TSpline.cxx:2480
 TSpline.cxx:2481
 TSpline.cxx:2482
 TSpline.cxx:2483
 TSpline.cxx:2484
 TSpline.cxx:2485
 TSpline.cxx:2486
 TSpline.cxx:2487
 TSpline.cxx:2488
 TSpline.cxx:2489
 TSpline.cxx:2490
 TSpline.cxx:2491
 TSpline.cxx:2492
 TSpline.cxx:2493
 TSpline.cxx:2494
 TSpline.cxx:2495
 TSpline.cxx:2496
 TSpline.cxx:2497
 TSpline.cxx:2498
 TSpline.cxx:2499
 TSpline.cxx:2500
 TSpline.cxx:2501
 TSpline.cxx:2502
 TSpline.cxx:2503
 TSpline.cxx:2504
 TSpline.cxx:2505
 TSpline.cxx:2506
 TSpline.cxx:2507
 TSpline.cxx:2508
 TSpline.cxx:2509
 TSpline.cxx:2510
 TSpline.cxx:2511
 TSpline.cxx:2512
 TSpline.cxx:2513
 TSpline.cxx:2514
 TSpline.cxx:2515
 TSpline.cxx:2516
 TSpline.cxx:2517
 TSpline.cxx:2518
 TSpline.cxx:2519
 TSpline.cxx:2520
 TSpline.cxx:2521
 TSpline.cxx:2522
 TSpline.cxx:2523
 TSpline.cxx:2524
 TSpline.cxx:2525
 TSpline.cxx:2526
 TSpline.cxx:2527
 TSpline.cxx:2528
 TSpline.cxx:2529
 TSpline.cxx:2530
 TSpline.cxx:2531
 TSpline.cxx:2532
 TSpline.cxx:2533
 TSpline.cxx:2534
 TSpline.cxx:2535
 TSpline.cxx:2536
 TSpline.cxx:2537
 TSpline.cxx:2538
 TSpline.cxx:2539
 TSpline.cxx:2540
 TSpline.cxx:2541
 TSpline.cxx:2542
 TSpline.cxx:2543
 TSpline.cxx:2544
 TSpline.cxx:2545