TSpline.cxx

Go to the documentation of this file.
00001 // @(#)root/hist:$Id: TSpline.cxx 36247 2010-10-10 10:16:42Z brun $
00002 // Author: Federico Carminati   28/02/2000
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 //////////////////////////////////////////////////////////////////////////
00013 //                                                                      //
00014 // TSpline                                                              //
00015 //                                                                      //
00016 // Base class for spline implementation containing the Draw/Paint       //
00017 // methods                                                              //
00018 //                                                                      //
00019 //////////////////////////////////////////////////////////////////////////
00020 
00021 #include "TROOT.h"
00022 #include "TSpline.h"
00023 #include "TVirtualPad.h"
00024 #include "TH1.h"
00025 #include "TF1.h"
00026 #include "TSystem.h"
00027 #include "Riostream.h"
00028 #include "TClass.h"
00029 #include "TMath.h"
00030 
00031 ClassImp(TSplinePoly)
00032 ClassImp(TSplinePoly3)
00033 ClassImp(TSplinePoly5)
00034 ClassImp(TSpline3)
00035 ClassImp(TSpline5)
00036 ClassImp(TSpline)
00037 
00038 
00039 //______________________________________________________________________________
00040 TSpline::TSpline(const TSpline &sp) :
00041   TNamed(sp),
00042   TAttLine(sp),
00043   TAttFill(sp),
00044   TAttMarker(sp),
00045   fDelta(sp.fDelta),
00046   fXmin(sp.fXmin),
00047   fXmax(sp.fXmax),
00048   fNp(sp.fNp),
00049   fKstep(sp.fKstep),
00050   fHistogram(0),
00051   fGraph(0),
00052   fNpx(sp.fNpx)
00053 {
00054    //copy constructor
00055 }
00056 
00057 
00058 //______________________________________________________________________________
00059 TSpline::~TSpline()
00060 {
00061    //destructor
00062    if(fHistogram) delete fHistogram;
00063    if(fGraph) delete fGraph;
00064 }
00065 
00066 
00067 //______________________________________________________________________________
00068 TSpline& TSpline::operator=(const TSpline &sp)
00069 {
00070    //assignment operator
00071    if(this!=&sp) {
00072       TNamed::operator=(sp);
00073       TAttLine::operator=(sp);
00074       TAttFill::operator=(sp);
00075       TAttMarker::operator=(sp);
00076       fDelta=sp.fDelta;
00077       fXmin=sp.fXmin;
00078       fXmax=sp.fXmax;
00079       fNp=sp.fNp;
00080       fKstep=sp.fKstep;
00081       fHistogram=0;
00082       fGraph=0;
00083       fNpx=sp.fNpx;
00084    }
00085    return *this;
00086 }
00087 
00088 
00089 //______________________________________________________________________________
00090 void TSpline::Draw(Option_t *option)
00091 {
00092    // Draw this function with its current attributes.
00093    //
00094    // Possible option values are:
00095    //   "SAME"  superimpose on top of existing picture
00096    //   "L"     connect all computed points with a straight line
00097    //   "C"     connect all computed points with a smooth curve.
00098    //   "P"     add a polymarker at each knot
00099    //
00100    // Note that the default value is "L". Therefore to draw on top
00101    // of an existing picture, specify option "LSAME"
00102 
00103    TString opt = option;
00104    opt.ToLower();
00105    if (gPad && !opt.Contains("same")) gPad->Clear();
00106 
00107    AppendPad(option);
00108 }
00109 
00110 
00111 //______________________________________________________________________________
00112 Int_t TSpline::DistancetoPrimitive(Int_t px, Int_t py)
00113 {
00114    // Compute distance from point px,py to a spline.
00115 
00116    if (!fHistogram) return 999;
00117    return fHistogram->DistancetoPrimitive(px, py);
00118 }
00119 
00120 
00121 //______________________________________________________________________________
00122 void TSpline::ExecuteEvent(Int_t event, Int_t px, Int_t py)
00123 {
00124    // Execute action corresponding to one event.
00125 
00126    if (!fHistogram) return;
00127    fHistogram->ExecuteEvent(event, px, py);
00128 }
00129 
00130 
00131 //______________________________________________________________________________
00132 void TSpline::Paint(Option_t *option)
00133 {
00134    // Paint this function with its current attributes.
00135 
00136    Int_t i;
00137    Double_t xv;
00138 
00139    TString opt = option;
00140    opt.ToLower();
00141    Double_t xmin, xmax, pmin, pmax;
00142    pmin = gPad->PadtoX(gPad->GetUxmin());
00143    pmax = gPad->PadtoX(gPad->GetUxmax());
00144    xmin = fXmin;
00145    xmax = fXmax;
00146    if (opt.Contains("same")) {
00147       if (xmax < pmin) return;  // Otto: completely outside
00148       if (xmin > pmax) return;
00149       if (xmin < pmin) xmin = pmin;
00150       if (xmax > pmax) xmax = pmax;
00151    } else {
00152       gPad->Clear();
00153    }
00154 
00155    //  Create a temporary histogram and fill each channel with the function value
00156    if (fHistogram)
00157       if ((!gPad->GetLogx() && fHistogram->TestBit(TH1::kLogX)) ||
00158           (gPad->GetLogx() && !fHistogram->TestBit(TH1::kLogX)))
00159          { delete fHistogram; fHistogram = 0;}
00160 
00161    if (fHistogram) {
00162       //if (xmin != fXmin || xmax != fXmax)
00163       fHistogram->GetXaxis()->SetLimits(xmin,xmax);
00164    } else {
00165    //      if logx, we must bin in logx and not in x !!!
00166    //      otherwise if several decades, one gets crazy results
00167       if (xmin > 0 && gPad->GetLogx()) {
00168          Double_t *xbins  = new Double_t[fNpx+1];
00169          Double_t xlogmin = TMath::Log10(xmin);
00170          Double_t xlogmax = TMath::Log10(xmax);
00171          Double_t dlogx   = (xlogmax-xlogmin)/((Double_t)fNpx);
00172          for (i=0;i<=fNpx;i++) {
00173             xbins[i] = gPad->PadtoX(xlogmin+ i*dlogx);
00174          }
00175          fHistogram = new TH1F("Spline",GetTitle(),fNpx,xbins);
00176          fHistogram->SetBit(TH1::kLogX);
00177          delete [] xbins;
00178       } else {
00179          fHistogram = new TH1F("Spline",GetTitle(),fNpx,xmin,xmax);
00180       }
00181       if (!fHistogram) return;
00182       fHistogram->SetDirectory(0);
00183    }
00184    for (i=1;i<=fNpx;i++) {
00185       xv = fHistogram->GetBinCenter(i);
00186       fHistogram->SetBinContent(i,this->Eval(xv));
00187    }
00188 
00189    // Copy Function attributes to histogram attributes
00190    fHistogram->SetBit(TH1::kNoStats);
00191    fHistogram->SetLineColor(GetLineColor());
00192    fHistogram->SetLineStyle(GetLineStyle());
00193    fHistogram->SetLineWidth(GetLineWidth());
00194    fHistogram->SetFillColor(GetFillColor());
00195    fHistogram->SetFillStyle(GetFillStyle());
00196    fHistogram->SetMarkerColor(GetMarkerColor());
00197    fHistogram->SetMarkerStyle(GetMarkerStyle());
00198    fHistogram->SetMarkerSize(GetMarkerSize());
00199 
00200    //  Draw the histogram
00201    //  but first strip off the 'p' option if any
00202    char *o = (char *) opt.Data();
00203    Int_t j=0;
00204    i=0;
00205    Bool_t graph=kFALSE;
00206    do
00207       if(o[i]=='p') graph=kTRUE ; else o[j++]=o[i];
00208    while(o[i++]);
00209    if (opt.Length() == 0 ) fHistogram->Paint("lf");
00210    else if (opt == "same") fHistogram->Paint("lfsame");
00211    else                    fHistogram->Paint(opt.Data());
00212 
00213    // Think about the graph, if demanded
00214    if(graph) {
00215       if(!fGraph) {
00216          Double_t *xx = new Double_t[fNp];
00217          Double_t *yy = new Double_t[fNp];
00218          for(i=0; i<fNp; ++i)
00219             GetKnot(i,xx[i],yy[i]);
00220          fGraph=new TGraph(fNp,xx,yy);
00221          delete [] xx;
00222          delete [] yy;
00223       }
00224       fGraph->SetMarkerColor(GetMarkerColor());
00225       fGraph->SetMarkerStyle(GetMarkerStyle());
00226       fGraph->SetMarkerSize(GetMarkerSize());
00227       fGraph->Paint("p");
00228    }
00229 }
00230 
00231 
00232 //______________________________________________________________________________
00233 void TSpline::Streamer(TBuffer &R__b)
00234 {
00235    // Stream an object of class TSpline.
00236 
00237    if (R__b.IsReading()) {
00238       UInt_t R__s, R__c;
00239       Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
00240       if (R__v > 1) {
00241          R__b.ReadClassBuffer(TSpline::Class(), this, R__v, R__s, R__c);
00242          return;
00243       }
00244       //====process old versions before automatic schema evolution
00245       TNamed::Streamer(R__b);
00246       TAttLine::Streamer(R__b);
00247       TAttFill::Streamer(R__b);
00248       TAttMarker::Streamer(R__b);
00249 
00250       fNp = 0;
00251       /*
00252       R__b >> fDelta;
00253       R__b >> fXmin;
00254       R__b >> fXmax;
00255       R__b >> fNp;
00256       R__b >> fKstep;
00257       R__b >> fHistogram;
00258       R__b >> fGraph;
00259       R__b >> fNpx;
00260       */
00261       R__b.CheckByteCount(R__s, R__c, TSpline::IsA());
00262       //====end of old versions
00263 
00264    } else {
00265       R__b.WriteClassBuffer(TSpline::Class(),this);
00266    }
00267 }
00268 
00269 
00270 //////////////////////////////////////////////////////////////////////////
00271 //                                                                      //
00272 // TSplinePoly                                                          //
00273 //                                                                      //
00274 // Base class for TSpline knot                                          //
00275 //                                                                      //
00276 //////////////////////////////////////////////////////////////////////////
00277 
00278 
00279 //______________________________________________________________________________
00280 TSplinePoly &TSplinePoly::operator=(TSplinePoly const &other)
00281 {
00282    //assignment operator
00283    if(this != &other) {
00284       TObject::operator=(other);
00285       CopyPoly(other);
00286    }
00287    return *this;
00288 }
00289 
00290 //______________________________________________________________________________
00291 void TSplinePoly::CopyPoly(TSplinePoly const &other)
00292 {
00293    //utility called by the copy constructors and = operator
00294    fX = other.fX;
00295    fY = other.fY;
00296 }
00297 
00298 //////////////////////////////////////////////////////////////////////////
00299 //                                                                      //
00300 // TSplinePoly3                                                         //
00301 //                                                                      //
00302 // Class for TSpline3 knot                                              //
00303 //                                                                      //
00304 //////////////////////////////////////////////////////////////////////////
00305 
00306 
00307 //______________________________________________________________________________
00308 TSplinePoly3 &TSplinePoly3::operator=(TSplinePoly3 const &other)
00309 {
00310    //assignment operator
00311    if(this != &other) {
00312       TSplinePoly::operator=(other);
00313       CopyPoly(other);
00314    }
00315    return *this;
00316 }
00317 
00318 //______________________________________________________________________________
00319 void TSplinePoly3::CopyPoly(TSplinePoly3 const &other)
00320 {
00321    //utility called by the copy constructors and = operator
00322    fB = other.fB;
00323    fC = other.fC;
00324    fD = other.fD;
00325 }
00326 
00327 //////////////////////////////////////////////////////////////////////////
00328 //                                                                      //
00329 // TSplinePoly5                                                         //
00330 //                                                                      //
00331 // Class for TSpline5 knot                                              //
00332 //                                                                      //
00333 //////////////////////////////////////////////////////////////////////////
00334 
00335 
00336 //______________________________________________________________________________
00337 TSplinePoly5 &TSplinePoly5::operator=(TSplinePoly5 const &other)
00338 {
00339    //assignment operator
00340    if(this != &other) {
00341       TSplinePoly::operator=(other);
00342       CopyPoly(other);
00343    }
00344    return *this;
00345 }
00346 
00347 //______________________________________________________________________________
00348 void TSplinePoly5::CopyPoly(TSplinePoly5 const &other)
00349 {
00350    //utility called by the copy constructors and = operator
00351    fB = other.fB;
00352    fC = other.fC;
00353    fD = other.fD;
00354    fE = other.fE;
00355    fF = other.fF;
00356 }
00357 
00358 
00359 
00360 //////////////////////////////////////////////////////////////////////////
00361 //                                                                      //
00362 // TSpline3                                                             //
00363 //                                                                      //
00364 // Class to create third splines to interpolate knots                   //
00365 // Arbitrary conditions can be introduced for first and second          //
00366 // derivatives at beginning and ending points                           //
00367 //                                                                      //
00368 //////////////////////////////////////////////////////////////////////////
00369 
00370 
00371 //______________________________________________________________________________
00372 TSpline3::TSpline3(const char *title,
00373                    Double_t x[], Double_t y[], Int_t n, const char *opt,
00374                    Double_t valbeg, Double_t valend) :
00375   TSpline(title,-1,x[0],x[n-1],n,kFALSE),
00376   fValBeg(valbeg), fValEnd(valend), fBegCond(0), fEndCond(0)
00377 {
00378    // Third spline creator given an array of
00379    // arbitrary knots in increasing abscissa order and
00380    // possibly end point conditions
00381 
00382    fName="Spline3";
00383 
00384    // Set endpoint conditions
00385    if(opt) SetCond(opt);
00386 
00387    // Create the plynomial terms and fill
00388    // them with node information
00389    fPoly = new TSplinePoly3[n];
00390    for (Int_t i=0; i<n; ++i) {
00391       fPoly[i].X() = x[i];
00392       fPoly[i].Y() = y[i];
00393    }
00394 
00395    // Build the spline coefficients
00396    BuildCoeff();
00397 }
00398 
00399 //______________________________________________________________________________
00400 TSpline3::TSpline3(const char *title,
00401                    Double_t xmin, Double_t xmax,
00402                    Double_t y[], Int_t n, const char *opt,
00403                    Double_t valbeg, Double_t valend) :
00404   TSpline(title,(xmax-xmin)/(n-1), xmin, xmax, n, kTRUE),
00405   fValBeg(valbeg), fValEnd(valend),
00406   fBegCond(0), fEndCond(0)
00407 {
00408    // Third spline creator given an array of
00409    // arbitrary function values on equidistant n abscissa
00410    // values from xmin to xmax and possibly end point conditions
00411 
00412    fName="Spline3";
00413 
00414    // Set endpoint conditions
00415    if(opt) SetCond(opt);
00416 
00417    // Create the plynomial terms and fill
00418    // them with node information
00419    fPoly = new TSplinePoly3[n];
00420    for (Int_t i=0; i<n; ++i) {
00421       fPoly[i].X() = fXmin+i*fDelta;
00422       fPoly[i].Y() = y[i];
00423    }
00424 
00425    // Build the spline coefficients
00426    BuildCoeff();
00427 }
00428 
00429 
00430 //______________________________________________________________________________
00431 TSpline3::TSpline3(const char *title,
00432                    Double_t x[], const TF1 *func, Int_t n, const char *opt,
00433                    Double_t valbeg, Double_t valend) :
00434   TSpline(title,-1, x[0], x[n-1], n, kFALSE),
00435   fValBeg(valbeg), fValEnd(valend),
00436   fBegCond(0), fEndCond(0)
00437 {
00438    // Third spline creator given an array of
00439    // arbitrary abscissas in increasing order and a function
00440    // to interpolate and possibly end point conditions
00441 
00442    fName="Spline3";
00443 
00444    // Set endpoint conditions
00445    if(opt) SetCond(opt);
00446 
00447    // Create the plynomial terms and fill
00448    // them with node information
00449    fPoly = new TSplinePoly3[n];
00450    for (Int_t i=0; i<n; ++i) {
00451       fPoly[i].X() = x[i];
00452       fPoly[i].Y() = ((TF1*)func)->Eval(x[i]);
00453    }
00454 
00455    // Build the spline coefficients
00456    BuildCoeff();
00457 }
00458 
00459 
00460 //______________________________________________________________________________
00461 TSpline3::TSpline3(const char *title,
00462                    Double_t xmin, Double_t xmax,
00463                    const TF1 *func, Int_t n, const char *opt,
00464                    Double_t valbeg, Double_t valend) :
00465   TSpline(title,(xmax-xmin)/(n-1), xmin, xmax, n, kTRUE),
00466   fValBeg(valbeg), fValEnd(valend),
00467   fBegCond(0), fEndCond(0)
00468 {
00469    // Third spline creator given a function to be
00470    // evaluated on n equidistand abscissa points between xmin
00471    // and xmax and possibly end point conditions
00472 
00473    fName="Spline3";
00474 
00475    // Set endpoint conditions
00476    if(opt) SetCond(opt);
00477 
00478    // Create the plynomial terms and fill
00479    // them with node information
00480    fPoly = new TSplinePoly3[n];
00481    //when func is null we return. In this case it is assumed that the spline
00482    //points will be given later via SetPoint and SetPointCoeff
00483    if (!func) {fKstep = kFALSE; fDelta = -1; return;}
00484    for (Int_t i=0; i<n; ++i) {
00485       Double_t x=fXmin+i*fDelta;
00486       fPoly[i].X() = x;
00487       fPoly[i].Y() = ((TF1*)func)->Eval(x);
00488    }
00489 
00490    // Build the spline coefficients
00491    BuildCoeff();
00492 }
00493 
00494 
00495 //______________________________________________________________________________
00496 TSpline3::TSpline3(const char *title,
00497                    const TGraph *g, const char *opt,
00498                    Double_t valbeg, Double_t valend) :
00499   TSpline(title,-1,0,0,g->GetN(),kFALSE),
00500   fValBeg(valbeg), fValEnd(valend),
00501   fBegCond(0), fEndCond(0)
00502 {
00503    // Third spline creator given a TGraph with
00504    // abscissa in increasing order and possibly end
00505    // point conditions
00506 
00507    fName="Spline3";
00508 
00509    // Set endpoint conditions
00510    if(opt) SetCond(opt);
00511 
00512    // Create the plynomial terms and fill
00513    // them with node information
00514    fPoly = new TSplinePoly3[fNp];
00515    for (Int_t i=0; i<fNp; ++i) {
00516       Double_t xx, yy;
00517       g->GetPoint(i,xx,yy);
00518       fPoly[i].X()=xx;
00519       fPoly[i].Y()=yy;
00520    }
00521    fXmin = fPoly[0].X();
00522    fXmax = fPoly[fNp-1].X();
00523 
00524    // Build the spline coefficients
00525    BuildCoeff();
00526 }
00527 
00528 
00529 //______________________________________________________________________________
00530 TSpline3::TSpline3(const TH1 *h, const char *opt,
00531                    Double_t valbeg, Double_t valend) :
00532   TSpline(h->GetTitle(),-1,0,0,h->GetNbinsX(),kFALSE),
00533   fValBeg(valbeg), fValEnd(valend),
00534   fBegCond(0), fEndCond(0)
00535 {
00536    // Third spline creator given a TH1 
00537 
00538    fName=h->GetName();
00539 
00540    // Set endpoint conditions
00541    if(opt) SetCond(opt);
00542 
00543    // Create the plynomial terms and fill
00544    // them with node information
00545    fPoly = new TSplinePoly3[fNp];
00546    for (Int_t i=0; i<fNp; ++i) {
00547       fPoly[i].X()=h->GetXaxis()->GetBinCenter(i+1);
00548       fPoly[i].Y()=h->GetBinContent(i+1);
00549    }
00550    fXmin = fPoly[0].X();
00551    fXmax = fPoly[fNp-1].X();
00552 
00553    // Build the spline coefficients
00554    BuildCoeff();
00555 }
00556 
00557 
00558 //______________________________________________________________________________
00559 TSpline3::TSpline3(const TSpline3& sp3) :
00560   TSpline(sp3),
00561   fPoly(0),
00562   fValBeg(sp3.fValBeg),
00563   fValEnd(sp3.fValEnd),
00564   fBegCond(sp3.fBegCond),
00565   fEndCond(sp3.fEndCond)
00566 {
00567    //copy constructor
00568    if (fNp > 0) fPoly = new TSplinePoly3[fNp];
00569    for (Int_t i=0; i<fNp; ++i) 
00570       fPoly[i] = sp3.fPoly[i];
00571 }
00572 
00573 
00574 //______________________________________________________________________________
00575 TSpline3& TSpline3::operator=(const TSpline3& sp3)
00576 {
00577    //assignment operator
00578    if(this!=&sp3) {
00579       TSpline::operator=(sp3);
00580       fPoly= 0;
00581       if (fNp > 0) fPoly = new TSplinePoly3[fNp];
00582       for (Int_t i=0; i<fNp; ++i) 
00583          fPoly[i] = sp3.fPoly[i];
00584       
00585       fValBeg=sp3.fValBeg;
00586       fValEnd=sp3.fValEnd;
00587       fBegCond=sp3.fBegCond;
00588       fEndCond=sp3.fEndCond;
00589    }
00590    return *this;
00591 }
00592 
00593 
00594 //______________________________________________________________________________
00595 void TSpline3::SetCond(const char *opt)
00596 {
00597    // Check the boundary conditions
00598 
00599    const char *b1 = strstr(opt,"b1");
00600    const char *e1 = strstr(opt,"e1");
00601    const char *b2 = strstr(opt,"b2");
00602    const char *e2 = strstr(opt,"e2");
00603    if (b1 && b2)
00604       Error("SetCond","Cannot specify first and second derivative at first point");
00605    if (e1 && e2)
00606       Error("SetCond","Cannot specify first and second derivative at last point");
00607    if (b1) fBegCond=1;
00608    else if (b2) fBegCond=2;
00609    if (e1) fEndCond=1;
00610    else if (e2) fEndCond=2;
00611 }
00612 
00613 
00614 //______________________________________________________________________________
00615 void TSpline3::Test()
00616 {
00617    // Test method for TSpline5
00618    //
00619    //   n          number of data points.
00620    //   m          2*m-1 is order of spline.
00621    //                 m = 2 always for third spline.
00622    //   nn,nm1,mm,
00623    //   mm1,i,k,
00624    //   j,jj       temporary integer variables.
00625    //   z,p        temporary double precision variables.
00626    //   x[n]       the sequence of knots.
00627    //   y[n]       the prescribed function values at the knots.
00628    //   a[200][4]  two dimensional array whose columns are
00629    //                 the computed spline coefficients
00630    //   diff[3]    maximum values of differences of values and
00631    //                 derivatives to right and left of knots.
00632    //   com[3]     maximum values of coefficients.
00633    //
00634    //
00635    //   test of TSpline3 with nonequidistant knots and
00636    //      equidistant knots follows.
00637 
00638    Double_t hx;
00639    Double_t diff[3];
00640    Double_t a[800], c[4];
00641    Int_t i, j, k, m, n;
00642    Double_t x[200], y[200], z;
00643    Int_t jj, mm;
00644    Int_t mm1, nm1;
00645    Double_t com[3];
00646    printf("1         TEST OF TSpline3 WITH NONEQUIDISTANT KNOTS\n");
00647    n = 5;
00648    x[0] = -3;
00649    x[1] = -1;
00650    x[2] = 0;
00651    x[3] = 3;
00652    x[4] = 4;
00653    y[0] = 7;
00654    y[1] = 11;
00655    y[2] = 26;
00656    y[3] = 56;
00657    y[4] = 29;
00658    m = 2;
00659    mm = m << 1;
00660    mm1 = mm-1;
00661    printf("\n-N = %3d    M =%2d\n",n,m);
00662    TSpline3 *spline = new TSpline3("Test",x,y,n);
00663    for (i = 0; i < n; ++i)
00664       spline->GetCoeff(i,hx, a[i],a[i+200],a[i+400],a[i+600]);
00665    delete spline;
00666    for (i = 0; i < mm1; ++i) diff[i] = com[i] = 0;
00667    for (k = 0; k < n; ++k) {
00668       for (i = 0; i < mm; ++i) c[i] = a[k+i*200];
00669       printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
00670       printf("%12.8f\n",x[k]);
00671       if (k == n-1) {
00672          printf("%16.8f\n",c[0]);
00673       } else {
00674          for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
00675          printf("\n");
00676          for (i = 0; i < mm1; ++i)
00677             if ((z=TMath::Abs(a[k+i*200])) > com[i]) com[i] = z;
00678          z = x[k+1]-x[k];
00679          for (i = 1; i < mm; ++i)
00680             for (jj = i; jj < mm; ++jj) {
00681                j = mm+i-jj;
00682                c[j-2] = c[j-1]*z+c[j-2];
00683             }
00684          for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
00685          printf("\n");
00686          for (i = 0; i < mm1; ++i)
00687             if (!(k >= n-2 && i != 0))
00688                if((z = TMath::Abs(c[i]-a[k+1+i*200]))
00689                   > diff[i]) diff[i] = z;
00690       }
00691    }
00692    printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
00693    for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
00694    printf("\n");
00695    printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
00696    if (TMath::Abs(c[0]) > com[0])
00697       com[0] = TMath::Abs(c[0]);
00698    for (i = 0; i < mm1; ++i) printf("%16.8f",com[i]);
00699    printf("\n");
00700    m = 2;
00701    for (n = 10; n <= 100; n += 10) {
00702       mm = m << 1;
00703       mm1 = mm-1;
00704       nm1 = n-1;
00705       for (i = 0; i < nm1; i += 2) {
00706          x[i] = i+1;
00707          x[i+1] = i+2;
00708          y[i] = 1;
00709          y[i+1] = 0;
00710       }
00711       if (n % 2 != 0) {
00712          x[n-1] = n;
00713          y[n-1] = 1;
00714       }
00715       printf("\n-N = %3d    M =%2d\n",n,m);
00716       spline = new TSpline3("Test",x,y,n);
00717       for (i = 0; i < n; ++i)
00718          spline->GetCoeff(i,hx,a[i],a[i+200],a[i+400],a[i+600]);
00719       delete spline;
00720       for (i = 0; i < mm1; ++i)
00721          diff[i] = com[i] = 0;
00722       for (k = 0; k < n; ++k) {
00723          for (i = 0; i < mm; ++i)
00724             c[i] = a[k+i*200];
00725          if (n < 11) {
00726             printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
00727             printf("%12.8f\n",x[k]);
00728             if (k == n-1) printf("%16.8f\n",c[0]);
00729          }
00730          if (k == n-1) break;
00731          if (n <= 10) {
00732             for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
00733             printf("\n");
00734          }
00735          for (i = 0; i < mm1; ++i)
00736          if ((z=TMath::Abs(a[k+i*200])) > com[i])
00737             com[i] = z;
00738          z = x[k+1]-x[k];
00739          for (i = 1; i < mm; ++i)
00740             for (jj = i; jj < mm; ++jj) {
00741                j = mm+i-jj;
00742                c[j-2] = c[j-1]*z+c[j-2];
00743             }
00744          if (n <= 10) {
00745             for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
00746             printf("\n");
00747          }
00748          for (i = 0; i < mm1; ++i)
00749          if (!(k >= n-2 && i != 0))
00750             if ((z = TMath::Abs(c[i]-a[k+1+i*200]))
00751                > diff[i]) diff[i] = z;
00752       }
00753       printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
00754       for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
00755       printf("\n");
00756       printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
00757       if (TMath::Abs(c[0]) > com[0])
00758          com[0] = TMath::Abs(c[0]);
00759       for (i = 0; i < mm1; ++i) printf("%16.8E",com[i]);
00760          printf("\n");
00761    }
00762 }
00763 
00764 
00765 //______________________________________________________________________________
00766 Int_t TSpline3::FindX(Double_t x) const
00767 {
00768    // Find X
00769 
00770    Int_t klow=0, khig=fNp-1;
00771    //
00772    // If out of boundaries, extrapolate
00773    // It may be badly wrong
00774    if(x<=fXmin) klow=0;
00775    else if(x>=fXmax) klow=khig;
00776    else {
00777       if(fKstep) {
00778          //
00779          // Equidistant knots, use histogramming
00780          klow = TMath::FloorNint((x-fXmin)/fDelta);
00781          // Correction for rounding errors
00782          if (x < fPoly[klow].X())
00783             klow = TMath::Max(klow-1,0);
00784          else if (klow < khig) {
00785             if (x > fPoly[klow+1].X()) ++klow;
00786          }
00787       } else {
00788          Int_t khalf;
00789          //
00790          // Non equidistant knots, binary search
00791          while(khig-klow>1)
00792             if(x>fPoly[khalf=(klow+khig)/2].X())
00793                klow=khalf;
00794             else
00795                khig=khalf;
00796          //
00797          // This could be removed, sanity check
00798          if(!(fPoly[klow].X()<=x && x<=fPoly[klow+1].X()))
00799             Error("Eval",
00800                   "Binary search failed x(%d) = %f < x= %f < x(%d) = %f\n",
00801                   klow,fPoly[klow].X(),x,klow+1,fPoly[klow+1].X());
00802       }
00803    }
00804    return klow;
00805 }
00806 
00807 
00808 //______________________________________________________________________________
00809 Double_t TSpline3::Eval(Double_t x) const
00810 {
00811    // Eval this spline at x
00812 
00813    Int_t klow=FindX(x);
00814    if (klow >= fNp-1) klow = fNp-2; //see: https://savannah.cern.ch/bugs/?71651
00815    return fPoly[klow].Eval(x);
00816 }
00817 
00818 
00819 //______________________________________________________________________________
00820 Double_t TSpline3::Derivative(Double_t x) const
00821 {
00822    // Derivative
00823 
00824    Int_t klow=FindX(x);
00825    if (klow >= fNp-1) klow = fNp-2; //see: https://savannah.cern.ch/bugs/?71651
00826    return fPoly[klow].Derivative(x);
00827 }
00828 
00829 
00830 //______________________________________________________________________________
00831 void TSpline3::SaveAs(const char *filename, Option_t * /*option*/) const
00832 {
00833    // write this spline as a C++ function that can be executed without ROOT
00834    // the name of the function is the name of the file up to the "." if any
00835 
00836    //open the file
00837    ofstream *f = new ofstream(filename,ios::out);
00838    if (f == 0 || gSystem->AccessPathName(filename,kWritePermission)) {
00839       Error("SaveAs","Cannot open file:%s\n",filename);
00840       return;
00841    }
00842 
00843    //write the function name and the spline constants
00844    char buffer[512];
00845    Int_t nch = strlen(filename);
00846    snprintf(buffer,512,"double %s",filename);
00847    char *dot = strstr(buffer,".");
00848    if (dot) *dot = 0;
00849    strlcat(buffer,"(double x) {\n",512);
00850    nch = strlen(buffer); f->write(buffer,nch);
00851    snprintf(buffer,512,"   const int fNp = %d, fKstep = %d;\n",fNp,fKstep);
00852    nch = strlen(buffer); f->write(buffer,nch);
00853    snprintf(buffer,512,"   const double fDelta = %g, fXmin = %g, fXmax = %g;\n",fDelta,fXmin,fXmax);
00854    nch = strlen(buffer); f->write(buffer,nch);
00855 
00856    //write the spline coefficients
00857    //array fX
00858    snprintf(buffer,512,"   const double fX[%d] = {",fNp);
00859    nch = strlen(buffer); f->write(buffer,nch);
00860    buffer[0] = 0;
00861    Int_t i;
00862    char numb[20];
00863    for (i=0;i<fNp;i++) {
00864       snprintf(numb,20," %g,",fPoly[i].X());
00865       nch = strlen(numb);
00866       if (i == fNp-1) numb[nch-1]=0;
00867       strlcat(buffer,numb,512);
00868       if (i%5 == 4 || i == fNp-1) {
00869          nch = strlen(buffer); f->write(buffer,nch);
00870          if (i != fNp-1) snprintf(buffer,512,"\n                       ");
00871       }
00872    }
00873    snprintf(buffer,512," };\n");
00874    nch = strlen(buffer); f->write(buffer,nch);
00875    //array fY
00876    snprintf(buffer,512,"   const double fY[%d] = {",fNp);
00877    nch = strlen(buffer); f->write(buffer,nch);
00878    buffer[0] = 0;
00879    for (i=0;i<fNp;i++) {
00880       snprintf(numb,20," %g,",fPoly[i].Y());
00881       nch = strlen(numb);
00882       if (i == fNp-1) numb[nch-1]=0;
00883       strlcat(buffer,numb,512);
00884       if (i%5 == 4 || i == fNp-1) {
00885          nch = strlen(buffer); f->write(buffer,nch);
00886          if (i != fNp-1) snprintf(buffer,512,"\n                       ");
00887       }
00888    }
00889    snprintf(buffer,512," };\n");
00890    nch = strlen(buffer); f->write(buffer,nch);
00891    //array fB
00892    snprintf(buffer,512,"   const double fB[%d] = {",fNp);
00893    nch = strlen(buffer); f->write(buffer,nch);
00894    buffer[0] = 0;
00895    for (i=0;i<fNp;i++) {
00896       snprintf(numb,20," %g,",fPoly[i].B());
00897       nch = strlen(numb);
00898       if (i == fNp-1) numb[nch-1]=0;
00899       strlcat(buffer,numb,512);
00900       if (i%5 == 4 || i == fNp-1) {
00901          nch = strlen(buffer); f->write(buffer,nch);
00902          if (i != fNp-1) snprintf(buffer,512,"\n                       ");
00903       }
00904    }
00905    snprintf(buffer,512," };\n");
00906    nch = strlen(buffer); f->write(buffer,nch);
00907    //array fC
00908    snprintf(buffer,512,"   const double fC[%d] = {",fNp);
00909    nch = strlen(buffer); f->write(buffer,nch);
00910    buffer[0] = 0;
00911    for (i=0;i<fNp;i++) {
00912       snprintf(numb,20," %g,",fPoly[i].C());
00913       nch = strlen(numb);
00914       if (i == fNp-1) numb[nch-1]=0;
00915       strlcat(buffer,numb,512);
00916       if (i%5 == 4 || i == fNp-1) {
00917          nch = strlen(buffer); f->write(buffer,nch);
00918          if (i != fNp-1) snprintf(buffer,512,"\n                       ");
00919       }
00920    }
00921    snprintf(buffer,512," };\n");
00922    nch = strlen(buffer); f->write(buffer,nch);
00923     //array fD
00924    snprintf(buffer,512,"   const double fD[%d] = {",fNp);
00925    nch = strlen(buffer); f->write(buffer,nch);
00926    buffer[0] = 0;
00927    for (i=0;i<fNp;i++) {
00928       snprintf(numb,20," %g,",fPoly[i].D());
00929       nch = strlen(numb);
00930       if (i == fNp-1) numb[nch-1]=0;
00931       strlcat(buffer,numb,512);
00932       if (i%5 == 4 || i == fNp-1) {
00933          nch = strlen(buffer); f->write(buffer,nch);
00934          if (i != fNp-1) snprintf(buffer,512,"\n                       ");
00935       }
00936    }
00937    snprintf(buffer,512," };\n");
00938    nch = strlen(buffer); f->write(buffer,nch);
00939 
00940    //generate code for the spline evaluation
00941    snprintf(buffer,512,"   int klow=0;\n");
00942    nch = strlen(buffer); f->write(buffer,nch);
00943 
00944    snprintf(buffer,512,"   // If out of boundaries, extrapolate. It may be badly wrong\n");
00945    snprintf(buffer,512,"   if(x<=fXmin) klow=0;\n");
00946    nch = strlen(buffer); f->write(buffer,nch);
00947    snprintf(buffer,512,"   else if(x>=fXmax) klow=fNp-1;\n");
00948    nch = strlen(buffer); f->write(buffer,nch);
00949    snprintf(buffer,512,"   else {\n");
00950    nch = strlen(buffer); f->write(buffer,nch);
00951    snprintf(buffer,512,"     if(fKstep) {\n");
00952    nch = strlen(buffer); f->write(buffer,nch);
00953 
00954    snprintf(buffer,512,"       // Equidistant knots, use histogramming\n");
00955    nch = strlen(buffer); f->write(buffer,nch);
00956    snprintf(buffer,512,"       klow = int((x-fXmin)/fDelta);\n");
00957    nch = strlen(buffer); f->write(buffer,nch);
00958    snprintf(buffer,512,"       if (klow < fNp-1) klow = fNp-1;\n");
00959    nch = strlen(buffer); f->write(buffer,nch);
00960    snprintf(buffer,512,"     } else {\n");
00961    nch = strlen(buffer); f->write(buffer,nch);
00962    snprintf(buffer,512,"       int khig=fNp-1, khalf;\n");
00963    nch = strlen(buffer); f->write(buffer,nch);
00964 
00965    snprintf(buffer,512,"       // Non equidistant knots, binary search\n");
00966    nch = strlen(buffer); f->write(buffer,nch);
00967    snprintf(buffer,512,"       while(khig-klow>1)\n");
00968    nch = strlen(buffer); f->write(buffer,nch);
00969    snprintf(buffer,512,"         if(x>fX[khalf=(klow+khig)/2]) klow=khalf;\n");
00970    nch = strlen(buffer); f->write(buffer,nch);
00971    snprintf(buffer,512,"         else khig=khalf;\n");
00972    nch = strlen(buffer); f->write(buffer,nch);
00973    snprintf(buffer,512,"     }\n");
00974    nch = strlen(buffer); f->write(buffer,nch);
00975    snprintf(buffer,512,"   }\n");
00976    nch = strlen(buffer); f->write(buffer,nch);
00977    snprintf(buffer,512,"   // Evaluate now\n");
00978    nch = strlen(buffer); f->write(buffer,nch);
00979    snprintf(buffer,512,"   double dx=x-fX[klow];\n");
00980    nch = strlen(buffer); f->write(buffer,nch);
00981    snprintf(buffer,512,"   return (fY[klow]+dx*(fB[klow]+dx*(fC[klow]+dx*fD[klow])));\n");
00982    nch = strlen(buffer); f->write(buffer,nch);
00983 
00984    //close file
00985    f->write("}\n",2);
00986 
00987    if (f) { f->close(); delete f;}
00988 }
00989 
00990 
00991 //______________________________________________________________________________
00992 void TSpline3::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
00993 {
00994     // Save primitive as a C++ statement(s) on output stream out
00995 
00996    char quote = '"';
00997    out<<"   "<<endl;
00998    if (gROOT->ClassSaved(TSpline3::Class())) {
00999       out<<"   ";
01000    } else {
01001       out<<"   TSpline3 *";
01002    }
01003    out<<"spline3 = new TSpline3("<<quote<<GetTitle()<<quote<<","
01004       <<fXmin<<","<<fXmax<<",(TF1*)0,"<<fNp<<","<<quote<<quote<<","
01005       <<fValBeg<<","<<fValEnd<<");"<<endl;
01006    out<<"   spline3->SetName("<<quote<<GetName()<<quote<<");"<<endl;
01007 
01008    SaveFillAttributes(out,"spline3",0,1001);
01009    SaveLineAttributes(out,"spline3",1,1,1);
01010    SaveMarkerAttributes(out,"spline3",1,1,1);
01011    if (fNpx != 100) out<<"   spline3->SetNpx("<<fNpx<<");"<<endl;
01012 
01013    for (Int_t i=0;i<fNp;i++) {
01014       out<<"   spline3->SetPoint("<<i<<","<<fPoly[i].X()<<","<<fPoly[i].Y()<<");"<<endl;
01015       out<<"   spline3->SetPointCoeff("<<i<<","<<fPoly[i].B()<<","<<fPoly[i].C()<<","<<fPoly[i].D()<<");"<<endl;
01016    }
01017    out<<"   spline3->Draw("<<quote<<option<<quote<<");"<<endl;
01018 }
01019 
01020 
01021 //______________________________________________________________________________
01022 void TSpline3::SetPoint(Int_t i, Double_t x, Double_t y)
01023 {
01024    //set point number i.
01025    
01026    if (i < 0 || i >= fNp) return;
01027    fPoly[i].X()= x;
01028    fPoly[i].Y()= y;
01029 }
01030 
01031 //______________________________________________________________________________
01032 void TSpline3::SetPointCoeff(Int_t i, Double_t b, Double_t c, Double_t d)
01033 {
01034    // set point coefficient number i
01035  
01036    if (i < 0 || i >= fNp) return;
01037    fPoly[i].B()= b;
01038    fPoly[i].C()= c;
01039    fPoly[i].D()= d;
01040 }
01041 
01042 //______________________________________________________________________________
01043 void TSpline3::BuildCoeff()
01044 {
01045    //      subroutine cubspl ( tau, c, n, ibcbeg, ibcend )
01046    //  from  * a practical guide to splines *  by c. de boor
01047    //     ************************  input  ***************************
01048    //     n = number of data points. assumed to be .ge. 2.
01049    //     (tau(i), c(1,i), i=1,...,n) = abscissae and ordinates of the
01050    //        data points. tau is assumed to be strictly increasing.
01051    //     ibcbeg, ibcend = boundary condition indicators, and
01052    //     c(2,1), c(2,n) = boundary condition information. specifically,
01053    //        ibcbeg = 0  means no boundary condition at tau(1) is given.
01054    //           in this case, the not-a-knot condition is used, i.e. the
01055    //           jump in the third derivative across tau(2) is forced to
01056    //           zero, thus the first and the second cubic polynomial pieces
01057    //           are made to coincide.)
01058    //        ibcbeg = 1  means that the slope at tau(1) is made to equal
01059    //           c(2,1), supplied by input.
01060    //        ibcbeg = 2  means that the second derivative at tau(1) is
01061    //           made to equal c(2,1), supplied by input.
01062    //        ibcend = 0, 1, or 2 has analogous meaning concerning the
01063    //           boundary condition at tau(n), with the additional infor-
01064    //           mation taken from c(2,n).
01065    //     ***********************  output  **************************
01066    //     c(j,i), j=1,...,4; i=1,...,l (= n-1) = the polynomial coefficients
01067    //        of the cubic interpolating spline with interior knots (or
01068    //        joints) tau(2), ..., tau(n-1). precisely, in the interval
01069    //        (tau(i), tau(i+1)), the spline f is given by
01070    //           f(x) = c(1,i)+h*(c(2,i)+h*(c(3,i)+h*c(4,i)/3.)/2.)
01071    //        where h = x - tau(i). the function program *ppvalu* may be
01072    //        used to evaluate f or its derivatives from tau,c, l = n-1,
01073    //        and k=4.
01074 
01075    Int_t i, j, l, m;
01076    Double_t   divdf1,divdf3,dtau,g=0;
01077    //***** a tridiagonal linear system for the unknown slopes s(i) of
01078    //  f  at tau(i), i=1,...,n, is generated and then solved by gauss elim-
01079    //  ination, with s(i) ending up in c(2,i), all i.
01080    //     c(3,.) and c(4,.) are used initially for temporary storage.
01081    l = fNp-1;
01082    // compute first differences of x sequence and store in C also,
01083    // compute first divided difference of data and store in D.
01084    for (m=1; m<fNp ; ++m) {
01085       fPoly[m].C() = fPoly[m].X() - fPoly[m-1].X();
01086       fPoly[m].D() = (fPoly[m].Y() - fPoly[m-1].Y())/fPoly[m].C();
01087    }
01088    // construct first equation from the boundary condition, of the form
01089    //             D[0]*s[0] + C[0]*s[1] = B[0]
01090    if(fBegCond==0) {
01091       if(fNp == 2) {
01092          //     no condition at left end and n = 2.
01093          fPoly[0].D() = 1.;
01094          fPoly[0].C() = 1.;
01095          fPoly[0].B() = 2.*fPoly[1].D();
01096       } else {
01097          //     not-a-knot condition at left end and n .gt. 2.
01098          fPoly[0].D() = fPoly[2].C();
01099          fPoly[0].C() = fPoly[1].C() + fPoly[2].C();
01100          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();
01101       }
01102    } else if (fBegCond==1) {
01103       //     slope prescribed at left end.
01104       fPoly[0].B() = fValBeg;
01105       fPoly[0].D() = 1.;
01106       fPoly[0].C() = 0.;
01107    } else if (fBegCond==2) {
01108       //     second derivative prescribed at left end.
01109       fPoly[0].D() = 2.;
01110       fPoly[0].C() = 1.;
01111       fPoly[0].B() = 3.*fPoly[1].D() - fPoly[1].C()/2.*fValBeg;
01112    }
01113    if(fNp > 2) {
01114       //  if there are interior knots, generate the corresp. equations and car-
01115       //  ry out the forward pass of gauss elimination, after which the m-th
01116       //  equation reads    D[m]*s[m] + C[m]*s[m+1] = B[m].
01117       for (m=1; m<l; ++m) {
01118          g = -fPoly[m+1].C()/fPoly[m-1].D();
01119          fPoly[m].B() = g*fPoly[m-1].B() + 3.*(fPoly[m].C()*fPoly[m+1].D()+fPoly[m+1].C()*fPoly[m].D());
01120          fPoly[m].D() = g*fPoly[m-1].C() + 2.*(fPoly[m].C() + fPoly[m+1].C());
01121       }
01122       // construct last equation from the second boundary condition, of the form
01123       //           (-g*D[n-2])*s[n-2] + D[n-1]*s[n-1] = B[n-1]
01124       //     if slope is prescribed at right end, one can go directly to back-
01125       //     substitution, since c array happens to be set up just right for it
01126       //     at this point.
01127       if(fEndCond == 0) {
01128          if (fNp > 3 || fBegCond != 0) {
01129             //     not-a-knot and n .ge. 3, and either n.gt.3 or  also not-a-knot at
01130             //     left end point.
01131             g = fPoly[fNp-2].C() + fPoly[fNp-1].C();
01132             fPoly[fNp-1].B() = ((fPoly[fNp-1].C()+2.*g)*fPoly[fNp-1].D()*fPoly[fNp-2].C()
01133                          + fPoly[fNp-1].C()*fPoly[fNp-1].C()*(fPoly[fNp-2].Y()-fPoly[fNp-3].Y())/fPoly[fNp-2].C())/g;
01134             g = -g/fPoly[fNp-2].D();
01135             fPoly[fNp-1].D() = fPoly[fNp-2].C();
01136          } else {
01137             //     either (n=3 and not-a-knot also at left) or (n=2 and not not-a-
01138             //     knot at left end point).
01139             fPoly[fNp-1].B() = 2.*fPoly[fNp-1].D();
01140             fPoly[fNp-1].D() = 1.;
01141             g = -1./fPoly[fNp-2].D();
01142          }
01143       } else if (fEndCond == 1) {
01144          fPoly[fNp-1].B() = fValEnd;
01145          goto L30;
01146       } else if (fEndCond == 2) {
01147          //     second derivative prescribed at right endpoint.
01148          fPoly[fNp-1].B() = 3.*fPoly[fNp-1].D() + fPoly[fNp-1].C()/2.*fValEnd;
01149          fPoly[fNp-1].D() = 2.;
01150          g = -1./fPoly[fNp-2].D();
01151       }
01152    } else {
01153       if(fEndCond == 0) {
01154          if (fBegCond > 0) {
01155             //     either (n=3 and not-a-knot also at left) or (n=2 and not not-a-
01156             //     knot at left end point).
01157             fPoly[fNp-1].B() = 2.*fPoly[fNp-1].D();
01158             fPoly[fNp-1].D() = 1.;
01159             g = -1./fPoly[fNp-2].D();
01160          } else {
01161             //     not-a-knot at right endpoint and at left endpoint and n = 2.
01162             fPoly[fNp-1].B() = fPoly[fNp-1].D();
01163             goto L30;
01164          }
01165       } else if(fEndCond == 1) {
01166          fPoly[fNp-1].B() = fValEnd;
01167          goto L30;
01168       } else if(fEndCond == 2) {
01169          //     second derivative prescribed at right endpoint.
01170          fPoly[fNp-1].B() = 3.*fPoly[fNp-1].D() + fPoly[fNp-1].C()/2.*fValEnd;
01171          fPoly[fNp-1].D() = 2.;
01172          g = -1./fPoly[fNp-2].D();
01173       }
01174    }
01175    // complete forward pass of gauss elimination.
01176    fPoly[fNp-1].D() = g*fPoly[fNp-2].C() + fPoly[fNp-1].D();
01177    fPoly[fNp-1].B() = (g*fPoly[fNp-2].B() + fPoly[fNp-1].B())/fPoly[fNp-1].D();
01178    // carry out back substitution
01179 L30: j = l-1;
01180    do {
01181       fPoly[j].B() = (fPoly[j].B() - fPoly[j].C()*fPoly[j+1].B())/fPoly[j].D();
01182       --j;
01183    }  while (j>=0);
01184    //****** generate cubic coefficients in each interval, i.e., the deriv.s
01185    //  at its left endpoint, from value and slope at its endpoints.
01186    for (i=1; i<fNp; ++i) {
01187       dtau = fPoly[i].C();
01188       divdf1 = (fPoly[i].Y() - fPoly[i-1].Y())/dtau;
01189       divdf3 = fPoly[i-1].B() + fPoly[i].B() - 2.*divdf1;
01190       fPoly[i-1].C() = (divdf1 - fPoly[i-1].B() - divdf3)/dtau;
01191       fPoly[i-1].D() = (divdf3/dtau)/dtau;
01192    }
01193 }
01194 
01195 
01196 //______________________________________________________________________________
01197 void TSpline3::Streamer(TBuffer &R__b)
01198 {
01199    // Stream an object of class TSpline3.
01200 
01201    if (R__b.IsReading()) {
01202       UInt_t R__s, R__c;
01203       Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
01204       if (R__v > 1) {
01205          R__b.ReadClassBuffer(TSpline3::Class(), this, R__v, R__s, R__c);
01206          return;
01207       }
01208       //====process old versions before automatic schema evolution
01209       TSpline::Streamer(R__b);
01210       if (fNp > 0) {
01211          fPoly = new TSplinePoly3[fNp];
01212          for(Int_t i=0; i<fNp; ++i) {
01213             fPoly[i].Streamer(R__b);
01214          }
01215       }
01216       //      R__b >> fPoly;
01217       R__b >> fValBeg;
01218       R__b >> fValEnd;
01219       R__b >> fBegCond;
01220       R__b >> fEndCond;
01221    } else {
01222       R__b.WriteClassBuffer(TSpline3::Class(),this);
01223    }
01224 }
01225 
01226 
01227 //////////////////////////////////////////////////////////////////////////
01228 //                                                                      //
01229 // TSpline5                                                             //
01230 //                                                                      //
01231 // Class to create quintic natural splines to interpolate knots         //
01232 // Arbitrary conditions can be introduced for first and second          //
01233 // derivatives using double knots (see BuildCoeff) for more on this.    //
01234 // Double knots are automatically introduced at ending points           //
01235 //                                                                      //
01236 //////////////////////////////////////////////////////////////////////////
01237 
01238 
01239 //______________________________________________________________________________
01240 TSpline5::TSpline5(const char *title,
01241                    Double_t x[], Double_t y[], Int_t n,
01242                    const char *opt, Double_t b1, Double_t e1,
01243                    Double_t b2, Double_t e2) :
01244   TSpline(title,-1, x[0], x[n-1], n, kFALSE)
01245 {
01246    // Quintic natural spline creator given an array of
01247    // arbitrary knots in increasing abscissa order and
01248    // possibly end point conditions
01249 
01250    Int_t beg, end;
01251    const char *cb1, *ce1, *cb2, *ce2;
01252    fName="Spline5";
01253 
01254    // Check endpoint conditions
01255    BoundaryConditions(opt,beg,end,cb1,ce1,cb2,ce2);
01256 
01257    // Create the plynomial terms and fill
01258    // them with node information
01259    fPoly = new TSplinePoly5[fNp];
01260    for (Int_t i=0; i<n; ++i) {
01261       fPoly[i+beg].X() = x[i];
01262       fPoly[i+beg].Y() = y[i];
01263    }
01264 
01265    // Set the double knots at boundaries
01266    SetBoundaries(b1,e1,b2,e2,cb1,ce1,cb2,ce2);
01267 
01268    // Build the spline coefficients
01269    BuildCoeff();
01270 }
01271 
01272 
01273 //______________________________________________________________________________
01274 TSpline5::TSpline5(const char *title,
01275                    Double_t xmin, Double_t xmax,
01276                    Double_t y[], Int_t n,
01277                    const char *opt, Double_t b1, Double_t e1,
01278                    Double_t b2, Double_t e2) :
01279   TSpline(title,(xmax-xmin)/(n-1), xmin, xmax, n, kTRUE)
01280 {
01281    // Quintic natural spline creator given an array of
01282    // arbitrary function values on equidistant n abscissa
01283    // values from xmin to xmax and possibly end point conditions
01284 
01285    Int_t beg, end;
01286    const char *cb1, *ce1, *cb2, *ce2;
01287    fName="Spline5";
01288 
01289    // Check endpoint conditions
01290    BoundaryConditions(opt,beg,end,cb1,ce1,cb2,ce2);
01291 
01292    // Create the plynomial terms and fill
01293    // them with node information
01294    fPoly = new TSplinePoly5[fNp];
01295    for (Int_t i=0; i<n; ++i) {
01296       fPoly[i+beg].X() = fXmin+i*fDelta;
01297       fPoly[i+beg].Y() = y[i];
01298    }
01299 
01300    // Set the double knots at boundaries
01301    SetBoundaries(b1,e1,b2,e2,cb1,ce1,cb2,ce2);
01302 
01303    // Build the spline coefficients
01304    BuildCoeff();
01305 }
01306 
01307 
01308 //______________________________________________________________________________
01309 TSpline5::TSpline5(const char *title,
01310                    Double_t x[], const TF1 *func, Int_t n,
01311                    const char *opt, Double_t b1, Double_t e1,
01312                    Double_t b2, Double_t e2) :
01313   TSpline(title,-1, x[0], x[n-1], n, kFALSE)
01314 {
01315    // Quintic natural spline creator given an array of
01316    // arbitrary abscissas in increasing order and a function
01317    // to interpolate and possibly end point conditions
01318 
01319    Int_t beg, end;
01320    const char *cb1, *ce1, *cb2, *ce2;
01321    fName="Spline5";
01322 
01323    // Check endpoint conditions
01324    BoundaryConditions(opt,beg,end,cb1,ce1,cb2,ce2);
01325 
01326    // Create the plynomial terms and fill
01327    // them with node information
01328    fPoly = new TSplinePoly5[fNp];
01329    for (Int_t i=0; i<n; i++) {
01330       fPoly[i+beg].X() = x[i];
01331       fPoly[i+beg].Y() = ((TF1*)func)->Eval(x[i]);
01332    }
01333 
01334    // Set the double knots at boundaries
01335    SetBoundaries(b1,e1,b2,e2,cb1,ce1,cb2,ce2);
01336 
01337    // Build the spline coefficients
01338    BuildCoeff();
01339 }
01340 
01341 
01342 //______________________________________________________________________________
01343 TSpline5::TSpline5(const char *title,
01344                    Double_t xmin, Double_t xmax,
01345                    const TF1 *func, Int_t n,
01346                    const char *opt, Double_t b1, Double_t e1,
01347                    Double_t b2, Double_t e2) :
01348   TSpline(title,(xmax-xmin)/(n-1), xmin, xmax, n, kTRUE)
01349 {
01350    // Quintic natural spline creator given a function to be
01351    // evaluated on n equidistand abscissa points between xmin
01352    // and xmax and possibly end point conditions
01353 
01354    Int_t beg, end;
01355    const char *cb1, *ce1, *cb2, *ce2;
01356    fName="Spline5";
01357 
01358    // Check endpoint conditions
01359    BoundaryConditions(opt,beg,end,cb1,ce1,cb2,ce2);
01360 
01361    // Create the plynomial terms and fill
01362    // them with node information
01363    fPoly = new TSplinePoly5[fNp];
01364    for (Int_t i=0; i<n; ++i) {
01365       Double_t x=fXmin+i*fDelta;
01366       fPoly[i+beg].X() = x;
01367       if (func) fPoly[i+beg].Y() = ((TF1*)func)->Eval(x);
01368    }
01369    if (!func) {fDelta = -1; fKstep = kFALSE;}
01370 
01371    // Set the double knots at boundaries
01372    SetBoundaries(b1,e1,b2,e2,cb1,ce1,cb2,ce2);
01373 
01374    // Build the spline coefficients
01375    if (func) BuildCoeff();
01376 }
01377 
01378 
01379 //______________________________________________________________________________
01380 TSpline5::TSpline5(const char *title,
01381                    const TGraph *g,
01382                    const char *opt, Double_t b1, Double_t e1,
01383                    Double_t b2, Double_t e2) :
01384   TSpline(title,-1,0,0,g->GetN(),kFALSE)
01385 {
01386    // Quintic natural spline creator given a TGraph with
01387    // abscissa in increasing order and possibly end
01388    // point conditions
01389 
01390    Int_t beg, end;
01391    const char *cb1, *ce1, *cb2, *ce2;
01392    fName="Spline5";
01393 
01394    // Check endpoint conditions
01395    BoundaryConditions(opt,beg,end,cb1,ce1,cb2,ce2);
01396 
01397    // Create the plynomial terms and fill
01398    // them with node information
01399    fPoly = new TSplinePoly5[fNp];
01400    for (Int_t i=0; i<fNp-beg; ++i) {
01401       Double_t xx, yy;
01402       g->GetPoint(i,xx,yy);
01403       fPoly[i+beg].X()=xx;
01404       fPoly[i+beg].Y()=yy;
01405    }
01406 
01407    // Set the double knots at boundaries
01408    SetBoundaries(b1,e1,b2,e2,cb1,ce1,cb2,ce2);
01409    fXmin = fPoly[0].X();
01410    fXmax = fPoly[fNp-1].X();
01411 
01412    // Build the spline coefficients
01413    BuildCoeff();
01414 }
01415 
01416 
01417 //______________________________________________________________________________
01418 TSpline5::TSpline5(const TH1 *h,
01419                    const char *opt, Double_t b1, Double_t e1,
01420                    Double_t b2, Double_t e2) :
01421   TSpline(h->GetTitle(),-1,0,0,h->GetNbinsX(),kFALSE)
01422 {
01423    // Quintic natural spline creator given a TH1
01424 
01425    Int_t beg, end;
01426    const char *cb1, *ce1, *cb2, *ce2;
01427    fName=h->GetName();
01428 
01429    // Check endpoint conditions
01430    BoundaryConditions(opt,beg,end,cb1,ce1,cb2,ce2);
01431 
01432    // Create the plynomial terms and fill
01433    // them with node information
01434    fPoly = new TSplinePoly5[fNp];
01435    for (Int_t i=0; i<fNp-beg; ++i) {
01436       fPoly[i+beg].X()=h->GetXaxis()->GetBinCenter(i+1);
01437       fPoly[i+beg].Y()=h->GetBinContent(i+1);
01438    }
01439 
01440    // Set the double knots at boundaries
01441    SetBoundaries(b1,e1,b2,e2,cb1,ce1,cb2,ce2);
01442    fXmin = fPoly[0].X();
01443    fXmax = fPoly[fNp-1].X();
01444 
01445    // Build the spline coefficients
01446    BuildCoeff();
01447 }
01448 
01449 
01450 //______________________________________________________________________________
01451 TSpline5::TSpline5(const TSpline5& sp5) :
01452   TSpline(sp5),
01453   fPoly(0)
01454 {
01455    //copy constructor
01456    if (fNp > 0) fPoly = new TSplinePoly5[fNp];
01457    for (Int_t i=0; i<fNp; ++i) {
01458       fPoly[i] = sp5.fPoly[i];
01459    }
01460 }
01461 
01462 
01463 //______________________________________________________________________________
01464 TSpline5& TSpline5::operator=(const TSpline5& sp5)
01465 {
01466    //assignment operator
01467    if(this!=&sp5) {
01468       TSpline::operator=(sp5);
01469       fPoly=0;
01470       if (fNp > 0) fPoly = new TSplinePoly5[fNp];
01471       for (Int_t i=0; i<fNp; ++i) {
01472          fPoly[i] = sp5.fPoly[i];
01473       }
01474    }
01475    return *this;
01476 }
01477 
01478 
01479 //______________________________________________________________________________
01480 void TSpline5::BoundaryConditions(const char *opt,Int_t &beg,Int_t &end,
01481                                   const char *&cb1,const char *&ce1,
01482                                   const char *&cb2,const char *&ce2)
01483 {
01484    // Check the boundary conditions and the
01485    // amount of extra double knots needed
01486 
01487    cb1=ce1=cb2=ce2=0;
01488    beg=end=0;
01489    if(opt) {
01490       cb1 = strstr(opt,"b1");
01491       ce1 = strstr(opt,"e1");
01492       cb2 = strstr(opt,"b2");
01493       ce2 = strstr(opt,"e2");
01494       if(cb2) {
01495          fNp=fNp+2;
01496          beg=2;
01497       } else if(cb1) {
01498          fNp=fNp+1;
01499          beg=1;
01500       }
01501       if(ce2) {
01502          fNp=fNp+2;
01503          end=2;
01504       } else if(ce1) {
01505          fNp=fNp+1;
01506          end=1;
01507       }
01508    }
01509 }
01510 
01511 
01512 //______________________________________________________________________________
01513 void TSpline5::SetBoundaries(Double_t b1, Double_t e1, Double_t b2, Double_t e2,
01514                              const char *cb1, const char *ce1, const char *cb2,
01515                              const char *ce2)
01516 {
01517    // Set the boundary conditions at double/triple knots
01518 
01519    if(cb2) {
01520 
01521       // Second derivative at the beginning
01522       fPoly[0].X() = fPoly[1].X() = fPoly[2].X();
01523       fPoly[0].Y() = fPoly[2].Y();
01524       fPoly[2].Y()=b2;
01525 
01526       // If first derivative not given, we take the finite
01527       // difference from first and second point... not recommended
01528       if(cb1)
01529          fPoly[1].Y()=b1;
01530       else
01531          fPoly[1].Y()=(fPoly[3].Y()-fPoly[0].Y())/(fPoly[3].X()-fPoly[2].X());
01532    } else if(cb1) {
01533 
01534       // First derivative at the end
01535       fPoly[0].X() = fPoly[1].X();
01536       fPoly[0].Y() = fPoly[1].Y();
01537       fPoly[1].Y()=b1;
01538    }
01539    if(ce2) {
01540 
01541       // Second derivative at the end
01542       fPoly[fNp-1].X() = fPoly[fNp-2].X() = fPoly[fNp-3].X();
01543       fPoly[fNp-1].Y()=e2;
01544 
01545       // If first derivative not given, we take the finite
01546       // difference from first and second point... not recommended
01547       if(ce1)
01548          fPoly[fNp-2].Y()=e1;
01549       else
01550          fPoly[fNp-2].Y()=
01551          (fPoly[fNp-3].Y()-fPoly[fNp-4].Y())
01552          /(fPoly[fNp-3].X()-fPoly[fNp-4].X());
01553    } else if(ce1) {
01554 
01555       // First derivative at the end
01556       fPoly[fNp-1].X() = fPoly[fNp-2].X();
01557       fPoly[fNp-1].Y()=e1;
01558    }
01559 }
01560 
01561 
01562 //______________________________________________________________________________
01563 Int_t TSpline5::FindX(Double_t x) const
01564 {
01565    // Find X
01566 
01567    Int_t klow=0;
01568 
01569    // If out of boundaries, extrapolate
01570    // It may be badly wrong
01571    if(x<=fXmin) klow=0;
01572    else if(x>=fXmax) klow=fNp-1;
01573    else {
01574       if(fKstep) {
01575 
01576          // Equidistant knots, use histogramming
01577          klow = TMath::Min(Int_t((x-fXmin)/fDelta),fNp-1);
01578       } else {
01579          Int_t khig=fNp-1, khalf;
01580 
01581          // Non equidistant knots, binary search
01582          while(khig-klow>1)
01583             if(x>fPoly[khalf=(klow+khig)/2].X())
01584                klow=khalf;
01585             else
01586                khig=khalf;
01587       }
01588 
01589       // This could be removed, sanity check
01590       if(!(fPoly[klow].X()<=x && x<=fPoly[klow+1].X()))
01591          Error("Eval",
01592                "Binary search failed x(%d) = %f < x(%d) = %f\n",
01593                 klow,fPoly[klow].X(),klow+1,fPoly[klow+1].X());
01594    }
01595    return klow;
01596 }
01597 
01598 
01599 //______________________________________________________________________________
01600 Double_t TSpline5::Eval(Double_t x) const
01601 {
01602    // Eval this spline at x
01603 
01604    Int_t klow=FindX(x);
01605    return fPoly[klow].Eval(x);
01606 }
01607 
01608 
01609 //______________________________________________________________________________
01610 Double_t TSpline5::Derivative(Double_t x) const
01611 {
01612    // Derivative
01613 
01614    Int_t klow=FindX(x);
01615    return fPoly[klow].Derivative(x);
01616 }
01617 
01618 
01619 //______________________________________________________________________________
01620 void TSpline5::SaveAs(const char *filename, Option_t * /*option*/) const
01621 {
01622    // write this spline as a C++ function that can be executed without ROOT
01623    // the name of the function is the name of the file up to the "." if any
01624 
01625    //open the file
01626    ofstream *f = new ofstream(filename,ios::out);
01627    if (f == 0 || gSystem->AccessPathName(filename,kWritePermission)) {
01628       Error("SaveAs","Cannot open file:%s\n",filename);
01629       return;
01630    }
01631 
01632    //write the function name and the spline constants
01633    char buffer[512];
01634    Int_t nch = strlen(filename);
01635    snprintf(buffer,512,"double %s",filename);
01636    char *dot = strstr(buffer,".");
01637    if (dot) *dot = 0;
01638    strlcat(buffer,"(double x) {\n",512);
01639    nch = strlen(buffer); f->write(buffer,nch);
01640    snprintf(buffer,512,"   const int fNp = %d, fKstep = %d;\n",fNp,fKstep);
01641    nch = strlen(buffer); f->write(buffer,nch);
01642    snprintf(buffer,512,"   const double fDelta = %g, fXmin = %g, fXmax = %g;\n",fDelta,fXmin,fXmax);
01643    nch = strlen(buffer); f->write(buffer,nch);
01644 
01645    //write the spline coefficients
01646    //array fX
01647    snprintf(buffer,512,"   const double fX[%d] = {",fNp);
01648    nch = strlen(buffer); f->write(buffer,nch);
01649    buffer[0] = 0;
01650    Int_t i;
01651    char numb[20];
01652    for (i=0;i<fNp;i++) {
01653       snprintf(numb,20," %g,",fPoly[i].X());
01654       nch = strlen(numb);
01655       if (i == fNp-1) numb[nch-1]=0;
01656       strlcat(buffer,numb,512);
01657       if (i%5 == 4 || i == fNp-1) {
01658          nch = strlen(buffer); f->write(buffer,nch);
01659          if (i != fNp-1) snprintf(buffer,512,"\n                       ");
01660       }
01661    }
01662    snprintf(buffer,512," };\n");
01663    nch = strlen(buffer); f->write(buffer,nch);
01664    //array fY
01665    snprintf(buffer,512,"   const double fY[%d] = {",fNp);
01666    nch = strlen(buffer); f->write(buffer,nch);
01667    buffer[0] = 0;
01668    for (i=0;i<fNp;i++) {
01669       snprintf(numb,20," %g,",fPoly[i].Y());
01670       nch = strlen(numb);
01671       if (i == fNp-1) numb[nch-1]=0;
01672       strlcat(buffer,numb,512);
01673       if (i%5 == 4 || i == fNp-1) {
01674          nch = strlen(buffer); f->write(buffer,nch);
01675          if (i != fNp-1) snprintf(buffer,512,"\n                       ");
01676       }
01677    }
01678    snprintf(buffer,512," };\n");
01679    nch = strlen(buffer); f->write(buffer,nch);
01680    //array fB
01681    snprintf(buffer,512,"   const double fB[%d] = {",fNp);
01682    nch = strlen(buffer); f->write(buffer,nch);
01683    buffer[0] = 0;
01684    for (i=0;i<fNp;i++) {
01685       snprintf(numb,20," %g,",fPoly[i].B());
01686       nch = strlen(numb);
01687       if (i == fNp-1) numb[nch-1]=0;
01688       strlcat(buffer,numb,512);
01689       if (i%5 == 4 || i == fNp-1) {
01690          nch = strlen(buffer); f->write(buffer,nch);
01691          if (i != fNp-1) snprintf(buffer,512,"\n                       ");
01692       }
01693    }
01694    snprintf(buffer,512," };\n");
01695    nch = strlen(buffer); f->write(buffer,nch);
01696    //array fC
01697    snprintf(buffer,512,"   const double fC[%d] = {",fNp);
01698    nch = strlen(buffer); f->write(buffer,nch);
01699    buffer[0] = 0;
01700    for (i=0;i<fNp;i++) {
01701       snprintf(numb,20," %g,",fPoly[i].C());
01702       nch = strlen(numb);
01703       if (i == fNp-1) numb[nch-1]=0;
01704       strlcat(buffer,numb,512);
01705       if (i%5 == 4 || i == fNp-1) {
01706          nch = strlen(buffer); f->write(buffer,nch);
01707          if (i != fNp-1) snprintf(buffer,512,"\n                       ");
01708       }
01709    }
01710    snprintf(buffer,512," };\n");
01711    nch = strlen(buffer); f->write(buffer,nch);
01712     //array fD
01713    snprintf(buffer,512,"   const double fD[%d] = {",fNp);
01714    nch = strlen(buffer); f->write(buffer,nch);
01715    buffer[0] = 0;
01716    for (i=0;i<fNp;i++) {
01717       snprintf(numb,20," %g,",fPoly[i].D());
01718       nch = strlen(numb);
01719       if (i == fNp-1) numb[nch-1]=0;
01720       strlcat(buffer,numb,512);
01721       if (i%5 == 4 || i == fNp-1) {
01722          nch = strlen(buffer); f->write(buffer,nch);
01723          if (i != fNp-1) snprintf(buffer,512,"\n                       ");
01724       }
01725    }
01726    snprintf(buffer,512," };\n");
01727    nch = strlen(buffer); f->write(buffer,nch);
01728     //array fE
01729    snprintf(buffer,512,"   const double fE[%d] = {",fNp);
01730    nch = strlen(buffer); f->write(buffer,nch);
01731    buffer[0] = 0;
01732    for (i=0;i<fNp;i++) {
01733       snprintf(numb,20," %g,",fPoly[i].E());
01734       nch = strlen(numb);
01735       if (i == fNp-1) numb[nch-1]=0;
01736       strlcat(buffer,numb,512);
01737       if (i%5 == 4 || i == fNp-1) {
01738          nch = strlen(buffer); f->write(buffer,nch);
01739          if (i != fNp-1) snprintf(buffer,512,"\n                       ");
01740       }
01741    }
01742    snprintf(buffer,512," };\n");
01743    nch = strlen(buffer); f->write(buffer,nch);
01744     //array fF
01745    snprintf(buffer,512,"   const double fF[%d] = {",fNp);
01746    nch = strlen(buffer); f->write(buffer,nch);
01747    buffer[0] = 0;
01748    for (i=0;i<fNp;i++) {
01749       snprintf(numb,20," %g,",fPoly[i].F());
01750       nch = strlen(numb);
01751       if (i == fNp-1) numb[nch-1]=0;
01752       strlcat(buffer,numb,512);
01753       if (i%5 == 4 || i == fNp-1) {
01754          nch = strlen(buffer); f->write(buffer,nch);
01755          if (i != fNp-1) snprintf(buffer,512,"\n                       ");
01756       }
01757    }
01758    snprintf(buffer,512," };\n");
01759    nch = strlen(buffer); f->write(buffer,nch);
01760 
01761    //generate code for the spline evaluation
01762    snprintf(buffer,512,"   int klow=0;\n");
01763    nch = strlen(buffer); f->write(buffer,nch);
01764 
01765    snprintf(buffer,512,"   // If out of boundaries, extrapolate. It may be badly wrong\n");
01766    snprintf(buffer,512,"   if(x<=fXmin) klow=0;\n");
01767    nch = strlen(buffer); f->write(buffer,nch);
01768    snprintf(buffer,512,"   else if(x>=fXmax) klow=fNp-1;\n");
01769    nch = strlen(buffer); f->write(buffer,nch);
01770    snprintf(buffer,512,"   else {\n");
01771    nch = strlen(buffer); f->write(buffer,nch);
01772    snprintf(buffer,512,"     if(fKstep) {\n");
01773    nch = strlen(buffer); f->write(buffer,nch);
01774 
01775    snprintf(buffer,512,"       // Equidistant knots, use histogramming\n");
01776    nch = strlen(buffer); f->write(buffer,nch);
01777    snprintf(buffer,512,"       klow = int((x-fXmin)/fDelta);\n");
01778    nch = strlen(buffer); f->write(buffer,nch);
01779    snprintf(buffer,512,"       if (klow < fNp-1) klow = fNp-1;\n");
01780    nch = strlen(buffer); f->write(buffer,nch);
01781    snprintf(buffer,512,"     } else {\n");
01782    nch = strlen(buffer); f->write(buffer,nch);
01783    snprintf(buffer,512,"       int khig=fNp-1, khalf;\n");
01784    nch = strlen(buffer); f->write(buffer,nch);
01785 
01786    snprintf(buffer,512,"       // Non equidistant knots, binary search\n");
01787    nch = strlen(buffer); f->write(buffer,nch);
01788    snprintf(buffer,512,"       while(khig-klow>1)\n");
01789    nch = strlen(buffer); f->write(buffer,nch);
01790    snprintf(buffer,512,"         if(x>fX[khalf=(klow+khig)/2]) klow=khalf;\n");
01791    nch = strlen(buffer); f->write(buffer,nch);
01792    snprintf(buffer,512,"         else khig=khalf;\n");
01793    nch = strlen(buffer); f->write(buffer,nch);
01794    snprintf(buffer,512,"     }\n");
01795    nch = strlen(buffer); f->write(buffer,nch);
01796    snprintf(buffer,512,"   }\n");
01797    nch = strlen(buffer); f->write(buffer,nch);
01798    snprintf(buffer,512,"   // Evaluate now\n");
01799    nch = strlen(buffer); f->write(buffer,nch);
01800    snprintf(buffer,512,"   double dx=x-fX[klow];\n");
01801    nch = strlen(buffer); f->write(buffer,nch);
01802    snprintf(buffer,512,"   return (fY[klow]+dx*(fB[klow]+dx*(fC[klow]+dx*(fD[klow]+dx*(fE[klow]+dx*fF[klow])))));\n");
01803    nch = strlen(buffer); f->write(buffer,nch);
01804 
01805    //close file
01806    f->write("}\n",2);
01807 
01808    if (f) { f->close(); delete f;}
01809 }
01810 
01811 
01812 //______________________________________________________________________________
01813 void TSpline5::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
01814 {
01815     // Save primitive as a C++ statement(s) on output stream out
01816 
01817    char quote = '"';
01818    out<<"   "<<endl;
01819    if (gROOT->ClassSaved(TSpline5::Class())) {
01820       out<<"   ";
01821    } else {
01822       out<<"   TSpline5 *";
01823    }
01824    Double_t b1 = fPoly[1].Y();
01825    Double_t e1 = fPoly[fNp-1].Y();
01826    Double_t b2 = fPoly[2].Y();
01827    Double_t e2 = fPoly[fNp-1].Y();
01828    out<<"spline5 = new TSpline5("<<quote<<GetTitle()<<quote<<","
01829       <<fXmin<<","<<fXmax<<",(TF1*)0,"<<fNp<<","<<quote<<quote<<","
01830       <<b1<<","<<e1<<","<<b2<<","<<e2<<");"<<endl;
01831    out<<"   spline5->SetName("<<quote<<GetName()<<quote<<");"<<endl;
01832 
01833    SaveFillAttributes(out,"spline5",0,1001);
01834    SaveLineAttributes(out,"spline5",1,1,1);
01835    SaveMarkerAttributes(out,"spline5",1,1,1);
01836    if (fNpx != 100) out<<"   spline5->SetNpx("<<fNpx<<");"<<endl;
01837 
01838    for (Int_t i=0;i<fNp;i++) {
01839       out<<"   spline5->SetPoint("<<i<<","<<fPoly[i].X()<<","<<fPoly[i].Y()<<");"<<endl;
01840       out<<"   spline5->SetPointCoeff("<<i<<","<<fPoly[i].B()<<","<<fPoly[i].C()<<","<<fPoly[i].D()<<","<<fPoly[i].E()<<","<<fPoly[i].F()<<");"<<endl;
01841    }
01842    out<<"   spline5->Draw("<<quote<<option<<quote<<");"<<endl;
01843 }
01844 
01845 
01846 //______________________________________________________________________________
01847 void TSpline5::SetPoint(Int_t i, Double_t x, Double_t y)
01848 {
01849    //set point number i.
01850    
01851    
01852    if (i < 0 || i >= fNp) return;
01853    fPoly[i].X()= x;
01854    fPoly[i].Y()= y;
01855 }
01856 
01857 //______________________________________________________________________________
01858 void TSpline5::SetPointCoeff(Int_t i, Double_t b, Double_t c, Double_t d,
01859                              Double_t e, Double_t f)
01860 {
01861    // set point coefficient number i
01862    
01863    if (i < 0 || i >= fNp) return;
01864    fPoly[i].B()= b;
01865    fPoly[i].C()= c;
01866    fPoly[i].D()= d;
01867    fPoly[i].E()= e;
01868    fPoly[i].F()= f;
01869 }
01870 
01871 
01872 
01873 //______________________________________________________________________________
01874 void TSpline5::BuildCoeff()
01875 {
01876    //     algorithm 600, collected algorithms from acm.
01877    //     algorithm appeared in acm-trans. math. software, vol.9, no. 2,
01878    //     jun., 1983, p. 258-259.
01879    //
01880    //     TSpline5 computes the coefficients of a quintic natural quintic spli
01881    //     s(x) with knots x(i) interpolating there to given function values:
01882    //               s(x(i)) = y(i)  for i = 1,2, ..., n.
01883    //     in each interval (x(i),x(i+1)) the spline function s(xx) is a
01884    //     polynomial of fifth degree:
01885    //     s(xx) = ((((f(i)*p+e(i))*p+d(i))*p+c(i))*p+b(i))*p+y(i)    (*)
01886    //           = ((((-f(i)*q+e(i+1))*q-d(i+1))*q+c(i+1))*q-b(i+1))*q+y(i+1)
01887    //     where  p = xx - x(i)  and  q = x(i+1) - xx.
01888    //     (note the first subscript in the second expression.)
01889    //     the different polynomials are pieced together so that s(x) and
01890    //     its derivatives up to s"" are continuous.
01891    //
01892    //        input:
01893    //
01894    //     n          number of data points, (at least three, i.e. n > 2)
01895    //     x(1:n)     the strictly increasing or decreasing sequence of
01896    //                knots.  the spacing must be such that the fifth power
01897    //                of x(i+1) - x(i) can be formed without overflow or
01898    //                underflow of exponents.
01899    //     y(1:n)     the prescribed function values at the knots.
01900    //
01901    //        output:
01902    //
01903    //     b,c,d,e,f  the computed spline coefficients as in (*).
01904    //         (1:n)  specifically
01905    //                b(i) = s'(x(i)), c(i) = s"(x(i))/2, d(i) = s"'(x(i))/6,
01906    //                e(i) = s""(x(i))/24,  f(i) = s""'(x(i))/120.
01907    //                f(n) is neither used nor altered.  the five arrays
01908    //                b,c,d,e,f must always be distinct.
01909    //
01910    //        option:
01911    //
01912    //     it is possible to specify values for the first and second
01913    //     derivatives of the spline function at arbitrarily many knots.
01914    //     this is done by relaxing the requirement that the sequence of
01915    //     knots be strictly increasing or decreasing.  specifically:
01916    //
01917    //     if x(j) = x(j+1) then s(x(j)) = y(j) and s'(x(j)) = y(j+1),
01918    //     if x(j) = x(j+1) = x(j+2) then in addition s"(x(j)) = y(j+2).
01919    //
01920    //     note that s""(x) is discontinuous at a double knot and, in
01921    //     addition, s"'(x) is discontinuous at a triple knot.  the
01922    //     subroutine assigns y(i) to y(i+1) in these cases and also to
01923    //     y(i+2) at a triple knot.  the representation (*) remains
01924    //     valid in each open interval (x(i),x(i+1)).  at a double knot,
01925    //     x(j) = x(j+1), the output coefficients have the following values:
01926    //       y(j) = s(x(j))          = y(j+1)
01927    //       b(j) = s'(x(j))         = b(j+1)
01928    //       c(j) = s"(x(j))/2       = c(j+1)
01929    //       d(j) = s"'(x(j))/6      = d(j+1)
01930    //       e(j) = s""(x(j)-0)/24     e(j+1) = s""(x(j)+0)/24
01931    //       f(j) = s""'(x(j)-0)/120   f(j+1) = s""'(x(j)+0)/120
01932    //     at a triple knot, x(j) = x(j+1) = x(j+2), the output
01933    //     coefficients have the following values:
01934    //       y(j) = s(x(j))         = y(j+1)    = y(j+2)
01935    //       b(j) = s'(x(j))        = b(j+1)    = b(j+2)
01936    //       c(j) = s"(x(j))/2      = c(j+1)    = c(j+2)
01937    //       d(j) = s"'((x(j)-0)/6    d(j+1) = 0  d(j+2) = s"'(x(j)+0)/6
01938    //       e(j) = s""(x(j)-0)/24    e(j+1) = 0  e(j+2) = s""(x(j)+0)/24
01939    //       f(j) = s""'(x(j)-0)/120  f(j+1) = 0  f(j+2) = s""'(x(j)+0)/120
01940 
01941    Int_t i, m;
01942    Double_t pqqr, p, q, r, s, t, u, v,
01943       b1, p2, p3, q2, q3, r2, pq, pr, qr;
01944 
01945    if (fNp <= 2) {
01946       return;
01947    }
01948 
01949    //     coefficients of a positive definite, pentadiagonal matrix,
01950    //     stored in D, E, F from 1 to n-3.
01951    m = fNp-2;
01952    q = fPoly[1].X()-fPoly[0].X();
01953    r = fPoly[2].X()-fPoly[1].X();
01954    q2 = q*q;
01955    r2 = r*r;
01956    qr = q+r;
01957    fPoly[0].D() = fPoly[0].E() = 0;
01958    if (q) fPoly[1].D() = q*6.*q2/(qr*qr);
01959    else fPoly[1].D() = 0;
01960 
01961    if (m > 1) {
01962       for (i = 1; i < m; ++i) {
01963          p = q;
01964          q = r;
01965          r = fPoly[i+2].X()-fPoly[i+1].X();
01966          p2 = q2;
01967          q2 = r2;
01968          r2 = r*r;
01969          pq = qr;
01970          qr = q+r;
01971          if (q) {
01972             q3 = q2*q;
01973             pr = p*r;
01974             pqqr = pq*qr;
01975             fPoly[i+1].D() = q3*6./(qr*qr);
01976             fPoly[i].D() += (q+q)*(pr*15.*pr+(p+r)*q
01977                                  *(pr* 20.+q2*7.)+q2*
01978                                  ((p2+r2)*8.+pr*21.+q2+q2))/(pqqr*pqqr);
01979             fPoly[i-1].D() += q3*6./(pq*pq);
01980             fPoly[i].E() = q2*(p*qr+pq*3.*(qr+r+r))/(pqqr*qr);
01981             fPoly[i-1].E() += q2*(r*pq+qr*3.*(pq+p+p))/(pqqr*pq);
01982             fPoly[i-1].F() = q3/pqqr;
01983          } else
01984             fPoly[i+1].D() = fPoly[i].E() = fPoly[i-1].F() = 0;
01985       }
01986    }
01987    if (r) fPoly[m-1].D() += r*6.*r2/(qr*qr);
01988 
01989    //     First and second order divided differences of the given function
01990    //     values, stored in b from 2 to n and in c from 3 to n
01991    //     respectively. care is taken of double and triple knots.
01992    for (i = 1; i < fNp; ++i) {
01993       if (fPoly[i].X() != fPoly[i-1].X()) {
01994          fPoly[i].B() =
01995             (fPoly[i].Y()-fPoly[i-1].Y())/(fPoly[i].X()-fPoly[i-1].X());
01996       } else {
01997          fPoly[i].B() = fPoly[i].Y();
01998          fPoly[i].Y() = fPoly[i-1].Y();
01999       }
02000    }
02001    for (i = 2; i < fNp; ++i) {
02002       if (fPoly[i].X() != fPoly[i-2].X()) {
02003          fPoly[i].C() =
02004             (fPoly[i].B()-fPoly[i-1].B())/(fPoly[i].X()-fPoly[i-2].X());
02005       } else {
02006          fPoly[i].C() = fPoly[i].B()*.5;
02007          fPoly[i].B() = fPoly[i-1].B();
02008       }
02009    }
02010 
02011    //     Solve the linear system with c(i+2) - c(i+1) as right-hand side. */
02012    if (m > 1) {
02013       p=fPoly[0].C()=fPoly[m-1].E()=fPoly[0].F()
02014          =fPoly[m-2].F()=fPoly[m-1].F()=0;
02015       fPoly[1].C() = fPoly[3].C()-fPoly[2].C();
02016       fPoly[1].D() = 1./fPoly[1].D();
02017 
02018       if (m > 2) {
02019          for (i = 2; i < m; ++i) {
02020             q = fPoly[i-1].D()*fPoly[i-1].E();
02021             fPoly[i].D() = 1./(fPoly[i].D()-p*fPoly[i-2].F()-q*fPoly[i-1].E());
02022             fPoly[i].E() -= q*fPoly[i-1].F();
02023             fPoly[i].C() = fPoly[i+2].C()-fPoly[i+1].C()-p*fPoly[i-2].C()
02024                            -q*fPoly[i-1].C();
02025             p = fPoly[i-1].D()*fPoly[i-1].F();
02026          }
02027       }
02028    }
02029 
02030    fPoly[fNp-2].C() = fPoly[fNp-1].C() = 0;
02031    if (fNp > 3)
02032       for (i=fNp-3; i > 0; --i)
02033          fPoly[i].C() = (fPoly[i].C()-fPoly[i].E()*fPoly[i+1].C()
02034                         -fPoly[i].F()*fPoly[i+2].C())*fPoly[i].D();
02035 
02036    //     Integrate the third derivative of s(x)
02037    m = fNp-1;
02038    q = fPoly[1].X()-fPoly[0].X();
02039    r = fPoly[2].X()-fPoly[1].X();
02040    b1 = fPoly[1].B();
02041    q3 = q*q*q;
02042    qr = q+r;
02043    if (qr) {
02044       v = fPoly[1].C()/qr;
02045       t = v;
02046    } else
02047       v = t = 0;
02048    if (q) fPoly[0].F() = v/q;
02049    else fPoly[0].F() = 0;
02050    for (i = 1; i < m; ++i) {
02051       p = q;
02052       q = r;
02053       if (i != m-1) r = fPoly[i+2].X()-fPoly[i+1].X();
02054       else r = 0;
02055       p3 = q3;
02056       q3 = q*q*q;
02057       pq = qr;
02058       qr = q+r;
02059       s = t;
02060       if (qr) t = (fPoly[i+1].C()-fPoly[i].C())/qr;
02061       else t = 0;
02062       u = v;
02063       v = t-s;
02064       if (pq) {
02065          fPoly[i].F() = fPoly[i-1].F();
02066          if (q) fPoly[i].F() = v/q;
02067          fPoly[i].E() = s*5.;
02068          fPoly[i].D() = (fPoly[i].C()-q*s)*10;
02069          fPoly[i].C() =
02070          fPoly[i].D()*(p-q)+(fPoly[i+1].B()-fPoly[i].B()+(u-fPoly[i].E())*
02071                             p3-(v+fPoly[i].E())*q3)/pq;
02072          fPoly[i].B() = (p*(fPoly[i+1].B()-v*q3)+q*(fPoly[i].B()-u*p3))/pq-p
02073          *q*(fPoly[i].D()+fPoly[i].E()*(q-p));
02074       } else {
02075          fPoly[i].C() = fPoly[i-1].C();
02076          fPoly[i].D() = fPoly[i].E() = fPoly[i].F() = 0;
02077       }
02078    }
02079 
02080    //     End points x(1) and x(n)
02081    p = fPoly[1].X()-fPoly[0].X();
02082    s = fPoly[0].F()*p*p*p;
02083    fPoly[0].E() = fPoly[0].D() = 0;
02084    fPoly[0].C() = fPoly[1].C()-s*10;
02085    fPoly[0].B() = b1-(fPoly[0].C()+s)*p;
02086 
02087    q = fPoly[fNp-1].X()-fPoly[fNp-2].X();
02088    t = fPoly[fNp-2].F()*q*q*q;
02089    fPoly[fNp-1].E() = fPoly[fNp-1].D() = 0;
02090    fPoly[fNp-1].C() = fPoly[fNp-2].C()+t*10;
02091    fPoly[fNp-1].B() += (fPoly[fNp-1].C()-t)*q;
02092 }
02093 
02094 
02095 //______________________________________________________________________________
02096 void TSpline5::Test()
02097 {
02098    // Test method for TSpline5
02099    //
02100    //
02101    //   n          number of data points.
02102    //   m          2*m-1 is order of spline.
02103    //                 m = 3 always for quintic spline.
02104    //   nn,nm1,mm,
02105    //   mm1,i,k,
02106    //   j,jj       temporary integer variables.
02107    //   z,p        temporary double precision variables.
02108    //   x[n]       the sequence of knots.
02109    //   y[n]       the prescribed function values at the knots.
02110    //   a[200][6]  two dimensional array whose columns are
02111    //                 the computed spline coefficients
02112    //   diff[5]    maximum values of differences of values and
02113    //                 derivatives to right and left of knots.
02114    //   com[5]     maximum values of coefficients.
02115    //
02116    //
02117    //   test of TSpline5 with nonequidistant knots and
02118    //      equidistant knots follows.
02119 
02120    Double_t hx;
02121    Double_t diff[5];
02122    Double_t a[1200], c[6];
02123    Int_t i, j, k, m, n;
02124    Double_t p, x[200], y[200], z;
02125    Int_t jj, mm, nn;
02126    Int_t mm1, nm1;
02127    Double_t com[5];
02128 
02129    printf("1         TEST OF TSpline5 WITH NONEQUIDISTANT KNOTS\n");
02130    n = 5;
02131    x[0] = -3;
02132    x[1] = -1;
02133    x[2] = 0;
02134    x[3] = 3;
02135    x[4] = 4;
02136    y[0] = 7;
02137    y[1] = 11;
02138    y[2] = 26;
02139    y[3] = 56;
02140    y[4] = 29;
02141    m = 3;
02142    mm = m << 1;
02143    mm1 = mm-1;
02144    printf("\n-N = %3d    M =%2d\n",n,m);
02145    TSpline5 *spline = new TSpline5("Test",x,y,n);
02146    for (i = 0; i < n; ++i)
02147       spline->GetCoeff(i,hx, a[i],a[i+200],a[i+400],
02148                        a[i+600],a[i+800],a[i+1000]);
02149    delete spline;
02150    for (i = 0; i < mm1; ++i) diff[i] = com[i] = 0;
02151    for (k = 0; k < n; ++k) {
02152       for (i = 0; i < mm; ++i) c[i] = a[k+i*200];
02153       printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
02154       printf("%12.8f\n",x[k]);
02155       if (k == n-1) {
02156          printf("%16.8f\n",c[0]);
02157       } else {
02158          for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
02159          printf("\n");
02160          for (i = 0; i < mm1; ++i)
02161             if ((z=TMath::Abs(a[k+i*200])) > com[i]) com[i] = z;
02162          z = x[k+1]-x[k];
02163          for (i = 1; i < mm; ++i)
02164             for (jj = i; jj < mm; ++jj) {
02165                j = mm+i-jj;
02166                c[j-2] = c[j-1]*z+c[j-2];
02167             }
02168          for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
02169          printf("\n");
02170          for (i = 0; i < mm1; ++i)
02171             if (!(k >= n-2 && i != 0))
02172                if((z = TMath::Abs(c[i]-a[k+1+i*200]))
02173                   > diff[i]) diff[i] = z;
02174       }
02175    }
02176    printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
02177    for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
02178    printf("\n");
02179    printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
02180    if (TMath::Abs(c[0]) > com[0])
02181       com[0] = TMath::Abs(c[0]);
02182    for (i = 0; i < mm1; ++i) printf("%16.8f",com[i]);
02183    printf("\n");
02184    m = 3;
02185    for (n = 10; n <= 100; n += 10) {
02186       mm = m << 1;
02187       mm1 = mm-1;
02188       nm1 = n-1;
02189       for (i = 0; i < nm1; i += 2) {
02190          x[i] = i+1;
02191          x[i+1] = i+2;
02192          y[i] = 1;
02193          y[i+1] = 0;
02194       }
02195       if (n % 2 != 0) {
02196          x[n-1] = n;
02197          y[n-1] = 1;
02198       }
02199       printf("\n-N = %3d    M =%2d\n",n,m);
02200       spline = new TSpline5("Test",x,y,n);
02201       for (i = 0; i < n; ++i)
02202          spline->GetCoeff(i,hx,a[i],a[i+200],a[i+400],
02203                           a[i+600],a[i+800],a[i+1000]);
02204       delete spline;
02205       for (i = 0; i < mm1; ++i)
02206          diff[i] = com[i] = 0;
02207       for (k = 0; k < n; ++k) {
02208          for (i = 0; i < mm; ++i)
02209             c[i] = a[k+i*200];
02210          if (n < 11) {
02211             printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
02212             printf("%12.8f\n",x[k]);
02213             if (k == n-1) printf("%16.8f\n",c[0]);
02214          }
02215          if (k == n-1) break;
02216          if (n <= 10) {
02217             for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
02218             printf("\n");
02219          }
02220          for (i = 0; i < mm1; ++i)
02221             if ((z=TMath::Abs(a[k+i*200])) > com[i])
02222                com[i] = z;
02223          z = x[k+1]-x[k];
02224          for (i = 1; i < mm; ++i)
02225             for (jj = i; jj < mm; ++jj) {
02226                j = mm+i-jj;
02227                c[j-2] = c[j-1]*z+c[j-2];
02228             }
02229          if (n <= 10) {
02230             for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
02231             printf("\n");
02232          }
02233          for (i = 0; i < mm1; ++i)
02234             if (!(k >= n-2 && i != 0))
02235                if ((z = TMath::Abs(c[i]-a[k+1+i*200]))
02236                   > diff[i]) diff[i] = z;
02237       }
02238       printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
02239       for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
02240       printf("\n");
02241       printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
02242       if (TMath::Abs(c[0]) > com[0])
02243          com[0] = TMath::Abs(c[0]);
02244       for (i = 0; i < mm1; ++i) printf("%16.8E",com[i]);
02245       printf("\n");
02246    }
02247 
02248    //     Test of TSpline5 with nonequidistant double knots follows
02249    printf("1  TEST OF TSpline5 WITH NONEQUIDISTANT DOUBLE KNOTS\n");
02250    n = 5;
02251    x[0] = -3;
02252    x[1] = -3;
02253    x[2] = -1;
02254    x[3] = -1;
02255    x[4] = 0;
02256    x[5] = 0;
02257    x[6] = 3;
02258    x[7] = 3;
02259    x[8] = 4;
02260    x[9] = 4;
02261    y[0] = 7;
02262    y[1] = 2;
02263    y[2] = 11;
02264    y[3] = 15;
02265    y[4] = 26;
02266    y[5] = 10;
02267    y[6] = 56;
02268    y[7] = -27;
02269    y[8] = 29;
02270    y[9] = -30;
02271    m = 3;
02272    nn = n << 1;
02273    mm = m << 1;
02274    mm1 = mm-1;
02275    printf("-N = %3d    M =%2d\n",n,m);
02276    spline = new TSpline5("Test",x,y,nn);
02277    for (i = 0; i < nn; ++i)
02278       spline->GetCoeff(i,hx,a[i],a[i+200],a[i+400],
02279                        a[i+600],a[i+800],a[i+1000]);
02280    delete spline;
02281    for (i = 0; i < mm1; ++i)
02282       diff[i] = com[i] = 0;
02283    for (k = 0; k < nn; ++k) {
02284       for (i = 0; i < mm; ++i)
02285          c[i] = a[k+i*200];
02286       printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
02287       printf("%12.8f\n",x[k]);
02288       if (k == nn-1) {
02289          printf("%16.8f\n",c[0]);
02290          break;
02291       }
02292       for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
02293       printf("\n");
02294       for (i = 0; i < mm1; ++i)
02295          if ((z=TMath::Abs(a[k+i*200])) > com[i]) com[i] = z;
02296       z = x[k+1]-x[k];
02297       for (i = 1; i < mm; ++i)
02298          for (jj = i; jj < mm; ++jj) {
02299             j = mm+i-jj;
02300             c[j-2] = c[j-1]*z+c[j-2];
02301          }
02302       for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
02303       printf("\n");
02304       for (i = 0; i < mm1; ++i)
02305          if (!(k >= nn-2 && i != 0))
02306             if ((z = TMath::Abs(c[i]-a[k+1+i*200]))
02307                > diff[i]) diff[i] = z;
02308    }
02309    printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
02310    for (i = 1; i <= mm1; ++i) {
02311       printf("%18.9E",diff[i-1]);
02312    }
02313    printf("\n");
02314    if (TMath::Abs(c[0]) > com[0])
02315       com[0] = TMath::Abs(c[0]);
02316    printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
02317    for (i = 0; i < mm1; ++i) printf("%16.8f",com[i]);
02318    printf("\n");
02319    m = 3;
02320    for (n = 10; n <= 100; n += 10) {
02321       nn = n << 1;
02322       mm = m << 1;
02323       mm1 = mm-1;
02324       p = 0;
02325       for (i = 0; i < n; ++i) {
02326          p += TMath::Abs(TMath::Sin(i+1));
02327          x[(i << 1)] = p;
02328          x[(i << 1)+1] = p;
02329          y[(i << 1)] = TMath::Cos(i+1)-.5;
02330          y[(i << 1)+1] = TMath::Cos((i << 1)+2)-.5;
02331       }
02332       printf("-N = %3d    M =%2d\n",n,m);
02333       spline = new TSpline5("Test",x,y,nn);
02334       for (i = 0; i < nn; ++i)
02335          spline->GetCoeff(i,hx,a[i],a[i+200],a[i+400],
02336                           a[i+600],a[i+800],a[i+1000]);
02337       delete spline;
02338       for (i = 0; i < mm1; ++i)
02339          diff[i] = com[i] = 0;
02340       for (k = 0; k < nn; ++k) {
02341          for (i = 0; i < mm; ++i)
02342             c[i] = a[k+i*200];
02343          if (n < 11) {
02344             printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
02345             printf("%12.8f\n",x[k]);
02346             if (k == nn-1) printf("%16.8f\n",c[0]);
02347          }
02348          if (k == nn-1) break;
02349          if (n <= 10) {
02350             for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
02351             printf("\n");
02352          }
02353          for (i = 0; i < mm1; ++i)
02354             if ((z=TMath::Abs(a[k+i*200])) > com[i]) com[i] = z;
02355          z = x[k+1]-x[k];
02356          for (i = 1; i < mm; ++i) {
02357             for (jj = i; jj < mm; ++jj) {
02358                j = mm+i-jj;
02359                c[j-2] = c[j-1]*z+c[j-2];
02360             }
02361          }
02362          if (n <= 10) {
02363             for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
02364             printf("\n");
02365          }
02366          for (i = 0; i < mm1; ++i)
02367             if (!(k >= nn-2 && i != 0))
02368                if ((z = TMath::Abs(c[i]-a[k+1+i*200]))
02369                   > diff[i]) diff[i] = z;
02370       }
02371       printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
02372       for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
02373       printf("\n");
02374       printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
02375       if (TMath::Abs(c[0]) > com[0])
02376          com[0] = TMath::Abs(c[0]);
02377       for (i = 0; i < mm1; ++i) printf("%18.9E",com[i]);
02378       printf("\n");
02379    }
02380 
02381    //     test of TSpline5 with nonequidistant knots, one double knot,
02382    //        one triple knot, follows.
02383    printf("1         TEST OF TSpline5 WITH NONEQUIDISTANT KNOTS,\n");
02384    printf("             ONE DOUBLE, ONE TRIPLE KNOT\n");
02385    n = 8;
02386    x[0] = -3;
02387    x[1] = -1;
02388    x[2] = -1;
02389    x[3] = 0;
02390    x[4] = 3;
02391    x[5] = 3;
02392    x[6] = 3;
02393    x[7] = 4;
02394    y[0] = 7;
02395    y[1] = 11;
02396    y[2] = 15;
02397    y[3] = 26;
02398    y[4] = 56;
02399    y[5] = -30;
02400    y[6] = -7;
02401    y[7] = 29;
02402    m = 3;
02403    mm = m << 1;
02404    mm1 = mm-1;
02405    printf("-N = %3d    M =%2d\n",n,m);
02406    spline=new TSpline5("Test",x,y,n);
02407    for (i = 0; i < n; ++i)
02408       spline->GetCoeff(i,hx,a[i],a[i+200],a[i+400],
02409                        a[i+600],a[i+800],a[i+1000]);
02410    delete spline;
02411    for (i = 0; i < mm1; ++i)
02412       diff[i] = com[i] = 0;
02413    for (k = 0; k < n; ++k) {
02414       for (i = 0; i < mm; ++i)
02415          c[i] = a[k+i*200];
02416       printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
02417       printf("%12.8f\n",x[k]);
02418       if (k == n-1) {
02419          printf("%16.8f\n",c[0]);
02420          break;
02421       }
02422       for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
02423       printf("\n");
02424       for (i = 0; i < mm1; ++i)
02425          if ((z=TMath::Abs(a[k+i*200])) > com[i]) com[i] = z;
02426       z = x[k+1]-x[k];
02427       for (i = 1; i < mm; ++i)
02428          for (jj = i; jj < mm; ++jj) {
02429             j = mm+i-jj;
02430             c[j-2] = c[j-1]*z+c[j-2];
02431          }
02432       for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
02433       printf("\n");
02434       for (i = 0; i < mm1; ++i)
02435          if (!(k >= n-2 && i != 0))
02436             if ((z = TMath::Abs(c[i]-a[k+1+i*200]))
02437                > diff[i]) diff[i] = z;
02438    }
02439    printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
02440    for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
02441    printf("\n");
02442    printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
02443    if (TMath::Abs(c[0]) > com[0])
02444       com[0] = TMath::Abs(c[0]);
02445    for (i = 0; i < mm1; ++i) printf("%16.8f",com[i]);
02446    printf("\n");
02447 
02448    //     Test of TSpline5 with nonequidistant knots, two double knots,
02449    //        one triple knot,follows.
02450    printf("1         TEST OF TSpline5 WITH NONEQUIDISTANT KNOTS,\n");
02451    printf("             TWO DOUBLE, ONE TRIPLE KNOT\n");
02452    n = 10;
02453    x[0] = 0;
02454    x[1] = 2;
02455    x[2] = 2;
02456    x[3] = 3;
02457    x[4] = 3;
02458    x[5] = 3;
02459    x[6] = 5;
02460    x[7] = 8;
02461    x[8] = 9;
02462    x[9] = 9;
02463    y[0] = 163;
02464    y[1] = 237;
02465    y[2] = -127;
02466    y[3] = 119;
02467    y[4] = -65;
02468    y[5] = 192;
02469    y[6] = 293;
02470    y[7] = 326;
02471    y[8] = 0;
02472    y[9] = -414;
02473    m = 3;
02474    mm = m << 1;
02475    mm1 = mm-1;
02476    printf("-N = %3d    M =%2d\n",n,m);
02477    spline = new TSpline5("Test",x,y,n);
02478    for (i = 0; i < n; ++i)
02479       spline->GetCoeff(i,hx,a[i],a[i+200],a[i+400],
02480                        a[i+600],a[i+800],a[i+1000]);
02481    delete spline;
02482    for (i = 0; i < mm1; ++i)
02483       diff[i] = com[i] = 0;
02484    for (k = 0; k < n; ++k) {
02485       for (i = 0; i < mm; ++i)
02486          c[i] = a[k+i*200];
02487       printf(" ---------------------------------------%3d --------------------------------------------\n",k+1);
02488       printf("%12.8f\n",x[k]);
02489       if (k == n-1) {
02490          printf("%16.8f\n",c[0]);
02491          break;
02492       }
02493       for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
02494       printf("\n");
02495       for (i = 0; i < mm1; ++i)
02496          if ((z=TMath::Abs(a[k+i*200])) > com[i]) com[i] = z;
02497       z = x[k+1]-x[k];
02498       for (i = 1; i < mm; ++i)
02499          for (jj = i; jj < mm; ++jj) {
02500             j = mm+i-jj;
02501             c[j-2] = c[j-1]*z+c[j-2];
02502          }
02503       for (i = 0; i < mm; ++i) printf("%16.8f",c[i]);
02504       printf("\n");
02505       for (i = 0; i < mm1; ++i)
02506          if (!(k >= n-2 && i != 0))
02507             if((z = TMath::Abs(c[i]-a[k+1+i*200]))
02508                > diff[i]) diff[i] = z;
02509    }
02510    printf("  MAXIMUM ABSOLUTE VALUES OF DIFFERENCES \n");
02511    for (i = 0; i < mm1; ++i) printf("%18.9E",diff[i]);
02512    printf("\n");
02513    printf("  MAXIMUM ABSOLUTE VALUES OF COEFFICIENTS \n");
02514    if (TMath::Abs(c[0]) > com[0])
02515       com[0] = TMath::Abs(c[0]);
02516    for (i = 0; i < mm1; ++i) printf("%16.8f",com[i]);
02517    printf("\n");
02518 }
02519 
02520 
02521 //______________________________________________________________________________
02522 void TSpline5::Streamer(TBuffer &R__b)
02523 {
02524    // Stream an object of class TSpline5.
02525 
02526    if (R__b.IsReading()) {
02527       UInt_t R__s, R__c;
02528       Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
02529       if (R__v > 1) {
02530          R__b.ReadClassBuffer(TSpline5::Class(), this, R__v, R__s, R__c);
02531          return;
02532       }
02533       //====process old versions before automatic schema evolution
02534       TSpline::Streamer(R__b);
02535       if (fNp > 0) {
02536          fPoly = new TSplinePoly5[fNp];
02537          for(Int_t i=0; i<fNp; ++i) {
02538             fPoly[i].Streamer(R__b);
02539          }
02540       }
02541       //      R__b >> fPoly;
02542    } else {
02543       R__b.WriteClassBuffer(TSpline5::Class(),this);
02544    }
02545 }

Generated on Tue Jul 5 14:24:22 2011 for ROOT_528-00b_version by  doxygen 1.5.1