36 ClassImp(RooStats::HistFactory::FlexibleInterpVar);
38 using namespace RooStats;
39 using namespace HistFactory;
44 FlexibleInterpVar::FlexibleInterpVar()
46 _paramIter = _paramList.createIterator() ;
56 FlexibleInterpVar::FlexibleInterpVar(
const char* name,
const char* title,
57 const RooArgList& paramList,
58 Double_t nominal, vector<double> low, vector<double> high) :
59 RooAbsReal(name, title),
60 _paramList(
"paramList",
"List of paramficients",this),
61 _nominal(nominal), _low(low), _high(high), _interpBoundary(1.)
64 _paramIter = _paramList.createIterator() ;
67 TIterator* paramIter = paramList.createIterator() ;
69 while((param = (RooAbsArg*)paramIter->Next())) {
70 if (!dynamic_cast<RooAbsReal*>(param)) {
71 coutE(InputArguments) <<
"FlexibleInterpVar::ctor(" << GetName() <<
") ERROR: paramficient " << param->GetName()
72 <<
" is not of type RooAbsReal" << endl ;
75 _paramList.add(*param) ;
76 _interpCode.push_back(0);
78 if (
int(_low.size() ) != _paramList.getSize() || _low.size() != _high.size()) {
79 coutE(InputArguments) <<
"FlexibleInterpVar::ctor(" << GetName() <<
") invalid input low/high vectors " << endl;
80 R__ASSERT(
int(_low.size() ) == _paramList.getSize());
81 R__ASSERT(_low.size() == _high.size());
91 FlexibleInterpVar::FlexibleInterpVar(
const char* name,
const char* title,
92 const RooArgList& paramList,
93 double nominal,
const RooArgList& low,
const RooArgList& high) :
94 RooAbsReal(name, title),
95 _paramList(
"paramList",
"List of paramficients",this),
96 _nominal(nominal), _interpBoundary(1.)
98 RooFIter lowIter = low.fwdIterator() ;
100 while ((val = (RooAbsReal*) lowIter.next())) {
101 _low.push_back(val->getVal()) ;
104 RooFIter highIter = high.fwdIterator() ;
105 while ((val = (RooAbsReal*) highIter.next())) {
106 _high.push_back(val->getVal()) ;
111 _paramIter = _paramList.createIterator() ;
114 TIterator* paramIter = paramList.createIterator() ;
116 while((param = (RooAbsArg*)paramIter->Next())) {
117 if (!dynamic_cast<RooAbsReal*>(param)) {
118 coutE(InputArguments) <<
"FlexibleInterpVar::ctor(" << GetName() <<
") ERROR: paramficient " << param->GetName()
119 <<
" is not of type RooAbsReal" << endl ;
122 _paramList.add(*param) ;
123 _interpCode.push_back(0);
125 if (
int(_low.size() ) != _paramList.getSize() || _low.size() != _high.size()) {
126 coutE(InputArguments) <<
"FlexibleInterpVar::ctor(" << GetName() <<
") invalid input low/high lists " << endl;
127 R__ASSERT(
int(_low.size() ) == _paramList.getSize());
128 R__ASSERT(_low.size() == _high.size());
141 FlexibleInterpVar::FlexibleInterpVar(
const char* name,
const char* title,
142 const RooArgList& paramList,
143 double nominal, vector<double> low, vector<double> high,
145 RooAbsReal(name, title),
146 _paramList(
"paramList",
"List of paramficients",this),
147 _nominal(nominal), _low(low), _high(high), _interpCode(code), _interpBoundary(1.)
150 _paramIter = _paramList.createIterator() ;
153 TIterator* paramIter = paramList.createIterator() ;
155 while((param = (RooAbsArg*)paramIter->Next())) {
156 if (!dynamic_cast<RooAbsReal*>(param)) {
157 coutE(InputArguments) <<
"FlexibleInterpVar::ctor(" << GetName() <<
") ERROR: paramficient " << param->GetName()
158 <<
" is not of type RooAbsReal" << endl ;
162 _paramList.add(*param) ;
164 if (
int(_low.size() ) != _paramList.getSize() || _low.size() != _high.size() || _low.size() != _interpCode.size()) {
165 coutE(InputArguments) <<
"FlexibleInterpVar::ctor(" << GetName() <<
") invalid input vectors " << endl;
166 R__ASSERT(
int(_low.size() ) == _paramList.getSize());
167 R__ASSERT(_low.size() == _high.size());
168 R__ASSERT(_low.size() == _interpCode.size());
178 FlexibleInterpVar::FlexibleInterpVar(
const char* name,
const char* title) :
179 RooAbsReal(name, title),
180 _paramList(
"paramList",
"List of coefficients",this),
181 _nominal(0), _interpBoundary(1.)
184 _paramIter = _paramList.createIterator() ;
190 FlexibleInterpVar::FlexibleInterpVar(
const FlexibleInterpVar& other,
const char* name) :
191 RooAbsReal(other, name),
192 _paramList(
"paramList",this,other._paramList),
193 _nominal(other._nominal), _low(other._low), _high(other._high), _interpCode(other._interpCode), _interpBoundary(other._interpBoundary)
198 _paramIter = _paramList.createIterator() ;
207 FlexibleInterpVar::~FlexibleInterpVar()
216 void FlexibleInterpVar::setInterpCode(RooAbsReal& param,
int code){
217 int index = _paramList.index(¶m);
219 coutE(InputArguments) <<
"FlexibleInterpVar::setInterpCode ERROR: " << param.GetName()
220 <<
" is not in list" << endl ;
222 coutW(InputArguments) <<
"FlexibleInterpVar::setInterpCode : " << param.GetName()
223 <<
" is now " << code << endl ;
224 _interpCode.at(index) = code;
233 void FlexibleInterpVar::setAllInterpCodes(
int code){
234 for(
unsigned int i=0; i<_interpCode.size(); ++i){
235 _interpCode.at(i) = code;
245 void FlexibleInterpVar::setNominal(Double_t newNominal){
246 coutW(InputArguments) <<
"FlexibleInterpVar::setNominal : nominal is now " << newNominal << endl ;
247 _nominal = newNominal;
256 void FlexibleInterpVar::setLow(RooAbsReal& param, Double_t newLow){
257 int index = _paramList.index(¶m);
259 coutE(InputArguments) <<
"FlexibleInterpVar::setLow ERROR: " << param.GetName()
260 <<
" is not in list" << endl ;
262 coutW(InputArguments) <<
"FlexibleInterpVar::setLow : " << param.GetName()
263 <<
" is now " << newLow << endl ;
264 _low.at(index) = newLow;
274 void FlexibleInterpVar::setHigh(RooAbsReal& param, Double_t newHigh){
275 int index = _paramList.index(¶m);
277 coutE(InputArguments) <<
"FlexibleInterpVar::setHigh ERROR: " << param.GetName()
278 <<
" is not in list" << endl ;
280 coutW(InputArguments) <<
"FlexibleInterpVar::setHigh : " << param.GetName()
281 <<
" is now " << newHigh << endl ;
282 _high.at(index) = newHigh;
291 void FlexibleInterpVar::printAllInterpCodes(){
292 for(
unsigned int i=0; i<_interpCode.size(); ++i){
293 coutI(InputArguments) <<
"interp code for " << _paramList.at(i)->GetName() <<
" = " << _interpCode.at(i) <<endl;
295 if( _low.at(i) <= 0.001 ) coutE(InputArguments) << GetName() <<
", " << _paramList.at(i)->GetName() <<
": low value = " << _low.at(i) << endl;
296 if( _high.at(i) <= 0.001 ) coutE(InputArguments) << GetName() <<
", " << _paramList.at(i)->GetName() <<
": high value = " << _high.at(i) << endl;
303 double FlexibleInterpVar::PolyInterpValue(
int i,
double x)
const {
306 double boundary = _interpBoundary;
308 double x0 = boundary;
317 unsigned int n = _low.size();
318 assert(n == _high.size() );
320 _polCoeff.resize(n*6) ;
322 for (
unsigned int j = 0; j < n ; j++) {
325 double * coeff = &_polCoeff[j * 6];
328 double pow_up = std::pow(_high[j]/_nominal, x0);
329 double pow_down = std::pow(_low[j]/_nominal, x0);
330 double logHi = std::log(_high[j]) ;
331 double logLo = std::log(_low[j] );
332 double pow_up_log = _high[j] <= 0.0 ? 0.0 : pow_up * logHi;
333 double pow_down_log = _low[j] <= 0.0 ? 0.0 : -pow_down * logLo;
334 double pow_up_log2 = _high[j] <= 0.0 ? 0.0 : pow_up_log * logHi;
335 double pow_down_log2= _low[j] <= 0.0 ? 0.0 : -pow_down_log* logLo;
337 double S0 = (pow_up+pow_down)/2;
338 double A0 = (pow_up-pow_down)/2;
339 double S1 = (pow_up_log+pow_down_log)/2;
340 double A1 = (pow_up_log-pow_down_log)/2;
341 double S2 = (pow_up_log2+pow_down_log2)/2;
342 double A2 = (pow_up_log2-pow_down_log2)/2;
347 coeff[0] = 1./(8*x0) *( 15*A0 - 7*x0*S1 + x0*x0*A2);
348 coeff[1] = 1./(8*x0*x0) *(-24 + 24*S0 - 9*x0*A1 + x0*x0*S2);
349 coeff[2] = 1./(4*pow(x0, 3))*( - 5*A0 + 5*x0*S1 - x0*x0*A2);
350 coeff[3] = 1./(4*pow(x0, 4))*( 12 - 12*S0 + 7*x0*A1 - x0*x0*S2);
351 coeff[4] = 1./(8*pow(x0, 5))*( + 3*A0 - 3*x0*S1 + x0*x0*A2);
352 coeff[5] = 1./(8*pow(x0, 6))*( -8 + 8*S0 - 5*x0*A1 + x0*x0*S2);
364 assert(
int(_polCoeff.size()) > i );
365 const double * coefficients = &_polCoeff.front() + 6*i;
367 double a = coefficients[0];
368 double b = coefficients[1];
369 double c = coefficients[2];
370 double d = coefficients[3];
371 double e = coefficients[4];
372 double f = coefficients[5];
376 double value = 1. + x * (a + x * ( b + x * ( c + x * ( d + x * ( e + x * f ) ) ) ) );
383 Double_t FlexibleInterpVar::evaluate()
const
385 Double_t total(_nominal) ;
386 _paramIter->Reset() ;
396 while((param=(RooAbsReal*)_paramIter->Next())) {
400 Int_t icode = _interpCode[i] ;
406 if(param->getVal()>0)
407 total += param->getVal()*(_high[i] - _nominal );
409 total += param->getVal()*(_nominal - _low[i]);
414 if(param->getVal()>=0)
415 total *= pow(_high[i]/_nominal, +param->getVal());
417 total *= pow(_low[i]/_nominal, -param->getVal());
422 double a = 0.5*(_high[i]+_low[i])-_nominal;
423 double b = 0.5*(_high[i]-_low[i]);
425 if(param->getVal()>1 ){
426 total += (2*a+b)*(param->getVal()-1)+_high[i]-_nominal;
427 }
else if(param->getVal()<-1 ) {
428 total += -1*(2*a-b)*(param->getVal()+1)+_low[i]-_nominal;
430 total += a*pow(param->getVal(),2) + b*param->getVal()+c;
436 double a = 0.5*(_high[i]+_low[i])-_nominal;
437 double b = 0.5*(_high[i]-_low[i]);
439 if(param->getVal()>1 ){
440 total += (2*a+b)*(param->getVal()-1)+_high[i]-_nominal;
441 }
else if(param->getVal()<-1 ) {
442 total += -1*(2*a-b)*(param->getVal()+1)+_low[i]-_nominal;
444 total += a*pow(param->getVal(),2) + b*param->getVal()+c;
450 double boundary = _interpBoundary;
451 double x = param->getVal();
456 total *= std::pow(_high[i]/_nominal, +param->getVal());
458 else if (x <= -boundary)
460 total *= std::pow(_low[i]/_nominal, -param->getVal());
464 total *= PolyInterpValue(i, x);
469 coutE(InputArguments) <<
"FlexibleInterpVar::evaluate ERROR: " << param->GetName()
470 <<
" with unknown interpolation code" << endl ;
477 total= TMath::Limits<double>::Min();
483 void FlexibleInterpVar::printMultiline(ostream& os, Int_t contents,
484 Bool_t verbose, TString indent)
const
486 RooAbsReal::printMultiline(os,contents,verbose,indent);
487 os << indent <<
"--- FlexibleInterpVar ---" << endl;
488 printFlexibleInterpVars(os);
491 void FlexibleInterpVar::printFlexibleInterpVars(ostream& os)
const
494 for (
int i=0;i<(int)_low.size();i++) {
495 RooAbsReal* param=(RooAbsReal*)_paramIter->Next();
496 os << setw(36) << param->GetName()<<
": "<<setw(7) << _low[i]<<
" "<<setw(7) << _high[i]