Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
FlexibleInterpVar.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id: cranmer $
2// Author: Kyle Cranmer, Akira Shibata
3// Author: Giovanni Petrucciani (UCSD) (log-interpolation)
4/*************************************************************************
5 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12////////////////////////////////////////////////////////////////////////////////
13
14/** \class RooStats::HistFactory::FlexibleInterpVar
15 * \ingroup HistFactory
16 */
17
18#include <RooMsgService.h>
19#include <RooTrace.h>
20
23
24#include <Riostream.h>
25#include <TMath.h>
26
27
28using namespace RooStats;
29using namespace HistFactory;
30
31
32////////////////////////////////////////////////////////////////////////////////
33/// Default constructor
34
39
40
41////////////////////////////////////////////////////////////////////////////////
42
43FlexibleInterpVar::FlexibleInterpVar(const char* name, const char* title,
44 const RooArgList& paramList,
45 double argNominal, std::vector<double> const& lowVec, std::vector<double> const& highVec) :
46 FlexibleInterpVar{name, title, paramList, argNominal, lowVec, highVec, std::vector<int>(lowVec.size(), 0)}
47{
48}
49
50
51////////////////////////////////////////////////////////////////////////////////
52
53FlexibleInterpVar::FlexibleInterpVar(const char* name, const char* title,
54 const RooArgList& paramList,
55 double argNominal, std::vector<double> const& lowVec, std::vector<double> const& highVec,
56 std::vector<int> const& codes) :
57 RooAbsReal(name, title),
58 _paramList("paramList","List of paramficients",this),
59 _nominal(argNominal), _low(lowVec), _high(highVec)
60{
61 for (auto param : paramList) {
62 if (!dynamic_cast<RooAbsReal*>(param)) {
63 coutE(InputArguments) << "FlexibleInterpVar::ctor(" << GetName() << ") ERROR: paramficient " << param->GetName()
64 << " is not of type RooAbsReal" << std::endl ;
65 // use R__ASSERT which remains also in release mode
66 R__ASSERT(0) ;
67 }
68 _paramList.add(*param) ;
69 }
70
71 _interpCode.resize(_paramList.size());
72 for (std::size_t i = 0; i < codes.size(); ++i) {
74 }
75
76 if (_low.size() != _paramList.size() || _low.size() != _high.size() || _low.size() != _interpCode.size()) {
77 coutE(InputArguments) << "FlexibleInterpVar::ctor(" << GetName() << ") invalid input std::vectors " << std::endl;
78 R__ASSERT(_low.size() == _paramList.size());
79 R__ASSERT(_low.size() == _high.size());
80 R__ASSERT(_low.size() == _interpCode.size());
81 }
82
84}
85
86////////////////////////////////////////////////////////////////////////////////
87/// Constructor of flat polynomial function
88
89FlexibleInterpVar::FlexibleInterpVar(const char* name, const char* title) :
90 RooAbsReal(name, title),
91 _paramList("paramList","List of coefficients",this)
92{
94}
95
96////////////////////////////////////////////////////////////////////////////////
97
100 _paramList("paramList",this,other._paramList),
101 _nominal(other._nominal), _low(other._low), _high(other._high), _interpCode(other._interpCode), _interpBoundary(other._interpBoundary)
102
103{
105}
106
107
108////////////////////////////////////////////////////////////////////////////////
109/// Destructor
110
115
117{
118 int index = _paramList.index(&param);
119 if (index < 0) {
120 coutE(InputArguments) << "FlexibleInterpVar::setInterpCode ERROR: " << param.GetName() << " is not in list"
121 << std::endl;
122 return;
123 }
125}
126
128{
129 for (std::size_t i = 0; i < _interpCode.size(); ++i) {
130 setInterpCodeForParam(i, code);
131 }
132}
133
135{
136 RooAbsArg const &param = _paramList[iParam];
137 if (code < 0 || code > 5) {
138 coutE(InputArguments) << "FlexibleInterpVar::setInterpCode ERROR: " << param.GetName()
139 << " with unknown interpolation code " << code << ", keeping current code "
140 << _interpCode[iParam] << std::endl;
141 return;
142 }
143 if (code == 3) {
144 // In the past, code 3 was equivalent to code 2, which confused users.
145 // Now, we just say that code 3 doesn't exist and default to code 2 in
146 // that case for backwards compatible behavior.
147 coutE(InputArguments) << "FlexibleInterpVar::setInterpCode ERROR: " << param.GetName()
148 << " with unknown interpolation code " << code << ", defaulting to code 2" << std::endl;
149 code = 2;
150 }
151 _interpCode.at(iParam) = code;
153}
154
155////////////////////////////////////////////////////////////////////////////////
156
158 coutW(InputArguments) << "FlexibleInterpVar::setNominal : nominal is now " << newNominal << std::endl ;
160
162}
163
164////////////////////////////////////////////////////////////////////////////////
165
167 int index = _paramList.index(&param);
168 if(index<0){
169 coutE(InputArguments) << "FlexibleInterpVar::setLow ERROR: " << param.GetName()
170 << " is not in list" << std::endl ;
171 } else {
172 coutW(InputArguments) << "FlexibleInterpVar::setLow : " << param.GetName()
173 << " is now " << newLow << std::endl ;
174 _low.at(index) = newLow;
175 }
176
178}
179
180////////////////////////////////////////////////////////////////////////////////
181
183 int index = _paramList.index(&param);
184 if(index<0){
185 coutE(InputArguments) << "FlexibleInterpVar::setHigh ERROR: " << param.GetName()
186 << " is not in list" << std::endl ;
187 } else {
188 coutW(InputArguments) << "FlexibleInterpVar::setHigh : " << param.GetName()
189 << " is now " << newHigh << std::endl ;
190 _high.at(index) = newHigh;
191 }
192
194}
195
196////////////////////////////////////////////////////////////////////////////////
197
199 for(unsigned int i=0; i<_interpCode.size(); ++i){
200 coutI(InputArguments) <<"interp code for " << _paramList.at(i)->GetName() << " = " << _interpCode.at(i) << std::endl;
201 // GHL: Adding suggestion by Swagato:
202 if( _low.at(i) <= 0.001 ) coutE(InputArguments) << GetName() << ", " << _paramList.at(i)->GetName() << ": low value = " << _low.at(i) << std::endl;
203 if( _high.at(i) <= 0.001 ) coutE(InputArguments) << GetName() << ", " << _paramList.at(i)->GetName() << ": high value = " << _high.at(i) << std::endl;
204 }
205
206}
207
208
209////////////////////////////////////////////////////////////////////////////////
210/// Calculate and return value of polynomial
211
213{
214 double total(_nominal);
215 for (std::size_t i = 0; i < _paramList.size(); ++i) {
216 int code = _interpCode[i];
217 // To get consistent codes with the PiecewiseInterpolation
218 if (code == 4) {
219 code = 5;
220 }
221 double paramVal = static_cast<const RooAbsReal *>(&_paramList[i])->getVal();
223 paramVal, total);
224 }
225
226 if (total <= 0) {
228 }
229
230 return total;
231}
232
234{
235 double total(_nominal);
236
237 for (std::size_t i = 0; i < _paramList.size(); ++i) {
238 int code = _interpCode[i];
239 // To get consistent codes with the PiecewiseInterpolation
240 if (code == 4) {
241 code = 5;
242 }
244 ctx.at(&_paramList[i])[0], total);
245 }
246
247 if (total <= 0) {
249 }
250
251 ctx.output()[0] = total;
252}
253
254void FlexibleInterpVar::printMultiline(std::ostream& os, Int_t contents,
255 bool verbose, TString indent) const
256{
257 RooAbsReal::printMultiline(os,contents,verbose,indent);
258 os << indent << "--- FlexibleInterpVar ---" << std::endl;
260}
261
263{
264 for (int i=0;i<(int)_low.size();i++) {
265 auto& param = static_cast<RooAbsReal&>(_paramList[i]);
266 os << std::setw(36) << param.GetName()<<": "<<std::setw(7) << _low[i]<<" "<<std::setw(7) << _high[i]
267 <<std::endl;
268 }
269}
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define coutI(a)
#define coutW(a)
#define coutE(a)
#define TRACE_DESTROY
Definition RooTrace.h:24
#define TRACE_CREATE
Definition RooTrace.h:23
static void indent(ostringstream &buf, int indent_level)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
static unsigned int total
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
char name[80]
Definition TGX11.cxx:110
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
void setValueDirty()
Mark the element dirty. This forces a re-evaluation when a value is requested.
Definition RooAbsArg.h:431
Int_t index(const RooAbsArg *arg) const
Returns index of given arg, or -1 if arg is not in the collection.
Storage_t::size_type size() const
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Structure printing.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
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...
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
std::span< double > output()
void doEval(RooFit::EvalContext &) const override
Base function for computing multiple values of a RooAbsReal.
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Interface for detailed printing of object.
void setInterpCode(RooAbsReal &param, int code)
void setLow(RooAbsReal &param, double newLow)
void setInterpCodeForParam(int iParam, int code)
void setHigh(RooAbsReal &param, double newHigh)
double evaluate() const override
Calculate and return value of polynomial.
virtual void printFlexibleInterpVars(std::ostream &os) const
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Basic string class.
Definition TString.h:139
double flexibleInterpSingle(unsigned int code, double low, double high, double boundary, double nominal, double paramVal, double res)
Definition MathFuncs.h:232
Namespace for the RooStats classes.
Definition CodegenImpl.h:58
static T Min()
Returns maximum representation for type T.
Definition TMath.h:929