00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027 #include "RooFit.h"
00028 #include "Riostream.h"
00029
00030 #include "RooGenFitStudy.h"
00031 #include "RooWorkspace.h"
00032 #include "RooMsgService.h"
00033 #include "RooDataSet.h"
00034 #include "RooAbsPdf.h"
00035 #include "RooRealVar.h"
00036 #include "RooGlobalFunc.h"
00037 #include "RooFitResult.h"
00038
00039
00040 using namespace std ;
00041
00042 ClassImp(RooGenFitStudy)
00043 ;
00044
00045
00046
00047 RooGenFitStudy::RooGenFitStudy(const char* name, const char* title) :
00048 RooAbsStudy(name?name:"RooGenFitStudy",title?title:"RooGenFitStudy"),
00049 _genPdf(0),
00050 _fitPdf(0),
00051 _genSpec(0),
00052 _nllVar(0),
00053 _ngenVar(0),
00054 _params(0),
00055 _initParams(0)
00056 {
00057
00058 }
00059
00060
00061
00062
00063 RooGenFitStudy::RooGenFitStudy(const RooGenFitStudy& other) :
00064 RooAbsStudy(other),
00065 _genPdfName(other._genPdfName),
00066 _genObsName(other._genObsName),
00067 _fitPdfName(other._fitPdfName),
00068 _fitObsName(other._fitObsName),
00069 _genPdf(0),
00070 _fitPdf(0),
00071 _genSpec(0),
00072 _nllVar(0),
00073 _ngenVar(0),
00074 _params(0),
00075 _initParams(0)
00076 {
00077
00078 TIterator* giter = other._genOpts.MakeIterator() ;
00079 TObject* o ;
00080 while((o=giter->Next())) {
00081 _genOpts.Add(o->Clone()) ;
00082 }
00083 delete giter ;
00084
00085 TIterator* fiter = other._fitOpts.MakeIterator() ;
00086 while((o=fiter->Next())) {
00087 _fitOpts.Add(o->Clone()) ;
00088 }
00089 delete fiter ;
00090
00091 }
00092
00093
00094
00095
00096 RooGenFitStudy::~RooGenFitStudy()
00097 {
00098 if (_params) delete _params ;
00099 }
00100
00101
00102
00103
00104 Bool_t RooGenFitStudy::attach(RooWorkspace& w)
00105 {
00106
00107 Bool_t ret = kFALSE ;
00108
00109 RooAbsPdf* pdf = w.pdf(_genPdfName.c_str()) ;
00110 if (pdf) {
00111 _genPdf = pdf ;
00112 } else {
00113 coutE(InputArguments) << "RooGenFitStudy(" << GetName() << ") ERROR: generator p.d.f named " << _genPdfName << " not found in workspace " << w.GetName() << endl ;
00114 ret = kTRUE ;
00115 }
00116
00117 _genObs.add(w.argSet(_genObsName.c_str())) ;
00118 if (_genObs.getSize()==0) {
00119 coutE(InputArguments) << "RooGenFitStudy(" << GetName() << ") ERROR: no generator observables defined" << endl ;
00120 ret = kTRUE ;
00121 }
00122
00123 pdf = w.pdf(_fitPdfName.c_str()) ;
00124 if (pdf) {
00125 _fitPdf = pdf ;
00126 } else {
00127 coutE(InputArguments) << "RooGenFitStudy(" << GetName() << ") ERROR: fitting p.d.f named " << _fitPdfName << " not found in workspace " << w.GetName() << endl ;
00128 ret = kTRUE ;
00129 }
00130
00131 _fitObs.add(w.argSet(_fitObsName.c_str())) ;
00132 if (_fitObs.getSize()==0) {
00133 coutE(InputArguments) << "RooGenFitStudy(" << GetName() << ") ERROR: no fitting observables defined" << endl ;
00134 ret = kTRUE ;
00135 }
00136
00137 return ret ;
00138 }
00139
00140
00141
00142
00143 void RooGenFitStudy::setGenConfig(const char* pdfName, const char* obsName, const RooCmdArg& arg1,const RooCmdArg& arg2,const RooCmdArg& arg3)
00144 {
00145 _genPdfName = pdfName ;
00146 _genObsName = obsName ;
00147 _genOpts.Add(arg1.Clone()) ;
00148 _genOpts.Add(arg2.Clone()) ;
00149 _genOpts.Add(arg3.Clone()) ;
00150 }
00151
00152
00153
00154
00155 void RooGenFitStudy::setFitConfig(const char* pdfName, const char* obsName, const RooCmdArg& arg1,const RooCmdArg& arg2,const RooCmdArg& arg3)
00156 {
00157 _fitPdfName = pdfName ;
00158 _fitObsName = obsName ;
00159 _fitOpts.Add(arg1.Clone()) ;
00160 _fitOpts.Add(arg2.Clone()) ;
00161 _fitOpts.Add(arg3.Clone()) ;
00162 }
00163
00164
00165
00166
00167 Bool_t RooGenFitStudy::initialize()
00168 {
00169
00170
00171 _nllVar = new RooRealVar("NLL","-log(Likelihood)",0) ;
00172 _ngenVar = new RooRealVar("ngen","number of generated events",0) ;
00173
00174 _params = _fitPdf->getParameters(_genObs) ;
00175 _initParams = (RooArgSet*) _params->snapshot() ;
00176 _params->add(*_nllVar) ;
00177 _params->add(*_ngenVar) ;
00178
00179 _genSpec = _genPdf->prepareMultiGen(_genObs,(RooCmdArg&)*_genOpts.At(0),(RooCmdArg&)*_genOpts.At(1),(RooCmdArg&)*_genOpts.At(2)) ;
00180
00181 registerSummaryOutput(*_params) ;
00182 return kFALSE ;
00183 }
00184
00185
00186
00187
00188 Bool_t RooGenFitStudy::execute()
00189 {
00190
00191 *_params = *_initParams ;
00192 RooDataSet* data = _genPdf->generate(*_genSpec) ;
00193 RooFitResult* fr = _fitPdf->fitTo(*data,RooFit::Save(kTRUE),(RooCmdArg&)*_fitOpts.At(0),(RooCmdArg&)*_fitOpts.At(1),(RooCmdArg&)*_fitOpts.At(2)) ;
00194
00195 if (fr->status()==0) {
00196 _ngenVar->setVal(data->sumEntries()) ;
00197 _nllVar->setVal(fr->minNll()) ;
00198 storeSummaryOutput(*_params) ;
00199 storeDetailedOutput(*fr) ;
00200 }
00201
00202 delete data ;
00203 return kFALSE ;
00204 }
00205
00206
00207
00208
00209 Bool_t RooGenFitStudy::finalize()
00210 {
00211
00212 delete _params ;
00213 delete _nllVar ;
00214 delete _ngenVar ;
00215 delete _initParams ;
00216 delete _genSpec ;
00217 _params = 0 ;
00218 _nllVar = 0 ;
00219 _ngenVar = 0 ;
00220 _initParams = 0 ;
00221 _genSpec = 0 ;
00222
00223
00224 return kFALSE ;
00225 }
00226
00227
00228
00229 void RooGenFitStudy::Print(Option_t* ) const
00230 {
00231 }
00232
00233