Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooFormulaVar.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17//////////////////////////////////////////////////////////////////////////////
18/// \class RooFormulaVar
19///
20/// A RooFormulaVar is a generic implementation of a real-valued object,
21/// which takes a RooArgList of servers and a C++ expression string defining how
22/// its value should be calculated from the given list of servers.
23/// RooFormulaVar uses a RooFormula object to perform the expression evaluation.
24///
25/// If RooAbsPdf objects are supplied to RooFormulaVar as servers, their
26/// raw (unnormalized) values will be evaluated. Use RooGenericPdf, which
27/// constructs generic PDF functions, to access their properly normalized
28/// values.
29///
30/// The string expression can be any valid TFormula expression referring to the
31/// listed servers either by name or by their ordinal list position. These three are
32/// equivalent:
33/// ```
34/// RooFormulaVar("gen", "x*y", RooArgList(x,y)) // reference by name
35/// RooFormulaVar("gen", "@0*@1", RooArgList(x,y)) // reference by ordinal with @
36/// RooFormulaVar("gen", "x[0]*x[1]", RooArgList(x,y)) // TFormula-builtin reference by ordinal
37/// ```
38/// Note that `x[i]` is an expression reserved for TFormula. All variable references
39/// are automatically converted to the TFormula-native format. If a variable with
40/// the name `x` is given, the RooFormula interprets `x[i]` as a list position,
41/// but `x` without brackets as the name of a RooFit object.
42///
43/// The last two versions, while slightly less readable, are more versatile because
44/// the names of the arguments are not hard coded.
45///
46
47
48#include "Riostream.h"
49
50#include "RooFormulaVar.h"
51#include "RooStreamParser.h"
52#include "RooNLLVar.h"
53#include "RooChi2Var.h"
54#include "RooMsgService.h"
55#include "RooTrace.h"
56#include "RooFormula.h"
57
58
59using namespace std;
60
62
64
66{
67 if(_formula) delete _formula;
68}
69
70////////////////////////////////////////////////////////////////////////////////
71/// Constructor with formula expression and list of input variables.
72/// \param[in] name Name of the formula.
73/// \param[in] title Title of the formula.
74/// \param[in] inFormula Expression to be evaluated.
75/// \param[in] dependents Variables that should be passed to the formula.
76/// \param[in] checkVariables Check that all variables from `dependents` are used in the expression.
77RooFormulaVar::RooFormulaVar(const char *name, const char *title, const char* inFormula, const RooArgList& dependents,
78 bool checkVariables) :
79 RooAbsReal(name,title),
80 _actualVars("actualVars","Variables used by formula expression",this),
81 _formExpr(inFormula)
82{
83 if (dependents.empty()) {
84 _value = traceEval(nullptr);
85 } else {
86 _formula = new RooFormula(GetName(), _formExpr, dependents, checkVariables);
87 _formExpr = _formula->formulaString().c_str();
89 }
90}
91
92
93
94////////////////////////////////////////////////////////////////////////////////
95/// Constructor with formula expression, title and list of input variables.
96/// \param[in] name Name of the formula.
97/// \param[in] title Formula expression. Will also be used as the title.
98/// \param[in] dependents Variables that should be passed to the formula.
99/// \param[in] checkVariables Check that all variables from `dependents` are used in the expression.
100RooFormulaVar::RooFormulaVar(const char *name, const char *title, const RooArgList& dependents,
101 bool checkVariables) :
102 RooAbsReal(name,title),
103 _actualVars("actualVars","Variables used by formula expression",this),
104 _formExpr(title)
105{
106 if (dependents.empty()) {
107 _value = traceEval(nullptr);
108 } else {
109 _formula = new RooFormula(GetName(), _formExpr, dependents, checkVariables);
110 _formExpr = _formula->formulaString().c_str();
112 }
113}
114
115
116
117////////////////////////////////////////////////////////////////////////////////
118/// Copy constructor
119
121 RooAbsReal(other, name),
122 _actualVars("actualVars",this,other._actualVars),
123 _formExpr(other._formExpr)
124{
125 if (other._formula && other._formula->ok()) {
126 _formula = new RooFormula(*other._formula);
127 _formExpr = _formula->formulaString().c_str();
128 }
129}
130
131
132////////////////////////////////////////////////////////////////////////////////
133/// Return reference to internal RooFormula object.
134/// If it doesn't exist, create it on the fly.
136{
137 if (!_formula) {
138 // After being read from file, the formula object might not exist, yet:
140 const_cast<TString&>(_formExpr) = _formula->formulaString().c_str();
141 }
142
143 return *_formula;
144}
145
146
147bool RooFormulaVar::ok() const { return getFormula().ok() ; }
148
149
151
152
153////////////////////////////////////////////////////////////////////////////////
154/// Calculate current value of object from internal formula
155
157{
158 return getFormula().eval(_actualVars.nset());
159}
160
161
162void RooFormulaVar::computeBatch(double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
163{
164 getFormula().computeBatch(output, nEvents, dataMap);
165}
166
167
168////////////////////////////////////////////////////////////////////////////////
169/// Propagate server change information to embedded RooFormula object
170
171bool RooFormulaVar::redirectServersHook(const RooAbsCollection& newServerList, bool mustReplaceAll, bool nameChange, bool isRecursive)
172{
173 bool error = getFormula().changeDependents(newServerList,mustReplaceAll,nameChange);
174
176 return error || RooAbsReal::redirectServersHook(newServerList, mustReplaceAll, nameChange, isRecursive);
177}
178
179
180
181////////////////////////////////////////////////////////////////////////////////
182/// Print info about this object to the specified stream.
183
184void RooFormulaVar::printMultiline(ostream& os, Int_t contents, bool verbose, TString indent) const
185{
186 RooAbsReal::printMultiline(os,contents,verbose,indent);
187 if(verbose) {
188 indent.Append(" ");
189 os << indent;
190 getFormula().printMultiline(os,contents,verbose,indent);
191 }
192}
193
194
195
196////////////////////////////////////////////////////////////////////////////////
197/// Add formula expression as meta argument in printing interface
198
199void RooFormulaVar::printMetaArgs(ostream& os) const
200{
201 os << "formula=\"" << _formExpr << "\" " ;
202}
203
204
205
206
207////////////////////////////////////////////////////////////////////////////////
208/// Read object contents from given stream
209
210bool RooFormulaVar::readFromStream(istream& /*is*/, bool /*compact*/, bool /*verbose*/)
211{
212 coutE(InputArguments) << "RooFormulaVar::readFromStream(" << GetName() << "): can't read" << endl ;
213 return true ;
214}
215
216
217
218////////////////////////////////////////////////////////////////////////////////
219/// Write object contents to given stream
220
221void RooFormulaVar::writeToStream(ostream& os, bool compact) const
222{
223 if (compact) {
224 cout << getVal() << endl ;
225 } else {
226 os << GetTitle() ;
227 }
228}
229
230
231
232////////////////////////////////////////////////////////////////////////////////
233/// Forward the plot sampling hint from the p.d.f. that defines the observable obs
234
235std::list<double>* RooFormulaVar::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
236{
237 for (const auto par : _actualVars) {
238 auto func = static_cast<const RooAbsReal*>(par);
239 list<double>* binb = nullptr;
240
241 if (func && (binb = func->binBoundaries(obs,xlo,xhi)) ) {
242 return binb;
243 }
244 }
245
246 return nullptr;
247}
248
249
250
251////////////////////////////////////////////////////////////////////////////////
252/// Forward the plot sampling hint from the p.d.f. that defines the observable obs
253
254std::list<double>* RooFormulaVar::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
255{
256 for (const auto par : _actualVars) {
257 auto func = dynamic_cast<const RooAbsReal*>(par);
258 list<double>* hint = nullptr;
259
260 if (func && (hint = func->plotSamplingHint(obs,xlo,xhi)) ) {
261 return hint;
262 }
263 }
264
265 return nullptr;
266}
267
268
269
270////////////////////////////////////////////////////////////////////////////////
271/// Return the default error level for MINUIT error analysis
272/// If the formula contains one or more RooNLLVars and
273/// no RooChi2Vars, return the defaultErrorLevel() of
274/// RooNLLVar. If the addition contains one ore more RooChi2Vars
275/// and no RooNLLVars, return the defaultErrorLevel() of
276/// RooChi2Var. If the addition contains neither or both
277/// issue a warning message and return a value of 1
278
280{
281 RooAbsReal* nllArg(nullptr) ;
282 RooAbsReal* chi2Arg(nullptr) ;
283
284 for (const auto arg : _actualVars) {
285 if (dynamic_cast<RooNLLVar*>(arg)) {
286 nllArg = (RooAbsReal*)arg ;
287 }
288 if (dynamic_cast<RooChi2Var*>(arg)) {
289 chi2Arg = (RooAbsReal*)arg ;
290 }
291 }
292
293 if (nllArg && !chi2Arg) {
294 coutI(Minimization) << "RooFormulaVar::defaultErrorLevel(" << GetName()
295 << ") Formula contains a RooNLLVar, using its error level" << endl ;
296 return nllArg->defaultErrorLevel() ;
297 } else if (chi2Arg && !nllArg) {
298 coutI(Minimization) << "RooFormulaVar::defaultErrorLevel(" << GetName()
299 << ") Formula contains a RooChi2Var, using its error level" << endl ;
300 return chi2Arg->defaultErrorLevel() ;
301 } else if (!nllArg && !chi2Arg) {
302 coutI(Minimization) << "RooFormulaVar::defaultErrorLevel(" << GetName() << ") WARNING: "
303 << "Formula contains neither RooNLLVar nor RooChi2Var server, using default level of 1.0" << endl ;
304 } else {
305 coutI(Minimization) << "RooFormulaVar::defaultErrorLevel(" << GetName() << ") WARNING: "
306 << "Formula contains BOTH RooNLLVar and RooChi2Var server, using default level of 1.0" << endl ;
307 }
308
309 return 1.0 ;
310}
311
312
313
314
#define coutI(a)
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:377
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition TGX11.cxx:110
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
const RooArgSet * nset() const
Definition RooAbsProxy.h:52
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
virtual double defaultErrorLevel() const
Definition RooAbsReal.h:248
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Structure printing.
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursiveStep) override
Function that is called at the end of redirectServers().
double _value
Cache for current value of object.
Definition RooAbsReal.h:543
double traceEval(const RooArgSet *set) const
Calculate current value of object, with error tracing wrapper.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooChi2Var implements a simple calculation from a binned dataset and a PDF.
Definition RooChi2Var.h:25
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
~RooFormulaVar() override
RooListProxy _actualVars
Actual parameters used by formula engine.
RooFormula & getFormula() const
Return reference to internal RooFormula object.
RooFormula * _formula
! Formula engine
void dumpFormula()
Dump the formula to stdout.
double defaultErrorLevel() const override
Return the default error level for MINUIT error analysis If the formula contains one or more RooNLLVa...
bool redirectServersHook(const RooAbsCollection &newServerList, bool mustReplaceAll, bool nameChange, bool isRecursive) override
Propagate server change information to embedded RooFormula object.
bool ok() const
const RooArgList & dependents() const
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
Forward the plot sampling hint from the p.d.f. that defines the observable obs.
bool readFromStream(std::istream &is, bool compact, bool verbose=false) override
Read object contents from given stream.
void computeBatch(double *output, size_t nEvents, RooFit::Detail::DataMap const &dataMap) const override
Base function for computing multiple values of a RooAbsReal.
TString _formExpr
Formula expression string.
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Print info about this object to the specified stream.
double evaluate() const override
Calculate current value of object from internal formula.
void writeToStream(std::ostream &os, bool compact) const override
Write object contents to given stream.
void printMetaArgs(std::ostream &os) const override
Add formula expression as meta argument in printing interface.
std::list< double > * plotSamplingHint(RooAbsRealLValue &, double, double) const override
Forward the plot sampling hint from the p.d.f. that defines the observable obs.
RooFormula internally uses ROOT's TFormula to compute user-defined expressions of RooAbsArgs.
Definition RooFormula.h:27
std::string formulaString() const
Definition RooFormula.h:71
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Printing interface.
RooArgSet actualDependents() const
Return list of arguments which are used in the formula.
Definition RooFormula.h:39
bool ok() const
Definition RooFormula.h:50
void computeBatch(double *output, size_t nEvents, RooFit::Detail::DataMap const &) const
void dump() const
DEBUG: Dump state information.
double eval(const RooArgSet *nset=nullptr) const
Evaluate all parameters/observables, and then evaluate formula.
bool changeDependents(const RooAbsCollection &newDeps, bool mustReplaceAll, bool nameChange)
Change used variables to those with the same name in given list.
Class RooNLLVar implements a -log(likelihood) calculation from a dataset and a PDF.
Definition RooNLLVar.h:25
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
Basic string class.
Definition TString.h:139
static void output()