#ifndef ROOT_TSpline
#define ROOT_TSpline
#ifndef ROOT_TH1
#include "TH1.h"
#endif
#ifndef ROOT_TGraph
#include "TGraph.h"
#endif
class TF1;
class TSpline : public TNamed, public TAttLine,
public TAttFill, public TAttMarker
{
protected:
Double_t fDelta;
Double_t fXmin;
Double_t fXmax;
Int_t fNp;
Bool_t fKstep;
TH1F *fHistogram;
TGraph *fGraph;
Int_t fNpx;
TSpline(const TSpline&);
TSpline& operator=(const TSpline&);
virtual void BuildCoeff()=0;
public:
TSpline() : fDelta(-1), fXmin(0), fXmax(0),
fNp(0), fKstep(kFALSE), fHistogram(0), fGraph(0), fNpx(100) {}
TSpline(const char *title, Double_t delta, Double_t xmin,
Double_t xmax, Int_t np, Bool_t step) :
TNamed("Spline",title), TAttFill(0,1),
fDelta(delta), fXmin(xmin),
fXmax(xmax), fNp(np), fKstep(step),
fHistogram(0), fGraph(0), fNpx(100) {}
virtual ~TSpline() {if(fHistogram) delete fHistogram; if(fGraph) delete fGraph;}
virtual void GetKnot(Int_t i, Double_t &x, Double_t &y) const =0;
virtual void Draw(Option_t *option="");
virtual Int_t GetNpx() const {return fNpx;}
virtual void Paint(Option_t *option="");
virtual Double_t Eval(Double_t x) const=0;
virtual void SaveAs(const char * ,Option_t * ) const {;}
void SetNpx(Int_t n) {fNpx=n;}
ClassDef (TSpline,2)
};
class TSplinePoly : public TObject
{
protected:
Double_t fX;
Double_t fY;
public:
TSplinePoly() :
fX(0), fY(0) {}
TSplinePoly(Double_t x, Double_t y) :
fX(x), fY(y) {}
Double_t &X() {return fX;}
Double_t &Y() {return fY;}
void GetKnot(Double_t &x, Double_t &y) const {x=fX; y=fY;}
virtual Double_t Eval(Double_t) const {return fY;}
ClassDef(TSplinePoly,1)
};
class TSplinePoly3 : public TSplinePoly
{
private:
Double_t fB;
Double_t fC;
Double_t fD;
public:
TSplinePoly3() :
fB(0), fC(0), fD(0) {}
TSplinePoly3(Double_t x, Double_t y, Double_t b, Double_t c, Double_t d) :
TSplinePoly(x,y), fB(b), fC(c), fD(d) {}
Double_t &B() {return fB;}
Double_t &C() {return fC;}
Double_t &D() {return fD;}
Double_t Eval(Double_t x) const {
Double_t dx=x-fX;
return (fY+dx*(fB+dx*(fC+dx*fD)));
}
Double_t Derivative(Double_t x) const {
Double_t dx=x-fX;
return (fB+2*fC*dx+3*fD*dx*dx);
}
ClassDef(TSplinePoly3,1)
};
class TSplinePoly5 : public TSplinePoly
{
private:
Double_t fB;
Double_t fC;
Double_t fD;
Double_t fE;
Double_t fF;
public:
TSplinePoly5() :
fB(0), fC(0), fD(0), fE(0), fF(0) {}
TSplinePoly5(Double_t x, Double_t y, Double_t b, Double_t c,
Double_t d, Double_t e, Double_t f) :
TSplinePoly(x,y), fB(b), fC(c), fD(d), fE(e), fF(f) {}
Double_t &B() {return fB;}
Double_t &C() {return fC;}
Double_t &D() {return fD;}
Double_t &E() {return fE;}
Double_t &F() {return fF;}
Double_t Eval(Double_t x) const {
Double_t dx=x-fX;
return (fY+dx*(fB+dx*(fC+dx*(fD+dx*(fE+dx*fF)))));
}
Double_t Derivative(Double_t x) const{
Double_t dx=x-fX;
return (fB+2*fC*dx+3*fD*dx*dx+4*fE*dx*dx*dx+5*fF*dx*dx*dx*dx);
}
ClassDef(TSplinePoly5,1)
};
class TSpline3 : public TSpline
{
private:
TSplinePoly3 *fPoly;
Double_t fValBeg;
Double_t fValEnd;
Int_t fBegCond;
Int_t fEndCond;
void BuildCoeff();
void SetCond(const char *opt);
protected:
TSpline3(const TSpline3&);
TSpline3& operator=(const TSpline3&);
public:
TSpline3() : fPoly(0), fValBeg(0), fValEnd(0),
fBegCond(-1), fEndCond(-1) {}
TSpline3(const char *title,
Double_t x[], Double_t y[], Int_t n, const char *opt=0,
Double_t valbeg=0, Double_t valend=0);
TSpline3(const char *title,
Double_t xmin, Double_t xmax,
Double_t y[], Int_t n, const char *opt=0,
Double_t valbeg=0, Double_t valend=0);
TSpline3(const char *title,
Double_t x[], const TF1 *func, Int_t n, const char *opt=0,
Double_t valbeg=0, Double_t valend=0);
TSpline3(const char *title,
Double_t xmin, Double_t xmax,
const TF1 *func, Int_t n, const char *opt=0,
Double_t valbeg=0, Double_t valend=0);
TSpline3(const char *title,
const TGraph *g, const char *opt=0,
Double_t valbeg=0, Double_t valend=0);
Int_t FindX(Double_t x) const;
Double_t Eval(Double_t x) const;
Double_t Derivative(Double_t x) const;
virtual ~TSpline3() {if (fPoly) delete [] fPoly;}
void GetCoeff(Int_t i, Double_t &x, Double_t &y, Double_t &b,
Double_t &c, Double_t &d) {x=fPoly[i].X();y=fPoly[i].Y();
b=fPoly[i].B();c=fPoly[i].C();d=fPoly[i].D();}
void GetKnot(Int_t i, Double_t &x, Double_t &y) const
{x=fPoly[i].X(); y=fPoly[i].Y();}
virtual void SaveAs(const char *filename,Option_t *option="") const;
static void Test();
ClassDef (TSpline3,2)
};
class TSpline5 : public TSpline
{
private:
TSplinePoly5 *fPoly;
void BuildCoeff();
void BoundaryConditions(const char *opt, Int_t &beg, Int_t &end,
const char *&cb1, const char *&ce1, const char *&cb2,
const char *&ce2);
void SetBoundaries(Double_t b1, Double_t e1, Double_t b2, Double_t e2,
const char *cb1, const char *ce1, const char *cb2,
const char *ce2);
protected:
TSpline5(const TSpline5&);
TSpline5& operator=(const TSpline5&);
public:
TSpline5() : fPoly(0) {}
TSpline5(const char *title,
Double_t x[], Double_t y[], Int_t n,
const char *opt=0, Double_t b1=0, Double_t e1=0,
Double_t b2=0, Double_t e2=0);
TSpline5(const char *title,
Double_t xmin, Double_t xmax,
Double_t y[], Int_t n,
const char *opt=0, Double_t b1=0, Double_t e1=0,
Double_t b2=0, Double_t e2=0);
TSpline5(const char *title,
Double_t x[], const TF1 *func, Int_t n,
const char *opt=0, Double_t b1=0, Double_t e1=0,
Double_t b2=0, Double_t e2=0);
TSpline5(const char *title,
Double_t xmin, Double_t xmax,
const TF1 *func, Int_t n,
const char *opt=0, Double_t b1=0, Double_t e1=0,
Double_t b2=0, Double_t e2=0);
TSpline5(const char *title,
const TGraph *g,
const char *opt=0, Double_t b1=0, Double_t e1=0,
Double_t b2=0, Double_t e2=0);
Int_t FindX(Double_t x) const;
Double_t Eval(Double_t x) const;
Double_t Derivative(Double_t x) const;
~TSpline5() {if (fPoly) delete [] fPoly;}
void GetCoeff(Int_t i, Double_t &x, Double_t &y, Double_t &b,
Double_t &c, Double_t &d, Double_t &e, Double_t &f)
{x=fPoly[i].X();y=fPoly[i].Y();b=fPoly[i].B();
c=fPoly[i].C();d=fPoly[i].D();
e=fPoly[i].E();f=fPoly[i].F();}
void GetKnot(Int_t i, Double_t &x, Double_t &y) const
{x=fPoly[i].X(); y=fPoly[i].Y();}
virtual void SaveAs(const char *filename,Option_t *option="") const;
static void Test();
ClassDef (TSpline5,2)
};
#endif
ROOT page - Class index - Class Hierarchy - Top of the page
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.