Logo ROOT  
Reference Guide
LinInterpVar.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id: cranmer $
2// Author: Kyle Cranmer, Akira Shibata
3/*************************************************************************
4 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11////////////////////////////////////////////////////////////////////////////////
12
13/** \class RooStats::HistFactory::LinInterpVar
14 * \ingroup HistFactory
15 * RooAbsReal that does piecewise-linear interpolations.
16 */
17
18#include "RooFit.h"
19
20#include "Riostream.h"
21#include <math.h>
22
23#include "RooAbsReal.h"
24#include "RooRealVar.h"
25#include "RooArgList.h"
26#include "RooMsgService.h"
27#include "TMath.h"
28
30
31using namespace std;
32
34
35using namespace RooStats;
36using namespace HistFactory;
37
38////////////////////////////////////////////////////////////////////////////////
39/// Default constructor
40
41LinInterpVar::LinInterpVar()
42{
44 _nominal = 0 ;
45}
46
47
48////////////////////////////////////////////////////////////////////////////////
49
50LinInterpVar::LinInterpVar(const char* name, const char* title,
51 const RooArgList& paramList,
52 double nominal, vector<double> low, vector<double> high) :
53 RooAbsReal(name, title),
54 _paramList("paramList","List of paramficients",this),
55 _nominal(nominal), _low(low), _high(high)
56{
58
59
60 TIterator* paramIter = paramList.createIterator() ;
61 RooAbsArg* param ;
62 while((param = (RooAbsArg*)paramIter->Next())) {
63 if (!dynamic_cast<RooAbsReal*>(param)) {
64 coutE(InputArguments) << "LinInterpVar::ctor(" << GetName() << ") ERROR: paramficient " << param->GetName()
65 << " is not of type RooAbsReal" << endl ;
66 assert(0) ;
67 }
68 _paramList.add(*param) ;
69 }
70 delete paramIter ;
71
72}
73
74////////////////////////////////////////////////////////////////////////////////
75/// Constructor of flat polynomial function
76
77LinInterpVar::LinInterpVar(const char* name, const char* title) :
78 RooAbsReal(name, title),
79 _paramList("paramList","List of coefficients",this),
80 _nominal(0)
81{
83}
84
85////////////////////////////////////////////////////////////////////////////////
86
87LinInterpVar::LinInterpVar(const LinInterpVar& other, const char* name) :
88 RooAbsReal(other, name),
89 _paramList("paramList",this,other._paramList),
90 _nominal(other._nominal), _low(other._low), _high(other._high)
91
92{
93 // Copy constructor
95
96}
97
98
99////////////////////////////////////////////////////////////////////////////////
100/// Destructor
101
103{
104 delete _paramIter ;
105}
106
107
108
109
110////////////////////////////////////////////////////////////////////////////////
111/// Calculate and return value of polynomial
112
114{
116 _paramIter->Reset() ;
117
118 RooAbsReal* param ;
119 //const RooArgSet* nset = _paramList.nset() ;
120 int i=0;
121
122 while((param=(RooAbsReal*)_paramIter->Next())) {
123 // param->Print("v");
124
125 if(param->getVal()>0)
126 sum += param->getVal()*(_high.at(i) - _nominal );
127 else
128 sum += param->getVal()*(_nominal - _low.at(i));
129
130 ++i;
131 }
132
133 if(sum<=0) {
134 sum=1E-9;
135 }
136
137 return sum;
138}
139
140
141
#define coutE(a)
Definition: RooMsgService.h:33
#define ClassImp(name)
Definition: Rtypes.h:361
char name[80]
Definition: TGX11.cxx:109
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:73
TIterator * createIterator(Bool_t dir=kIterForward) const
TIterator-style iteration over contained elements.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:60
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:90
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
RooAbsReal that does piecewise-linear interpolations.
Definition: LinInterpVar.h:24
Double_t evaluate() const
do not persist
LinInterpVar()
Default constructor.
Iterator abstract base class.
Definition: TIterator.h:30
virtual void Reset()=0
virtual TObject * Next()=0
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
@ HistFactory
Definition: RooGlobalFunc.h:69
@ InputArguments
Definition: RooGlobalFunc.h:68
Namespace for the RooStats classes.
Definition: Asimov.h:19
constexpr Double_t E()
Base of natural log:
Definition: TMath.h:97
static long int sum(long int i)
Definition: Factory.cxx:2275