MaCh3  2.4.2
Reference Guide
SplineStructs.h
Go to the documentation of this file.
1 #pragma once
2 
3 // MaCh3 includes
6 
7 #include <cmath>
8 
9 #pragma GCC diagnostic push
10 #pragma GCC diagnostic ignored "-Wuseless-cast"
11 #pragma GCC diagnostic ignored "-Wfloat-conversion"
12 #pragma GCC diagnostic ignored "-Wnull-dereference"
13 
18 
19 // *******************
22 // *******************
25  nPts = -999;
26  CurrSegment = 0;
27  splineParsPointer = nullptr;
28  }
29 
31  virtual ~FastSplineInfo() = default;
32 
35 
37  std::vector<M3::float_t> xPts;
38 
41 
43  const double* splineParsPointer;
44 };
45 
46 
47 // ***************************************************************************
49 inline void ApplyKnotWeightCap(TGraph* xsecgraph, const int splineParsIndex, ParameterHandlerGeneric* XsecCov) {
50 // ***************************************************************************
51  if(xsecgraph == nullptr){
52  MACH3LOG_ERROR("hmmm looks like you're trying to apply capping for spline parameter {} but it hasn't been set in xsecgraph yet", XsecCov->GetParFancyName(splineParsIndex));
53  throw MaCh3Exception(__FILE__ , __LINE__ );
54  }
55 
56  // EM: cap the weights of the knots if specified in the config
57  if(
58  (XsecCov->GetParSplineKnotUpperBound(splineParsIndex) != M3::DefSplineKnotUpBound)
59  || (XsecCov->GetParSplineKnotLowerBound(splineParsIndex) != M3::DefSplineKnotLowBound))
60  {
61  for(int knotId = 0; knotId < xsecgraph->GetN(); knotId++){
62  double x,y;
63 
64  // EM: get the x and y of the point. Also double check that the requested point was valid just to be super safe
65  if( xsecgraph->GetPoint(knotId, x, y) == -1) {
66  MACH3LOG_ERROR("Invalid knot requested: {}", knotId);
67  throw MaCh3Exception(__FILE__ , __LINE__ );
68  }
69 
70  y = std::min(y, XsecCov->GetParSplineKnotUpperBound(splineParsIndex));
71  y = std::max(y, XsecCov->GetParSplineKnotLowerBound(splineParsIndex));
72 
73  xsecgraph->SetPoint(knotId, x, y);
74  }
75 
76  // EM: check if our cap made the spline flat, if so we set to 1 knot to avoid problems later on
77  bool isFlat = true;
78  for(int knotId = 0; knotId < xsecgraph->GetN(); knotId++){
79  double x,y;
80 
81  // EM: get the x and y of the point. Also double check that the requested point was valid just to be super safe
82  if( xsecgraph->GetPoint(knotId, x, y) == -1) {
83  MACH3LOG_ERROR("Invalid knot requested: {}", knotId);
84  throw MaCh3Exception(__FILE__ , __LINE__ );
85  }
86  if(std::abs(y - 1.0) > 1e-5) isFlat = false;
87  }
88 
89  if(isFlat){
90  xsecgraph->Set(1);
91  xsecgraph->SetPoint(0, 0.0, 1.0);
92  }
93  }
94 }
95 
96 // ***************************************************************************
98 inline void ApplyKnotWeightCapTSpline3(TSpline3* &Spline, const int splineParsIndex, ParameterHandlerGeneric* XsecCov) {
99 // ***************************************************************************
100  if(Spline == nullptr) {
101  MACH3LOG_ERROR("hmmm looks like you're trying to apply capping for spline parameter {} but it hasn't been set in Spline yet", XsecCov->GetParFancyName(splineParsIndex));
102  throw MaCh3Exception(__FILE__ , __LINE__ );
103  }
104 
105  std::string oldName = Spline->GetName();
106  // EM: cap the weights of the knots if specified in the config
107  if((XsecCov->GetParSplineKnotUpperBound(splineParsIndex) != M3::DefSplineKnotUpBound) || (XsecCov->GetParSplineKnotLowerBound(splineParsIndex) != M3::DefSplineKnotLowBound))
108  {
109  const int NValues = Spline->GetNp();
110  std::vector<double> XVals(NValues);
111  std::vector<double> YVals(NValues);
112  for(int knotId = 0; knotId < NValues; knotId++){
113  double x,y;
114  // Extract X and Y
115  Spline->GetKnot(knotId, x, y);
116 
117  y = std::min(y, XsecCov->GetParSplineKnotUpperBound(splineParsIndex));
118  y = std::max(y, XsecCov->GetParSplineKnotLowerBound(splineParsIndex));
119 
120  XVals[knotId] = x;
121  YVals[knotId] = y;
122  }
123  delete Spline;
124  // If we capped we have to make new TSpline3 to recalculate coefficients
125  Spline = new TSpline3(oldName.c_str(), XVals.data(), YVals.data(), NValues);
126  }
127 }
128 
129 // ************************
132 // ************************
133 public:
139  virtual double Eval(const double var)=0;
141  virtual void Print()=0;
143  virtual M3::int_t GetNp()=0;
144 };
145 
146 // ************************
149 // ************************
150 public:
153  length = 0;
154  Par = nullptr;
155  }
156 
158  virtual ~TF1_red() {
159  if (Par != NULL) {
160  delete[] Par;
161  Par = nullptr;
162  }
163  }
164 
167  length = nSize;
168  for (int i = 0; i < length; ++i) {
169  Par[i] = Array[i];
170  }
171  }
172 
174  TF1_red(TF1* &Function) : TResponseFunction_red() {
175  Par = nullptr;
176  SetFunc(Function);
177  }
178 
180  inline void SetFunc(TF1* &Func) {
181  length = M3::int_t(Func->GetNpar());
182  if (Par != nullptr) delete[] Par;
183  Par = new M3::float_t[length];
184  for (int i = 0; i < length; ++i) {
185  Par[i] = M3::float_t(Func->GetParameter(i));
186  }
187  delete Func;
188  Func = nullptr;
189  }
190 
192  inline double Eval(const double var) override {
193  return Par[1]+Par[0]*var;
194 
195  /* FIXME in future we might introduce more TF1
196  //If we have 5 parameters we're using a fifth order polynomial
197  if (Type == kFifthOrderPolynomial) {
198  return 1+Par[0]*var+Par[1]*var*var+Par[2]*var*var*var+Par[3]*var*var*var*var+Par[4]*var*var*var*var*var;
199  } else if (Type == kTwoLinears) {
200  return (var <= 0)*(Par[2]+Par[0]*var)+(var > 0)*(Par[2]+Par[1]*var);
201  } else if (Type == kLinear) {
202  return (Par[1]+Par[0]*var);
203  } else if (Type == kPseudoHeaviside) {
204  return (var <= 0)*(1+Par[0]*var) + (1 >= var)*(var > 0)*(1+Par[1]*var) + (var > 1)*(Par[3]+Par[2]*var);
205  }else {
206  MACH3LOG_ERROR(" Class only knows about 5th order polynomial, two superposed linear functions, linear function, or pseudo Heaviside.");
207  MACH3LOG_ERROR(" You have tried something else than this, which remains unimplemented.");
208  MACH3LOG_ERROR("{}: {}", __FILE__, __LINE__);
209  throw MaCh3Exception(__FILE__ , __LINE__ );
210  }
211  */
212  }
213 
215  inline void SetParameter(M3::int_t Parameter, M3::float_t Value) {
216  Par[Parameter] = Value;
217  }
218 
220  double GetParameter(M3::int_t Parameter) {
221  if (Parameter > length) {
222  MACH3LOG_ERROR("You requested parameter number {} but length is {} parameters", Parameter, length);
223  throw MaCh3Exception(__FILE__ , __LINE__ );
224  return -999.999;
225  }
226  return Par[Parameter];
227  }
228 
230  inline void SetSize(M3::int_t nSpline) {
231  length = nSpline;
232  Par = new M3::float_t[length];
233  }
235  inline int GetSize() { return length; }
237  inline void Print() override {
238  MACH3LOG_INFO("Printing TF1_red:");
239  MACH3LOG_INFO(" Length = {}", length);
240  for (int i = 0; i < length; i++) {
241  MACH3LOG_INFO(" Coeff {} = {}", i, Par[i]);
242  }
243  }
244 
246  inline TF1* ConstructTF1(const std::string& function, const int xmin, const int xmax) {
247  TF1 *func = new TF1("TF1", function.c_str(), xmin, xmax);
248  for(int i = 0; i < length; ++i) {
249  func->SetParameter(i, Par[i]);
250  }
251  return func;
252  }
253 
255  inline M3::int_t GetNp() override { return length; }
256 
257 private:
261 };
262 
263 // ************************
266 // ************************
267 public:
270  nPoints = 0;
271  Par = nullptr;
272  XPos = nullptr;
273  YResp = nullptr;
274  }
275 
277  TSpline3_red(TSpline3* &spline, SplineInterpolation InterPolation = kTSpline3) : TResponseFunction_red() {
278  Par = nullptr;
279  XPos = nullptr;
280  YResp = nullptr;
281  SetFunc(spline, InterPolation);
282  }
283 
286  nPoints = N;
287  // Save the parameters for each knot
288  Par = new M3::float_t*[nPoints];
289  // Save the positions of the knots
290  XPos = new M3::float_t[nPoints];
291  // Save the y response at each knot
292  YResp = new M3::float_t[nPoints];
293  for(int j = 0; j < N; ++j){
294  Par[j] = new M3::float_t[3];
295  Par[j][0] = P[j][0];
296  Par[j][1] = P[j][1];
297  Par[j][2] = P[j][2];
298  XPos[j] = X[j];
299  YResp[j] = Y[j];
300 
301  if((Par[j][0] == -999) | (Par[j][1] ==-999) | (Par[j][2] ==-999) | (XPos[j] ==-999) | (YResp[j] ==-999)){
302  MACH3LOG_ERROR("******************* Bad parameter values when constructing TSpline3_red *********************");
303  MACH3LOG_ERROR("Passed values (i, x, y, b, c, d): {}, {}, {}, {}, {}, {}", j, X[j], Y[j], P[j][0], P[j][1], P[j][2]);
304  MACH3LOG_ERROR("Set values (i, x, y, b, c, d): {}, {}, {}, {}, {}, {}", j, XPos[j], YResp[j], Par[j][0], Par[j][1], Par[j][2]);
305  MACH3LOG_ERROR("*********************************************************************************************");
306  }
307  }
308  }
311  inline void SetFunc(TSpline3* &spline, SplineInterpolation InterPolation = kTSpline3) {
312  nPoints = M3::int_t(spline->GetNp());
313  if (Par != nullptr) {
314  for (int i = 0; i < nPoints; ++i) {
315  delete[] Par[i];
316  Par[i] = nullptr;
317  }
318  delete[] Par;
319  Par = nullptr;
320  }
321  if (XPos != nullptr) delete[] XPos;
322  if (YResp != nullptr) delete[] YResp;
323  // Save the parameters for each knot
324  Par = new M3::float_t*[nPoints];
325  // Save the positions of the knots
326  XPos = new M3::float_t[nPoints];
327  // Save the y response at each knot
328  YResp = new M3::float_t[nPoints];
329 
330  //KS: Default TSpline3 ROOT implementation
331  if(InterPolation == kTSpline3)
332  {
333  for (int i = 0; i < nPoints; ++i) {
334  // 3 is the size of the TSpline3 coefficients
335  Par[i] = new M3::float_t[3];
336  double x = -999.99, y = -999.99, b = -999.99, c = -999.99, d = -999.99;
337  spline->GetCoeff(i, x, y, b, c, d);
338  XPos[i] = M3::float_t(x);
339  YResp[i] = M3::float_t(y);
340  Par[i][0] = M3::float_t(b);
341  Par[i][1] = M3::float_t(c);
342  Par[i][2] = M3::float_t(d);
343  }
344  }
345  //CW: Reduce to use linear spline interpolation for certain parameters
346  // Not the most elegant way: use TSpline3 object but set coefficients to zero and recalculate spline points; the smart way (but more human intensive) would be to save memory here and simply not store the zeros at all
347  // Get which parameters should be linear from the fit manager
348  // Convert the spline number to global xsec parameter
349  // Loop over the splines points
350  // KS: kLinearFunc should be used with TF1, this is just as safety
351  else if(InterPolation == kLinear || InterPolation == kLinearFunc)
352  {
353  for (int k = 0; k < nPoints; ++k) {
354  Par[k] = new M3::float_t[3];
355 
356  Double_t x1, y1, b1, c1, d1, x2, y2, b2, c2, d2 = 0;
357  spline->GetCoeff(k, x1, y1, b1, c1, d1);
358 
359  double tempb = 0;
360  if (k == nPoints - 1) {
361  tempb = Par[k-1][0];
362  } else {
363  spline->GetCoeff(k + 1, x2, y2, b2, c2, d2);
364  tempb = (y2-y1)/(x2-x1);
365  }
366  XPos[k] = M3::float_t(x1);
367  YResp[k] = M3::float_t(y1);
368  Par[k][0] = M3::float_t(tempb); // linear slope
369  Par[k][1] = M3::float_t(0);
370  Par[k][2] = M3::float_t(0);
371  }
372  }
373  //EM: Akima spline is similar to regular cubic spline but is allowed to be discontinuous in 2nd derivative and coefficients in any segment
374  // only depend on th 2 nearest points on either side
375  else if(InterPolation == kAkima)
376  {
377  // get the knot values for the spline
378  for (int i = 0; i < nPoints; ++i) {
379  // 3 is the size of the TSpline3 coefficients
380  Par[i] = new M3::float_t[3];
381 
382  double x = -999.99, y = -999.99;
383  spline->GetKnot(i, x, y);
384 
385  XPos[i] = M3::float_t(x);
386  YResp[i] = M3::float_t(y);
387  }
388 
389  M3::float_t* mvals = new M3::float_t[nPoints + 3];
390  M3::float_t* svals = new M3::float_t[nPoints + 1];
391 
392  for (int i = -2; i <= nPoints; ++i) {
393  // if segment is first or last or 2nd to first or last, needs to be dealt with slightly differently;
394  // need to estimate the values for additinal points which would lie outside of the spline
395  if(i ==-2){
396  mvals[i+2] = M3::float_t(3.0 * (YResp[1] - YResp[0]) / (XPos[1] - XPos[0]) - 2.0*(YResp[2] - YResp[1]) / (XPos[2] - XPos[1]));
397  }
398  else if(i==-1){
399  mvals[i+2] = M3::float_t(2.0 * (YResp[1] - YResp[0]) / (XPos[1] - XPos[0]) - (YResp[2] - YResp[1]) / (XPos[2] - XPos[1]));
400  }
401  else if(i==nPoints){
402  mvals[i+2] = M3::float_t(3.0 * (YResp[nPoints-1] - YResp[nPoints-2]) / (XPos[nPoints-1] - XPos[nPoints-2]) - 2.0*(YResp[nPoints-2] - YResp[nPoints-3]) / (XPos[nPoints-2] - XPos[nPoints-3]));
403  }
404  else if(i == nPoints - 1){
405  mvals[i+2] = M3::float_t(2.0 * (YResp[nPoints-1] - YResp[nPoints-2]) / (XPos[nPoints-1] - XPos[nPoints-2]) - (YResp[nPoints-2] - YResp[nPoints-3]) / (XPos[nPoints-2] - XPos[nPoints-3]));
406  }
407  //standard internal segment
408  else{
409  mvals[i+2] = (YResp[i+1] - YResp[i])/ (XPos[i+1] - XPos[i]);
410  }
411  }
412 
413  for(int i =2; i<=nPoints+2; i++){
414  if (std::abs(mvals[i+1] - mvals[i]) + std::abs(mvals[i-1] - mvals[i-2]) != 0.0){
415  svals[i-2] = (std::abs(mvals[i+1] - mvals[i]) * mvals[i-1] + std::abs(mvals[i-1] - mvals[i-2]) *mvals[i]) / (std::abs(mvals[i+1] - mvals[i]) + std::abs(mvals[i-1] - mvals[i-2]));
416  }
417  else{svals[i-2] = mvals[i];}
418  }
419 
420  // calculate the coefficients for the spline
421  for(int i = 0; i <nPoints; i++){
422  M3::float_t b, c, d = M3::float_t(-999.999);
423 
424  b = svals[i];
425  c = M3::float_t(3.0* (YResp[i+1] - YResp[i]) / (XPos[i+1] - XPos[i]) -2.0 *svals[i] - svals[i +1]) /(XPos[i+1] - XPos[i]);
426  d = M3::float_t((svals[i + 1] +svals[i]) - 2.0*(YResp[i+1] - YResp[i]) / (XPos[i+1] - XPos[i])) / ((XPos[i+1] - XPos[i]) * (XPos[i+1] - XPos[i]));
427 
428  Par[i][0] = b;
429  Par[i][1] = c;
430  Par[i][2] = d;
431  }
432 
433  // check the input spline for linear segments, if there are any then overwrite the calculated coefficients
434  // this will pretty much only ever be the case if they are set to be linear in samplePDFND i.e. the user wants it to be linear
435  for(int i = 0; i <nPoints-1; i++){
436  double x = -999.99, y = -999.99, b = -999.99, c = -999.99, d = -999.99;
437  spline->GetCoeff(i, x, y, b, c, d);
438 
439  if((c == 0.0 && d == 0.0)){
440  Par[i][0] = M3::float_t(b);
441  Par[i][1] = M3::float_t(0.0);
442  Par[i][2] = M3::float_t(0.0);
443  }
444  }
445  delete[] mvals;
446  delete[] svals;
447  }
448  //EM: Monotone spline is similar to regular cubic spline but enforce the condition that the interpolated value at any point
449  // must be between its two nearest knots, DOES NOT make the entire spline monotonic, only the segments
450  else if(InterPolation == kMonotonic)
451  {
452  // values of the secants at each point (for calculating monotone spline)
453  M3::float_t * Secants = new M3::float_t[nPoints -1];
454  // values of the tangens at each point (for calculating monotone spline)
455  M3::float_t * Tangents = new M3::float_t[nPoints];
456 
457  // get the knot values for the spline
458  for (int i = 0; i < nPoints; ++i) {
459  // 3 is the size of the TSpline3 coefficients
460  Par[i] = new M3::float_t[3];
461 
462  double x = -999.99, y = -999.99;
463  spline->GetKnot(i, x, y);
464 
465  XPos[i] = M3::float_t(x);
466  YResp[i] = M3::float_t(y);
467 
468  Tangents[i] = 0.0;
469  }
470 
471  // deal with the case of two points (just do linear interpolation between them)
472  if (nPoints == 2){
473  Par[0][0] = (YResp[1] - YResp[0]) / (XPos[1] - XPos[0]);
474  Par[0][1] = 0.0;
475  Par[0][2] = 0.0;
476  // extra "virtual" segment at end to make Par array shape fit with knot arrays shapes
477  Par[1][1] = 0.0;
478  Par[1][2] = 0.0;
479 
480  return;
481  } // if nPoints !=2 do full monotonic spline treatment:
482  else
483  {
484  // first pass over knots to calculate the secants
485  for (int i = 0; i < nPoints-1; ++i) {
486  Secants[i] = (YResp[i+1] - YResp[i]) / (XPos[i+1] - XPos[i]);
487  MACH3LOG_TRACE("Secant {}: {}", i, Secants[i]);
488  }
489 
490  Tangents[0] = Secants[0];
491  Tangents[nPoints-1] = Secants[nPoints -2];
492 
493  M3::float_t alpha;
494  M3::float_t beta;
495 
496  // second pass over knots to calculate tangents
497  for (int i = 1; i < nPoints-1; ++i) {
498  if ((Secants[i-1] >= 0.0 && Secants[i] >= 0.0) | (Secants[i-1] < 0.0 && Secants[i] < 0.0)){ //check for same sign
499  Tangents[i] = M3::float_t((Secants[i-1] + Secants[i]) /2.0);
500  }
501  }
502 
503  // third pass over knots to rescale tangents
504  for (int i = 0; i < nPoints-1; ++i) {
505  if (Secants[i] == 0.0){
506  Tangents[i] = 0.0;
507  Tangents[i+1] = 0.0;
508  }
509 
510  else{
511  alpha = Tangents[i] / Secants[i];
512  beta = Tangents[i+1] / Secants[i];
513 
514  if (alpha <0.0){
515  Tangents[i] = 0.0;
516  }
517  if (beta < 0.0){
518  Tangents[i+1] = 0.0;
519  }
520 
521  if (alpha * alpha + beta * beta >9.0){
522  M3::float_t tau = M3::float_t(3.0 / std::sqrt(alpha * alpha + beta * beta));
523  Tangents[i] = tau * alpha * Secants[i];
524  Tangents[i+1] = tau * beta * Secants[i];
525  }
526  }
527  } // finished rescaling tangents
528  // fourth pass over knots to calculate the coefficients for the spline
529  M3::float_t dx;
530  for(int i = 0; i <nPoints-1; i++){
531  M3::float_t b, c, d = M3::float_t(-999.999);
532  dx = XPos[i+1] - XPos[i];
533 
534  b = Tangents[i] * dx;
535  c = M3::float_t(3.0* (YResp[i+1] - YResp[i]) -2.0 *dx * Tangents[i] - dx * Tangents[i +1]);
536  d = M3::float_t(2.0* (YResp[i] - YResp[i+1]) + dx * (Tangents[i] + Tangents[i+1]));
537 
538  Par[i][0] = b / dx;
539  Par[i][1] = c / (dx * dx);
540  Par[i][2] = d / (dx * dx * dx);
541 
542  if((Par[i][0] == -999) | (Par[i][1] == -999) | (Par[i][2] ==-999) | (Par[i][0] == -999.999) | (Par[i][1] == -999.999) | (Par[i][2] ==-999.999)){
543  MACH3LOG_INFO("Bad spline parameters for segment {}: (b, c, d) = {}, {}, {}. This will cause problems with GPU.",
544  i, Par[i][0], Par[i][1], Par[i][2]);
545  }
546  MACH3LOG_TRACE("b: {}", b);
547  MACH3LOG_TRACE("dx: {}, x_0: {}, x_1: {}", dx, XPos[i], XPos[i+1]);
548  MACH3LOG_TRACE(" y_0: {}, y_1: {}", YResp[i], YResp[i+1]);
549  }
550 
551  // include params for final "segment" outside of the spline so that par array fits with x and y arrays,
552  // should never actually get used but if not set then the GPU code gets very angry
553  Par[nPoints-1][0] = 0.0;
554  Par[nPoints-1][1] = 0.0;
555  Par[nPoints-1][2] = 0.0;
556 
557  // check the input spline for linear segments, if there are any then overwrite the calculated coefficients
558  // this will pretty much only ever be the case if they are set to be linear in samplePDFND i.e. the user wants it to be linear
559  for(int i = 0; i <nPoints-1; i++){
560  double x = -999.99, y = -999.99, b = -999.99, c = -999.99, d = -999.99;
561  spline->GetCoeff(i, x, y, b, c, d);
562 
563  if((c == 0.0 && d == 0.0)){
564  Par[i][0] = M3::float_t(b);
565  Par[i][1] = 0.0;
566  Par[i][2] = 0.0;
567  }
568  }
569  delete[] Secants;
570  delete[] Tangents;
571  } // end of if(nPoints !=2)
572  }
573  else if (InterPolation == kKochanekBartels)
574  {
575  // Allocate memory for tangents and coefficients
576  M3::float_t* Tangents = new M3::float_t[nPoints];
577  for (int i = 0; i < nPoints; ++i)
578  {
579  Par[i] = new M3::float_t[3];
580  double x = -999.99, y = -999.99;
581  spline->GetKnot(i, x, y);
582  XPos[i] = M3::float_t(x);
583  YResp[i] = M3::float_t(y);
584  Tangents[i] = 0.0;
585  }
586 
587  // KS: Setting all parameters to 0 gives a Catmull-Rom spline.
588  // Tension (T):
589  // - T = 0.0: Default smooth interpolation (Catmull-Rom).
590  // - T = 1.0: Maximum tension; curve becomes linear between knots (no curvature).
591  // - T = -1.0: Overshooting; creates loops or "bouncy" effects around knots.
592  constexpr M3::float_t T = 0.0;
593 
594  // Continuity (C):
595  // - C = 0.0: Default smooth transition (continuous first derivative).
596  // - C = 1.0: Sharp corner/crease at knot (discontinuous first derivative).
597  // - C = -1.0: Inverted curve; creates loops or kinks at knot.
598  constexpr M3::float_t C = 0.0;
599 
600  // Bias (B):
601  // - B = 0.0: Default symmetric curve around knot.
602  // - B = 1.0: Curve is "pulled" toward the next knot (leading effect).
603  // - B = -1.0: Curve is "pulled" toward the previous knot (lagging effect).
604  constexpr M3::float_t B = 0.0;
605 
606  // Calculate tangents for internal knots
607  for (int i = 1; i < nPoints - 1; ++i)
608  {
609  M3::float_t d0 = (YResp[i] - YResp[i-1]) / (XPos[i] - XPos[i-1]);
610  M3::float_t d1 = (YResp[i+1] - YResp[i]) / (XPos[i+1] - XPos[i]);
611 
612  M3::float_t term1 = (1.0 - T) * (1.0 + B) * (1.0 + C) * 0.5 * d0;
613  M3::float_t term2 = (1.0 - T) * (1.0 - B) * (1.0 - C) * 0.5 * d1;
614 
615  Tangents[i] = term1 + term2;
616  }
617 
618  // Boundary conditions (simple choice: secant slopes)
619  Tangents[0] = (YResp[1] - YResp[0]) / (XPos[1] - XPos[0]);
620  Tangents[nPoints-1]= (YResp[nPoints-1] - YResp[nPoints-2]) / (XPos[nPoints-1] - XPos[nPoints-2]);
621 
622  // Compute cubic coefficients for each segment
623  for (int i = 0; i < nPoints - 1; ++i)
624  {
625  M3::float_t dx = XPos[i+1] - XPos[i];
626  M3::float_t dy = YResp[i+1] - YResp[i];
627  M3::float_t m0 = Tangents[i];
628  M3::float_t m1 = Tangents[i+1];
629 
630  Par[i][0] = m0; // b
631  Par[i][1] = (3*dy/(dx*dx)) - (2*m0 + m1)/dx; // c
632  Par[i][2] = (m0 + m1 - 2*dy/dx) / (dx*dx); // d
633 
634  MACH3LOG_TRACE("KB segment {}: dx={}, dy={}, b={}, c={}, d={}",
635  i, dx, dy, Par[i][0], Par[i][1], Par[i][2]);
636  }
637  // Last "virtual" segment
638  Par[nPoints-1][0] = 0.0;
639  Par[nPoints-1][1] = 0.0;
640  Par[nPoints-1][2] = 0.0;
641 
642  delete[] Tangents;
643  }
644  else
645  {
646  MACH3LOG_ERROR("Unsupported interpolation type {}", static_cast<int>(InterPolation));
647  throw MaCh3Exception(__FILE__ , __LINE__ );
648  }
649 
650  delete spline;
651  spline = nullptr;
652  }
653 
655  virtual ~TSpline3_red() {
656  if(Par != nullptr) {
657  for (int i = 0; i < nPoints; ++i) {
658  if (Par[i] != nullptr) {
659  delete[] Par[i];
660  }
661  }
662  delete[] Par;
663  }
664  if(XPos != nullptr) delete[] XPos;
665  if(YResp != nullptr) delete[] YResp;
666  Par = nullptr;
667  XPos = YResp = nullptr;
668  }
669 
673  inline int FindX(double x) {
674  // The segment we're interested in (klow in ROOT code)
675  int segment = 0;
676  int kHigh = nPoints-1;
677  // If the variation is below the lowest saved spline point
678  if (x <= XPos[0]){
679  segment = 0;
680  // If the variation is above the highest saved spline point
681  } else if (x >= XPos[nPoints-1]) {
682  // Yes, the -2 is indeed correct, see TSpline.cxx:814 and //see: https://savannah.cern.ch/bugs/?71651
683  segment = kHigh;
684  // If the variation is between the maximum and minimum, perform a binary search
685  } else {
686  // The top point we've got
687  int kHalf = 0;
688  // While there is still a difference in the points (we haven't yet found the segment)
689  // This is a binary search, incrementing segment and decrementing kHalf until we've found the segment
690  while (kHigh - segment > 1) {
691  // Increment the half-step
692  kHalf = (segment + kHigh)/2;
693  // If our variation is above the kHalf, set the segment to kHalf
694  if (x > XPos[kHalf]) {
695  segment = kHalf;
696  // Else move kHigh down
697  } else {
698  kHigh = kHalf;
699  }
700  } // End the while: we've now done our binary search
701  } // End the else: we've now found our point
702  if (segment >= nPoints-1 && nPoints > 1) segment = nPoints-2;
703  return segment;
704  }
705 
707  inline double Eval(double var) override {
708  // Get the segment for this variation
709  int segment = FindX(var);
710  // The get the coefficients for this variation
711  M3::float_t x = M3::float_t(-999.99), y = M3::float_t(-999.99), b = M3::float_t(-999.99), c = M3::float_t(-999.99), d = M3::float_t(-999.99);
712  GetCoeff(segment, x, y, b, c, d);
713  double dx = var - x;
714  // Evaluate the third order polynomial
715  double weight = y+dx*(b+dx*(c+d*dx));
716  return weight;
717  }
718 
720  inline M3::int_t GetNp() override { return nPoints; }
721  // Get the ith knot's x and y position
722  inline void GetKnot(int i, M3::float_t &xtmp, M3::float_t &ytmp) {
723  xtmp = XPos[i];
724  ytmp = YResp[i];
725  }
726 
728  inline void GetCoeff(int segment, M3::float_t &x, M3::float_t &y, M3::float_t &b, M3::float_t &c, M3::float_t &d) {
729  b = Par[segment][0];
730  c = Par[segment][1];
731  d = Par[segment][2];
732  x = XPos[segment];
733  y = YResp[segment];
734  }
735 
737  inline TSpline3* ConstructTSpline3() {
738  // KS: Sadly ROOT only accepts double...
739  #ifdef _LOW_MEMORY_STRUCTS_
740  std::vector<Double_t> xPosDoubles(nPoints);
741  std::vector<Double_t> yPosDoubles(nPoints);
742  for (Int_t i = 0; i < nPoints; ++i) {
743  xPosDoubles[i] = static_cast<Double_t>(XPos[i]); // Convert float to double
744  yPosDoubles[i] = static_cast<Double_t>(YResp[i]); // Convert float to double
745  }
746  TSpline3 *spline = new TSpline3("Spline", xPosDoubles.data(), yPosDoubles.data(), static_cast<int>(nPoints));
747  #else
748  TSpline3 *spline = new TSpline3("Spline", XPos, YResp, nPoints);
749  #endif
750  for (Int_t i = 0; i < nPoints; ++i) {
751  spline->SetPointCoeff(i, Par[i][0], Par[i][1], Par[i][2]);
752  }
753 
754  return spline;
755  }
756 
758  inline void Print() override {
759  MACH3LOG_INFO("Printing TSpline_red:");
760  MACH3LOG_INFO(" Nknots = {}", nPoints);
761  for (int i = 0; i < nPoints; ++i) {
762  MACH3LOG_INFO(" i = {} x = {} y = {} b = {} c = {} d = {}",
763  i, XPos[i], YResp[i], Par[i][0], Par[i][1], Par[i][2]);
764  }
765  }
766 
767  protected: //changed to protected from private so can be accessed by derived classes
776 };
777 
778 // *****************************************
781 inline bool isFlat(TSpline3_red* &spl) {
782 // *****************************************
783  int Np = spl->GetNp();
784  M3::float_t x, y, b, c, d;
785  // Go through spline segment parameters,
786  // Get y values for each spline knot,
787  // Every knot must evaluate to 1.0 to create a flat spline
788  for(int i = 0; i < Np; i++) {
789  spl->GetCoeff(i, x, y, b, c, d);
790  if (y != 1) {
791  return false;
792  }
793  }
794  return true;
795 }
796 
797 // *********************************
800 inline std::vector<std::vector<TSpline3_red*> > ReduceTSpline3(std::vector<std::vector<TSpline3*> > &MasterSpline) {
801 // *********************************
802  std::vector<std::vector<TSpline3*> >::iterator OuterIt;
803  std::vector<TSpline3*>::iterator InnerIt;
804 
805  // The return vector
806  std::vector<std::vector<TSpline3_red*> > ReducedVector;
807  ReducedVector.reserve(MasterSpline.size());
808 
809  // Loop over each parameter
810  int OuterCounter = 0;
811  for (OuterIt = MasterSpline.begin(); OuterIt != MasterSpline.end(); ++OuterIt, ++OuterCounter) {
812  // Make the temp vector
813  std::vector<TSpline3_red*> TempVector;
814  TempVector.reserve(OuterIt->size());
815  int InnerCounter = 0;
816  // Loop over each TSpline3 pointer
817  for (InnerIt = OuterIt->begin(); InnerIt != OuterIt->end(); ++InnerIt, ++InnerCounter) {
818  // Here's our delicious TSpline3 object
819  TSpline3 *spline = (*InnerIt);
820  // Now make the reduced TSpline3 pointer
821  TSpline3_red *red = nullptr;
822  if (spline != NULL) {
823  red = new TSpline3_red(spline);
824  (*InnerIt) = spline;
825  }
826  // Push back onto new vector
827  TempVector.push_back(red);
828  } // End inner for loop
829  ReducedVector.push_back(TempVector);
830  } // End outer for loop
831  // Now have the reduced vector
832  return ReducedVector;
833 }
834 
835 // *********************************
838 inline std::vector<std::vector<TF1_red*> > ReduceTF1(std::vector<std::vector<TF1*> > &MasterSpline) {
839 // *********************************
840  std::vector<std::vector<TF1*> >::iterator OuterIt;
841  std::vector<TF1*>::iterator InnerIt;
842 
843  // The return vector
844  std::vector<std::vector<TF1_red*> > ReducedVector;
845  ReducedVector.reserve(MasterSpline.size());
846 
847  // Loop over each parameter
848  int OuterCounter = 0;
849  for (OuterIt = MasterSpline.begin(); OuterIt != MasterSpline.end(); ++OuterIt, ++OuterCounter) {
850  // Make the temp vector
851  std::vector<TF1_red*> TempVector;
852  TempVector.reserve(OuterIt->size());
853  int InnerCounter = 0;
854  // Loop over each TSpline3 pointer
855  for (InnerIt = OuterIt->begin(); InnerIt != OuterIt->end(); ++InnerIt, ++InnerCounter) {
856  // Here's our delicious TSpline3 object
857  TF1* spline = (*InnerIt);
858  // Now make the reduced TSpline3 pointer (which deleted TSpline3)
859  TF1_red* red = nullptr;
860  if (spline != NULL) {
861  red = new TF1_red(spline);
862  (*InnerIt) = spline;
863  }
864  // Push back onto new vector
865  TempVector.push_back(red);
866  } // End inner for loop
867  ReducedVector.push_back(TempVector);
868  } // End outer for loop
869  // Now have the reduced vector
870  return ReducedVector;
871 }
872 
873 // *********************************
880  const RespFuncType SplineRespFuncType,
881  const SplineInterpolation SplineInterpolationType,
882  const std::string& Title) {
883 // *********************************
884  TResponseFunction_red* RespFunc = nullptr;
885 
886  if (graph && graph->GetN() > 1)
887  {
888  if(SplineRespFuncType == kTSpline3_red)
889  {
890  // Here's the TSpline3
891  TSpline3* spline = nullptr;
892  TSpline3_red *spline_red = nullptr;
893 
894  // Create the TSpline3* from the TGraph* and build the coefficients
895  spline = new TSpline3(Title.c_str(), graph);
896  spline->SetNameTitle(Title.c_str(), Title.c_str());
897 
898  // Make the reduced TSpline3 format and delete the old spline
899  spline_red = new TSpline3_red(spline, SplineInterpolationType);
900 
901  RespFunc = spline_red;
902  }
903  else if(SplineRespFuncType == kTF1_red)
904  {
905  // The TF1 object we build from fitting the TGraph
906  TF1 *Fitter = nullptr;
907  TF1_red *tf1_red = nullptr;
908 
909  if(graph->GetN() != 2) {
910  MACH3LOG_ERROR("Trying to make TF1 from more than 2 knots. Knots = {}", graph->GetN());
911  MACH3LOG_ERROR("Currently support only linear with 2 knots :(");
912  throw MaCh3Exception(__FILE__ , __LINE__ );
913  }
914 
915  // Try simple linear function
916  Fitter = new TF1(Title.c_str(), "([1]+[0]*x)", graph->GetX()[0], graph->GetX()[graph->GetN()-1]);
917  //CW: For 2p2h shape C and O we can't fit a polynomial: try a linear combination of two linear functions around 0
918  //Fitter = new TF1(Title.c_str(), "(x<=0)*(1+[0]*x)+(x>0)*([1]*x+1)", graph->GetX()[0], graph->GetX()[graph->GetN()-1]);
919  // Fit 5hd order polynomial for all other parameters
920  //Fitter = new TF1(Title.c_str(), "1+[0]*x+[1]*x*x+[2]*x*x*x+[3]*x*x*x*x+[4]*x*x*x*x*x", graph->GetX()[0], graph->GetX()[graph->GetN()-1]);
921  //Pseudo Heaviside for Pauli Blocking
922  //Fitter = new TF1(Title.c_str(), "(x <= 0)*(1+[0]*x) + (1 >= x)*(x > 0)*(1+[1]*x) + (x > 1)*([3]+[2]*x)", graph->GetX()[0], graph->GetX()[graph->GetN()-1]);
923 
924  // Fit the TF1 to the graph
925  graph->Fit(Fitter, "Q0");
926  // Make the reduced TF1 if we want
927  tf1_red = new TF1_red(Fitter);
928 
929  RespFunc = tf1_red;
930  }
931  else
932  {
933  MACH3LOG_ERROR("Unsupported response function type");
934  throw MaCh3Exception(__FILE__ , __LINE__ );
935  }
936  }
937  else
938  {
939  RespFunc = nullptr;
940  }
941  return RespFunc;
942 }
943 
944 #pragma GCC diagnostic pop
945 
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
#define MACH3LOG_TRACE
Definition: MaCh3Logger.h:33
Definitions of generic parameter structs and utility templates for MaCh3.
SplineInterpolation
Make an enum of the spline interpolation type.
@ kTSpline3
Default TSpline3 interpolation.
@ kMonotonic
EM: DOES NOT make the entire spline monotonic, only the segments.
@ kKochanekBartels
KS: Kochanek-Bartels spline: allows local control of tension, continuity, and bias.
@ kLinear
Linear interpolation between knots.
@ kLinearFunc
Liner interpolation using TF1 not spline.
@ kAkima
EM: Akima spline iis allowed to be discontinuous in 2nd derivative and coefficients in any segment.
RespFuncType
Make an enum of the spline interpolation type.
@ kTF1_red
Uses TF1_red for interpolation.
@ kTSpline3_red
Uses TSpline3_red for interpolation.
TResponseFunction_red * CreateResponseFunction(TGraph *&graph, const RespFuncType SplineRespFuncType, const SplineInterpolation SplineInterpolationType, const std::string &Title)
KS: Create Response Function using TGraph.
void ApplyKnotWeightCap(TGraph *xsecgraph, const int splineParsIndex, ParameterHandlerGeneric *XsecCov)
EM: Apply capping to knot weight for specified spline parameter. param graph needs to have been set i...
Definition: SplineStructs.h:49
void ApplyKnotWeightCapTSpline3(TSpline3 *&Spline, const int splineParsIndex, ParameterHandlerGeneric *XsecCov)
EM: Apply capping to knot weight for specified spline parameter. param graph needs to have been set i...
Definition: SplineStructs.h:98
std::vector< std::vector< TF1_red * > > ReduceTF1(std::vector< std::vector< TF1 * > > &MasterSpline)
CW: Reduced the TF1 to TF1_red.
std::vector< std::vector< TSpline3_red * > > ReduceTSpline3(std::vector< std::vector< TSpline3 * > > &MasterSpline)
CW: Reduced the TSpline3 to TSpline3_red.
bool isFlat(TSpline3_red *&spl)
CW: Helper function used in the constructor, tests to see if the spline is flat.
Custom exception class used throughout MaCh3.
std::string GetParFancyName(const int i) const
Get fancy name of the Parameter.
Class responsible for handling of systematic error parameters with different types defined in the con...
double GetParSplineKnotUpperBound(const int i) const
EM: value at which we cap spline knot weight.
double GetParSplineKnotLowerBound(const int i) const
EM: value at which we cap spline knot weight.
CW: A reduced TF1 class only. Only saves parameters for each TF1 and how many parameters each paramet...
void SetSize(M3::int_t nSpline)
Set the size.
M3::int_t length
void SetFunc(TF1 *&Func)
Set the function.
TF1_red()
Empty constructor.
void Print() override
Print detailed info.
M3::int_t GetNp() override
DL: Get number of points.
TF1 * ConstructTF1(const std::string &function, const int xmin, const int xmax)
KS: Make a TF1 from the reduced TF1.
TF1_red(TF1 *&Function)
The TF1 constructor with deep copy.
M3::float_t * Par
The parameters.
virtual ~TF1_red()
Empty destructor.
double Eval(const double var) override
Evaluate a variation.
TF1_red(M3::int_t nSize, M3::float_t *Array)
The useful constructor with deep copy.
double GetParameter(M3::int_t Parameter)
Get a parameter value.
void SetParameter(M3::int_t Parameter, M3::float_t Value)
Set a parameter to a value.
int GetSize()
Get the size.
KS: A reduced ResponseFunction Generic function used for evaluating weight.
virtual ~TResponseFunction_red()
Empty destructor.
TResponseFunction_red()
Empty constructor.
virtual M3::int_t GetNp()=0
DL: Get number of points.
virtual double Eval(const double var)=0
Evaluate a variation.
virtual void Print()=0
KS: Printer.
CW: Reduced TSpline3 class.
double Eval(double var) override
CW: Evaluate the weight from a variation.
void GetKnot(int i, M3::float_t &xtmp, M3::float_t &ytmp)
int FindX(double x)
Find the segment relevant to this variation in x.
virtual ~TSpline3_red()
Empty destructor.
M3::float_t ** Par
Always uses a third order polynomial, so hard-code the number of coefficients in implementation.
TSpline3_red(M3::float_t *X, M3::float_t *Y, M3::int_t N, M3::float_t **P)
constructor taking parameters
M3::int_t nPoints
Number of points/knot in TSpline3.
TSpline3_red(TSpline3 *&spline, SplineInterpolation InterPolation=kTSpline3)
The constructor that takes a TSpline3 pointer and copies in to memory.
M3::int_t GetNp() override
CW: Get the number of points.
TSpline3 * ConstructTSpline3()
CW: Make a TSpline3 from the reduced splines.
void GetCoeff(int segment, M3::float_t &x, M3::float_t &y, M3::float_t &b, M3::float_t &c, M3::float_t &d)
CW: Get the coefficient of a given segment.
void Print() override
Print detailed info.
void SetFunc(TSpline3 *&spline, SplineInterpolation InterPolation=kTSpline3)
Set the function .
M3::float_t * XPos
Positions of each x for each knot.
TSpline3_red()
Empty constructor.
M3::float_t * YResp
y-value for each knot
constexpr static const double DefSplineKnotUpBound
Default value for spline knot capping, default mean not capping is being applied.
Definition: Core.h:86
double float_t
Definition: Core.h:37
constexpr static const double DefSplineKnotLowBound
Default value for spline knot capping, default mean not capping is being applied.
Definition: Core.h:88
int int_t
Definition: Core.h:38
CW: Add a struct to hold info about the splinified xsec parameters and help with FindSplineSegment.
Definition: SplineStructs.h:21
virtual ~FastSplineInfo()=default
Destructor.
std::vector< M3::float_t > xPts
Array of the knots positions.
Definition: SplineStructs.h:37
const double * splineParsPointer
Array of the knots positions.
Definition: SplineStructs.h:43
M3::int_t CurrSegment
Array of what segment of spline we're currently interested in. Gets updated once per MCMC iteration.
Definition: SplineStructs.h:40
FastSplineInfo()
Constructor.
Definition: SplineStructs.h:24
M3::int_t nPts
Number of points in spline.
Definition: SplineStructs.h:34