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
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044 #include "RooFit.h"
00045 #include "Riostream.h"
00046
00047 #include "RooMCStudy.h"
00048 #include "RooAbsMCStudyModule.h"
00049
00050 #include "RooGenContext.h"
00051 #include "RooAbsPdf.h"
00052 #include "RooDataSet.h"
00053 #include "RooDataHist.h"
00054 #include "RooRealVar.h"
00055 #include "RooFitResult.h"
00056 #include "RooErrorVar.h"
00057 #include "RooFormulaVar.h"
00058 #include "RooArgList.h"
00059 #include "RooPlot.h"
00060 #include "RooGenericPdf.h"
00061 #include "RooRandom.h"
00062 #include "RooCmdConfig.h"
00063 #include "RooGlobalFunc.h"
00064 #include "RooPullVar.h"
00065 #include "RooMsgService.h"
00066 #include "RooProdPdf.h"
00067
00068 using namespace std ;
00069
00070 ClassImp(RooMCStudy)
00071 ;
00072
00073
00074
00075 RooMCStudy::RooMCStudy(const RooAbsPdf& model, const RooArgSet& observables,
00076 const RooCmdArg& arg1, const RooCmdArg& arg2,
00077 const RooCmdArg& arg3,const RooCmdArg& arg4,const RooCmdArg& arg5,
00078 const RooCmdArg& arg6,const RooCmdArg& arg7,const RooCmdArg& arg8) : TNamed("mcstudy","mcstudy")
00079
00080 {
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112 RooLinkedList cmdList;
00113 cmdList.Add(const_cast<RooCmdArg*>(&arg1)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg2)) ;
00114 cmdList.Add(const_cast<RooCmdArg*>(&arg3)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg4)) ;
00115 cmdList.Add(const_cast<RooCmdArg*>(&arg5)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg6)) ;
00116 cmdList.Add(const_cast<RooCmdArg*>(&arg7)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg8)) ;
00117
00118
00119 RooCmdConfig pc(Form("RooMCStudy::RooMCStudy(%s)",model.GetName())) ;
00120
00121 pc.defineObject("fitModel","FitModel",0,0) ;
00122 pc.defineObject("condObs","ProjectedDependents",0,0) ;
00123 pc.defineObject("protoData","PrototypeData",0,0) ;
00124 pc.defineSet("cPars","Constrain",0,0) ;
00125 pc.defineSet("extCons","ExternalConstraints",0,0) ;
00126 pc.defineInt("silence","Silence",0,0) ;
00127 pc.defineInt("randProtoData","PrototypeData",0,0) ;
00128 pc.defineInt("verboseGen","Verbose",0,0) ;
00129 pc.defineInt("extendedGen","Extended",0,0) ;
00130 pc.defineInt("binGenData","Binned",0,0) ;
00131 pc.defineString("fitOpts","FitOptions",0,"") ;
00132 pc.defineInt("dummy","FitOptArgs",0,0) ;
00133 pc.defineMutex("FitOptions","FitOptArgs") ;
00134 pc.defineMutex("Constrain","FitOptions") ;
00135 pc.defineMutex("ExternalConstraints","FitOptions") ;
00136
00137
00138 pc.process(cmdList) ;
00139 if (!pc.ok(kTRUE)) {
00140
00141 throw std::string("RooMCStudy::RooMCStudy() Error in parsing arguments passed to contructor") ;
00142 return ;
00143 }
00144
00145
00146 if (pc.hasProcessed("FitOptArgs")) {
00147 RooCmdArg* fitOptArg = static_cast<RooCmdArg*>(cmdList.FindObject("FitOptArgs")) ;
00148 for (Int_t i=0 ; i<fitOptArg->subArgs().GetSize() ;i++) {
00149 _fitOptList.Add(new RooCmdArg(static_cast<RooCmdArg&>(*fitOptArg->subArgs().At(i)))) ;
00150 }
00151 }
00152
00153
00154 _silence = pc.getInt("silence") ;
00155 _verboseGen = pc.getInt("verboseGen") ;
00156 _extendedGen = pc.getInt("extendedGen") ;
00157 _binGenData = pc.getInt("binGenData") ;
00158 _randProto = pc.getInt("randProtoData") ;
00159
00160
00161 const RooArgSet* cParsTmp = pc.getSet("cPars") ;
00162 const RooArgSet* extCons = pc.getSet("extCons") ;
00163
00164 RooArgSet* cPars = new RooArgSet ;
00165 if (cParsTmp) {
00166 cPars->add(*cParsTmp) ;
00167 }
00168
00169
00170 if (cPars) {
00171 _fitOptList.Add(RooFit::Constrain(*cPars).Clone()) ;
00172 }
00173 if (extCons) {
00174 _fitOptList.Add(RooFit::ExternalConstraints(*extCons).Clone()) ;
00175 }
00176
00177
00178 RooArgSet allConstraints ;
00179 RooArgSet consPars ;
00180 if (cPars) {
00181 RooArgSet* constraints = model.getConstraints(observables,*cPars,kTRUE) ;
00182 if (constraints) {
00183 allConstraints.add(*constraints) ;
00184 delete constraints ;
00185 }
00186 }
00187
00188
00189 if (allConstraints.getSize()>0) {
00190 _constrPdf = new RooProdPdf("mcs_constr_prod","RooMCStudy constraints product",allConstraints) ;
00191
00192 if (cPars) {
00193 consPars.add(*cPars) ;
00194 } else {
00195 RooArgSet* params = model.getParameters(observables) ;
00196 RooArgSet* cparams = _constrPdf->getObservables(*params) ;
00197 consPars.add(*cparams) ;
00198 delete params ;
00199 delete cparams ;
00200 }
00201 _constrGenContext = _constrPdf->genContext(consPars,0,0,_verboseGen) ;
00202
00203 _perExptGenParams = kTRUE ;
00204
00205 coutI(Generation) << "RooMCStudy::RooMCStudy: INFO have pdf with constraints, will generate paramaters from constraint pdf for each experiment" << endl ;
00206
00207
00208 } else {
00209 _constrPdf = 0 ;
00210 _constrGenContext=0 ;
00211
00212 _perExptGenParams = kFALSE ;
00213 }
00214
00215
00216
00217 _genModel = const_cast<RooAbsPdf*>(&model) ;
00218 _genSample = 0 ;
00219 RooAbsPdf* fitModel = static_cast<RooAbsPdf*>(pc.getObject("fitModel",0)) ;
00220 _fitModel = fitModel ? fitModel : _genModel ;
00221
00222
00223 _genProtoData = static_cast<RooDataSet*>(pc.getObject("protoData",0)) ;
00224 if (pc.getObject("condObs",0)) {
00225 _projDeps.add(static_cast<RooArgSet&>(*pc.getObject("condObs",0))) ;
00226 }
00227
00228 _dependents.add(observables) ;
00229
00230 _allDependents.add(_dependents) ;
00231 _fitOptions = pc.getString("fitOpts") ;
00232 _canAddFitResults = kTRUE ;
00233
00234 if (_extendedGen && _genProtoData && !_randProto) {
00235 oocoutW(_fitModel,Generation) << "RooMCStudy::RooMCStudy: WARNING Using generator option 'e' (Poisson distribution of #events) together " << endl
00236 << " with a prototype dataset implies incomplete sampling or oversampling of proto data." << endl
00237 << " Use option \"r\" to randomize prototype dataset order and thus to randomize" << endl
00238 << " the set of over/undersampled prototype events for each generation cycle." << endl ;
00239 }
00240
00241 _genParams = _genModel->getParameters(&_dependents) ;
00242 if (!_binGenData) {
00243 _genContext = _genModel->genContext(_dependents,_genProtoData,0,_verboseGen) ;
00244 _genContext->attach(*_genParams) ;
00245 } else {
00246 _genContext = 0 ;
00247 }
00248
00249 _genInitParams = (RooArgSet*) _genParams->snapshot(kFALSE) ;
00250
00251
00252 _fitParams = _fitModel->getParameters(&_dependents) ;
00253 _fitInitParams = (RooArgSet*) _fitParams->snapshot(kTRUE) ;
00254
00255 _nExpGen = _extendedGen ? _genModel->expectedEvents(&_dependents) : 0 ;
00256
00257
00258 _nllVar = new RooRealVar("NLL","-log(Likelihood)",0) ;
00259
00260
00261 _ngenVar = new RooRealVar("ngen","number of generated events",0) ;
00262
00263
00264 RooArgSet tmp2(*_fitParams) ;
00265 tmp2.add(*_nllVar) ;
00266 tmp2.add(*_ngenVar) ;
00267
00268
00269 tmp2.setAttribAll("StoreError",kTRUE) ;
00270 tmp2.setAttribAll("StoreAsymError",kTRUE) ;
00271 TString fpdName ;
00272 if (_fitModel==_genModel) {
00273 fpdName = Form("fitParData_%s",_fitModel->GetName()) ;
00274 } else {
00275 fpdName= Form("fitParData_%s_%s",_fitModel->GetName(),_genModel->GetName()) ;
00276 }
00277
00278 _fitParData = new RooDataSet(fpdName.Data(),"Fit Parameters DataSet",tmp2) ;
00279 tmp2.setAttribAll("StoreError",kFALSE) ;
00280 tmp2.setAttribAll("StoreAsymError",kFALSE) ;
00281
00282 if (_perExptGenParams) {
00283 _genParData = new RooDataSet("genParData","Generated Parameters dataset",*_genParams) ;
00284 } else {
00285 _genParData = 0 ;
00286 }
00287
00288
00289 if (_genProtoData) {
00290 _allDependents.add(*_genProtoData->get(),kTRUE) ;
00291 }
00292
00293
00294 list<RooAbsMCStudyModule*>::iterator iter ;
00295 for (iter=_modList.begin() ; iter!= _modList.end() ; ++iter) {
00296 Bool_t ok = (*iter)->doInitializeInstance(*this) ;
00297 if (!ok) {
00298 oocoutE(_fitModel,Generation) << "RooMCStudy::ctor: removing study module " << (*iter)->GetName() << " from analysis chain because initialization failed" << endl ;
00299 iter = _modList.erase(iter) ;
00300 }
00301 }
00302
00303 }
00304
00305
00306
00307 RooMCStudy::RooMCStudy(const RooAbsPdf& genModel, const RooAbsPdf& fitModel,
00308 const RooArgSet& dependents, const char* genOptions,
00309 const char* fitOptions, const RooDataSet* genProtoData,
00310 const RooArgSet& projDeps) :
00311 TNamed("mcstudy","mcstudy"),
00312 _genModel((RooAbsPdf*)&genModel),
00313 _genProtoData(genProtoData),
00314 _projDeps(projDeps),
00315 _constrPdf(0),
00316 _constrGenContext(0),
00317 _dependents(dependents),
00318 _allDependents(dependents),
00319 _fitModel((RooAbsPdf*)&fitModel),
00320 _nllVar(0),
00321 _ngenVar(0),
00322 _genParData(0),
00323 _fitOptions(fitOptions),
00324 _canAddFitResults(kTRUE),
00325 _perExptGenParams(0),
00326 _silence(kFALSE)
00327 {
00328
00329
00330
00331
00332
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344
00345 TString genOpt(genOptions) ;
00346 genOpt.ToLower() ;
00347 _verboseGen = genOpt.Contains("v") ;
00348 _extendedGen = genOpt.Contains("e") ;
00349 _binGenData = genOpt.Contains("b") ;
00350 _randProto = genOpt.Contains("r") ;
00351
00352 if (_extendedGen && genProtoData && !_randProto) {
00353 oocoutE(_fitModel,Generation) << "RooMCStudy::RooMCStudy: WARNING Using generator option 'e' (Poisson distribution of #events) together " << endl
00354 << " with a prototype dataset implies incomplete sampling or oversampling of proto data." << endl
00355 << " Use option \"r\" to randomize prototype dataset order and thus to randomize" << endl
00356 << " the set of over/undersampled prototype events for each generation cycle." << endl ;
00357 }
00358
00359 if (!_binGenData) {
00360 _genContext = genModel.genContext(dependents,genProtoData,0,_verboseGen) ;
00361 } else {
00362 _genContext = 0 ;
00363 }
00364 _genParams = _genModel->getParameters(&_dependents) ;
00365 _genSample = 0 ;
00366 RooArgSet* tmp = genModel.getParameters(&dependents) ;
00367 _genInitParams = (RooArgSet*) tmp->snapshot(kFALSE) ;
00368 delete tmp ;
00369
00370
00371 _fitParams = fitModel.getParameters(&dependents) ;
00372 _fitInitParams = (RooArgSet*) _fitParams->snapshot(kTRUE) ;
00373
00374 _nExpGen = _extendedGen ? genModel.expectedEvents(&dependents) : 0 ;
00375
00376
00377 _nllVar = new RooRealVar("NLL","-log(Likelihood)",0) ;
00378
00379
00380 _ngenVar = new RooRealVar("ngen","number of generated events",0) ;
00381
00382
00383 RooArgSet tmp2(*_fitParams) ;
00384 tmp2.add(*_nllVar) ;
00385 tmp2.add(*_ngenVar) ;
00386
00387
00388 tmp2.setAttribAll("StoreError",kTRUE) ;
00389 tmp2.setAttribAll("StoreAsymError",kTRUE) ;
00390 _fitParData = new RooDataSet("fitParData","Fit Parameters DataSet",tmp2) ;
00391 tmp2.setAttribAll("StoreError",kFALSE) ;
00392 tmp2.setAttribAll("StoreAsymError",kFALSE) ;
00393
00394
00395 if (genProtoData) {
00396 _allDependents.add(*genProtoData->get(),kTRUE) ;
00397 }
00398
00399
00400 list<RooAbsMCStudyModule*>::iterator iter ;
00401 for (iter=_modList.begin() ; iter!= _modList.end() ; ++iter) {
00402 Bool_t ok = (*iter)->doInitializeInstance(*this) ;
00403 if (!ok) {
00404 oocoutE(_fitModel,Generation) << "RooMCStudy::ctor: removing study module " << (*iter)->GetName() << " from analysis chain because initialization failed" << endl ;
00405 iter = _modList.erase(iter) ;
00406 }
00407 }
00408
00409 }
00410
00411
00412
00413
00414 RooMCStudy::~RooMCStudy()
00415 {
00416
00417
00418 _genDataList.Delete() ;
00419 _fitResList.Delete() ;
00420 _fitOptList.Delete() ;
00421 delete _ngenVar ;
00422 delete _fitParData ;
00423 delete _genParData ;
00424 delete _fitInitParams ;
00425 delete _fitParams ;
00426 delete _genInitParams ;
00427 delete _genParams ;
00428 delete _genContext ;
00429 delete _nllVar ;
00430 delete _constrPdf ;
00431 delete _constrGenContext ;
00432 }
00433
00434
00435
00436
00437 void RooMCStudy::addModule(RooAbsMCStudyModule& module)
00438 {
00439
00440
00441
00442 module.doInitializeInstance(*this) ;
00443 _modList.push_back(&module) ;
00444 }
00445
00446
00447
00448
00449 Bool_t RooMCStudy::run(Bool_t doGenerate, Bool_t DoFit, Int_t nSamples, Int_t nEvtPerSample, Bool_t keepGenData, const char* asciiFilePat)
00450 {
00451
00452
00453
00454
00455
00456
00457
00458
00459
00460
00461
00462
00463 RooFit::MsgLevel oldLevel(RooFit::FATAL) ;
00464 if (_silence) {
00465 oldLevel = RooMsgService::instance().globalKillBelow() ;
00466 RooMsgService::instance().setGlobalKillBelow(RooFit::PROGRESS) ;
00467 }
00468
00469 list<RooAbsMCStudyModule*>::iterator iter ;
00470 for (iter=_modList.begin() ; iter!= _modList.end() ; ++iter) {
00471 (*iter)->initializeRun(nSamples) ;
00472 }
00473
00474 Int_t prescale = nSamples>100 ? Int_t(nSamples/100) : 1 ;
00475
00476 while(nSamples--) {
00477
00478 if (nSamples%prescale==0) {
00479 oocoutP(_fitModel,Generation) << "RooMCStudy::run: " ;
00480 if (doGenerate) ooccoutI(_fitModel,Generation) << "Generating " ;
00481 if (doGenerate && DoFit) ooccoutI(_fitModel,Generation) << "and " ;
00482 if (DoFit) ooccoutI(_fitModel,Generation) << "fitting " ;
00483 ooccoutP(_fitModel,Generation) << "sample " << nSamples << endl ;
00484 }
00485
00486 _genSample = 0;
00487 Bool_t existingData = kFALSE ;
00488 if (doGenerate) {
00489
00490 Int_t nEvt(nEvtPerSample) ;
00491
00492
00493 *_genParams = *_genInitParams ;
00494
00495
00496 if (_constrPdf) {
00497 RooDataSet* tmp = _constrGenContext->generate(1) ;
00498 *_genParams = *tmp->get() ;
00499 delete tmp ;
00500 }
00501
00502
00503 if (_genParData) {
00504 _genParData->add(*_genParams) ;
00505 }
00506
00507
00508 list<RooAbsMCStudyModule*>::iterator iter2 ;
00509 for (iter2=_modList.begin() ; iter2!= _modList.end() ; ++iter2) {
00510 (*iter2)->processBeforeGen(nSamples) ;
00511 }
00512
00513 if (_binGenData) {
00514
00515
00516 if (_extendedGen) {
00517 _nExpGen = _genModel->expectedEvents(&_dependents) ;
00518 nEvt = RooRandom::randomGenerator()->Poisson(nEvtPerSample==0?_nExpGen:nEvtPerSample) ;
00519 }
00520
00521
00522 _genSample = _genModel->generateBinned(_dependents,nEvt) ;
00523
00524 } else {
00525
00526
00527 if (_extendedGen) {
00528 _nExpGen = _genModel->expectedEvents(&_dependents) ;
00529 nEvt = RooRandom::randomGenerator()->Poisson(nEvtPerSample==0?_nExpGen:nEvtPerSample) ;
00530 }
00531
00532
00533 if (_randProto && _genProtoData && _genProtoData->numEntries()!=nEvt) {
00534 oocoutI(_fitModel,Generation) << "RooMCStudy: (Re)randomizing event order in prototype dataset (Nevt=" << nEvt << ")" << endl ;
00535 Int_t* newOrder = _genModel->randomizeProtoOrder(_genProtoData->numEntries(),nEvt) ;
00536 _genContext->setProtoDataOrder(newOrder) ;
00537 delete[] newOrder ;
00538 }
00539
00540
00541 if (nEvt>0) {
00542 _genSample = _genContext->generate(nEvt) ;
00543 } else {
00544
00545 _genSample = new RooDataSet("emptySample","emptySample",_dependents) ;
00546 }
00547 }
00548
00549
00550
00551 } else if (asciiFilePat) {
00552
00553
00554 char asciiFile[1024] ;
00555 snprintf(asciiFile,1024,asciiFilePat,nSamples) ;
00556 RooArgList depList(_allDependents) ;
00557 _genSample = RooDataSet::read(asciiFile,depList,"q") ;
00558
00559 } else {
00560
00561
00562 _genSample = (RooDataSet*) _genDataList.At(nSamples) ;
00563 existingData = kTRUE ;
00564 if (!_genSample) {
00565 oocoutW(_fitModel,Generation) << "RooMCStudy::run: WARNING: Sample #" << nSamples << " not loaded, skipping" << endl ;
00566 continue ;
00567 }
00568 }
00569
00570
00571 _ngenVar->setVal(_genSample->sumEntries()) ;
00572
00573
00574 list<RooAbsMCStudyModule*>::iterator iter3 ;
00575 for (iter3=_modList.begin() ; iter3!= _modList.end() ; ++iter3) {
00576 (*iter3)->processBetweenGenAndFit(nSamples) ;
00577 }
00578
00579 if (DoFit) fitSample(_genSample) ;
00580
00581
00582 for (iter3=_modList.begin() ; iter3!= _modList.end() ; ++iter3) {
00583 (*iter3)->processAfterFit(nSamples) ;
00584 }
00585
00586
00587 if (doGenerate && asciiFilePat && *asciiFilePat) {
00588 char asciiFile[1024] ;
00589 snprintf(asciiFile,1024,asciiFilePat,nSamples) ;
00590 RooDataSet* unbinnedData = dynamic_cast<RooDataSet*>(_genSample) ;
00591 if (unbinnedData) {
00592 unbinnedData->write(asciiFile) ;
00593 } else {
00594 coutE(InputArguments) << "RooMCStudy::run(" << GetName() << ") ERROR: ASCII writing of binned datasets is not supported" << endl ;
00595 }
00596 }
00597
00598
00599 if (!existingData) {
00600 if (keepGenData) {
00601 _genDataList.Add(_genSample) ;
00602 } else {
00603 delete _genSample ;
00604 }
00605 }
00606 }
00607
00608 for (iter=_modList.begin() ; iter!= _modList.end() ; ++iter) {
00609 RooDataSet* auxData = (*iter)->finalizeRun() ;
00610 if (auxData) {
00611 _fitParData->merge(auxData) ;
00612 }
00613 }
00614
00615 _canAddFitResults = kFALSE ;
00616
00617 if (_genParData) {
00618 const RooArgSet* genPars = _genParData->get() ;
00619 TIterator* iter2 = genPars->createIterator() ;
00620 RooAbsArg* arg ;
00621 while((arg=(RooAbsArg*)iter2->Next())) {
00622 _genParData->changeObservableName(arg->GetName(),Form("%s_gen",arg->GetName())) ;
00623 }
00624 delete iter2 ;
00625
00626 _fitParData->merge(_genParData) ;
00627 }
00628
00629 if (DoFit) calcPulls() ;
00630
00631 if (_silence) {
00632 RooMsgService::instance().setGlobalKillBelow(oldLevel) ;
00633 }
00634
00635 return kFALSE ;
00636 }
00637
00638
00639
00640
00641
00642
00643
00644 Bool_t RooMCStudy::generateAndFit(Int_t nSamples, Int_t nEvtPerSample, Bool_t keepGenData, const char* asciiFilePat)
00645 {
00646
00647
00648
00649
00650
00651
00652
00653
00654
00655
00656 _fitResList.Delete() ;
00657 _genDataList.Delete() ;
00658 _fitParData->reset() ;
00659
00660 return run(kTRUE,kTRUE,nSamples,nEvtPerSample,keepGenData,asciiFilePat) ;
00661 }
00662
00663
00664
00665
00666 Bool_t RooMCStudy::generate(Int_t nSamples, Int_t nEvtPerSample, Bool_t keepGenData, const char* asciiFilePat)
00667 {
00668
00669
00670
00671
00672
00673
00674
00675
00676
00677
00678 _genDataList.Delete() ;
00679
00680 return run(kTRUE,kFALSE,nSamples,nEvtPerSample,keepGenData,asciiFilePat) ;
00681 }
00682
00683
00684
00685
00686 Bool_t RooMCStudy::fit(Int_t nSamples, const char* asciiFilePat)
00687 {
00688
00689
00690
00691
00692
00693
00694
00695 _fitResList.Delete() ;
00696 _fitParData->reset() ;
00697
00698 return run(kFALSE,kTRUE,nSamples,0,kFALSE,asciiFilePat) ;
00699 }
00700
00701
00702
00703
00704 Bool_t RooMCStudy::fit(Int_t nSamples, TList& dataSetList)
00705 {
00706
00707
00708
00709
00710 _fitResList.Delete() ;
00711 _genDataList.Delete() ;
00712 _fitParData->reset() ;
00713
00714
00715 TIterator* iter = dataSetList.MakeIterator() ;
00716 RooAbsData* gset ;
00717 while((gset=(RooAbsData*)iter->Next())) {
00718 _genDataList.Add(gset) ;
00719 }
00720 delete iter ;
00721
00722 return run(kFALSE,kTRUE,nSamples,0,kTRUE,0) ;
00723 }
00724
00725
00726
00727
00728 void RooMCStudy::resetFitParams()
00729 {
00730
00731
00732
00733 *_fitParams = *_fitInitParams ;
00734 }
00735
00736
00737
00738
00739 RooFitResult* RooMCStudy::doFit(RooAbsData* genSample)
00740 {
00741
00742
00743
00744 TString fitOpt2(_fitOptions) ; fitOpt2.Append("r") ;
00745 if (_silence) {
00746 fitOpt2.Append("b") ;
00747 }
00748
00749
00750 RooAbsData* data ;
00751 if (_binGenData) {
00752 RooArgSet* depList = _fitModel->getObservables(genSample) ;
00753 data = new RooDataHist(genSample->GetName(),genSample->GetTitle(),*depList,*genSample) ;
00754 delete depList ;
00755 } else {
00756 data = genSample ;
00757 }
00758
00759 RooFitResult* fr ;
00760 if (_fitOptList.GetSize()==0) {
00761 if (_projDeps.getSize()>0) {
00762 fr = (RooFitResult*) _fitModel->fitTo(*data,RooFit::ConditionalObservables(_projDeps),RooFit::FitOptions(fitOpt2)) ;
00763 } else {
00764 fr = (RooFitResult*) _fitModel->fitTo(*data,RooFit::FitOptions(fitOpt2)) ;
00765 }
00766 } else {
00767 RooCmdArg save = RooFit::Save() ;
00768 RooCmdArg condo = RooFit::ConditionalObservables(_projDeps) ;
00769 RooCmdArg plevel = RooFit::PrintLevel(-1) ;
00770 RooLinkedList fitOptList(_fitOptList) ;
00771 fitOptList.Add(&save) ;
00772 if (_projDeps.getSize()>0) {
00773 fitOptList.Add(&condo) ;
00774 }
00775 if (_silence) {
00776 fitOptList.Add(&plevel) ;
00777 }
00778 fr = (RooFitResult*) _fitModel->fitTo(*data,fitOptList) ;
00779 }
00780
00781 if (_binGenData) delete data ;
00782
00783 return fr ;
00784 }
00785
00786
00787
00788
00789 RooFitResult* RooMCStudy::refit(RooAbsData* genSample)
00790 {
00791
00792
00793
00794 if (!genSample) {
00795 genSample = _genSample ;
00796 }
00797
00798 RooFitResult* fr(0) ;
00799 if (genSample->sumEntries()>0) {
00800 fr = doFit(genSample) ;
00801 }
00802
00803 return fr ;
00804 }
00805
00806
00807
00808
00809 Bool_t RooMCStudy::fitSample(RooAbsData* genSample)
00810 {
00811
00812
00813
00814
00815
00816
00817
00818
00819
00820
00821
00822 resetFitParams() ;
00823
00824
00825 Bool_t ok ;
00826 RooFitResult* fr(0) ;
00827 if (genSample->sumEntries()>0) {
00828 fr = doFit(genSample) ;
00829 ok = (fr->status()==0) ;
00830 } else {
00831 ok = kFALSE ;
00832 }
00833
00834
00835 if (ok) {
00836 _nllVar->setVal(fr->minNll()) ;
00837 RooArgSet tmp(*_fitParams) ;
00838 tmp.add(*_nllVar) ;
00839 tmp.add(*_ngenVar) ;
00840 _fitParData->add(tmp) ;
00841 }
00842
00843
00844 Bool_t userSaveRequest = kFALSE ;
00845 if (_fitOptList.GetSize()>0) {
00846 if (_fitOptList.FindObject("Save")) userSaveRequest = kTRUE ;
00847 } else {
00848 if (_fitOptions.Contains("r")) userSaveRequest = kTRUE ;
00849 }
00850
00851 if (userSaveRequest) {
00852 _fitResList.Add(fr) ;
00853 } else {
00854 delete fr ;
00855 }
00856
00857 return !ok ;
00858 }
00859
00860
00861
00862
00863 Bool_t RooMCStudy::addFitResult(const RooFitResult& fr)
00864 {
00865
00866
00867
00868
00869
00870
00871
00872
00873
00874 if (!_canAddFitResults) {
00875 oocoutE(_fitModel,InputArguments) << "RooMCStudy::addFitResult: ERROR cannot add fit results in current state" << endl ;
00876 return kTRUE ;
00877 }
00878
00879
00880 *_fitParams = RooArgSet(fr.floatParsFinal()) ;
00881
00882
00883 Bool_t ok = (fr.status()==0) ;
00884 if (ok) {
00885 _nllVar->setVal(fr.minNll()) ;
00886 RooArgSet tmp(*_fitParams) ;
00887 tmp.add(*_nllVar) ;
00888 tmp.add(*_ngenVar) ;
00889 _fitParData->add(tmp) ;
00890 }
00891
00892
00893 if (_fitOptions.Contains("r")) {
00894 _fitResList.Add((TObject*)&fr) ;
00895 }
00896
00897 return kFALSE ;
00898 }
00899
00900
00901
00902
00903 void RooMCStudy::calcPulls()
00904 {
00905
00906
00907
00908 TIterator* iter = _fitParams->createIterator() ;
00909 RooRealVar* par ;
00910 while((par=(RooRealVar*)iter->Next())) {
00911
00912 RooErrorVar* err = par->errorVar() ;
00913 _fitParData->addColumn(*err) ;
00914 delete err ;
00915
00916 TString name(par->GetName()), title(par->GetTitle()) ;
00917 name.Append("pull") ;
00918 title.Append(" Pull") ;
00919
00920
00921 RooAbsReal* genParOrig = (RooAbsReal*) _fitParData->get()->find(Form("%s_gen",par->GetName())) ;
00922 if (genParOrig && _perExptGenParams) {
00923
00924 RooPullVar pull(name,title,*par,*genParOrig) ;
00925 _fitParData->addColumn(pull,kFALSE) ;
00926
00927 } else {
00928
00929 genParOrig = (RooAbsReal*)_genInitParams->find(par->GetName()) ;
00930
00931 if (genParOrig) {
00932 RooAbsReal* genPar = (RooAbsReal*) genParOrig->Clone("truth") ;
00933 RooPullVar pull(name,title,*par,*genPar) ;
00934
00935 _fitParData->addColumn(pull,kFALSE) ;
00936 delete genPar ;
00937
00938 }
00939
00940 }
00941
00942 }
00943 delete iter ;
00944
00945 }
00946
00947
00948
00949
00950
00951 const RooDataSet& RooMCStudy::fitParDataSet()
00952 {
00953
00954
00955
00956
00957 if (_canAddFitResults) {
00958 calcPulls() ;
00959 _canAddFitResults = kFALSE ;
00960 }
00961
00962 return *_fitParData ;
00963 }
00964
00965
00966
00967
00968 const RooArgSet* RooMCStudy::fitParams(Int_t sampleNum) const
00969 {
00970
00971
00972
00973
00974
00975
00976
00977
00978 if (sampleNum<0 || sampleNum>=_fitParData->numEntries()) {
00979 oocoutE(_fitModel,InputArguments) << "RooMCStudy::fitParams: ERROR, invalid sample number: " << sampleNum << endl ;
00980 return 0 ;
00981 }
00982
00983 return _fitParData->get(sampleNum) ;
00984 }
00985
00986
00987
00988
00989 const RooFitResult* RooMCStudy::fitResult(Int_t sampleNum) const
00990 {
00991
00992
00993
00994 if (sampleNum<0 || sampleNum>=_fitResList.GetSize()) {
00995 oocoutE(_fitModel,InputArguments) << "RooMCStudy::fitResult: ERROR, invalid sample number: " << sampleNum << endl ;
00996 return 0 ;
00997 }
00998
00999
01000 const RooFitResult* fr = (RooFitResult*) _fitResList.At(sampleNum) ;
01001 if (fr) {
01002 return fr ;
01003 } else {
01004 oocoutE(_fitModel,InputArguments) << "RooMCStudy::fitResult: ERROR, no fit result saved for sample "
01005 << sampleNum << ", did you use the 'r; fit option?" << endl ;
01006 }
01007 return 0 ;
01008 }
01009
01010
01011
01012
01013 const RooDataSet* RooMCStudy::genData(Int_t sampleNum) const
01014 {
01015
01016
01017
01018
01019 if (_genDataList.GetSize()==0) {
01020 oocoutE(_fitModel,InputArguments) << "RooMCStudy::genData() ERROR, generated data was not saved" << endl ;
01021 return 0 ;
01022 }
01023
01024
01025 if (sampleNum<0 || sampleNum>=_genDataList.GetSize()) {
01026 oocoutE(_fitModel,InputArguments) << "RooMCStudy::genData() ERROR, invalid sample number: " << sampleNum << endl ;
01027 return 0 ;
01028 }
01029
01030 return (RooDataSet*) _genDataList.At(sampleNum) ;
01031 }
01032
01033
01034
01035
01036 RooPlot* RooMCStudy::plotParamOn(RooPlot* frame, const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4,
01037 const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
01038 {
01039
01040
01041
01042
01043
01044
01045
01046
01047 _fitParData->plotOn(frame,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8) ;
01048 return frame ;
01049 }
01050
01051
01052
01053
01054 RooPlot* RooMCStudy::plotParam(const char* paramName, const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4,
01055 const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
01056 {
01057
01058
01059
01060
01061
01062
01063
01064
01065
01066
01067
01068
01069
01070 RooRealVar* param = static_cast<RooRealVar*>(_fitParData->get()->find(paramName)) ;
01071 if (!param) {
01072 oocoutE(_fitModel,InputArguments) << "RooMCStudy::plotParam: ERROR: no parameter defined with name " << paramName << endl ;
01073 return 0 ;
01074 }
01075
01076
01077 return plotParam(*param,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8) ;
01078 }
01079
01080
01081
01082
01083 RooPlot* RooMCStudy::plotParam(const RooRealVar& param, const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4,
01084 const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
01085 {
01086
01087
01088
01089
01090
01091
01092
01093
01094
01095
01096
01097
01098 RooLinkedList cmdList;
01099 cmdList.Add(const_cast<RooCmdArg*>(&arg1)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg2)) ;
01100 cmdList.Add(const_cast<RooCmdArg*>(&arg3)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg4)) ;
01101 cmdList.Add(const_cast<RooCmdArg*>(&arg5)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg6)) ;
01102 cmdList.Add(const_cast<RooCmdArg*>(&arg7)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg8)) ;
01103
01104 RooPlot* frame = makeFrameAndPlotCmd(param, cmdList) ;
01105 if (frame) {
01106 _fitParData->plotOn(frame, cmdList) ;
01107 }
01108
01109 return frame ;
01110 }
01111
01112
01113
01114
01115 RooPlot* RooMCStudy::plotNLL(const RooCmdArg& arg1, const RooCmdArg& arg2,
01116 const RooCmdArg& arg3, const RooCmdArg& arg4,
01117 const RooCmdArg& arg5, const RooCmdArg& arg6,
01118 const RooCmdArg& arg7, const RooCmdArg& arg8)
01119 {
01120
01121
01122
01123
01124
01125
01126
01127
01128
01129
01130
01131 return plotParam(*_nllVar,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8) ;
01132 }
01133
01134
01135
01136
01137 RooPlot* RooMCStudy::plotError(const RooRealVar& param, const RooCmdArg& arg1, const RooCmdArg& arg2,
01138 const RooCmdArg& arg3, const RooCmdArg& arg4,
01139 const RooCmdArg& arg5, const RooCmdArg& arg6,
01140 const RooCmdArg& arg7, const RooCmdArg& arg8)
01141 {
01142
01143
01144
01145
01146
01147
01148
01149
01150
01151
01152
01153 if (_canAddFitResults) {
01154 calcPulls() ;
01155 _canAddFitResults=kFALSE ;
01156 }
01157
01158 RooErrorVar* evar = param.errorVar() ;
01159 RooRealVar* evar_rrv = static_cast<RooRealVar*>(evar->createFundamental()) ;
01160 RooPlot* frame = plotParam(*evar_rrv,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8) ;
01161 delete evar_rrv ;
01162 delete evar ;
01163 return frame ;
01164 }
01165
01166
01167
01168
01169 RooPlot* RooMCStudy::plotPull(const RooRealVar& param, const RooCmdArg& arg1, const RooCmdArg& arg2,
01170 const RooCmdArg& arg3, const RooCmdArg& arg4,
01171 const RooCmdArg& arg5, const RooCmdArg& arg6,
01172 const RooCmdArg& arg7, const RooCmdArg& arg8)
01173 {
01174
01175
01176
01177
01178
01179
01180
01181
01182
01183
01184
01185
01186
01187
01188 RooLinkedList cmdList;
01189 cmdList.Add(const_cast<RooCmdArg*>(&arg1)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg2)) ;
01190 cmdList.Add(const_cast<RooCmdArg*>(&arg3)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg4)) ;
01191 cmdList.Add(const_cast<RooCmdArg*>(&arg5)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg6)) ;
01192 cmdList.Add(const_cast<RooCmdArg*>(&arg7)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg8)) ;
01193
01194 TString name(param.GetName()), title(param.GetTitle()) ;
01195 name.Append("pull") ; title.Append(" Pull") ;
01196 RooRealVar pvar(name,title,-100,100) ;
01197 pvar.setBins(100) ;
01198
01199
01200 RooPlot* frame = makeFrameAndPlotCmd(pvar, cmdList, kTRUE) ;
01201 if (frame) {
01202
01203
01204 RooCmdConfig pc(Form("RooMCStudy::plotPull(%s)",_genModel->GetName())) ;
01205 pc.defineInt("fitGauss","FitGauss",0,0) ;
01206 pc.allowUndefined() ;
01207 pc.process(cmdList) ;
01208 Bool_t fitGauss=pc.getInt("fitGauss") ;
01209
01210
01211 pc.stripCmdList(cmdList,"FitGauss") ;
01212 _fitParData->plotOn(frame,cmdList) ;
01213
01214
01215 if (fitGauss) {
01216 RooRealVar pullMean("pullMean","Mean of pull",0,-10,10) ;
01217 RooRealVar pullSigma("pullSigma","Width of pull",1,0.1,5) ;
01218 RooGenericPdf pullGauss("pullGauss","Gaussian of pull",
01219 "exp(-0.5*(@0-@1)*(@0-@1)/(@2*@2))",
01220 RooArgSet(pvar,pullMean,pullSigma)) ;
01221 pullGauss.fitTo(*_fitParData,RooFit::Minos(0),RooFit::PrintLevel(-1)) ;
01222 pullGauss.plotOn(frame) ;
01223 pullGauss.paramOn(frame,_fitParData) ;
01224 }
01225 }
01226 return frame ; ;
01227 }
01228
01229
01230
01231
01232 RooPlot* RooMCStudy::makeFrameAndPlotCmd(const RooRealVar& param, RooLinkedList& cmdList, Bool_t symRange) const
01233 {
01234
01235
01236
01237
01238
01239 RooCmdConfig pc(Form("RooMCStudy::plotParam(%s)",_genModel->GetName())) ;
01240 pc.defineInt("nbins","Bins",0,0) ;
01241 pc.defineDouble("xlo","Range",0,0) ;
01242 pc.defineDouble("xhi","Range",1,0) ;
01243 pc.defineInt("dummy","FrameArgs",0,0) ;
01244 pc.defineMutex("Bins","FrameArgs") ;
01245 pc.defineMutex("Range","FrameArgs") ;
01246
01247
01248 pc.allowUndefined() ;
01249 pc.process(cmdList) ;
01250 if (!pc.ok(kTRUE)) {
01251 return 0 ;
01252 }
01253
01254
01255 Int_t nbins = pc.getInt("nbins") ;
01256 Double_t xlo = pc.getDouble("xlo") ;
01257 Double_t xhi = pc.getDouble("xhi") ;
01258 RooPlot* frame ;
01259
01260 if (pc.hasProcessed("FrameArgs")) {
01261
01262 RooCmdArg* frameArg = static_cast<RooCmdArg*>(cmdList.FindObject("FrameArgs")) ;
01263 frame = param.frame(frameArg->subArgs()) ;
01264 } else {
01265
01266 RooCmdArg bins = RooFit::Bins(nbins) ;
01267 RooCmdArg range = RooFit::Range(xlo,xhi) ;
01268 RooCmdArg autor = symRange ? RooFit::AutoSymRange(*_fitParData,0.2) : RooFit::AutoRange(*_fitParData,0.2) ;
01269 RooLinkedList frameCmdList ;
01270
01271 if (pc.hasProcessed("Bins")) frameCmdList.Add(&bins) ;
01272 if (pc.hasProcessed("Range")) {
01273 frameCmdList.Add(&range) ;
01274 } else {
01275 frameCmdList.Add(&autor) ;
01276 }
01277 frame = param.frame(frameCmdList) ;
01278 }
01279
01280
01281 pc.stripCmdList(cmdList,"FrameArgs,Bins,Range") ;
01282
01283 return frame ;
01284 }
01285
01286
01287
01288
01289 RooPlot* RooMCStudy::plotNLL(Double_t lo, Double_t hi, Int_t nBins)
01290 {
01291
01292
01293
01294 RooPlot* frame = _nllVar->frame(lo,hi,nBins) ;
01295
01296 _fitParData->plotOn(frame) ;
01297 return frame ;
01298 }
01299
01300
01301
01302
01303 RooPlot* RooMCStudy::plotError(const RooRealVar& param, Double_t lo, Double_t hi, Int_t nbins)
01304 {
01305
01306
01307
01308 if (_canAddFitResults) {
01309 calcPulls() ;
01310 _canAddFitResults=kFALSE ;
01311 }
01312
01313 RooErrorVar* evar = param.errorVar() ;
01314 RooPlot* frame = evar->frame(lo,hi,nbins) ;
01315 _fitParData->plotOn(frame) ;
01316
01317 delete evar ;
01318 return frame ;
01319 }
01320
01321
01322
01323
01324 RooPlot* RooMCStudy::plotPull(const RooRealVar& param, Double_t lo, Double_t hi, Int_t nbins, Bool_t fitGauss)
01325 {
01326
01327
01328
01329
01330
01331
01332 if (_canAddFitResults) {
01333 calcPulls() ;
01334 _canAddFitResults=kFALSE ;
01335 }
01336
01337
01338 TString name(param.GetName()), title(param.GetTitle()) ;
01339 name.Append("pull") ; title.Append(" Pull") ;
01340 RooRealVar pvar(name,title,lo,hi) ;
01341 pvar.setBins(nbins) ;
01342
01343 RooPlot* frame = pvar.frame() ;
01344 _fitParData->plotOn(frame) ;
01345
01346 if (fitGauss) {
01347 RooRealVar pullMean("pullMean","Mean of pull",0,lo,hi) ;
01348 RooRealVar pullSigma("pullSigma","Width of pull",1,0,5) ;
01349 RooGenericPdf pullGauss("pullGauss","Gaussian of pull",
01350 "exp(-0.5*(@0-@1)*(@0-@1)/(@2*@2))",
01351 RooArgSet(pvar,pullMean,pullSigma)) ;
01352 pullGauss.fitTo(*_fitParData,"mh") ;
01353 pullGauss.plotOn(frame) ;
01354 pullGauss.paramOn(frame,_fitParData) ;
01355 }
01356
01357 return frame ;
01358 }
01359
01360
01361