ROOT  6.06/09
Reference Guide
WrapperRooPdf.h
Go to the documentation of this file.
1 // wrapper class for a RooPdf
2 
3 #ifndef ROOT_WrapperRooPdf
4 #define ROOT_WrapperRooPdf
5 
6 #include <RooAbsPdf.h>
7 #include <RooRealVar.h>
8 #include <RooArgSet.h>
9 #include <RooGaussian.h>
10 #include <TF1.h>
11 
12 #include <Math/IParamFunction.h>
13 
14 #include <cassert>
15 
17 
18 public:
19 
20  /**
21  for pdf with only 1D observables using as default the name x
22  */
23  WrapperRooPdf(RooAbsPdf * pdf, const std::string xvar = "x", bool norm = true) :
24  fNorm(norm),
25  fPdf(pdf),
26  fX(0),
27  fParams(0)
28  {
29  assert(fPdf != 0);
30 
31  RooArgSet *vars = fPdf->getVariables();
32  RooAbsArg * arg = vars->find(xvar.c_str()); // code should abort if not found
33  if (!arg) std::cout <<"Error - observable " << xvar << "is not in the list of pdf variables" << std::endl;
34  assert(arg != 0);
35  RooArgSet obsList(*arg);
36  //arg.setDirtyInhibit(true); // do have faster setter of values
37  fX = fPdf->getObservables(obsList);
38  fParams = fPdf->getParameters(obsList);
39  assert(fX!=0);
40  assert(fParams!=0);
41  delete vars;
42 #ifdef DEBUG
43  fX->Print("v");
44  fParams->Print("v");
45 #endif
46  }
47 
48 
49 
50 
51  /**
52  for pdf with multi-dim observables specifying observables in the RooArgSet
53  */
54  WrapperRooPdf(RooAbsPdf * pdf, const RooArgSet & obsList, bool norm = true ) :
55  fNorm(norm),
56  fPdf(pdf),
57  fX(0),
58  fParams(0)
59  {
60  assert(fPdf != 0);
61 
62  fX = fPdf->getObservables(obsList);
63  fParams = fPdf->getParameters(obsList);
64  assert(fX!=0);
65  assert(fParams!=0);
66 #ifdef DEBUG
67  fX->Print("v");
68  fParams->Print("v");
69 #endif
70 // // iterate on fX
71 // TIterator* itr = fX->createIterator() ;
72 // RooAbsArg* arg = 0;
73 // while( ( arg = dynamic_cast<RooAbsArg*>(itr->Next() ) ) ) {
74 // assert(arg != 0);
75 // arg->setDirtyInhibit(true); // for having faster setter later in DoEval
76 // }
77 
78  }
79 
80 
82  // need to delete observables and parameter list
83  if (fX) delete fX;
84  if (fParams) delete fParams;
85  }
86 
87  /**
88  clone the function
89  */
90 #ifndef _WIN32
92 #else
94 #endif
95  * Clone() const {
96  // copy the pdf function pointer
97  return new WrapperRooPdf(fPdf, *fX, fNorm);
98  }
99 
100  unsigned int NPar() const {
101  return fParams->getSize();
102  }
103  unsigned int NDim() const {
104  return fX->getSize();
105  }
106  const double * Parameters() const {
107  if (fParamValues.size() != NPar() )
108  fParamValues.resize(NPar() );
109 
110  // iterate on parameters and set values
111  TIterator* itr = fParams->createIterator() ;
112  std::vector<double>::iterator vpitr = fParamValues.begin();
113 
114  RooRealVar* var = 0;
115  while( ( var = dynamic_cast<RooRealVar*>(itr->Next() ) ) ) {
116  assert(var != 0);
117  *vpitr++ = var->getVal();
118  }
119  return &fParamValues.front();
120  }
121 
122  std::string ParameterName(unsigned int i) const {
123  // iterate on parameters and set values
124  TIterator* itr = fParams->createIterator() ;
125  RooRealVar* var = 0;
126  unsigned int index = 0;
127  while( ( var = dynamic_cast<RooRealVar*>(itr->Next() ) ) ) {
128  assert(var != 0);
129  if (index == i) return std::string(var->GetName() );
130  index++;
131  }
132  return "not_found";
133  }
134 
135 
136  /**
137  set parameters. Order of parameter is the one defined by the RooPdf and must be checked by user
138  */
139 
140  void SetParameters(const double * p) {
141  DoSetParameters(p);
142  }
143 
144 // double operator() (double *x, double * p = 0) {
145 // if (p != 0) SetParameters(p);
146 // // iterate on observables
147 // TIterator* itr = fX->createIterator() ;
148 // RooRealVar* var = 0;
149 // while( ( var = dynamic_cast<RooRealVar*>(itr->Next() ) ) ) {
150 // assert(var != 0);
151 // var->setVal(*x++);
152 // }
153 // // debug
154 // //fX->Print("v");
155 
156 // if (fNorm)
157 // return fPdf->getVal(fX);
158 // else
159 // return fPdf->getVal(); // get unnormalized value
160 // }
161 
162 
163 private:
164 
165  double DoEvalPar(const double * x, const double * p) const {
166 
167  // should maybe be optimized ???
168  DoSetParameters(p);
169 
170  // iterate on observables
171  TIterator* itr = fX->createIterator() ;
172  RooRealVar* var = 0;
173  while( ( var = dynamic_cast<RooRealVar*>(itr->Next() ) ) ) {
174  assert(var != 0);
175 #ifndef _WIN32
176  var->setDirtyInhibit(true);
177 #endif
178  var->setVal(*x++);
179  }
180  // debug
181  //fX->Print("v");
182 
183  if (fNorm)
184  return fPdf->getVal(fX);
185  else
186  return fPdf->getVal(); // get unnormalized value
187 
188  }
189 
190 
191  void DoSetParameters(const double * p) const {
192  // iterate on parameters and set values
193  TIterator* itr = fParams->createIterator() ;
194  RooRealVar* var = 0;
195  while( ( var = dynamic_cast<RooRealVar*>(itr->Next() ) ) ) {
196  assert(var != 0);
197  var->setVal(*p++);
198  }
199  // debug
200  //fParams->Print("v");
201  }
202 
203 
204  bool fNorm;
205  mutable RooAbsPdf * fPdf;
206  mutable RooArgSet * fX;
207  mutable RooArgSet * fParams;
208  mutable std::vector<double> fParamValues;
209 
210 
211 };
212 
213 
214 
215 #endif
double DoEvalPar(const double *x, const double *p) const
Implementation of the evaluation function using the x values and the parameters.
const double * Parameters() const
Access the parameter values.
unsigned int NPar() const
Return the number of Parameters.
RooAbsPdf * fPdf
#define assert(cond)
Definition: unittest.h:542
void DoSetParameters(const double *p) const
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Definition: RooAbsArg.h:194
std::string ParameterName(unsigned int i) const
Return the name of the i-th parameter (starting from zero) Overwrite if want to avoid the default nam...
Iterator abstract base class.
Definition: TIterator.h:32
Double_t x[n]
Definition: legend1.C:17
WrapperRooPdf(RooAbsPdf *pdf, const std::string xvar="x", bool norm=true)
for pdf with only 1D observables using as default the name x
Definition: WrapperRooPdf.h:23
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:2081
void SetParameters(const double *p)
set parameters.
std::vector< double > fParamValues
TIterator * createIterator(Bool_t dir=kIterForward) const
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
RooArgSet * getParameters(const RooAbsData *data, Bool_t stripDisconnected=kTRUE) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
Definition: RooAbsArg.cxx:559
RooArgSet * fParams
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Definition: RooRealVar.cxx:202
unsigned int NDim() const
Retrieve the dimension of the function.
RooAbsArg * find(const char *name) const
Find object with given name in list.
IParamFunction interface (abstract class) describing multi-dimensional parameteric functions It is a ...
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
WrapperRooPdf(RooAbsPdf *pdf, const RooArgSet &obsList, bool norm=true)
for pdf with multi-dim observables specifying observables in the RooArgSet
Definition: WrapperRooPdf.h:54
static void setDirtyInhibit(Bool_t flag)
Control global dirty inhibit mode.
Definition: RooAbsArg.cxx:238
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
virtual TObject * Next()=0
Int_t getSize() const
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
WrapperRooPdf * Clone() const
clone the function
Definition: WrapperRooPdf.h:95
Documentation for the abstract class IBaseFunctionMultiDim.
Definition: IFunction.h:63
RooArgSet * fX
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40