15 #include "libxml/parser.h"
16 #include "libxml/xmlmemory.h"
21 #include <TSQLServer.h>
34 using std::setprecision;
38 using namespace genie;
60 string fmt = (is_xml) ?
"XML" :
"ASCII";
63 <<
"Constructing spline from data in " << fmt <<
" file: " << filename;
76 LOG(
"Spline",
pDEBUG) <<
"Constructing spline from data in a TNtuple";
85 LOG(
"Spline",
pDEBUG) <<
"Constructing spline from data in a TTree";
94 LOG(
"Spline",
pDEBUG) <<
"Constructing spline from data in a MySQL server";
104 <<
"Constructing spline from the arrays passed to the ctor";
114 <<
"Constructing spline from the arrays passed to the ctor";
116 double * dblx =
new double[nentries];
117 double * dbly =
new double[nentries];
119 for(
int i = 0; i < nentries; i++) {
120 dblx[i] = double( x[i] );
121 dbly[i] = double( y[i] );
132 TObject(), fInterpolator(0)
134 LOG(
"Spline",
pDEBUG) <<
"Spline copy constructor";
140 TObject(), fInterpolator(0)
143 <<
"Constructing spline from the input TSpline3 object";
157 LOG(
"Spline",
pDEBUG) <<
"Retrieving data from file: " << filename;
159 xmlDocPtr xml_doc = xmlParseFile(filename.c_str());
163 <<
"XML file could not be parsed! [filename: " << filename <<
"]";
167 xmlNodePtr xmlCur = xmlDocGetRootElement(xml_doc);
171 <<
"XML doc. has null root element! [filename: " << filename <<
"]";
174 if( xmlStrcmp(xmlCur->name, (
const xmlChar *)
"spline") ) {
176 <<
"XML doc. has invalid root element! [filename: " << filename <<
"]";
184 int nknots = atoi( snkn.c_str() );
187 <<
"Parsing XML spline: " << name <<
", nknots = " << nknots;
190 double * vx =
new double[nknots];
191 double * vy =
new double[nknots];
193 xmlNodePtr xmlSplChild = xmlCur->xmlChildrenNode;
196 while (xmlSplChild != NULL) {
198 <<
"Got <spline> children node: " << xmlSplChild->name;
201 if( (!xmlStrcmp(xmlSplChild->name, (
const xmlChar *)
"knot")) ) {
203 xmlNodePtr xmlKnotChild = xmlSplChild->xmlChildrenNode;
206 while (xmlKnotChild != NULL) {
208 <<
"Got <knot> children node: " << xmlKnotChild->name;
211 const xmlChar * tag = xmlKnotChild->name;
212 bool is_xtag = ! xmlStrcmp(tag,(
const xmlChar *) xtag.c_str());
213 bool is_ytag = ! xmlStrcmp(tag,(
const xmlChar *) ytag.c_str());
214 if (is_xtag || is_ytag) {
215 xmlNodePtr xmlValTagChild = xmlKnotChild->xmlChildrenNode;
217 xmlNodeListGetString(xml_doc, xmlValTagChild, 1));
219 if (is_xtag) vx[iknot] = atof(val.c_str());
220 if (is_ytag) vy[iknot] = atof(val.c_str());
222 xmlFree(xmlValTagChild);
223 LOG(
"Spline",
pDEBUG) <<
"tag: " << tag <<
", value: " << val;
226 xmlKnotChild = xmlKnotChild->next;
228 xmlFree(xmlKnotChild);
231 xmlSplChild = xmlSplChild->next;
233 xmlFree(xmlSplChild);
244 LOG(
"Spline",
pDEBUG) <<
"Retrieving data from file: " << filename;
246 TNtupleD nt(
"ntuple",
"",
"x:y");
247 nt.ReadFile(filename.c_str());
255 if(!nt)
return false;
257 TTree * tree =
dynamic_cast<TTree *
> (nt);
263 LOG(
"Spline",
pDEBUG) <<
"Retrieving data from tree: " << tree->GetName();
265 if(!cut.size()) tree->Draw(var.c_str(),
"",
"GOFF");
266 else tree->Draw(var.c_str(), cut.c_str(),
"GOFF");
268 TH2F * hst = (TH2F*)gROOT->FindObject(
"htemp");
269 if(hst) { hst->SetDirectory(0);
delete hst; }
276 int n = tree->GetSelectedRows();
278 int * idx =
new int[n];
279 double * x =
new double[n];
280 double * y =
new double[n];
282 TMath::Sort(n,tree->GetV1(),idx,
false);
284 for(
int i=0; i<n; i++) {
286 x[i] = (tree->GetV1())[ii];
287 y[i] = (tree->GetV2())[ii];
301 LOG(
"Spline",
pDEBUG) <<
"Retrieving data from data-base: ";
307 int nentries = nknots;
309 double * x =
new double[nentries];
310 double * y =
new double[nentries];
311 double cx = 0, cy = 0;
313 for(
int i = 0; i < nentries; i++) {
314 spline.GetKnot(i, cx, cy);
329 LOG(
"Spline",
pWARN) <<
"Spline has not been built yet!";
338 LOG(
"Spline",
pWARN) <<
"Spline has not been built yet!";
349 LOG(
"Spline",
pWARN) <<
"Spline has not been built yet!";
359 bool is_in_range = (
fXMin <= x && x <=
fXMax);
365 LOG(
"Spline",
pDEBUG) <<
"Evaluating spline at point x = " << x;
366 assert(!TMath::IsNaN(x));
378 LOG(
"Spline",
pDEBUG) <<
"Point is between non-zero knots";
389 LOG(
"Spline",
pDEBUG) <<
"Point is between zero knots";
394 <<
"Point has zero" << (is0n ?
" left " :
" right ") <<
"knot";
395 double xpknot=0, ypknot=0, xnknot=0, ynknot=0;
398 if(is0n) y = ypknot * (x-xnknot)/(xpknot-xnknot);
399 else y = ynknot * (x-xnknot)/(xpknot-xnknot);
405 <<
" is not within spline range [" <<
fXMin <<
", " <<
fXMax <<
"]";
409 LOG(
"Spline",
pINFO) <<
"Negative y (" << y <<
")";
410 LOG(
"Spline",
pINFO) <<
"x = " << x;
414 LOG(
"Spline",
pDEBUG) <<
"Spline(x = " << x <<
") = " << y;
420 string filename,
string xtag,
string ytag,
string name)
const
422 ofstream outxml(filename.c_str());
423 if(!outxml.is_open()) {
424 LOG(
"Spline",
pERROR) <<
"Couldn't create file = " << filename;
427 string spline_name = (name.size()>0 ? name :
fName);
429 this->
SaveAsXml(outxml, xtag, ytag, spline_name);
436 ofstream & ofs,
string xtag,
string ytag,
string name)
const
438 string spline_name = (name.size()>0 ? name :
fName);
441 int nknots = this->
NKnots();
442 string padding =
" ";
443 ofs << padding <<
"<spline name=\"" << spline_name
444 <<
"\" nknots=\"" << nknots <<
"\">" << endl;
448 for(
int iknot = 0; iknot < nknots; iknot++)
452 ofs << std::fixed << setprecision(5);
454 <<
" <" << xtag <<
"> " << setfill(
' ')
455 << setw(10) << x <<
" </" << xtag <<
">";
456 ofs << std::scientific << setprecision(10);
457 ofs <<
" <" << ytag <<
"> " << setfill(
' ')
458 << setw(10) << y <<
" </" << ytag <<
">"
459 <<
" </knot>" << endl;
461 ofs << padding <<
"</spline>" << endl;
466 ofstream outtxt(filename.c_str());
467 if(!outtxt.is_open()) {
468 LOG(
"Spline",
pERROR) <<
"Couldn't create file = " << filename;
471 int nknots = this->
NKnots();
472 outtxt << nknots << endl;
475 for(
int iknot = 0; iknot < nknots; iknot++) {
478 sprintf(line,format.c_str(),x,y);
479 outtxt << line << endl;
487 string spline_name = (name.size()>0 ? name :
fName);
489 string opt = ( (recreate) ?
"RECREATE" :
"UPDATE" );
491 TFile f(filename.c_str(), opt.c_str());
497 int np,
bool scale_with_x,
bool in_log,
double fx,
double fy)
const
502 np = TMath::Max(np,2);
504 bool use_log = in_log && (
fXMin>0) && (
fXMax>0);
507 xmin = TMath::Log10(xmin);
508 xmax = TMath::Log10(xmax);
511 double dx = (xmax - xmin) / (np-1);
513 double * x =
new double[np];
514 double * y =
new double[np];
516 for(
int i=0; i<np; i++) {
517 x[i] = ( (use_log) ? TMath::Power(10, xmin+i*dx) : xmin + i*dx );
521 if (scale_with_x) y[i] /= x[i];
528 TGraph * graph =
new TGraph(np, x, y);
535 double x,
double & xknot,
double & yknot, Option_t * opt)
const
539 bool pos = (option.find(
"+") != string::npos);
540 bool neg = (option.find(
"-") != string::npos);
542 if(!pos && !neg)
return;
546 double xp=0, yp=0, xn=0, yn=0;
550 bool p = (TMath::Abs(x-xp) < TMath::Abs(x-xn));
553 if(p) { xknot = xp; yknot = yp; }
554 else { xknot = xn; yknot = yn; }
556 if(pos) { xknot = xp; yknot = yp; }
557 if(neg) { xknot = xn; yknot = yn; }
563 double xknot = 0, yknot = 0;
571 int nknots = this->
NKnots();
572 double xmin = this->
XMin();
573 double xmax = this->
XMax();
576 stream <<
"** Spline: " << this->
Name() << endl;
577 stream <<
"Has " << nknots
578 <<
" knots in the [" << xmin <<
", " << xmax <<
"] range" << endl;
580 for(
int i=0; i<nknots; i++) {
582 stream <<
"-- knot : " << i
583 <<
" -> (x = " << x <<
", y = " << y <<
")" << endl;
590 double xmin = this->
XMin();
591 double xmax = this->
XMax();
594 LOG(
"Spline",
pERROR) <<
"** Can't add splines. X-range mismatch!";
598 int nknots = this->
NKnots();
599 double * x =
new double[nknots];
600 double * y =
new double[nknots];
602 for(
int i=0; i<nknots; i++) {
615 double xmin = this->
XMin();
616 double xmax = this->
XMax();
619 LOG(
"Spline",
pERROR) <<
"** Can't multiply splines. X-range mismatch!";
623 int nknots = this->
NKnots();
624 double * x =
new double[nknots];
625 double * y =
new double[nknots];
627 for(
int i=0; i<nknots; i++) {
640 double xmin = this->
XMin();
641 double xmax = this->
XMax();
644 LOG(
"Spline",
pERROR) <<
"** Can't divide splines. X-range mismatch!";
648 int nknots = this->
NKnots();
649 double * x =
new double[nknots];
650 double * y =
new double[nknots];
652 for(
int i=0; i<nknots; i++) {
654 double denom = c * spl.
Evaluate(x[i]);
655 bool denom_is_zero = TMath::Abs(denom) < DBL_EPSILON;
657 LOG(
"Spline",
pERROR) <<
"** Refusing to divide spline knot by 0";
672 int nknots = this->
NKnots();
673 double * x =
new double[nknots];
674 double * y =
new double[nknots];
676 for(
int i=0; i<nknots; i++) {
688 int nknots = this->
NKnots();
689 double * x =
new double[nknots];
690 double * y =
new double[nknots];
692 for(
int i=0; i<nknots; i++) {
704 bool a_is_zero = TMath::Abs(a) < DBL_EPSILON;
706 LOG(
"Spline",
pERROR) <<
"** Refusing to divide spline by 0";
709 int nknots = this->
NKnots();
710 double * x =
new double[nknots];
711 double * y =
new double[nknots];
713 for(
int i=0; i<nknots; i++) {
725 LOG(
"Spline",
pDEBUG) <<
"Initializing spline...";
727 fName =
"genie-spline";
739 LOG(
"Spline",
pDEBUG) <<
"...done initializing spline";
752 LOG(
"Spline",
pDEBUG) <<
"Building spline...";
755 double xmax = x[nentries-1];
760 fYMax = y[ TMath::LocMax(nentries, y) ];
766 LOG(
"Spline",
pDEBUG) <<
"...done building spline";
775 ROOT::Math::Interpolation::Type gsltype;
796 gsltype = ROOT::Math::Interpolation::kLINEAR;
800 gsltype = ROOT::Math::Interpolation::kPOLYNOMIAL;
804 gsltype = ROOT::Math::Interpolation::kCSPLINE;
808 gsltype = ROOT::Math::Interpolation::kCSPLINE_PERIODIC;
812 gsltype = ROOT::Math::Interpolation::kAKIMA;
816 gsltype = ROOT::Math::Interpolation::kAKIMA_PERIODIC;
822 <<
"Unknown interpolator type. Setting it to default [TSpline3].";
TGraph * GetAsTGraph(int np=500, bool xscaling=false, bool inlog=false, double fx=1., double fy=1.) const
ROOT::Math::Interpolator * fGSLInterpolator
bool AreEqual(double x1, double x2)
void SetType(string type)
string TrimSpaces(xmlChar *xmls)
A numeric analysis tool class for interpolating 1-D functions.
void Print(ostream &stream) const
bool LoadFromAsciiFile(string filename)
void BuildSpline(int nentries, double x[], double y[])
double Evaluate(double x) const
bool ClosestKnotValueIsZero(double x, Option_t *opt="-+") const
bool IsWithinValidRange(double x) const
TSpline3 * GetAsTSpline(void) const
bool LoadFromXmlFile(string filename, string xtag, string ytag)
double GetKnotX(int iknot) const
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
void SaveAsROOT(string filename, string name="", bool recreate=false) const
void Add(const Spline &spl, double c=1)
void Multiply(const Spline &spl, double c=1)
void GetKnot(int iknot, double &x, double &y) const
string TrimSpaces(string input)
double GetKnotY(int iknot) const
bool LoadFromTSpline3(const TSpline3 &spline, int nknots)
void SaveAsText(string filename, string format="%10.6f\t%10.6f") const
void SaveAsXml(string filename, string xtag, string ytag, string name="") const
TSpline5 * fInterpolator5
bool LoadFromNtuple(TNtupleD *nt, string xy, string cut="")
bool LoadFromDBase(TSQLServer *db, string query)
ostream & operator<<(ostream &stream, const AlgConfigPool &config_pool)
bool LoadFromTree(TTree *tr, string xy, string cut="")
string ToUpper(string input)
void FindClosestKnot(double x, double &xknot, double &yknot, Option_t *opt="-+") const
string GetAttribute(xmlNodePtr xml_cur, string attr_name)
void Divide(const Spline &spl, double c=1)