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 "RooFit.h"
19
20#include "Riostream.h"
21#include <math.h>
22#include "TMath.h"
23
24#include "RooAbsReal.h"
25#include "RooRealVar.h"
26#include "RooArgList.h"
27#include "RooMsgService.h"
28#include "RooTrace.h"
29
31
32using namespace std;
33
35
36using namespace RooStats;
37using namespace HistFactory;
38
39////////////////////////////////////////////////////////////////////////////////
40/// Default constructor
41
43{
44 _nominal = 0;
48}
49
50
51////////////////////////////////////////////////////////////////////////////////
52
53FlexibleInterpVar::FlexibleInterpVar(const char* name, const char* title,
54 const RooArgList& paramList,
55 Double_t argNominal, vector<double> lowVec, vector<double> highVec) :
56 RooAbsReal(name, title),
57 _paramList("paramList","List of paramficients",this),
58 _nominal(argNominal), _low(lowVec), _high(highVec), _interpBoundary(1.)
59{
61
62 for (auto param : paramList) {
63 if (!dynamic_cast<RooAbsReal*>(param)) {
64 coutE(InputArguments) << "FlexibleInterpVar::ctor(" << GetName() << ") ERROR: paramficient " << param->GetName()
65 << " is not of type RooAbsReal" << endl ;
66 R__ASSERT(0) ;
67 }
68 _paramList.add(*param) ;
69 _interpCode.push_back(0); // default code
70 }
71 if (int(_low.size() ) != _paramList.getSize() || _low.size() != _high.size()) {
72 coutE(InputArguments) << "FlexibleInterpVar::ctor(" << GetName() << ") invalid input low/high vectors " << endl;
73 R__ASSERT(int(_low.size() ) == _paramList.getSize());
74 R__ASSERT(_low.size() == _high.size());
75 }
76
78}
79
80////////////////////////////////////////////////////////////////////////////////
81
82FlexibleInterpVar::FlexibleInterpVar(const char* name, const char* title,
83 const RooArgList& paramList,
84 double argNominal, const RooArgList& lowList, const RooArgList& highList) :
85 RooAbsReal(name, title),
86 _paramList("paramList","List of paramficients",this),
87 _nominal(argNominal), _interpBoundary(1.)
88{
89 RooFIter lowIter = lowList.fwdIterator() ;
90 RooAbsReal* val ;
91 while ((val = (RooAbsReal*) lowIter.next())) {
92 _low.push_back(val->getVal()) ;
93 }
94
95 RooFIter highIter = highList.fwdIterator() ;
96 while ((val = (RooAbsReal*) highIter.next())) {
97 _high.push_back(val->getVal()) ;
98 }
99
100
101 _logInit = kFALSE ;
102
103
104 for (auto param : paramList) {
105 if (!dynamic_cast<RooAbsReal*>(param)) {
106 coutE(InputArguments) << "FlexibleInterpVar::ctor(" << GetName() << ") ERROR: paramficient " << param->GetName()
107 << " is not of type RooAbsReal" << endl ;
108 R__ASSERT(0) ;
109 }
110 _paramList.add(*param) ;
111 _interpCode.push_back(0); // default code
112 }
113 if (int(_low.size() ) != _paramList.getSize() || _low.size() != _high.size()) {
114 coutE(InputArguments) << "FlexibleInterpVar::ctor(" << GetName() << ") invalid input low/high lists " << endl;
115 R__ASSERT(int(_low.size() ) == _paramList.getSize());
116 R__ASSERT(_low.size() == _high.size());
117 }
118
120}
121
122
123
124
125////////////////////////////////////////////////////////////////////////////////
126
127FlexibleInterpVar::FlexibleInterpVar(const char* name, const char* title,
128 const RooArgList& paramList,
129 double argNominal, vector<double> lowVec, vector<double> highVec,
130 vector<int> code) :
131 RooAbsReal(name, title),
132 _paramList("paramList","List of paramficients",this),
133 _nominal(argNominal), _low(lowVec), _high(highVec), _interpCode(code), _interpBoundary(1.)
134{
135 _logInit = kFALSE ;
136
137 for (auto param : paramList) {
138 if (!dynamic_cast<RooAbsReal*>(param)) {
139 coutE(InputArguments) << "FlexibleInterpVar::ctor(" << GetName() << ") ERROR: paramficient " << param->GetName()
140 << " is not of type RooAbsReal" << endl ;
141 // use R__ASSERT which remains also in release mode
142 R__ASSERT(0) ;
143 }
144 _paramList.add(*param) ;
145 }
146
147 if (int(_low.size() ) != _paramList.getSize() || _low.size() != _high.size() || _low.size() != _interpCode.size()) {
148 coutE(InputArguments) << "FlexibleInterpVar::ctor(" << GetName() << ") invalid input vectors " << endl;
149 R__ASSERT(int(_low.size() ) == _paramList.getSize());
150 R__ASSERT(_low.size() == _high.size());
151 R__ASSERT(_low.size() == _interpCode.size());
152 }
153
155}
156
157////////////////////////////////////////////////////////////////////////////////
158/// Constructor of flat polynomial function
159
160FlexibleInterpVar::FlexibleInterpVar(const char* name, const char* title) :
161 RooAbsReal(name, title),
162 _paramList("paramList","List of coefficients",this),
163 _nominal(0), _interpBoundary(1.)
164{
165 _logInit = kFALSE ;
167}
168
169////////////////////////////////////////////////////////////////////////////////
170
172 RooAbsReal(other, name),
173 _paramList("paramList",this,other._paramList),
174 _nominal(other._nominal), _low(other._low), _high(other._high), _interpCode(other._interpCode), _interpBoundary(other._interpBoundary)
175
176{
177 // Copy constructor
178 _logInit = kFALSE ;
180
181}
182
183
184////////////////////////////////////////////////////////////////////////////////
185/// Destructor
186
188{
190}
191
192
193////////////////////////////////////////////////////////////////////////////////
194
196 int index = _paramList.index(&param);
197 if(index<0){
198 coutE(InputArguments) << "FlexibleInterpVar::setInterpCode ERROR: " << param.GetName()
199 << " is not in list" << endl ;
200 } else {
201 coutW(InputArguments) << "FlexibleInterpVar::setInterpCode : " << param.GetName()
202 << " is now " << code << endl ;
203 _interpCode.at(index) = code;
204 }
205 // GHL: Adding suggestion by Swagato:
206 _logInit = kFALSE ;
208}
209
210////////////////////////////////////////////////////////////////////////////////
211
213 for(unsigned int i=0; i<_interpCode.size(); ++i){
214 _interpCode.at(i) = code;
215 }
216 // GHL: Adding suggestion by Swagato:
217 _logInit = kFALSE ;
219
220}
221
222////////////////////////////////////////////////////////////////////////////////
223
225 coutW(InputArguments) << "FlexibleInterpVar::setNominal : nominal is now " << newNominal << endl ;
226 _nominal = newNominal;
227
228 _logInit = kFALSE ;
229
231}
232
233////////////////////////////////////////////////////////////////////////////////
234
236 int index = _paramList.index(&param);
237 if(index<0){
238 coutE(InputArguments) << "FlexibleInterpVar::setLow ERROR: " << param.GetName()
239 << " is not in list" << endl ;
240 } else {
241 coutW(InputArguments) << "FlexibleInterpVar::setLow : " << param.GetName()
242 << " is now " << newLow << endl ;
243 _low.at(index) = newLow;
244 }
245
246 _logInit = kFALSE ;
247
249}
250
251////////////////////////////////////////////////////////////////////////////////
252
254 int index = _paramList.index(&param);
255 if(index<0){
256 coutE(InputArguments) << "FlexibleInterpVar::setHigh ERROR: " << param.GetName()
257 << " is not in list" << endl ;
258 } else {
259 coutW(InputArguments) << "FlexibleInterpVar::setHigh : " << param.GetName()
260 << " is now " << newHigh << endl ;
261 _high.at(index) = newHigh;
262 }
263
264 _logInit = kFALSE ;
266}
267
268////////////////////////////////////////////////////////////////////////////////
269
271 for(unsigned int i=0; i<_interpCode.size(); ++i){
272 coutI(InputArguments) <<"interp code for " << _paramList.at(i)->GetName() << " = " << _interpCode.at(i) <<endl;
273 // GHL: Adding suggestion by Swagato:
274 if( _low.at(i) <= 0.001 ) coutE(InputArguments) << GetName() << ", " << _paramList.at(i)->GetName() << ": low value = " << _low.at(i) << endl;
275 if( _high.at(i) <= 0.001 ) coutE(InputArguments) << GetName() << ", " << _paramList.at(i)->GetName() << ": high value = " << _high.at(i) << endl;
276 }
277
278}
279
280////////////////////////////////////////////////////////////////////////////////
281
282double FlexibleInterpVar::PolyInterpValue(int i, double x) const {
283 // code for polynomial interpolation used when interpCode=4
284
285 double boundary = _interpBoundary;
286
287 double x0 = boundary;
288
289
290 // cache the polynomial coefficient values
291 // which do not depend on x but on the boundaries values
292 if (!_logInit) {
293
295
296 unsigned int n = _low.size();
297 assert(n == _high.size() );
298
299 _polCoeff.resize(n*6) ;
300
301 for (unsigned int j = 0; j < n ; j++) {
302
303 // location of the 6 coefficient for the j-th variable
304 double * coeff = &_polCoeff[j * 6];
305
306 // GHL: Swagato's suggestions
307 double pow_up = std::pow(_high[j]/_nominal, x0);
308 double pow_down = std::pow(_low[j]/_nominal, x0);
309 double logHi = std::log(_high[j]) ;
310 double logLo = std::log(_low[j] );
311 double pow_up_log = _high[j] <= 0.0 ? 0.0 : pow_up * logHi;
312 double pow_down_log = _low[j] <= 0.0 ? 0.0 : -pow_down * logLo;
313 double pow_up_log2 = _high[j] <= 0.0 ? 0.0 : pow_up_log * logHi;
314 double pow_down_log2= _low[j] <= 0.0 ? 0.0 : -pow_down_log* logLo;
315
316 double S0 = (pow_up+pow_down)/2;
317 double A0 = (pow_up-pow_down)/2;
318 double S1 = (pow_up_log+pow_down_log)/2;
319 double A1 = (pow_up_log-pow_down_log)/2;
320 double S2 = (pow_up_log2+pow_down_log2)/2;
321 double A2 = (pow_up_log2-pow_down_log2)/2;
322
323 //fcns+der+2nd_der are eq at bd
324
325 // cache coefficient of the polynomial
326 coeff[0] = 1./(8*x0) *( 15*A0 - 7*x0*S1 + x0*x0*A2);
327 coeff[1] = 1./(8*x0*x0) *(-24 + 24*S0 - 9*x0*A1 + x0*x0*S2);
328 coeff[2] = 1./(4*pow(x0, 3))*( - 5*A0 + 5*x0*S1 - x0*x0*A2);
329 coeff[3] = 1./(4*pow(x0, 4))*( 12 - 12*S0 + 7*x0*A1 - x0*x0*S2);
330 coeff[4] = 1./(8*pow(x0, 5))*( + 3*A0 - 3*x0*S1 + x0*x0*A2);
331 coeff[5] = 1./(8*pow(x0, 6))*( -8 + 8*S0 - 5*x0*A1 + x0*x0*S2);
332
333 }
334
335 }
336
337 // GHL: Swagato's suggestions
338 // if( _low[i] == 0 ) _low[i] = 0.0001;
339 // if( _high[i] == 0 ) _high[i] = 0.0001;
340
341 // get pointer to location of coefficients in the vector
342
343 assert(int(_polCoeff.size()) > i );
344 const double * coefficients = &_polCoeff.front() + 6*i;
345
346 double a = coefficients[0];
347 double b = coefficients[1];
348 double c = coefficients[2];
349 double d = coefficients[3];
350 double e = coefficients[4];
351 double f = coefficients[5];
352
353
354 // evaluate the 6-th degree polynomial using Horner's method
355 double value = 1. + x * (a + x * ( b + x * ( c + x * ( d + x * ( e + x * f ) ) ) ) );
356 return value;
357}
358
359////////////////////////////////////////////////////////////////////////////////
360/// Const getters
361
363double FlexibleInterpVar::nominal() const { return _nominal; }
364const std::vector<double>& FlexibleInterpVar::low() const { return _low; }
365const std::vector<double>& FlexibleInterpVar::high() const { return _high; }
366
367////////////////////////////////////////////////////////////////////////////////
368/// Calculate and return value of polynomial
369
371{
373 int i=0;
374
375 for (auto arg : _paramList) {
376 auto param = static_cast<const RooAbsReal*>(arg);
377
378 Int_t icode = _interpCode[i] ;
379
380 switch(icode) {
381
382 case 0: {
383 // piece-wise linear
384 if(param->getVal()>0)
385 total += param->getVal()*(_high[i] - _nominal );
386 else
387 total += param->getVal()*(_nominal - _low[i]);
388 break ;
389 }
390 case 1: {
391 // pice-wise log
392 if(param->getVal()>=0)
393 total *= pow(_high[i]/_nominal, +param->getVal());
394 else
395 total *= pow(_low[i]/_nominal, -param->getVal());
396 break ;
397 }
398 case 2: {
399 // parabolic with linear
400 double a = 0.5*(_high[i]+_low[i])-_nominal;
401 double b = 0.5*(_high[i]-_low[i]);
402 double c = 0;
403 if(param->getVal()>1 ){
404 total += (2*a+b)*(param->getVal()-1)+_high[i]-_nominal;
405 } else if(param->getVal()<-1 ) {
406 total += -1*(2*a-b)*(param->getVal()+1)+_low[i]-_nominal;
407 } else {
408 total += a*pow(param->getVal(),2) + b*param->getVal()+c;
409 }
410 break ;
411 }
412 case 3: {
413 //parabolic version of log-normal
414 double a = 0.5*(_high[i]+_low[i])-_nominal;
415 double b = 0.5*(_high[i]-_low[i]);
416 double c = 0;
417 if(param->getVal()>1 ){
418 total += (2*a+b)*(param->getVal()-1)+_high[i]-_nominal;
419 } else if(param->getVal()<-1 ) {
420 total += -1*(2*a-b)*(param->getVal()+1)+_low[i]-_nominal;
421 } else {
422 total += a*pow(param->getVal(),2) + b*param->getVal()+c;
423 }
424 break ;
425 }
426
427 case 4: {
428 double boundary = _interpBoundary;
429 double x = param->getVal();
430 //std::cout << icode << " param " << param->GetName() << " " << param->getVal() << " boundary " << boundary << std::endl;
431
432 if(x >= boundary)
433 {
434 total *= std::pow(_high[i]/_nominal, +param->getVal());
435 }
436 else if (x <= -boundary)
437 {
438 total *= std::pow(_low[i]/_nominal, -param->getVal());
439 }
440 else if (x != 0)
441 {
442 total *= PolyInterpValue(i, x);
443 }
444 break ;
445 }
446 default: {
447 coutE(InputArguments) << "FlexibleInterpVar::evaluate ERROR: " << param->GetName()
448 << " with unknown interpolation code" << endl ;
449 }
450 }
451 ++i;
452 }
453
454 if(total<=0) {
456 }
457
458 return total;
459}
460
461void FlexibleInterpVar::printMultiline(ostream& os, Int_t contents,
462 Bool_t verbose, TString indent) const
463{
464 RooAbsReal::printMultiline(os,contents,verbose,indent);
465 os << indent << "--- FlexibleInterpVar ---" << endl;
467}
468
470{
471 for (int i=0;i<(int)_low.size();i++) {
472 auto& param = static_cast<RooAbsReal&>(_paramList[i]);
473 os << setw(36) << param.GetName()<<": "<<setw(7) << _low[i]<<" "<<setw(7) << _high[i]
474 <<endl;
475 }
476}
477
478
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define S0(x)
Definition RSha256.hxx:88
#define S1(x)
Definition RSha256.hxx:89
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
#define e(i)
Definition RSha256.hxx:103
#define coutI(a)
#define coutW(a)
#define coutE(a)
#define TRACE_DESTROY
Definition RooTrace.h:24
#define TRACE_CREATE
Definition RooTrace.h:23
const Bool_t kFALSE
Definition RtypesCore.h:101
const Bool_t kTRUE
Definition RtypesCore.h:100
#define ClassImp(name)
Definition Rtypes.h:364
static void indent(ostringstream &buf, int indent_level)
#define R__ASSERT(e)
Definition TError.h:118
static unsigned int total
char name[80]
Definition TGX11.cxx:110
void setValueDirty()
Mark the element dirty. This forces a re-evaluation when a value is requested.
Definition RooAbsArg.h:505
Int_t getSize() const
Int_t index(const RooAbsArg *arg) const
Returns index of given arg, or -1 if arg is not in the collection.
RooFIter fwdIterator() const
One-time forward iterator.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Structure printing.
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
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
A one-time forward iterator working on RooLinkedList or RooAbsCollection.
RooAbsArg * next()
Return next element or nullptr if at end.
RooListProxy is the concrete proxy for RooArgList objects.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Reimplementation of standard RooArgList::add()
std::vector< double > _polCoeff
! cached polynomial coefficients
Double_t evaluate() const
Calculate and return value of polynomial.
void setInterpCode(RooAbsReal &param, int code)
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Structure printing.
const std::vector< double > & high() const
void setLow(RooAbsReal &param, Double_t newLow)
const std::vector< double > & low() const
Bool_t _logInit
! flag used for caching polynomial coefficients
double PolyInterpValue(int i, double x) const
virtual void printFlexibleInterpVars(std::ostream &os) const
const RooListProxy & variables() const
Const getters.
void setHigh(RooAbsReal &param, Double_t newHigh)
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Basic string class.
Definition TString.h:136
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
Namespace for the RooStats classes.
Definition Asimov.h:19
static T Min()
Returns maximum representation for type T.
Definition TMath.h:871