Public Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
genie::Spline Class Reference

A numeric analysis tool class for interpolating 1-D functions. More...

#include <Spline.h>

Inheritance diagram for genie::Spline:

Public Member Functions

 Spline ()
 
 Spline (string filename, string xtag="", string ytag="", bool is_xml=false)
 
 Spline (TNtupleD *ntuple, string xy, string cut="")
 
 Spline (TTree *tree, string xy, string cut="")
 
 Spline (TSQLServer *db, string query)
 
 Spline (int nentries, double x[], double y[])
 
 Spline (int nentries, float x[], float y[])
 
 Spline (const Spline &spline)
 
 Spline (const TSpline3 &spline, int nknots)
 
virtual ~Spline ()
 
bool LoadFromXmlFile (string filename, string xtag, string ytag)
 
bool LoadFromAsciiFile (string filename)
 
bool LoadFromNtuple (TNtupleD *nt, string xy, string cut="")
 
bool LoadFromTree (TTree *tr, string xy, string cut="")
 
bool LoadFromDBase (TSQLServer *db, string query)
 
bool LoadFromTSpline3 (const TSpline3 &spline, int nknots)
 
int NKnots (void) const
 
void GetKnot (int iknot, double &x, double &y) const
 
double GetKnotX (int iknot) const
 
double GetKnotY (int iknot) const
 
double XMin (void) const
 
double XMax (void) const
 
double YMax (void) const
 
double Evaluate (double x) const
 
bool IsWithinValidRange (double x) const
 
void SetName (string name)
 
string Name (void) const
 
void YCanBeNegative (bool tf)
 
void SaveAsXml (string filename, string xtag, string ytag, string name="") const
 
void SaveAsXml (ofstream &str, string xtag, string ytag, string name="") const
 
void SaveAsText (string filename, string format="%10.6f\t%10.6f") const
 
void SaveAsROOT (string filename, string name="", bool recreate=false) const
 
TGraph * GetAsTGraph (int np=500, bool xscaling=false, bool inlog=false, double fx=1., double fy=1.) const
 
TSpline3 * GetAsTSpline (void) const
 
void FindClosestKnot (double x, double &xknot, double &yknot, Option_t *opt="-+") const
 
bool ClosestKnotValueIsZero (double x, Option_t *opt="-+") const
 
void Add (const Spline &spl, double c=1)
 
void Multiply (const Spline &spl, double c=1)
 
void Divide (const Spline &spl, double c=1)
 
void Add (double a)
 
void Multiply (double a)
 
void Divide (double a)
 
void Print (ostream &stream) const
 

Private Member Functions

void InitSpline (void)
 
void ResetSpline (void)
 
void BuildSpline (int nentries, double x[], double y[])
 

Private Attributes

string fName
 
int fNKnots
 
double fXMin
 
double fXMax
 
double fYMax
 
TSpline3 * fInterpolator
 
bool fYCanBeNegative
 

Friends

ostream & operator<< (ostream &stream, const Spline &spl)
 

Detailed Description

A numeric analysis tool class for interpolating 1-D functions.

Uses ROOT's TSpline3 for the actual interpolation and can retrieve function (x,y(x)) pairs from an XML file, a flat ascii file, a TNtuple, a TTree or an SQL database.

Author
Costas Andreopoulos <constantinos.andreopoulos cern.ch> University of Liverpool & STFC Rutherford Appleton Laboratory

May 04, 2004

Copyright (c) 2003-2020, The GENIE Collaboration For the full text of the license visit http://copyright.genie-mc.org

Definition at line 46 of file Spline.h.

Constructor & Destructor Documentation

Spline::Spline ( )

Definition at line 52 of file Spline.cxx.

53 {
54  this->InitSpline();
55 }
void InitSpline(void)
Definition: Spline.cxx:716
Spline::Spline ( string  filename,
string  xtag = "",
string  ytag = "",
bool  is_xml = false 
)

Definition at line 57 of file Spline.cxx.

57  :
58 TObject()
59 {
60  string fmt = (is_xml) ? "XML" : "ASCII";
61 
62  LOG("Spline", pDEBUG)
63  << "Constructing spline from data in " << fmt << " file: " << filename;
64 
65  this->InitSpline();
66 
67  if(is_xml)
68  this->LoadFromXmlFile(filename, xtag, ytag);
69  else
70  this->LoadFromAsciiFile(filename);
71 }
bool LoadFromAsciiFile(string filename)
Definition: Spline.cxx:240
string filename
Definition: train.py:213
bool LoadFromXmlFile(string filename, string xtag, string ytag)
Definition: Spline.cxx:153
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void InitSpline(void)
Definition: Spline.cxx:716
#define pDEBUG
Definition: Messenger.h:63
Spline::Spline ( TNtupleD *  ntuple,
string  xy,
string  cut = "" 
)

Definition at line 73 of file Spline.cxx.

73  :
74 TObject()
75 {
76  LOG("Spline", pDEBUG) << "Constructing spline from data in a TNtuple";
77 
78  this->InitSpline();
79  this->LoadFromNtuple(ntuple, var, cut);
80 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool LoadFromNtuple(TNtupleD *nt, string xy, string cut="")
Definition: Spline.cxx:251
int var
Definition: 018_def.c:9
void InitSpline(void)
Definition: Spline.cxx:716
#define pDEBUG
Definition: Messenger.h:63
Spline::Spline ( TTree *  tree,
string  xy,
string  cut = "" 
)

Definition at line 82 of file Spline.cxx.

82  :
83 TObject()
84 {
85  LOG("Spline", pDEBUG) << "Constructing spline from data in a TTree";
86 
87  this->InitSpline();
88  this->LoadFromTree(tree, var, cut);
89 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
int var
Definition: 018_def.c:9
void InitSpline(void)
Definition: Spline.cxx:716
bool LoadFromTree(TTree *tr, string xy, string cut="")
Definition: Spline.cxx:259
#define pDEBUG
Definition: Messenger.h:63
Spline::Spline ( TSQLServer *  db,
string  query 
)

Definition at line 91 of file Spline.cxx.

91  :
92 TObject()
93 {
94  LOG("Spline", pDEBUG) << "Constructing spline from data in a MySQL server";
95 
96  this->InitSpline();
97  this->LoadFromDBase(db, query);
98 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool LoadFromDBase(TSQLServer *db, string query)
Definition: Spline.cxx:297
void InitSpline(void)
Definition: Spline.cxx:716
query_result< Args... > query(sqlite3 *db, std::string const &ddl)
Definition: select.h:75
#define pDEBUG
Definition: Messenger.h:63
Spline::Spline ( int  nentries,
double  x[],
double  y[] 
)

Definition at line 100 of file Spline.cxx.

100  :
101 TObject()
102 {
103  LOG("Spline", pDEBUG)
104  << "Constructing spline from the arrays passed to the ctor";
105 
106  this->InitSpline();
107  this->BuildSpline(nentries, x, y);
108 }
void BuildSpline(int nentries, double x[], double y[])
Definition: Spline.cxx:738
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void InitSpline(void)
Definition: Spline.cxx:716
list x
Definition: train.py:276
#define pDEBUG
Definition: Messenger.h:63
Spline::Spline ( int  nentries,
float  x[],
float  y[] 
)

Definition at line 110 of file Spline.cxx.

110  :
111 TObject()
112 {
113  LOG("Spline", pDEBUG)
114  << "Constructing spline from the arrays passed to the ctor";
115 
116  double * dblx = new double[nentries];
117  double * dbly = new double[nentries];
118 
119  for(int i = 0; i < nentries; i++) {
120  dblx[i] = double( x[i] );
121  dbly[i] = double( y[i] );
122  }
123 
124  this->InitSpline();
125  this->BuildSpline(nentries, dblx, dbly);
126 
127  delete [] dblx;
128  delete [] dbly;
129 }
void BuildSpline(int nentries, double x[], double y[])
Definition: Spline.cxx:738
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void InitSpline(void)
Definition: Spline.cxx:716
list x
Definition: train.py:276
#define pDEBUG
Definition: Messenger.h:63
Spline::Spline ( const Spline spline)

Definition at line 131 of file Spline.cxx.

131  :
132  TObject(), fInterpolator(0)
133 {
134  LOG("Spline", pDEBUG) << "Spline copy constructor";
135 
136  this->LoadFromTSpline3( *spline.GetAsTSpline(), spline.NKnots() );
137 }
TSpline3 * GetAsTSpline(void) const
Definition: Spline.h:96
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool LoadFromTSpline3(const TSpline3 &spline, int nknots)
Definition: Spline.cxx:303
TSpline3 * fInterpolator
Definition: Spline.h:129
int NKnots(void) const
Definition: Spline.h:72
#define pDEBUG
Definition: Messenger.h:63
Spline::Spline ( const TSpline3 &  spline,
int  nknots 
)

Definition at line 139 of file Spline.cxx.

139  :
140  TObject(), fInterpolator(0)
141 {
142  LOG("Spline", pDEBUG)
143  << "Constructing spline from the input TSpline3 object";
144 
145  this->LoadFromTSpline3( spline, nknots );
146 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool LoadFromTSpline3(const TSpline3 &spline, int nknots)
Definition: Spline.cxx:303
TSpline3 * fInterpolator
Definition: Spline.h:129
#define pDEBUG
Definition: Messenger.h:63
Spline::~Spline ( )
virtual

Definition at line 148 of file Spline.cxx.

149 {
150  if(fInterpolator) delete fInterpolator;
151 }
TSpline3 * fInterpolator
Definition: Spline.h:129

Member Function Documentation

void Spline::Add ( const Spline spl,
double  c = 1 
)

Definition at line 580 of file Spline.cxx.

581 {
582  // continue only if the input spline is defined at a wider x-range
583  double xmin = this->XMin();
584  double xmax = this->XMax();
585  bool ok = spl.IsWithinValidRange(xmin) && spl.IsWithinValidRange(xmax);
586  if(!ok) {
587  LOG("Spline", pERROR) << "** Can't add splines. X-range mismatch!";
588  return;
589  }
590 
591  int nknots = this->NKnots();
592  double * x = new double[nknots];
593  double * y = new double[nknots];
594 
595  for(int i=0; i<nknots; i++) {
596  this->GetKnot(i,x[i],y[i]);
597  y[i] += (c * spl.Evaluate(x[i]));
598  }
599  this->ResetSpline();
600  this->BuildSpline(nknots,x,y);
601  delete [] x;
602  delete [] y;
603 }
#define pERROR
Definition: Messenger.h:59
void BuildSpline(int nentries, double x[], double y[])
Definition: Spline.cxx:738
double Evaluate(double x) const
Definition: Spline.cxx:361
bool IsWithinValidRange(double x) const
Definition: Spline.cxx:355
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void ResetSpline(void)
Definition: Spline.cxx:732
void GetKnot(int iknot, double &x, double &y) const
Definition: Spline.cxx:324
double XMax(void) const
Definition: Spline.h:77
int NKnots(void) const
Definition: Spline.h:72
list x
Definition: train.py:276
double XMin(void) const
Definition: Spline.h:76
void Spline::Add ( double  a)

Definition at line 663 of file Spline.cxx.

664 {
665  int nknots = this->NKnots();
666  double * x = new double[nknots];
667  double * y = new double[nknots];
668 
669  for(int i=0; i<nknots; i++) {
670  this->GetKnot(i,x[i],y[i]);
671  y[i]+=a;
672  }
673  this->ResetSpline();
674  this->BuildSpline(nknots,x,y);
675  delete [] x;
676  delete [] y;
677 }
void BuildSpline(int nentries, double x[], double y[])
Definition: Spline.cxx:738
void ResetSpline(void)
Definition: Spline.cxx:732
const double a
void GetKnot(int iknot, double &x, double &y) const
Definition: Spline.cxx:324
int NKnots(void) const
Definition: Spline.h:72
list x
Definition: train.py:276
void Spline::BuildSpline ( int  nentries,
double  x[],
double  y[] 
)
private

Definition at line 738 of file Spline.cxx.

739 {
740  LOG("Spline", pDEBUG) << "Building spline...";
741 
742  double xmin = x[0]; // minimum x in spline
743  double xmax = x[nentries-1]; // maximum x in spline
744 
745  fNKnots = nentries;
746  fXMin = xmin;
747  fXMax = xmax;
748  fYMax = y[ TMath::LocMax(nentries, y) ]; // maximum y in spline
749 
750  if(fInterpolator) delete fInterpolator;
751 
752  fInterpolator = new TSpline3("spl3", x, y, nentries, "0");
753 
754  LOG("Spline", pDEBUG) << "...done building spline";
755 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
double fXMax
Definition: Spline.h:127
double fYMax
Definition: Spline.h:128
int fNKnots
Definition: Spline.h:125
TSpline3 * fInterpolator
Definition: Spline.h:129
list x
Definition: train.py:276
double fXMin
Definition: Spline.h:126
#define pDEBUG
Definition: Messenger.h:63
bool Spline::ClosestKnotValueIsZero ( double  x,
Option_t *  opt = "-+" 
) const

Definition at line 554 of file Spline.cxx.

555 {
556  double xknot = 0, yknot = 0;
557  this->FindClosestKnot(x, xknot, yknot, opt);
558  if(utils::math::AreEqual(yknot,0)) return true;
559  return false;
560 }
bool AreEqual(double x1, double x2)
Definition: MathUtils.cxx:236
opt
Definition: train.py:196
list x
Definition: train.py:276
void FindClosestKnot(double x, double &xknot, double &yknot, Option_t *opt="-+") const
Definition: Spline.cxx:527
void Spline::Divide ( const Spline spl,
double  c = 1 
)

Definition at line 630 of file Spline.cxx.

631 {
632  // continue only if the input spline is defined at a wider x-range
633  double xmin = this->XMin();
634  double xmax = this->XMax();
635  bool ok = spl.IsWithinValidRange(xmin) && spl.IsWithinValidRange(xmax);
636  if(!ok) {
637  LOG("Spline", pERROR) << "** Can't divide splines. X-range mismatch!";
638  return;
639  }
640 
641  int nknots = this->NKnots();
642  double * x = new double[nknots];
643  double * y = new double[nknots];
644 
645  for(int i=0; i<nknots; i++) {
646  this->GetKnot(i,x[i],y[i]);
647  double denom = c * spl.Evaluate(x[i]);
648  bool denom_is_zero = TMath::Abs(denom) < DBL_EPSILON;
649  if(denom_is_zero) {
650  LOG("Spline", pERROR) << "** Refusing to divide spline knot by 0";
651  delete [] x;
652  delete [] y;
653  return;
654  }
655  y[i] /= denom;
656  }
657  this->ResetSpline();
658  this->BuildSpline(nknots,x,y);
659  delete [] x;
660  delete [] y;
661 }
#define pERROR
Definition: Messenger.h:59
void BuildSpline(int nentries, double x[], double y[])
Definition: Spline.cxx:738
double Evaluate(double x) const
Definition: Spline.cxx:361
bool IsWithinValidRange(double x) const
Definition: Spline.cxx:355
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void ResetSpline(void)
Definition: Spline.cxx:732
void GetKnot(int iknot, double &x, double &y) const
Definition: Spline.cxx:324
double XMax(void) const
Definition: Spline.h:77
int NKnots(void) const
Definition: Spline.h:72
list x
Definition: train.py:276
double XMin(void) const
Definition: Spline.h:76
void Spline::Divide ( double  a)

Definition at line 695 of file Spline.cxx.

696 {
697  bool a_is_zero = TMath::Abs(a) < DBL_EPSILON;
698  if(a_is_zero==0) {
699  LOG("Spline", pERROR) << "** Refusing to divide spline by 0";
700  return;
701  }
702  int nknots = this->NKnots();
703  double * x = new double[nknots];
704  double * y = new double[nknots];
705 
706  for(int i=0; i<nknots; i++) {
707  this->GetKnot(i,x[i],y[i]);
708  y[i]/=a;
709  }
710  this->ResetSpline();
711  this->BuildSpline(nknots,x,y);
712  delete [] x;
713  delete [] y;
714 }
#define pERROR
Definition: Messenger.h:59
void BuildSpline(int nentries, double x[], double y[])
Definition: Spline.cxx:738
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void ResetSpline(void)
Definition: Spline.cxx:732
const double a
void GetKnot(int iknot, double &x, double &y) const
Definition: Spline.cxx:324
int NKnots(void) const
Definition: Spline.h:72
list x
Definition: train.py:276
double Spline::Evaluate ( double  x) const

Definition at line 361 of file Spline.cxx.

362 {
363  LOG("Spline", pDEBUG) << "Evaluating spline at point x = " << x;
364  assert(!TMath::IsNaN(x));
365 
366  double y = 0;
367  if( this->IsWithinValidRange(x) ) {
368 
369  // we can interpolate within the range of spline knots - be careful with
370  // strange cubic spline behaviour when close to knots with y=0
371  bool is0p = this->ClosestKnotValueIsZero(x, "+");
372  bool is0n = this->ClosestKnotValueIsZero(x, "-");
373 
374  if(!is0p && !is0n) {
375  // both knots (on the left and right are non-zero) - just interpolate
376  LOG("Spline", pDEBUG) << "Point is between non-zero knots";
377  y = fInterpolator->Eval(x);
378  } else {
379  // at least one of the neighboring knots has y=0
380  if(is0p && is0n) {
381  // both neighboring knots have y=0
382  LOG("Spline", pDEBUG) << "Point is between zero knots";
383  y=0;
384  } else {
385  // just 1 neighboring knot has y=0 - do a linear interpolation
386  LOG("Spline", pDEBUG)
387  << "Point has zero" << (is0n ? " left " : " right ") << "knot";
388  double xpknot=0, ypknot=0, xnknot=0, ynknot=0;
389  this->FindClosestKnot(x, xnknot, ynknot, "-");
390  this->FindClosestKnot(x, xpknot, ypknot, "+");
391  if(is0n) y = ypknot * (x-xnknot)/(xpknot-xnknot);
392  else y = ynknot * (x-xnknot)/(xpknot-xnknot);
393  }
394  }
395 
396  } else {
397  LOG("Spline", pDEBUG) << "x = " << x
398  << " is not within spline range [" << fXMin << ", " << fXMax << "]";
399  }
400 
401  if(y<0 && !fYCanBeNegative) {
402  LOG("Spline", pINFO) << "Negative y (" << y << ")";
403  LOG("Spline", pINFO) << "x = " << x;
404  LOG("Spline", pINFO) << "spline range [" << fXMin << ", " << fXMax << "]";
405  }
406 
407  LOG("Spline", pDEBUG) << "Spline(x = " << x << ") = " << y;
408 
409  return y;
410 }
bool fYCanBeNegative
Definition: Spline.h:130
bool ClosestKnotValueIsZero(double x, Option_t *opt="-+") const
Definition: Spline.cxx:554
bool IsWithinValidRange(double x) const
Definition: Spline.cxx:355
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
double fXMax
Definition: Spline.h:127
#define pINFO
Definition: Messenger.h:62
TSpline3 * fInterpolator
Definition: Spline.h:129
list x
Definition: train.py:276
double fXMin
Definition: Spline.h:126
void FindClosestKnot(double x, double &xknot, double &yknot, Option_t *opt="-+") const
Definition: Spline.cxx:527
#define pDEBUG
Definition: Messenger.h:63
void Spline::FindClosestKnot ( double  x,
double &  xknot,
double &  yknot,
Option_t *  opt = "-+" 
) const

Definition at line 527 of file Spline.cxx.

529 {
530  string option(opt);
531 
532  bool pos = (option.find("+") != string::npos);
533  bool neg = (option.find("-") != string::npos);
534 
535  if(!pos && !neg) return;
536 
537  int iknot = fInterpolator->FindX(x);
538 
539  double xp=0, yp=0, xn=0, yn=0;
540  fInterpolator->GetKnot(iknot, xn,yn);
541  fInterpolator->GetKnot(iknot+1,xp,yp);
542 
543  bool p = (TMath::Abs(x-xp) < TMath::Abs(x-xn));
544 
545  if(pos&&neg) {
546  if(p) { xknot = xp; yknot = yp; }
547  else { xknot = xn; yknot = yn; }
548  } else {
549  if(pos) { xknot = xp; yknot = yp; }
550  if(neg) { xknot = xn; yknot = yn; }
551  }
552 }
opt
Definition: train.py:196
p
Definition: test.py:223
TSpline3 * fInterpolator
Definition: Spline.h:129
list x
Definition: train.py:276
TGraph * Spline::GetAsTGraph ( int  np = 500,
bool  xscaling = false,
bool  inlog = false,
double  fx = 1.,
double  fy = 1. 
) const

Definition at line 489 of file Spline.cxx.

491 {
492  double xmin = fXMin;
493  double xmax = fXMax;
494 
495  np = TMath::Max(np,2);
496 
497  bool use_log = in_log && (fXMin>0) && (fXMax>0);
498 
499  if(use_log) {
500  xmin = TMath::Log10(xmin);
501  xmax = TMath::Log10(xmax);
502  }
503 
504  double dx = (xmax - xmin) / (np-1);
505 
506  double * x = new double[np];
507  double * y = new double[np];
508 
509  for(int i=0; i<np; i++) {
510  x[i] = ( (use_log) ? TMath::Power(10, xmin+i*dx) : xmin + i*dx );
511  y[i] = this->Evaluate( x[i] );
512 
513  // scale with x if needed
514  if (scale_with_x) y[i] /= x[i];
515 
516  // apply x,y scaling
517  y[i] *= fy;
518  x[i] *= fx;
519  }
520 
521  TGraph * graph = new TGraph(np, x, y);
522  delete[] x;
523  delete[] y;
524  return graph;
525 }
def graph(desc, maker=maker)
Definition: apa.py:294
double Evaluate(double x) const
Definition: Spline.cxx:361
double fXMax
Definition: Spline.h:127
list x
Definition: train.py:276
double fXMin
Definition: Spline.h:126
TSpline3* genie::Spline::GetAsTSpline ( void  ) const
inline

Definition at line 96 of file Spline.h.

96 { return fInterpolator; }
TSpline3 * fInterpolator
Definition: Spline.h:129
void Spline::GetKnot ( int  iknot,
double &  x,
double &  y 
) const

Definition at line 324 of file Spline.cxx.

325 {
326  if(!fInterpolator) {
327  LOG("Spline", pWARN) << "Spline has not been built yet!";
328  return;
329  }
330  fInterpolator->GetKnot(iknot,x,y);
331 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pWARN
Definition: Messenger.h:60
TSpline3 * fInterpolator
Definition: Spline.h:129
list x
Definition: train.py:276
double Spline::GetKnotX ( int  iknot) const

Definition at line 333 of file Spline.cxx.

334 {
335  if(!fInterpolator) {
336  LOG("Spline", pWARN) << "Spline has not been built yet!";
337  return 0;
338  }
339  double x,y;
340  fInterpolator->GetKnot(iknot,x,y);
341  return x;
342 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pWARN
Definition: Messenger.h:60
TSpline3 * fInterpolator
Definition: Spline.h:129
list x
Definition: train.py:276
double Spline::GetKnotY ( int  iknot) const

Definition at line 344 of file Spline.cxx.

345 {
346  if(!fInterpolator) {
347  LOG("Spline", pWARN) << "Spline has not been built yet!";
348  return 0;
349  }
350  double x,y;
351  fInterpolator->GetKnot(iknot,x,y);
352  return y;
353 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pWARN
Definition: Messenger.h:60
TSpline3 * fInterpolator
Definition: Spline.h:129
list x
Definition: train.py:276
void Spline::InitSpline ( void  )
private

Definition at line 716 of file Spline.cxx.

717 {
718  LOG("Spline", pDEBUG) << "Initializing spline...";
719 
720  fName = "genie-spline";
721  fXMin = 0.0;
722  fXMax = 0.0;
723  fYMax = 0.0;
724 
725  fInterpolator = 0;
726 
727  fYCanBeNegative = false;
728 
729  LOG("Spline", pDEBUG) << "...done initializing spline";
730 }
bool fYCanBeNegative
Definition: Spline.h:130
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
double fXMax
Definition: Spline.h:127
double fYMax
Definition: Spline.h:128
string fName
Definition: Spline.h:124
TSpline3 * fInterpolator
Definition: Spline.h:129
double fXMin
Definition: Spline.h:126
#define pDEBUG
Definition: Messenger.h:63
bool Spline::IsWithinValidRange ( double  x) const

Definition at line 355 of file Spline.cxx.

356 {
357  bool is_in_range = (fXMin <= x && x <= fXMax);
358  return is_in_range;
359 }
double fXMax
Definition: Spline.h:127
list x
Definition: train.py:276
double fXMin
Definition: Spline.h:126
bool Spline::LoadFromAsciiFile ( string  filename)

Definition at line 240 of file Spline.cxx.

241 {
242  LOG("Spline", pDEBUG) << "Retrieving data from file: " << filename;
243 
244  TNtupleD nt("ntuple","","x:y");
245  nt.ReadFile(filename.c_str());
246 
247  this->LoadFromNtuple(&nt, "x:y");
248  return true;
249 }
string filename
Definition: train.py:213
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
bool LoadFromNtuple(TNtupleD *nt, string xy, string cut="")
Definition: Spline.cxx:251
#define pDEBUG
Definition: Messenger.h:63
bool Spline::LoadFromDBase ( TSQLServer *  db,
string  query 
)

Definition at line 297 of file Spline.cxx.

298 {
299  LOG("Spline", pDEBUG) << "Retrieving data from data-base: ";
300  return false;
301 }
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pDEBUG
Definition: Messenger.h:63
bool Spline::LoadFromNtuple ( TNtupleD *  nt,
string  xy,
string  cut = "" 
)

Definition at line 251 of file Spline.cxx.

252 {
253  if(!nt) return false;
254 
255  TTree * tree = dynamic_cast<TTree *> (nt);
256  return this->LoadFromTree(tree,var,cut);
257 }
int var
Definition: 018_def.c:9
bool LoadFromTree(TTree *tr, string xy, string cut="")
Definition: Spline.cxx:259
bool Spline::LoadFromTree ( TTree *  tr,
string  xy,
string  cut = "" 
)

Definition at line 259 of file Spline.cxx.

260 {
261  LOG("Spline", pDEBUG) << "Retrieving data from tree: " << tree->GetName();
262 
263  if(!cut.size()) tree->Draw(var.c_str(), "", "GOFF");
264  else tree->Draw(var.c_str(), cut.c_str(), "GOFF");
265 
266  TH2F * hst = (TH2F*)gROOT->FindObject("htemp");
267  if(hst) { hst->SetDirectory(0); delete hst; }
268 
269  // Now, take into account that the data retrieved from the ntuple would
270  // not be sorted in x and the resulting spline will be bogus...
271  // Sort the x array, use the x-sorting to re-arrange the y array and only
272  // then build the spline..
273 
274  int n = tree->GetSelectedRows();
275 
276  int * idx = new int[n];
277  double * x = new double[n];
278  double * y = new double[n];
279 
280  TMath::Sort(n,tree->GetV1(),idx,false);
281 
282  for(int i=0; i<n; i++) {
283  int ii = idx[i];
284  x[i] = (tree->GetV1())[ii];
285  y[i] = (tree->GetV2())[ii];
286  }
287 
288  this->BuildSpline(n,x,y);
289 
290  delete [] idx;
291  delete [] x;
292  delete [] y;
293 
294  return true;
295 }
void BuildSpline(int nentries, double x[], double y[])
Definition: Spline.cxx:738
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
std::void_t< T > n
int var
Definition: 018_def.c:9
list x
Definition: train.py:276
#define pDEBUG
Definition: Messenger.h:63
bool Spline::LoadFromTSpline3 ( const TSpline3 &  spline,
int  nknots 
)

Definition at line 303 of file Spline.cxx.

304 {
305  int nentries = nknots;
306 
307  double * x = new double[nentries];
308  double * y = new double[nentries];
309  double cx = 0, cy = 0;
310 
311  for(int i = 0; i < nentries; i++) {
312  spline.GetKnot(i, cx, cy);
313  x[i] = cx;
314  y[i] = cy;
315  }
316  this->BuildSpline(nentries, x, y);
317 
318  delete [] x;
319  delete [] y;
320 
321  return true;
322 }
void BuildSpline(int nentries, double x[], double y[])
Definition: Spline.cxx:738
list x
Definition: train.py:276
bool Spline::LoadFromXmlFile ( string  filename,
string  xtag,
string  ytag 
)

Definition at line 153 of file Spline.cxx.

154 {
155  LOG("Spline", pDEBUG) << "Retrieving data from file: " << filename;
156 
157  xmlDocPtr xml_doc = xmlParseFile(filename.c_str());
158 
159  if(xml_doc==NULL) {
160  LOG("Spline", pERROR)
161  << "XML file could not be parsed! [filename: " << filename << "]";
162  return false;
163  }
164 
165  xmlNodePtr xmlCur = xmlDocGetRootElement(xml_doc);
166 
167  if(xmlCur==NULL) {
168  LOG("Spline", pERROR)
169  << "XML doc. has null root element! [filename: " << filename << "]";
170  return false;
171  }
172  if( xmlStrcmp(xmlCur->name, (const xmlChar *) "spline") ) {
173  LOG("Spline", pERROR)
174  << "XML doc. has invalid root element! [filename: " << filename << "]";
175  return false;
176  }
177 
178  string name = utils::str::TrimSpaces(
179  utils::xml::GetAttribute(xmlCur, "name"));
180  string snkn = utils::str::TrimSpaces(
181  utils::xml::GetAttribute(xmlCur, "nknots"));
182  int nknots = atoi( snkn.c_str() );
183 
184  LOG("Spline", pINFO)
185  << "Parsing XML spline: " << name << ", nknots = " << nknots;
186 
187  int iknot = 0;
188  double * vx = new double[nknots];
189  double * vy = new double[nknots];
190 
191  xmlNodePtr xmlSplChild = xmlCur->xmlChildrenNode; // <knots>'s
192 
193  // loop over all xml tree nodes that are children of the <spline> node
194  while (xmlSplChild != NULL) {
195  LOG("Spline", pDEBUG)
196  << "Got <spline> children node: " << xmlSplChild->name;
197 
198  // enter everytime you find a <knot> tag
199  if( (!xmlStrcmp(xmlSplChild->name, (const xmlChar *) "knot")) ) {
200 
201  xmlNodePtr xmlKnotChild = xmlSplChild->xmlChildrenNode;
202 
203  // loop over all xml tree nodes that are children of this <knot>
204  while (xmlKnotChild != NULL) {
205  LOG("Spline", pDEBUG)
206  << "Got <knot> children node: " << xmlKnotChild->name;
207 
208  // enter everytime you find a <E> or a <xsec> tag
209  const xmlChar * tag = xmlKnotChild->name;
210  bool is_xtag = ! xmlStrcmp(tag,(const xmlChar *) xtag.c_str());
211  bool is_ytag = ! xmlStrcmp(tag,(const xmlChar *) ytag.c_str());
212  if (is_xtag || is_ytag) {
213  xmlNodePtr xmlValTagChild = xmlKnotChild->xmlChildrenNode;
214  string val = utils::xml::TrimSpaces(
215  xmlNodeListGetString(xml_doc, xmlValTagChild, 1));
216 
217  if (is_xtag) vx[iknot] = atof(val.c_str());
218  if (is_ytag) vy[iknot] = atof(val.c_str());
219 
220  xmlFree(xmlValTagChild);
221  LOG("Spline", pDEBUG) << "tag: " << tag << ", value: " << val;
222  }//if current tag is <E>,<xsec>
223 
224  xmlKnotChild = xmlKnotChild->next;
225  }//[end of] loop over tags within <knot>...</knot> tags
226  xmlFree(xmlKnotChild);
227  iknot++;
228  } // if <knot>
229  xmlSplChild = xmlSplChild->next;
230  } //[end of] loop over tags within <spline>...</spline> tags
231  xmlFree(xmlSplChild);
232 
233  this->BuildSpline(nknots, vx, vy);
234  delete [] vx;
235  delete [] vy;
236 
237  return true;
238 }
static QCString name
Definition: declinfo.cpp:673
#define pERROR
Definition: Messenger.h:59
string TrimSpaces(xmlChar *xmls)
void BuildSpline(int nentries, double x[], double y[])
Definition: Spline.cxx:738
string filename
Definition: train.py:213
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
#define pINFO
Definition: Messenger.h:62
string TrimSpaces(string input)
Definition: StringUtils.cxx:18
string GetAttribute(xmlNodePtr xml_cur, string attr_name)
#define pDEBUG
Definition: Messenger.h:63
void Spline::Multiply ( const Spline spl,
double  c = 1 
)

Definition at line 605 of file Spline.cxx.

606 {
607  // continue only if the input spline is defined at a wider x-range
608  double xmin = this->XMin();
609  double xmax = this->XMax();
610  bool ok = spl.IsWithinValidRange(xmin) && spl.IsWithinValidRange(xmax);
611  if(!ok) {
612  LOG("Spline", pERROR) << "** Can't multiply splines. X-range mismatch!";
613  return;
614  }
615 
616  int nknots = this->NKnots();
617  double * x = new double[nknots];
618  double * y = new double[nknots];
619 
620  for(int i=0; i<nknots; i++) {
621  this->GetKnot(i,x[i],y[i]);
622  y[i] *= (c * spl.Evaluate(x[i]));
623  }
624  this->ResetSpline();
625  this->BuildSpline(nknots,x,y);
626  delete [] x;
627  delete [] y;
628 }
#define pERROR
Definition: Messenger.h:59
void BuildSpline(int nentries, double x[], double y[])
Definition: Spline.cxx:738
double Evaluate(double x) const
Definition: Spline.cxx:361
bool IsWithinValidRange(double x) const
Definition: Spline.cxx:355
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void ResetSpline(void)
Definition: Spline.cxx:732
void GetKnot(int iknot, double &x, double &y) const
Definition: Spline.cxx:324
double XMax(void) const
Definition: Spline.h:77
int NKnots(void) const
Definition: Spline.h:72
list x
Definition: train.py:276
double XMin(void) const
Definition: Spline.h:76
void Spline::Multiply ( double  a)

Definition at line 679 of file Spline.cxx.

680 {
681  int nknots = this->NKnots();
682  double * x = new double[nknots];
683  double * y = new double[nknots];
684 
685  for(int i=0; i<nknots; i++) {
686  this->GetKnot(i,x[i],y[i]);
687  y[i]*=a;
688  }
689  this->ResetSpline();
690  this->BuildSpline(nknots,x,y);
691  delete [] x;
692  delete [] y;
693 }
void BuildSpline(int nentries, double x[], double y[])
Definition: Spline.cxx:738
void ResetSpline(void)
Definition: Spline.cxx:732
const double a
void GetKnot(int iknot, double &x, double &y) const
Definition: Spline.cxx:324
int NKnots(void) const
Definition: Spline.h:72
list x
Definition: train.py:276
string genie::Spline::Name ( void  ) const
inline

Definition at line 83 of file Spline.h.

83 { return fName; }
string fName
Definition: Spline.h:124
int genie::Spline::NKnots ( void  ) const
inline

Definition at line 72 of file Spline.h.

72 {return fNKnots;}
int fNKnots
Definition: Spline.h:125
void Spline::Print ( ostream &  stream) const

Definition at line 562 of file Spline.cxx.

563 {
564  int nknots = this->NKnots();
565  double xmin = this->XMin();
566  double xmax = this->XMax();
567 
568  stream << endl;
569  stream << "** Spline: " << this->Name() << endl;
570  stream << "Has " << nknots
571  << " knots in the [" << xmin << ", " << xmax << "] range" << endl;
572  double x,y;
573  for(int i=0; i<nknots; i++) {
574  this->GetKnot(i,x,y);
575  stream << "-- knot : " << i
576  << " -> (x = " << x << ", y = " << y << ")" << endl;
577  }
578 }
void GetKnot(int iknot, double &x, double &y) const
Definition: Spline.cxx:324
double XMax(void) const
Definition: Spline.h:77
string Name(void) const
Definition: Spline.h:83
int NKnots(void) const
Definition: Spline.h:72
list x
Definition: train.py:276
double XMin(void) const
Definition: Spline.h:76
QTextStream & endl(QTextStream &s)
void Spline::ResetSpline ( void  )
private

Definition at line 732 of file Spline.cxx.

733 {
734  if(fInterpolator) delete fInterpolator;
735  this->InitSpline();
736 }
TSpline3 * fInterpolator
Definition: Spline.h:129
void InitSpline(void)
Definition: Spline.cxx:716
void Spline::SaveAsROOT ( string  filename,
string  name = "",
bool  recreate = false 
) const

Definition at line 478 of file Spline.cxx.

479 {
480  string spline_name = (name.size()>0 ? name : fName);
481 
482  string opt = ( (recreate) ? "RECREATE" : "UPDATE" );
483 
484  TFile f(filename.c_str(), opt.c_str());
485  if(fInterpolator) fInterpolator->Write(spline_name.c_str());
486  f.Close();
487 }
static QCString name
Definition: declinfo.cpp:673
opt
Definition: train.py:196
uint size() const
Definition: qcstring.h:201
string filename
Definition: train.py:213
string fName
Definition: Spline.h:124
TSpline3 * fInterpolator
Definition: Spline.h:129
void Spline::SaveAsText ( string  filename,
string  format = "%10.6f\t%10.6f" 
) const

Definition at line 457 of file Spline.cxx.

458 {
459  ofstream outtxt(filename.c_str());
460  if(!outtxt.is_open()) {
461  LOG("Spline", pERROR) << "Couldn't create file = " << filename;
462  return;
463  }
464  int nknots = this->NKnots();
465  outtxt << nknots << endl;
466 
467  double x=0, y=0;
468  for(int iknot = 0; iknot < nknots; iknot++) {
469  fInterpolator->GetKnot(iknot, x, y);
470  char line[1024];
471  sprintf(line,format.c_str(),x,y);
472  outtxt << line << endl;
473  }
474  outtxt << endl;
475  outtxt.close();
476 }
#define pERROR
Definition: Messenger.h:59
static bool format(QChar::Decomposition tag, QString &str, int index, int len)
Definition: qstring.cpp:11496
string filename
Definition: train.py:213
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
TSpline3 * fInterpolator
Definition: Spline.h:129
void line(double t, double *p, double &x, double &y, double &z)
int NKnots(void) const
Definition: Spline.h:72
list x
Definition: train.py:276
QTextStream & endl(QTextStream &s)
void Spline::SaveAsXml ( string  filename,
string  xtag,
string  ytag,
string  name = "" 
) const

Definition at line 412 of file Spline.cxx.

414 {
415  ofstream outxml(filename.c_str());
416  if(!outxml.is_open()) {
417  LOG("Spline", pERROR) << "Couldn't create file = " << filename;
418  return;
419  }
420  string spline_name = (name.size()>0 ? name : fName);
421 
422  this->SaveAsXml(outxml, xtag, ytag, spline_name);
423 
424  outxml << endl;
425  outxml.close();
426 }
static QCString name
Definition: declinfo.cpp:673
#define pERROR
Definition: Messenger.h:59
uint size() const
Definition: qcstring.h:201
string filename
Definition: train.py:213
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
string fName
Definition: Spline.h:124
void SaveAsXml(string filename, string xtag, string ytag, string name="") const
Definition: Spline.cxx:412
QTextStream & endl(QTextStream &s)
void Spline::SaveAsXml ( ofstream &  str,
string  xtag,
string  ytag,
string  name = "" 
) const

Definition at line 428 of file Spline.cxx.

430 {
431  string spline_name = (name.size()>0 ? name : fName);
432 
433  // create a spline tag with the number of knots as an attribute
434  int nknots = this->NKnots();
435  string padding = " ";
436  ofs << padding << "<spline name=\"" << spline_name
437  << "\" nknots=\"" << nknots << "\">" << endl;
438 
439  // start printing the knots
440  double x=0, y=0;
441  for(int iknot = 0; iknot < nknots; iknot++)
442  {
443  fInterpolator->GetKnot(iknot, x, y);
444 
445  ofs << std::fixed << setprecision(5);
446  ofs << "\t<knot>"
447  << " <" << xtag << "> " << setfill(' ')
448  << setw(10) << x << " </" << xtag << ">";
449  ofs << std::scientific << setprecision(10);
450  ofs << " <" << ytag << "> " << setfill(' ')
451  << setw(10) << y << " </" << ytag << ">"
452  << " </knot>" << endl;
453  }
454  ofs << padding << "</spline>" << endl;
455 }
static QCString name
Definition: declinfo.cpp:673
uint size() const
Definition: qcstring.h:201
Q_EXPORT QTSManip setprecision(int p)
Definition: qtextstream.h:343
Q_EXPORT QTSManip setw(int w)
Definition: qtextstream.h:331
string fName
Definition: Spline.h:124
TSpline3 * fInterpolator
Definition: Spline.h:129
int NKnots(void) const
Definition: Spline.h:72
list x
Definition: train.py:276
Q_EXPORT QTSManip setfill(int f)
Definition: qtextstream.h:337
QTextStream & endl(QTextStream &s)
void genie::Spline::SetName ( string  name)
inline

Definition at line 82 of file Spline.h.

82 { fName = name; }
static QCString name
Definition: declinfo.cpp:673
string fName
Definition: Spline.h:124
double genie::Spline::XMax ( void  ) const
inline

Definition at line 77 of file Spline.h.

77 {return fXMax; }
double fXMax
Definition: Spline.h:127
double genie::Spline::XMin ( void  ) const
inline

Definition at line 76 of file Spline.h.

76 {return fXMin; }
double fXMin
Definition: Spline.h:126
void genie::Spline::YCanBeNegative ( bool  tf)
inline

Definition at line 85 of file Spline.h.

85 { fYCanBeNegative = tf; }
bool fYCanBeNegative
Definition: Spline.h:130
double genie::Spline::YMax ( void  ) const
inline

Definition at line 78 of file Spline.h.

78 {return fYMax; }
double fYMax
Definition: Spline.h:128

Friends And Related Function Documentation

ostream& operator<< ( ostream &  stream,
const Spline spl 
)
friend

Member Data Documentation

TSpline3* genie::Spline::fInterpolator
private

Definition at line 129 of file Spline.h.

string genie::Spline::fName
private

Definition at line 124 of file Spline.h.

int genie::Spline::fNKnots
private

Definition at line 125 of file Spline.h.

double genie::Spline::fXMax
private

Definition at line 127 of file Spline.h.

double genie::Spline::fXMin
private

Definition at line 126 of file Spline.h.

bool genie::Spline::fYCanBeNegative
private

Definition at line 130 of file Spline.h.

double genie::Spline::fYMax
private

Definition at line 128 of file Spline.h.


The documentation for this class was generated from the following files: