GENIEGenerator
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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.\n
11 
12  Update May 15, 2022 IK:
13  Adding as extra interpolators TSpline5 and
14  ROOT::Math::GSLInterpolator (LINEAR, POLYNOMIAL, CSPLINE, CSPLINE_PERIODIC,
15  AKIMA, AKIMA_PERIODIC)
16 
17 
18 \ref [1] GENIE docdb 297
19 
20 \author Costas Andreopoulos <c.andreopoulos \at cern.ch>
21  University of Liverpool \n
22  Igor Kakorin <kakorin@jinr.ru>
23  Joint Institute for Nuclear Research
24 
25 \created May 04, 2004
26 
27 \cpright Copyright (c) 2003-2024, The GENIE Collaboration
28  For the full text of the license visit http://copyright.genie-mc.org
29 */
30 //____________________________________________________________________________
31 
32 #ifndef _SPLINE_H_
33 #define _SPLINE_H_
34 
35 #include <vector>
36 #include <string>
37 #include <fstream>
38 #include <ostream>
39 
40 #include <TObject.h>
41 #include <TSpline.h>
42 #include <Math/Interpolator.h>
43 
44 class TNtupleD;
45 class TTree;
46 class TSQLServer;
47 class TGraph;
48 
49 using std::string;
50 using std::ostream;
51 using std::ofstream;
52 
53 namespace genie {
54 
55 class Spline;
56 ostream & operator << (ostream & stream, const Spline & spl);
57 
58 class Spline : public TObject {
59 
60 public:
61  using TObject::Print; // suppress clang 'hides overloaded virtual function [-Woverloaded-virtual]' warnings
62 
63  // ctors & dtor
64  Spline();
65  Spline(string filename, string xtag="", string ytag="", bool is_xml = false);
66  Spline(TNtupleD * ntuple, string xy, string cut="");
67  Spline(TTree * tree, string xy, string cut="");
68  Spline(TSQLServer * db, string query);
69  Spline(int nentries, double x[], double y[]);
70  Spline(int nentries, float x[], float y[]);
71  Spline(const Spline & spline);
72  Spline(const TSpline3 & spline, int nknots);
73  virtual ~Spline();
74 
75  // Load the Spline from XML, flat ASCII, ROOT ntuple/tree/tspline3, or SQL DB
76  bool LoadFromXmlFile (string filename, string xtag, string ytag);
77  bool LoadFromAsciiFile (string filename);
78  bool LoadFromNtuple (TNtupleD * nt, string xy, string cut = "");
79  bool LoadFromTree (TTree * tr, string xy, string cut = "");
80  bool LoadFromDBase (TSQLServer * db, string query);
81  bool LoadFromTSpline3 (const TSpline3 & spline, int nknots);
82 
83  // Get xmin,xmax,nknots, check x variable against valid range and evaluate spline
84  int NKnots (void) const {return fNKnots;}
85  void GetKnot (int iknot, double & x, double & y) const;
86  double GetKnotX (int iknot) const;
87  double GetKnotY (int iknot) const;
88  double XMin (void) const {return fXMin; }
89  double XMax (void) const {return fXMax; }
90  double YMax (void) const {return fYMax; }
91  double Evaluate (double x) const;
92  bool IsWithinValidRange (double x) const;
93 
94  void SetName (string name) { fName = name; }
95  string Name (void) const { return fName; }
96 
97  void YCanBeNegative(bool tf) { fYCanBeNegative = tf; }
98 
99  // Save the Spline in XML, flat ASCII or ROOT format
100  void SaveAsXml (string filename, string xtag, string ytag, string name="") const;
101  void SaveAsXml (ofstream & str, string xtag, string ytag, string name="") const;
102  void SaveAsText(string filename, string format="%10.6f\t%10.6f") const;
103  void SaveAsROOT(string filename, string name="", bool recreate=false) const;
104 
105  // Export Spline as TGraph or TSpline3
106  TGraph * GetAsTGraph (int np = 500, bool xscaling = false,
107  bool inlog=false, double fx=1., double fy=1.) const;
108  TSpline3 * GetAsTSpline (void) const { return fInterpolator; }
109 
110  // Knot manipulation methods in additions to the TSpline3 ones
111  void FindClosestKnot(double x, double & xknot, double & yknot, Option_t * opt="-+") const;
112  bool ClosestKnotValueIsZero(double x, Option_t * opt="-+") const;
113 
114  // Common mathematical operations applied simultaneously on all spline knots
115  void Add (const Spline & spl, double c=1);
116  void Multiply (const Spline & spl, double c=1);
117  void Divide (const Spline & spl, double c=1);
118  void Add (double a);
119  void Multiply (double a);
120  void Divide (double a);
121  void SetType (string type);
122  string GetType (void) { return fInterpolatorType; }
123 
124  // Print knots
125  void Print(ostream & stream) const;
126 
127  // Overloaded operators
128  friend ostream & operator << (ostream & stream, const Spline & spl);
129 
130 private:
131 
132  // Initialize and build spline
133  void InitSpline (void);
134  void ResetSpline (void);
135  void BuildSpline (int nentries, double x[], double y[]);
136 
137  // Private data members
138  string fName;
139  int fNKnots;
140  double fXMin;
141  double fXMax;
142  double fYMax;
143  TSpline3 * fInterpolator;
145  TSpline5 * fInterpolator5;
146  ROOT::Math::Interpolator * fGSLInterpolator;
148 
149 ClassDef(Spline,2)
150 };
151 
152 }
153 
154 #endif
TGraph * GetAsTGraph(int np=500, bool xscaling=false, bool inlog=false, double fx=1., double fy=1.) const
Definition: Spline.cxx:496
ROOT::Math::Interpolator * fGSLInterpolator
Definition: Spline.h:146
void SetType(string type)
Definition: Spline.cxx:769
A numeric analysis tool class for interpolating 1-D functions.
Definition: Spline.h:58
void Print(ostream &stream) const
Definition: Spline.cxx:569
bool LoadFromAsciiFile(string filename)
Definition: Spline.cxx:242
void BuildSpline(int nentries, double x[], double y[])
Definition: Spline.cxx:750
bool fYCanBeNegative
Definition: Spline.h:144
double Evaluate(double x) const
Definition: Spline.cxx:363
bool ClosestKnotValueIsZero(double x, Option_t *opt="-+") const
Definition: Spline.cxx:561
bool IsWithinValidRange(double x) const
Definition: Spline.cxx:357
double YMax(void) const
Definition: Spline.h:90
TSpline3 * GetAsTSpline(void) const
Definition: Spline.h:108
bool LoadFromXmlFile(string filename, string xtag, string ytag)
Definition: Spline.cxx:155
double GetKnotX(int iknot) const
Definition: Spline.cxx:335
string GetType(void)
Definition: Spline.h:122
string fInterpolatorType
Definition: Spline.h:147
void ResetSpline(void)
Definition: Spline.cxx:742
const double a
double fXMax
Definition: Spline.h:141
void SaveAsROOT(string filename, string name="", bool recreate=false) const
Definition: Spline.cxx:485
void Add(const Spline &spl, double c=1)
Definition: Spline.cxx:587
void Multiply(const Spline &spl, double c=1)
Definition: Spline.cxx:612
void GetKnot(int iknot, double &x, double &y) const
Definition: Spline.cxx:326
double fYMax
Definition: Spline.h:142
double GetKnotY(int iknot) const
Definition: Spline.cxx:346
double XMax(void) const
Definition: Spline.h:89
string fName
Definition: Spline.h:138
bool LoadFromTSpline3(const TSpline3 &spline, int nknots)
Definition: Spline.cxx:305
void SaveAsText(string filename, string format="%10.6f\t%10.6f") const
Definition: Spline.cxx:464
void SaveAsXml(string filename, string xtag, string ytag, string name="") const
Definition: Spline.cxx:419
int fNKnots
Definition: Spline.h:139
TSpline5 * fInterpolator5
Definition: Spline.h:145
TSpline3 * fInterpolator
Definition: Spline.h:143
void SetName(string name)
Definition: Spline.h:94
bool LoadFromNtuple(TNtupleD *nt, string xy, string cut="")
Definition: Spline.cxx:253
bool LoadFromDBase(TSQLServer *db, string query)
Definition: Spline.cxx:299
ostream & operator<<(ostream &stream, const AlgConfigPool &config_pool)
string Name(void) const
Definition: Spline.h:95
void YCanBeNegative(bool tf)
Definition: Spline.h:97
int NKnots(void) const
Definition: Spline.h:84
void InitSpline(void)
Definition: Spline.cxx:723
bool LoadFromTree(TTree *tr, string xy, string cut="")
Definition: Spline.cxx:261
virtual ~Spline()
Definition: Spline.cxx:148
double fXMin
Definition: Spline.h:140
void FindClosestKnot(double x, double &xknot, double &yknot, Option_t *opt="-+") const
Definition: Spline.cxx:534
double XMin(void) const
Definition: Spline.h:88
friend ostream & operator<<(ostream &stream, const Spline &spl)
Definition: Spline.cxx:45
void Divide(const Spline &spl, double c=1)
Definition: Spline.cxx:637