Spline.h
Go to the documentation of this file.
1 //____________________________________________________________________________
2 /*!
3 
4 \class genie::Spline
5 
6 \brief A numeric analysis tool class for interpolating 1-D functions.
7 
8  Uses ROOT's TSpline3 for the actual interpolation and can retrieve
9  function (x,y(x)) pairs from an XML file, a flat ascii file, a
10  TNtuple, a TTree or an SQL database.
11 
12 \author Costas Andreopoulos <constantinos.andreopoulos \at cern.ch>
13  University of Liverpool & STFC Rutherford Appleton Laboratory
14 
15 \created May 04, 2004
16 
17 \cpright Copyright (c) 2003-2020, The GENIE Collaboration
18  For the full text of the license visit http://copyright.genie-mc.org
19 */
20 //____________________________________________________________________________
21 
22 #ifndef _SPLINE_H_
23 #define _SPLINE_H_
24 
25 #include <string>
26 #include <fstream>
27 #include <ostream>
28 
29 #include <TObject.h>
30 #include <TSpline.h>
31 
32 class TNtupleD;
33 class TTree;
34 class TSQLServer;
35 class TGraph;
36 
37 using std::string;
38 using std::ostream;
39 using std::ofstream;
40 
41 namespace genie {
42 
43 class Spline;
44 ostream & operator << (ostream & stream, const Spline & spl);
45 
46 class Spline : public TObject {
47 
48 public:
49  using TObject::Print; // suppress clang 'hides overloaded virtual function [-Woverloaded-virtual]' warnings
50 
51  // ctors & dtor
52  Spline();
53  Spline(string filename, string xtag="", string ytag="", bool is_xml = false);
54  Spline(TNtupleD * ntuple, string xy, string cut="");
55  Spline(TTree * tree, string xy, string cut="");
56  Spline(TSQLServer * db, string query);
57  Spline(int nentries, double x[], double y[]);
58  Spline(int nentries, float x[], float y[]);
59  Spline(const Spline & spline);
60  Spline(const TSpline3 & spline, int nknots);
61  virtual ~Spline();
62 
63  // Load the Spline from XML, flat ASCII, ROOT ntuple/tree/tspline3, or SQL DB
64  bool LoadFromXmlFile (string filename, string xtag, string ytag);
65  bool LoadFromAsciiFile (string filename);
66  bool LoadFromNtuple (TNtupleD * nt, string xy, string cut = "");
67  bool LoadFromTree (TTree * tr, string xy, string cut = "");
68  bool LoadFromDBase (TSQLServer * db, string query);
69  bool LoadFromTSpline3 (const TSpline3 & spline, int nknots);
70 
71  // Get xmin,xmax,nknots, check x variable against valid range and evaluate spline
72  int NKnots (void) const {return fNKnots;}
73  void GetKnot (int iknot, double & x, double & y) const;
74  double GetKnotX (int iknot) const;
75  double GetKnotY (int iknot) const;
76  double XMin (void) const {return fXMin; }
77  double XMax (void) const {return fXMax; }
78  double YMax (void) const {return fYMax; }
79  double Evaluate (double x) const;
80  bool IsWithinValidRange (double x) const;
81 
82  void SetName (string name) { fName = name; }
83  string Name (void) const { return fName; }
84 
85  void YCanBeNegative(bool tf) { fYCanBeNegative = tf; }
86 
87  // Save the Spline in XML, flat ASCII or ROOT format
88  void SaveAsXml (string filename, string xtag, string ytag, string name="") const;
89  void SaveAsXml (ofstream & str, string xtag, string ytag, string name="") const;
90  void SaveAsText(string filename, string format="%10.6f\t%10.6f") const;
91  void SaveAsROOT(string filename, string name="", bool recreate=false) const;
92 
93  // Export Spline as TGraph or TSpline3
94  TGraph * GetAsTGraph (int np = 500, bool xscaling = false,
95  bool inlog=false, double fx=1., double fy=1.) const;
96  TSpline3 * GetAsTSpline (void) const { return fInterpolator; }
97 
98  // Knot manipulation methods in additions to the TSpline3 ones
99  void FindClosestKnot(double x, double & xknot, double & yknot, Option_t * opt="-+") const;
100  bool ClosestKnotValueIsZero(double x, Option_t * opt="-+") const;
101 
102  // Common mathematical operations applied simultaneously on all spline knots
103  void Add (const Spline & spl, double c=1);
104  void Multiply (const Spline & spl, double c=1);
105  void Divide (const Spline & spl, double c=1);
106  void Add (double a);
107  void Multiply (double a);
108  void Divide (double a);
109 
110  // Print knots
111  void Print(ostream & stream) const;
112 
113  // Overloaded operators
114  friend ostream & operator << (ostream & stream, const Spline & spl);
115 
116 private:
117 
118  // Initialize and build spline
119  void InitSpline (void);
120  void ResetSpline (void);
121  void BuildSpline (int nentries, double x[], double y[]);
122 
123  // Private data members
124  string fName;
125  int fNKnots;
126  double fXMin;
127  double fXMax;
128  double fYMax;
129  TSpline3 * fInterpolator;
131 
132 ClassDef(Spline,1)
133 };
134 
135 }
136 
137 #endif
static QCString name
Definition: declinfo.cpp:673
TGraph * GetAsTGraph(int np=500, bool xscaling=false, bool inlog=false, double fx=1., double fy=1.) const
Definition: Spline.cxx:489
THE MAIN GENIE PROJECT NAMESPACE
Definition: AlgCmp.h:25
std::string string
Definition: nybbler.cc:12
A numeric analysis tool class for interpolating 1-D functions.
Definition: Spline.h:46
opt
Definition: train.py:196
static bool format(QChar::Decomposition tag, QString &str, int index, int len)
Definition: qstring.cpp:11496
void Print(ostream &stream) const
Definition: Spline.cxx:562
bool LoadFromAsciiFile(string filename)
Definition: Spline.cxx:240
void BuildSpline(int nentries, double x[], double y[])
Definition: Spline.cxx:738
bool fYCanBeNegative
Definition: Spline.h:130
Definition: tf_graph.h:23
double Evaluate(double x) const
Definition: Spline.cxx:361
string filename
Definition: train.py:213
bool ClosestKnotValueIsZero(double x, Option_t *opt="-+") const
Definition: Spline.cxx:554
bool IsWithinValidRange(double x) const
Definition: Spline.cxx:355
double YMax(void) const
Definition: Spline.h:78
TSpline3 * GetAsTSpline(void) const
Definition: Spline.h:96
bool LoadFromXmlFile(string filename, string xtag, string ytag)
Definition: Spline.cxx:153
double GetKnotX(int iknot) const
Definition: Spline.cxx:333
void ResetSpline(void)
Definition: Spline.cxx:732
const double a
double fXMax
Definition: Spline.h:127
void SaveAsROOT(string filename, string name="", bool recreate=false) const
Definition: Spline.cxx:478
void Add(const Spline &spl, double c=1)
Definition: Spline.cxx:580
void Multiply(const Spline &spl, double c=1)
Definition: Spline.cxx:605
void GetKnot(int iknot, double &x, double &y) const
Definition: Spline.cxx:324
double fYMax
Definition: Spline.h:128
double GetKnotY(int iknot) const
Definition: Spline.cxx:344
double XMax(void) const
Definition: Spline.h:77
string fName
Definition: Spline.h:124
bool LoadFromTSpline3(const TSpline3 &spline, int nknots)
Definition: Spline.cxx:303
void SaveAsText(string filename, string format="%10.6f\t%10.6f") const
Definition: Spline.cxx:457
void SaveAsXml(string filename, string xtag, string ytag, string name="") const
Definition: Spline.cxx:412
int fNKnots
Definition: Spline.h:125
TSpline3 * fInterpolator
Definition: Spline.h:129
void SetName(string name)
Definition: Spline.h:82
bool LoadFromNtuple(TNtupleD *nt, string xy, string cut="")
Definition: Spline.cxx:251
bool LoadFromDBase(TSQLServer *db, string query)
Definition: Spline.cxx:297
ostream & operator<<(ostream &stream, const AlgConfigPool &config_pool)
string Name(void) const
Definition: Spline.h:83
void YCanBeNegative(bool tf)
Definition: Spline.h:85
int NKnots(void) const
Definition: Spline.h:72
void InitSpline(void)
Definition: Spline.cxx:716
bool LoadFromTree(TTree *tr, string xy, string cut="")
Definition: Spline.cxx:259
list x
Definition: train.py:276
virtual ~Spline()
Definition: Spline.cxx:148
double fXMin
Definition: Spline.h:126
void FindClosestKnot(double x, double &xknot, double &yknot, Option_t *opt="-+") const
Definition: Spline.cxx:527
double XMin(void) const
Definition: Spline.h:76
friend ostream & operator<<(ostream &stream, const Spline &spl)
query_result< Args... > query(sqlite3 *db, std::string const &ddl)
Definition: select.h:75
static QCString str
void Divide(const Spline &spl, double c=1)
Definition: Spline.cxx:630