00001
00002
00003 #ifndef ROOT_WrapperRooPdf
00004 #define ROOT_WrapperRooPdf
00005
00006 #include <RooAbsPdf.h>
00007 #include <RooRealVar.h>
00008 #include <RooArgSet.h>
00009 #include <RooGaussian.h>
00010 #include <TF1.h>
00011
00012 #include <Math/IParamFunction.h>
00013
00014 #include <cassert>
00015
00016 class WrapperRooPdf : public ROOT::Math::IParamMultiFunction {
00017
00018 public:
00019
00020
00021
00022
00023 WrapperRooPdf(RooAbsPdf * pdf, const std::string xvar = "x", bool norm = true) :
00024 fNorm(norm),
00025 fPdf(pdf),
00026 fX(0),
00027 fParams(0)
00028 {
00029 assert(fPdf != 0);
00030
00031 RooArgSet *vars = fPdf->getVariables();
00032 RooAbsArg & arg = (*vars)[xvar.c_str()];
00033 RooArgSet obsList(arg);
00034
00035 fX = fPdf->getObservables(obsList);
00036 fParams = fPdf->getParameters(obsList);
00037 assert(fX!=0);
00038 assert(fParams!=0);
00039 fX->Print("v");
00040 fParams->Print("v");
00041 }
00042
00043
00044
00045
00046
00047
00048
00049 WrapperRooPdf(RooAbsPdf * pdf, const RooArgSet & obsList, bool norm = true ) :
00050 fNorm(norm),
00051 fPdf(pdf),
00052 fX(0),
00053 fParams(0)
00054 {
00055 assert(fPdf != 0);
00056
00057 fX = fPdf->getObservables(obsList);
00058 fParams = fPdf->getParameters(obsList);
00059 assert(fX!=0);
00060 assert(fParams!=0);
00061 fX->Print("v");
00062 fParams->Print("v");
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072 }
00073
00074
00075 ~WrapperRooPdf() {
00076
00077 if (fX) delete fX;
00078 if (fParams) delete fParams;
00079 }
00080
00081
00082
00083
00084 #ifndef _WIN32
00085 WrapperRooPdf
00086 #else
00087 ROOT::Math::IMultiGenFunction
00088 #endif
00089 * Clone() const {
00090
00091 return new WrapperRooPdf(fPdf, *fX, fNorm);
00092 }
00093
00094 unsigned int NPar() const {
00095 return fParams->getSize();
00096 }
00097 unsigned int NDim() const {
00098 return fX->getSize();
00099 }
00100 const double * Parameters() const {
00101 if (fParamValues.size() != NPar() )
00102 fParamValues.resize(NPar() );
00103
00104
00105 TIterator* itr = fParams->createIterator() ;
00106 std::vector<double>::iterator vpitr = fParamValues.begin();
00107
00108 RooRealVar* var = 0;
00109 while( ( var = dynamic_cast<RooRealVar*>(itr->Next() ) ) ) {
00110 assert(var != 0);
00111 *vpitr++ = var->getVal();
00112 }
00113 return &fParamValues.front();
00114 }
00115
00116 std::string ParameterName(unsigned int i) const {
00117
00118 TIterator* itr = fParams->createIterator() ;
00119 RooRealVar* var = 0;
00120 unsigned int index = 0;
00121 while( ( var = dynamic_cast<RooRealVar*>(itr->Next() ) ) ) {
00122 assert(var != 0);
00123 if (index == i) return std::string(var->GetName() );
00124 index++;
00125 }
00126 return "not_found";
00127 }
00128
00129
00130
00131
00132
00133
00134 void SetParameters(const double * p) {
00135 DoSetParameters(p);
00136 }
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157 private:
00158
00159 double DoEvalPar(const double * x, const double * p) const {
00160
00161
00162 DoSetParameters(p);
00163
00164
00165 TIterator* itr = fX->createIterator() ;
00166 RooRealVar* var = 0;
00167 while( ( var = dynamic_cast<RooRealVar*>(itr->Next() ) ) ) {
00168 assert(var != 0);
00169 #ifndef _WIN32
00170 var->setDirtyInhibit(true);
00171 #endif
00172 var->setVal(*x++);
00173 }
00174
00175
00176
00177 if (fNorm)
00178 return fPdf->getVal(fX);
00179 else
00180 return fPdf->getVal();
00181
00182 }
00183
00184
00185 void DoSetParameters(const double * p) const {
00186
00187 TIterator* itr = fParams->createIterator() ;
00188 RooRealVar* var = 0;
00189 while( ( var = dynamic_cast<RooRealVar*>(itr->Next() ) ) ) {
00190 assert(var != 0);
00191 var->setVal(*p++);
00192 }
00193
00194
00195 }
00196
00197
00198 bool fNorm;
00199 mutable RooAbsPdf * fPdf;
00200 mutable RooArgSet * fX;
00201 mutable RooArgSet * fParams;
00202 mutable std::vector<double> fParamValues;
00203
00204
00205 };
00206
00207
00208
00209 #endif