GENIEGenerator
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Spline.cxx
Go to the documentation of this file.
1 //____________________________________________________________________________
2 /*
3  Copyright (c) 2003-2024, The GENIE Collaboration
4  For the full text of the license visit http://copyright.genie-mc.org
5 
6  Costas Andreopoulos <c.andreopoulos \at cern.ch>
7  University of Liverpool
8 */
9 //____________________________________________________________________________
10 
11 #include <cassert>
12 #include <iomanip>
13 #include <cfloat>
14 
15 #include "libxml/parser.h"
16 #include "libxml/xmlmemory.h"
17 
18 #include <TFile.h>
19 #include <TNtupleD.h>
20 #include <TTree.h>
21 #include <TSQLServer.h>
22 #include <TGraph.h>
23 #include <TMath.h>
24 #include <TH2F.h>
25 #include <TROOT.h>
26 
31 
32 using std::endl;
33 using std::setw;
34 using std::setprecision;
35 using std::setfill;
36 using std::ios;
37 
38 using namespace genie;
39 
41 
42 //___________________________________________________________________________
43 namespace genie
44 {
45  ostream & operator << (ostream & stream, const Spline & spl)
46  {
47  spl.Print(stream);
48  return stream;
49  }
50 }
51 //___________________________________________________________________________
53 {
54  this->InitSpline();
55 }
56 //___________________________________________________________________________
57 Spline::Spline(string filename, string xtag, string ytag, bool is_xml) :
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 }
72 //___________________________________________________________________________
73 Spline::Spline(TNtupleD * ntuple, string var, string cut) :
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 }
81 //___________________________________________________________________________
82 Spline::Spline(TTree * tree, string var, string cut) :
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 }
90 //___________________________________________________________________________
91 Spline::Spline(TSQLServer * db, string query) :
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 }
99 //___________________________________________________________________________
100 Spline::Spline(int nentries, double x[], double y[]) :
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 }
109 //___________________________________________________________________________
110 Spline::Spline(int nentries, float x[], float y[]) :
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 }
130 //___________________________________________________________________________
131 Spline::Spline(const Spline & spline) :
132  TObject(), fInterpolator(0)
133 {
134  LOG("Spline", pDEBUG) << "Spline copy constructor";
135  this->InitSpline();
136  this->LoadFromTSpline3( *spline.GetAsTSpline(), spline.NKnots() );
137 }
138 //___________________________________________________________________________
139 Spline::Spline(const TSpline3 & spline, int nknots) :
140  TObject(), fInterpolator(0)
141 {
142  LOG("Spline", pDEBUG)
143  << "Constructing spline from the input TSpline3 object";
144  this->InitSpline();
145  this->LoadFromTSpline3( spline, nknots );
146 }
147 //___________________________________________________________________________
149 {
150  if(fInterpolator) delete fInterpolator;
151  if(fInterpolator5) delete fInterpolator5;
153 }
154 //___________________________________________________________________________
155 bool Spline::LoadFromXmlFile(string filename, string xtag, string ytag)
156 {
157  LOG("Spline", pDEBUG) << "Retrieving data from file: " << filename;
158 
159  xmlDocPtr xml_doc = xmlParseFile(filename.c_str());
160 
161  if(xml_doc==NULL) {
162  LOG("Spline", pERROR)
163  << "XML file could not be parsed! [filename: " << filename << "]";
164  return false;
165  }
166 
167  xmlNodePtr xmlCur = xmlDocGetRootElement(xml_doc);
168 
169  if(xmlCur==NULL) {
170  LOG("Spline", pERROR)
171  << "XML doc. has null root element! [filename: " << filename << "]";
172  return false;
173  }
174  if( xmlStrcmp(xmlCur->name, (const xmlChar *) "spline") ) {
175  LOG("Spline", pERROR)
176  << "XML doc. has invalid root element! [filename: " << filename << "]";
177  return false;
178  }
179 
180  string name = utils::str::TrimSpaces(
181  utils::xml::GetAttribute(xmlCur, "name"));
182  string snkn = utils::str::TrimSpaces(
183  utils::xml::GetAttribute(xmlCur, "nknots"));
184  int nknots = atoi( snkn.c_str() );
185 
186  LOG("Spline", pINFO)
187  << "Parsing XML spline: " << name << ", nknots = " << nknots;
188 
189  int iknot = 0;
190  double * vx = new double[nknots];
191  double * vy = new double[nknots];
192 
193  xmlNodePtr xmlSplChild = xmlCur->xmlChildrenNode; // <knots>'s
194 
195  // loop over all xml tree nodes that are children of the <spline> node
196  while (xmlSplChild != NULL) {
197  LOG("Spline", pDEBUG)
198  << "Got <spline> children node: " << xmlSplChild->name;
199 
200  // enter everytime you find a <knot> tag
201  if( (!xmlStrcmp(xmlSplChild->name, (const xmlChar *) "knot")) ) {
202 
203  xmlNodePtr xmlKnotChild = xmlSplChild->xmlChildrenNode;
204 
205  // loop over all xml tree nodes that are children of this <knot>
206  while (xmlKnotChild != NULL) {
207  LOG("Spline", pDEBUG)
208  << "Got <knot> children node: " << xmlKnotChild->name;
209 
210  // enter everytime you find a <E> or a <xsec> tag
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;
216  string val = utils::xml::TrimSpaces(
217  xmlNodeListGetString(xml_doc, xmlValTagChild, 1));
218 
219  if (is_xtag) vx[iknot] = atof(val.c_str());
220  if (is_ytag) vy[iknot] = atof(val.c_str());
221 
222  xmlFree(xmlValTagChild);
223  LOG("Spline", pDEBUG) << "tag: " << tag << ", value: " << val;
224  }//if current tag is <E>,<xsec>
225 
226  xmlKnotChild = xmlKnotChild->next;
227  }//[end of] loop over tags within <knot>...</knot> tags
228  xmlFree(xmlKnotChild);
229  iknot++;
230  } // if <knot>
231  xmlSplChild = xmlSplChild->next;
232  } //[end of] loop over tags within <spline>...</spline> tags
233  xmlFree(xmlSplChild);
234 
235  this->BuildSpline(nknots, vx, vy);
236  delete [] vx;
237  delete [] vy;
238 
239  return true;
240 }
241 //___________________________________________________________________________
242 bool Spline::LoadFromAsciiFile(string filename)
243 {
244  LOG("Spline", pDEBUG) << "Retrieving data from file: " << filename;
245 
246  TNtupleD nt("ntuple","","x:y");
247  nt.ReadFile(filename.c_str());
248 
249  this->LoadFromNtuple(&nt, "x:y");
250  return true;
251 }
252 //___________________________________________________________________________
253 bool Spline::LoadFromNtuple(TNtupleD * nt, string var, string cut)
254 {
255  if(!nt) return false;
256 
257  TTree * tree = dynamic_cast<TTree *> (nt);
258  return this->LoadFromTree(tree,var,cut);
259 }
260 //___________________________________________________________________________
261 bool Spline::LoadFromTree(TTree * tree, string var, string cut)
262 {
263  LOG("Spline", pDEBUG) << "Retrieving data from tree: " << tree->GetName();
264 
265  if(!cut.size()) tree->Draw(var.c_str(), "", "GOFF");
266  else tree->Draw(var.c_str(), cut.c_str(), "GOFF");
267 
268  TH2F * hst = (TH2F*)gROOT->FindObject("htemp");
269  if(hst) { hst->SetDirectory(0); delete hst; }
270 
271  // Now, take into account that the data retrieved from the ntuple would
272  // not be sorted in x and the resulting spline will be bogus...
273  // Sort the x array, use the x-sorting to re-arrange the y array and only
274  // then build the spline..
275 
276  int n = tree->GetSelectedRows();
277 
278  int * idx = new int[n];
279  double * x = new double[n];
280  double * y = new double[n];
281 
282  TMath::Sort(n,tree->GetV1(),idx,false);
283 
284  for(int i=0; i<n; i++) {
285  int ii = idx[i];
286  x[i] = (tree->GetV1())[ii];
287  y[i] = (tree->GetV2())[ii];
288  }
289 
290  this->BuildSpline(n,x,y);
291 
292  delete [] idx;
293  delete [] x;
294  delete [] y;
295 
296  return true;
297 }
298 //___________________________________________________________________________
299 bool Spline::LoadFromDBase(TSQLServer * /*db*/, string /*query*/)
300 {
301  LOG("Spline", pDEBUG) << "Retrieving data from data-base: ";
302  return false;
303 }
304 //___________________________________________________________________________
305 bool Spline::LoadFromTSpline3(const TSpline3 & spline, int nknots)
306 {
307  int nentries = nknots;
308 
309  double * x = new double[nentries];
310  double * y = new double[nentries];
311  double cx = 0, cy = 0;
312 
313  for(int i = 0; i < nentries; i++) {
314  spline.GetKnot(i, cx, cy);
315  x[i] = cx;
316  y[i] = cy;
317  }
318  this->BuildSpline(nentries, x, y);
319 
320  delete [] x;
321  delete [] y;
322 
323  return true;
324 }
325 //___________________________________________________________________________
326 void Spline::GetKnot(int iknot, double & x, double & y) const
327 {
328  if(!fInterpolator) {
329  LOG("Spline", pWARN) << "Spline has not been built yet!";
330  return;
331  }
332  fInterpolator->GetKnot(iknot,x,y);
333 }
334 //___________________________________________________________________________
335 double Spline::GetKnotX(int iknot) const
336 {
337  if(!fInterpolator) {
338  LOG("Spline", pWARN) << "Spline has not been built yet!";
339  return 0;
340  }
341  double x,y;
342  fInterpolator->GetKnot(iknot,x,y);
343  return x;
344 }
345 //___________________________________________________________________________
346 double Spline::GetKnotY(int iknot) const
347 {
348  if(!fInterpolator) {
349  LOG("Spline", pWARN) << "Spline has not been built yet!";
350  return 0;
351  }
352  double x,y;
353  fInterpolator->GetKnot(iknot,x,y);
354  return y;
355 }
356 //___________________________________________________________________________
357 bool Spline::IsWithinValidRange(double x) const
358 {
359  bool is_in_range = (fXMin <= x && x <= fXMax);
360  return is_in_range;
361 }
362 //___________________________________________________________________________
363 double Spline::Evaluate(double x) const
364 {
365  LOG("Spline", pDEBUG) << "Evaluating spline at point x = " << x;
366  assert(!TMath::IsNaN(x));
367 
368  double y = 0;
369  if( this->IsWithinValidRange(x) ) {
370 
371  // we can interpolate within the range of spline knots - be careful with
372  // strange cubic spline behaviour when close to knots with y=0
373  bool is0p = this->ClosestKnotValueIsZero(x, "+");
374  bool is0n = this->ClosestKnotValueIsZero(x, "-");
375 
376  if(!is0p && !is0n) {
377  // both knots (on the left and right are non-zero) - just interpolate
378  LOG("Spline", pDEBUG) << "Point is between non-zero knots";
379  if (fInterpolatorType == "TSpline3")
380  y = fInterpolator->Eval(x);
381  else if (fInterpolatorType == "TSpline5")
382  y = fInterpolator5->Eval(x);
383  else
384  y = fGSLInterpolator->Eval(x);
385  } else {
386  // at least one of the neighboring knots has y=0
387  if(is0p && is0n) {
388  // both neighboring knots have y=0
389  LOG("Spline", pDEBUG) << "Point is between zero knots";
390  y=0;
391  } else {
392  // just 1 neighboring knot has y=0 - do a linear interpolation
393  LOG("Spline", pDEBUG)
394  << "Point has zero" << (is0n ? " left " : " right ") << "knot";
395  double xpknot=0, ypknot=0, xnknot=0, ynknot=0;
396  this->FindClosestKnot(x, xnknot, ynknot, "-");
397  this->FindClosestKnot(x, xpknot, ypknot, "+");
398  if(is0n) y = ypknot * (x-xnknot)/(xpknot-xnknot);
399  else y = ynknot * (x-xnknot)/(xpknot-xnknot);
400  }
401  }
402 
403  } else {
404  LOG("Spline", pDEBUG) << "x = " << x
405  << " is not within spline range [" << fXMin << ", " << fXMax << "]";
406  }
407 
408  if(y<0 && !fYCanBeNegative) {
409  LOG("Spline", pINFO) << "Negative y (" << y << ")";
410  LOG("Spline", pINFO) << "x = " << x;
411  LOG("Spline", pINFO) << "spline range [" << fXMin << ", " << fXMax << "]";
412  }
413 
414  LOG("Spline", pDEBUG) << "Spline(x = " << x << ") = " << y;
415 
416  return y;
417 }
418 //___________________________________________________________________________
420  string filename, string xtag, string ytag, string name) const
421 {
422  ofstream outxml(filename.c_str());
423  if(!outxml.is_open()) {
424  LOG("Spline", pERROR) << "Couldn't create file = " << filename;
425  return;
426  }
427  string spline_name = (name.size()>0 ? name : fName);
428 
429  this->SaveAsXml(outxml, xtag, ytag, spline_name);
430 
431  outxml << endl;
432  outxml.close();
433 }
434 //___________________________________________________________________________
436  ofstream & ofs, string xtag, string ytag, string name) const
437 {
438  string spline_name = (name.size()>0 ? name : fName);
439 
440  // create a spline tag with the number of knots as an attribute
441  int nknots = this->NKnots();
442  string padding = " ";
443  ofs << padding << "<spline name=\"" << spline_name
444  << "\" nknots=\"" << nknots << "\">" << endl;
445 
446  // start printing the knots
447  double x=0, y=0;
448  for(int iknot = 0; iknot < nknots; iknot++)
449  {
450  fInterpolator->GetKnot(iknot, x, y);
451 
452  ofs << std::fixed << setprecision(5);
453  ofs << "\t<knot>"
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;
460  }
461  ofs << padding << "</spline>" << endl;
462 }
463 //___________________________________________________________________________
464 void Spline::SaveAsText(string filename, string format) const
465 {
466  ofstream outtxt(filename.c_str());
467  if(!outtxt.is_open()) {
468  LOG("Spline", pERROR) << "Couldn't create file = " << filename;
469  return;
470  }
471  int nknots = this->NKnots();
472  outtxt << nknots << endl;
473 
474  double x=0, y=0;
475  for(int iknot = 0; iknot < nknots; iknot++) {
476  fInterpolator->GetKnot(iknot, x, y);
477  char line[1024];
478  sprintf(line,format.c_str(),x,y);
479  outtxt << line << endl;
480  }
481  outtxt << endl;
482  outtxt.close();
483 }
484 //___________________________________________________________________________
485 void Spline::SaveAsROOT(string filename, string name, bool recreate) const
486 {
487  string spline_name = (name.size()>0 ? name : fName);
488 
489  string opt = ( (recreate) ? "RECREATE" : "UPDATE" );
490 
491  TFile f(filename.c_str(), opt.c_str());
492  if(fInterpolator) fInterpolator->Write(spline_name.c_str());
493  f.Close();
494 }
495 //___________________________________________________________________________
497  int np, bool scale_with_x, bool in_log, double fx, double fy) const
498 {
499  double xmin = fXMin;
500  double xmax = fXMax;
501 
502  np = TMath::Max(np,2);
503 
504  bool use_log = in_log && (fXMin>0) && (fXMax>0);
505 
506  if(use_log) {
507  xmin = TMath::Log10(xmin);
508  xmax = TMath::Log10(xmax);
509  }
510 
511  double dx = (xmax - xmin) / (np-1);
512 
513  double * x = new double[np];
514  double * y = new double[np];
515 
516  for(int i=0; i<np; i++) {
517  x[i] = ( (use_log) ? TMath::Power(10, xmin+i*dx) : xmin + i*dx );
518  y[i] = this->Evaluate( x[i] );
519 
520  // scale with x if needed
521  if (scale_with_x) y[i] /= x[i];
522 
523  // apply x,y scaling
524  y[i] *= fy;
525  x[i] *= fx;
526  }
527 
528  TGraph * graph = new TGraph(np, x, y);
529  delete[] x;
530  delete[] y;
531  return graph;
532 }
533 //___________________________________________________________________________
535  double x, double & xknot, double & yknot, Option_t * opt) const
536 {
537  string option(opt);
538 
539  bool pos = (option.find("+") != string::npos);
540  bool neg = (option.find("-") != string::npos);
541 
542  if(!pos && !neg) return;
543 
544  int iknot = fInterpolator->FindX(x);
545 
546  double xp=0, yp=0, xn=0, yn=0;
547  fInterpolator->GetKnot(iknot, xn,yn);
548  fInterpolator->GetKnot(iknot+1,xp,yp);
549 
550  bool p = (TMath::Abs(x-xp) < TMath::Abs(x-xn));
551 
552  if(pos&&neg) {
553  if(p) { xknot = xp; yknot = yp; }
554  else { xknot = xn; yknot = yn; }
555  } else {
556  if(pos) { xknot = xp; yknot = yp; }
557  if(neg) { xknot = xn; yknot = yn; }
558  }
559 }
560 //___________________________________________________________________________
561 bool Spline::ClosestKnotValueIsZero(double x, Option_t * opt) const
562 {
563  double xknot = 0, yknot = 0;
564  this->FindClosestKnot(x, xknot, yknot, opt);
565  if(utils::math::AreEqual(yknot,0)) return true;
566  return false;
567 }
568 //___________________________________________________________________________
569 void Spline::Print(ostream & stream) const
570 {
571  int nknots = this->NKnots();
572  double xmin = this->XMin();
573  double xmax = this->XMax();
574 
575  stream << endl;
576  stream << "** Spline: " << this->Name() << endl;
577  stream << "Has " << nknots
578  << " knots in the [" << xmin << ", " << xmax << "] range" << endl;
579  double x,y;
580  for(int i=0; i<nknots; i++) {
581  this->GetKnot(i,x,y);
582  stream << "-- knot : " << i
583  << " -> (x = " << x << ", y = " << y << ")" << endl;
584  }
585 }
586 //___________________________________________________________________________
587 void Spline::Add(const Spline & spl, double c)
588 {
589  // continue only if the input spline is defined at a wider x-range
590  double xmin = this->XMin();
591  double xmax = this->XMax();
592  bool ok = spl.IsWithinValidRange(xmin) && spl.IsWithinValidRange(xmax);
593  if(!ok) {
594  LOG("Spline", pERROR) << "** Can't add splines. X-range mismatch!";
595  return;
596  }
597 
598  int nknots = this->NKnots();
599  double * x = new double[nknots];
600  double * y = new double[nknots];
601 
602  for(int i=0; i<nknots; i++) {
603  this->GetKnot(i,x[i],y[i]);
604  y[i] += (c * spl.Evaluate(x[i]));
605  }
606  this->ResetSpline();
607  this->BuildSpline(nknots,x,y);
608  delete [] x;
609  delete [] y;
610 }
611 //___________________________________________________________________________
612 void Spline::Multiply(const Spline & spl, double c)
613 {
614  // continue only if the input spline is defined at a wider x-range
615  double xmin = this->XMin();
616  double xmax = this->XMax();
617  bool ok = spl.IsWithinValidRange(xmin) && spl.IsWithinValidRange(xmax);
618  if(!ok) {
619  LOG("Spline", pERROR) << "** Can't multiply splines. X-range mismatch!";
620  return;
621  }
622 
623  int nknots = this->NKnots();
624  double * x = new double[nknots];
625  double * y = new double[nknots];
626 
627  for(int i=0; i<nknots; i++) {
628  this->GetKnot(i,x[i],y[i]);
629  y[i] *= (c * spl.Evaluate(x[i]));
630  }
631  this->ResetSpline();
632  this->BuildSpline(nknots,x,y);
633  delete [] x;
634  delete [] y;
635 }
636 //___________________________________________________________________________
637 void Spline::Divide(const Spline & spl, double c)
638 {
639  // continue only if the input spline is defined at a wider x-range
640  double xmin = this->XMin();
641  double xmax = this->XMax();
642  bool ok = spl.IsWithinValidRange(xmin) && spl.IsWithinValidRange(xmax);
643  if(!ok) {
644  LOG("Spline", pERROR) << "** Can't divide splines. X-range mismatch!";
645  return;
646  }
647 
648  int nknots = this->NKnots();
649  double * x = new double[nknots];
650  double * y = new double[nknots];
651 
652  for(int i=0; i<nknots; i++) {
653  this->GetKnot(i,x[i],y[i]);
654  double denom = c * spl.Evaluate(x[i]);
655  bool denom_is_zero = TMath::Abs(denom) < DBL_EPSILON;
656  if(denom_is_zero) {
657  LOG("Spline", pERROR) << "** Refusing to divide spline knot by 0";
658  delete [] x;
659  delete [] y;
660  return;
661  }
662  y[i] /= denom;
663  }
664  this->ResetSpline();
665  this->BuildSpline(nknots,x,y);
666  delete [] x;
667  delete [] y;
668 }
669 //___________________________________________________________________________
670 void Spline::Add(double a)
671 {
672  int nknots = this->NKnots();
673  double * x = new double[nknots];
674  double * y = new double[nknots];
675 
676  for(int i=0; i<nknots; i++) {
677  this->GetKnot(i,x[i],y[i]);
678  y[i]+=a;
679  }
680  this->ResetSpline();
681  this->BuildSpline(nknots,x,y);
682  delete [] x;
683  delete [] y;
684 }
685 //___________________________________________________________________________
686 void Spline::Multiply(double a)
687 {
688  int nknots = this->NKnots();
689  double * x = new double[nknots];
690  double * y = new double[nknots];
691 
692  for(int i=0; i<nknots; i++) {
693  this->GetKnot(i,x[i],y[i]);
694  y[i]*=a;
695  }
696  this->ResetSpline();
697  this->BuildSpline(nknots,x,y);
698  delete [] x;
699  delete [] y;
700 }
701 //___________________________________________________________________________
702 void Spline::Divide(double a)
703 {
704  bool a_is_zero = TMath::Abs(a) < DBL_EPSILON;
705  if(a_is_zero==0) {
706  LOG("Spline", pERROR) << "** Refusing to divide spline by 0";
707  return;
708  }
709  int nknots = this->NKnots();
710  double * x = new double[nknots];
711  double * y = new double[nknots];
712 
713  for(int i=0; i<nknots; i++) {
714  this->GetKnot(i,x[i],y[i]);
715  y[i]/=a;
716  }
717  this->ResetSpline();
718  this->BuildSpline(nknots,x,y);
719  delete [] x;
720  delete [] y;
721 }
722 //___________________________________________________________________________
724 {
725  LOG("Spline", pDEBUG) << "Initializing spline...";
726 
727  fName = "genie-spline";
728  fXMin = 0.0;
729  fXMax = 0.0;
730  fYMax = 0.0;
731 
732  fInterpolator = 0;
733  fInterpolator5 = 0;
734  fGSLInterpolator = 0;
735  fInterpolatorType = "TSpline3";
736 
737  fYCanBeNegative = false;
738 
739  LOG("Spline", pDEBUG) << "...done initializing spline";
740 }
741 //___________________________________________________________________________
743 {
744  if(fInterpolator) delete fInterpolator;
745  if(fInterpolator5) delete fInterpolator5;
747  this->InitSpline();
748 }
749 //___________________________________________________________________________
750 void Spline::BuildSpline(int nentries, double x[], double y[])
751 {
752  LOG("Spline", pDEBUG) << "Building spline...";
753 
754  double xmin = x[0]; // minimum x in spline
755  double xmax = x[nentries-1]; // maximum x in spline
756 
757  fNKnots = nentries;
758  fXMin = xmin;
759  fXMax = xmax;
760  fYMax = y[ TMath::LocMax(nentries, y) ]; // maximum y in spline
761 
762  if(fInterpolator) delete fInterpolator;
763 
764  fInterpolator = new TSpline3("spl3", x, y, nentries, "0");
765 
766  LOG("Spline", pDEBUG) << "...done building spline";
767 }
768 //___________________________________________________________________________
769 void Spline::SetType(string type)
770 {
771  if(!fInterpolator) return;
772 
774 
775  ROOT::Math::Interpolation::Type gsltype;
776 
777  if ( fInterpolatorType == "TSPLINE3" || fInterpolatorType == "" )
778  {
779  fInterpolatorType = "TSpline3";
780  return;
781  }
782  else if ( fInterpolatorType == "TSPLINE5" )
783  {
784  fInterpolatorType = "TSpline5";
785  if(fInterpolator5) delete fInterpolator5;
786  double x[fNKnots], y[fNKnots];
787  for (int i=0; i<fNKnots; i++)
788  {
789  fInterpolator->GetKnot(i, x[i], y[i]);
790  }
791  fInterpolator5 = new TSpline5("spl5", x, y, fNKnots, "0");
792  return;
793  }
794  else if ( fInterpolatorType == "LINEAR" )
795  {
796  gsltype = ROOT::Math::Interpolation::kLINEAR;
797  }
798  else if ( fInterpolatorType == "POLYNOMIAL" )
799  {
800  gsltype = ROOT::Math::Interpolation::kPOLYNOMIAL;
801  }
802  else if ( fInterpolatorType == "CSPLINE" )
803  {
804  gsltype = ROOT::Math::Interpolation::kCSPLINE;
805  }
806  else if ( fInterpolatorType == "CSPLINE_PERIODIC" )
807  {
808  gsltype = ROOT::Math::Interpolation::kCSPLINE_PERIODIC;
809  }
810  else if ( fInterpolatorType == "AKIMA" )
811  {
812  gsltype = ROOT::Math::Interpolation::kAKIMA;
813  }
814  else if ( fInterpolatorType == "AKIMA_PERIODIC" )
815  {
816  gsltype = ROOT::Math::Interpolation::kAKIMA_PERIODIC;
817  }
818  else
819  {
820  fInterpolatorType = "TSpline3";
821  LOG("Spline", pWARN)
822  << "Unknown interpolator type. Setting it to default [TSpline3].";
823  return;
824  }
825 
827  vector<double> x(fNKnots);
828  vector<double> y(fNKnots);
829  for (int i=0; i<fNKnots; i++)
830  {
831  fInterpolator->GetKnot(i, x[i], y[i]);
832  }
833  fGSLInterpolator = new ROOT::Math::Interpolator(x, y, gsltype);
834 
835 
836 }
837 //___________________________________________________________________________
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
#define pERROR
Definition: Messenger.h:59
bool AreEqual(double x1, double x2)
Definition: MathUtils.cxx:236
void SetType(string type)
Definition: Spline.cxx:769
string TrimSpaces(xmlChar *xmls)
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
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 fInterpolatorType
Definition: Spline.h:147
#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:742
const double a
double fXMax
Definition: Spline.h:141
void SaveAsROOT(string filename, string name="", bool recreate=false) const
Definition: Spline.cxx:485
#define pINFO
Definition: Messenger.h:62
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
#define pWARN
Definition: Messenger.h:60
double fYMax
Definition: Spline.h:142
string TrimSpaces(string input)
Definition: StringUtils.cxx:18
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
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
ClassImp(CacheBranchFx)
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
string ToUpper(string input)
Definition: StringUtils.cxx:92
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
string GetAttribute(xmlNodePtr xml_cur, string attr_name)
void Divide(const Spline &spl, double c=1)
Definition: Spline.cxx:637
#define pDEBUG
Definition: Messenger.h:63