Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooSpline.cxx
Go to the documentation of this file.
1/*
2 * Project: RooFit
3 * Author:
4 * Ruggero Turra <ruggero.turra@cern.ch>, 2016
5 *
6 * Copyright (c) 2023, CERN
7 *
8 * Redistribution and use in source and binary forms,
9 * with or without modification, are permitted according to the terms
10 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
11 */
12
13#include <RooSpline.h>
14
15#include <RooMsgService.h>
16
17#include <TGraph.h>
18
19#include <sstream>
20#include <string>
21#include <vector>
22
23/** \class RooSpline
24 \ingroup Roofit
25 \brief A RooFit class for creating spline functions.
26
27This class provides the functionality to create spline functions in RooFit,
28, using ROOT TSpline. It supports splines of order 3 or 5. It also support
29interpolation in the log-space (x or y), for example
30exp(spline({x0}, {log y0})), useful when you have something (as xsections)
31that is more similar to exponentials than polynomials.
32
33Usage example:
34~~~ {.cpp}
35RooRealVar x{"x", "x", 0, 5};
36
37std::vector<double> x0{1., 2, 3};
38std::vector<double> y0{10., 20, 50};
39
40RooSpline spline{"myspline", "my spline", x, x0, y0};
41
42auto frame = x.frame();
43spline.plotOn(frame);
44frame->Draw();
45~~~
46*/
47
49
50/// Constructor for creating a spline from a TGraph.
51/// \param[in] name The name of the spline.
52/// \param[in] title The title of the spline.
53/// \param[in] x The independent variable.
54/// \param[in] gr The input TGraph containing the data points.
55/// \param[in] order The order of the spline (3 or 5).
56/// \param[in] logx If true, the x values are logarithmically scaled before spline creation.
57/// \param[in] logy If true, the y values are logarithmically scaled before spline creation.
58RooSpline::RooSpline(const char *name, const char *title, RooAbsReal &x, const TGraph &gr, int order, bool logy,
59 bool logx)
60 : RooSpline(name, title, x, {gr.GetX(), gr.GetX() + gr.GetN()}, {gr.GetY(), gr.GetY() + gr.GetN()}, order, logx,
61 logy)
62{
63}
64
65/// Constructor for creating a spline from raw data.
66/// \param[in] name The name of the spline.
67/// \param[in] title The title of the spline.
68/// \param[in] x The independent variable.
69/// \param[in] x0 The array of x values for the spline points.
70/// \param[in] y0 The array of y values for the spline points.
71/// \param[in] order The order of the spline (3 or 5).
72/// \param[in] logx If true, the x values are logarithmically scaled before spline creation.
73/// \param[in] logy If true, the y values are logarithmically scaled before spline creation.
74RooSpline::RooSpline(const char *name, const char *title, RooAbsReal &x, std::span<const double> x0,
75 std::span<const double> y0, int order, bool logx, bool logy)
76 : RooAbsReal{name, title}, _x{"x", "x", this, x}, _logx{logx}, _logy{logy}
77{
78 const std::string title_spline = std::string(title) + "_spline";
79 if (x0.size() != y0.size()) {
80 std::stringstream errMsg;
81 errMsg << "RooSpline::ctor(" << GetName() << ") ERROR: size of x and y are not equal";
82 coutE(InputArguments) << errMsg.str() << std::endl;
83 throw std::invalid_argument(errMsg.str());
84 }
85
86 // To do the logarithm inplace if necessary.
87 std::vector<double> x0Copy;
88 x0Copy.assign(x0.begin(), x0.end());
89 std::vector<double> y0Copy;
90 y0Copy.assign(y0.begin(), y0.end());
91
92 if (_logx) {
93 for (auto &xRef : x0Copy) {
94 xRef = std::log(xRef);
95 }
96 }
97 if (_logy) {
98 for (auto &yRef : y0Copy) {
99 yRef = std::log(yRef);
100 }
101 }
102
103 if (order == 3) {
104 _spline = std::make_unique<TSpline3>(title_spline.c_str(), &x0Copy[0], &y0Copy[0], x0.size());
105 } else if (order == 5) {
106 _spline = std::make_unique<TSpline5>(title_spline.c_str(), &x0Copy[0], &y0Copy[0], x0.size());
107 } else {
108 std::stringstream errMsg;
109 errMsg << "supported orders are 3 or 5";
110 coutE(InputArguments) << errMsg.str() << std::endl;
111 throw std::invalid_argument(errMsg.str());
112 }
113}
114
115/// Copy constructor.
116/// \param[in] other The RooSpline object to copy from.
117/// \param[in] name The name of the new RooSpline object (optional).
118RooSpline::RooSpline(const RooSpline &other, const char *name)
119 : RooAbsReal(other, name),
120 _spline(static_cast<TSpline *>(other._spline->Clone())),
121 _x("x", this, other._x),
122 _logx(other._logx),
123 _logy(other._logy)
124{
125}
126
127/// Evaluate the spline function at the current point.
129{
130 const double x_val = (!_logx) ? _x : std::exp(_x);
131 return (!_logy) ? _spline->Eval(x_val) : std::exp(_spline->Eval(x_val));
132}
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:377
char name[80]
Definition TGX11.cxx:110
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
A RooFit class for creating spline functions.
Definition RooSpline.h:27
bool _logx
Flag indicating logarithmic scaling of x values.
Definition RooSpline.h:46
double evaluate() const override
Evaluate the spline function at the current point.
RooSpline()=default
RooRealProxy _x
The independent variable.
Definition RooSpline.h:45
std::unique_ptr< TSpline > _spline
The spline object.
Definition RooSpline.h:44
bool _logy
Flag indicating logarithmic scaling of y values.
Definition RooSpline.h:47
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
Double_t * GetY() const
Definition TGraph.h:138
Int_t GetN() const
Definition TGraph.h:130
Double_t * GetX() const
Definition TGraph.h:137
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Base class for spline implementation containing the Draw/Paint methods.
Definition TSpline.h:31
Double_t x[n]
Definition legend1.C:17
TGraphErrors * gr
Definition legend1.C:25