Logo ROOT   6.18/05
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 "Riostream.h"
22#include <math.h>
23#include "TMath.h"
24
25#include "RooAbsReal.h"
26#include "RooRealVar.h"
27#include "RooArgList.h"
28#include "RooMsgService.h"
29#include "TMath.h"
30
32
33using namespace std;
34
36
37using namespace RooStats;
38using namespace HistFactory;
39
40////////////////////////////////////////////////////////////////////////////////
41/// Default constructor
42
43LinInterpVar::LinInterpVar()
44{
46 _nominal = 0 ;
47}
48
49
50////////////////////////////////////////////////////////////////////////////////
51
52LinInterpVar::LinInterpVar(const char* name, const char* title,
53 const RooArgList& paramList,
54 double nominal, vector<double> low, vector<double> high) :
55 RooAbsReal(name, title),
56 _paramList("paramList","List of paramficients",this),
57 _nominal(nominal), _low(low), _high(high)
58{
60
61
62 TIterator* paramIter = paramList.createIterator() ;
63 RooAbsArg* param ;
64 while((param = (RooAbsArg*)paramIter->Next())) {
65 if (!dynamic_cast<RooAbsReal*>(param)) {
66 coutE(InputArguments) << "LinInterpVar::ctor(" << GetName() << ") ERROR: paramficient " << param->GetName()
67 << " is not of type RooAbsReal" << endl ;
68 assert(0) ;
69 }
70 _paramList.add(*param) ;
71 }
72 delete paramIter ;
73
74}
75
76////////////////////////////////////////////////////////////////////////////////
77/// Constructor of flat polynomial function
78
79LinInterpVar::LinInterpVar(const char* name, const char* title) :
80 RooAbsReal(name, title),
81 _paramList("paramList","List of coefficients",this),
82 _nominal(0)
83{
85}
86
87////////////////////////////////////////////////////////////////////////////////
88
89LinInterpVar::LinInterpVar(const LinInterpVar& other, const char* name) :
90 RooAbsReal(other, name),
91 _paramList("paramList",this,other._paramList),
92 _nominal(other._nominal), _low(other._low), _high(other._high)
93
94{
95 // Copy constructor
97
98}
99
100
101////////////////////////////////////////////////////////////////////////////////
102/// Destructor
103
105{
106 delete _paramIter ;
107}
108
109
110
111
112////////////////////////////////////////////////////////////////////////////////
113/// Calculate and return value of polynomial
114
116{
118 _paramIter->Reset() ;
119
120 RooAbsReal* param ;
121 //const RooArgSet* nset = _paramList.nset() ;
122 int i=0;
123
124 while((param=(RooAbsReal*)_paramIter->Next())) {
125 // param->Print("v");
126
127 if(param->getVal()>0)
128 sum += param->getVal()*(_high.at(i) - _nominal );
129 else
130 sum += param->getVal()*(_nominal - _low.at(i));
131
132 ++i;
133 }
134
135 if(sum<=0) {
136 sum=1E-9;
137 }
138
139 return sum;
140}
141
142
143
#define coutE(a)
Definition: RooMsgService.h:34
double Double_t
Definition: RtypesCore.h:55
#define ClassImp(name)
Definition: Rtypes.h:365
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:70
TIterator * createIterator(Bool_t dir=kIterForward) const R__SUGGEST_ALTERNATIVE("begin()
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:53
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:81
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
@ InputArguments
Definition: RooGlobalFunc.h:58
Namespace for the RooStats classes.
Definition: Asimov.h:20
constexpr Double_t E()
Base of natural log:
Definition: TMath.h:97
static long int sum(long int i)
Definition: Factory.cxx:2258