// RooMCStudy is a help class to facilitate Monte Carlo studies
// such as 'goodness-of-fit' studies, that involve fitting a PDF
// to multiple toy Monte Carlo sets generated from the same PDF
// or another PDF.
// <p>
// Given a fit PDF and a generator PDF, RooMCStudy can produce
// large numbers of toyMC samples and/or fit these samples
// and acculumate the final parameters of each fit in a dataset.
// <p>
// Additional plotting routines simplify the task of plotting
// the distribution of the minimized likelihood, each parameters fitted value,
// fitted error and pull distribution.
// <p>
// Class RooMCStudy provides the option to insert add-in modules
// that modify the generate and fit cycle and allow to perform
// extra steps in the cycle. Output of these modules can be stored
// alongside the fit results in the aggregate results dataset.
// These study modules should derive from classs RooAbsMCStudyModel
//
// END_HTML
#include "RooFit.h"
#include "Riostream.h"
#include "RooMCStudy.h"
#include "RooAbsMCStudyModule.h"
#include "RooGenContext.h"
#include "RooAbsPdf.h"
#include "RooDataSet.h"
#include "RooDataHist.h"
#include "RooRealVar.h"
#include "RooFitResult.h"
#include "RooErrorVar.h"
#include "RooFormulaVar.h"
#include "RooArgList.h"
#include "RooPlot.h"
#include "RooGenericPdf.h"
#include "RooRandom.h"
#include "RooCmdConfig.h"
#include "RooGlobalFunc.h"
#include "RooPullVar.h"
#include "RooMsgService.h"
#include "RooProdPdf.h"
using namespace std ;
ClassImp(RooMCStudy)
;
RooMCStudy::RooMCStudy(const RooAbsPdf& model, const RooArgSet& observables,
RooCmdArg arg1, RooCmdArg arg2,
RooCmdArg arg3,RooCmdArg arg4,RooCmdArg arg5,
RooCmdArg arg6,RooCmdArg arg7,RooCmdArg arg8) : TNamed("mcstudy","mcstudy")
{
RooLinkedList cmdList;
cmdList.Add(const_cast<RooCmdArg*>(&arg1)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg2)) ;
cmdList.Add(const_cast<RooCmdArg*>(&arg3)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg4)) ;
cmdList.Add(const_cast<RooCmdArg*>(&arg5)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg6)) ;
cmdList.Add(const_cast<RooCmdArg*>(&arg7)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg8)) ;
RooCmdConfig pc(Form("RooMCStudy::RooMCStudy(%s)",model.GetName())) ;
pc.defineObject("fitModel","FitModel",0,0) ;
pc.defineObject("condObs","ProjectedDependents",0,0) ;
pc.defineObject("protoData","PrototypeData",0,0) ;
pc.defineSet("cPars","Constrain",0,0) ;
pc.defineSet("extCons","ExternalConstraints",0,0) ;
pc.defineInt("silence","Silence",0,0) ;
pc.defineInt("randProtoData","PrototypeData",0,0) ;
pc.defineInt("verboseGen","Verbose",0,0) ;
pc.defineInt("extendedGen","Extended",0,0) ;
pc.defineInt("binGenData","Binned",0,0) ;
pc.defineString("fitOpts","FitOptions",0,"") ;
pc.defineInt("dummy","FitOptArgs",0,0) ;
pc.defineMutex("FitOptions","FitOptArgs") ;
pc.defineMutex("Constrain","FitOptions") ;
pc.defineMutex("ExternalConstraints","FitOptions") ;
pc.process(cmdList) ;
if (!pc.ok(kTRUE)) {
return ;
}
if (pc.hasProcessed("FitOptArgs")) {
RooCmdArg* fitOptArg = static_cast<RooCmdArg*>(cmdList.FindObject("FitOptArgs")) ;
for (Int_t i=0 ; i<fitOptArg->subArgs().GetSize() ;i++) {
_fitOptList.Add(new RooCmdArg(static_cast<RooCmdArg&>(*fitOptArg->subArgs().At(i)))) ;
}
}
_silence = pc.getInt("silence") ;
_verboseGen = pc.getInt("verboseGen") ;
_extendedGen = pc.getInt("extendedGen") ;
_binGenData = pc.getInt("binGenData") ;
_randProto = pc.getInt("randProtoData") ;
const RooArgSet* cPars = pc.getSet("cPars") ;
const RooArgSet* extCons = pc.getSet("extCons") ;
if (cPars) {
_fitOptList.Add(RooFit::Constrain(*cPars).Clone()) ;
}
if (extCons) {
_fitOptList.Add(RooFit::ExternalConstraints(*extCons).Clone()) ;
}
RooArgSet allConstraints ;
RooArgSet consPars ;
if (cPars) {
RooArgSet* constraints = model.getConstraints(observables,*cPars,kTRUE) ;
allConstraints.add(*constraints) ;
delete constraints ;
}
if (allConstraints.getSize()>0) {
_constrPdf = new RooProdPdf("mcs_constr_prod","RooMCStudy constraints product",allConstraints) ;
if (cPars) {
consPars.add(*cPars) ;
} else {
RooArgSet* params = model.getParameters(observables) ;
RooArgSet* cparams = _constrPdf->getObservables(*params) ;
consPars.add(*cparams) ;
delete params ;
delete cparams ;
}
_constrGenContext = _constrPdf->genContext(consPars,0,0,_verboseGen) ;
_perExptGenParams = kTRUE ;
} else {
_constrPdf = 0 ;
_constrGenContext=0 ;
_perExptGenParams = kFALSE ;
}
_genModel = const_cast<RooAbsPdf*>(&model) ;
RooAbsPdf* fitModel = static_cast<RooAbsPdf*>(pc.getObject("fitModel",0)) ;
_fitModel = fitModel ? fitModel : _genModel ;
_genProtoData = static_cast<RooDataSet*>(pc.getObject("protoData",0)) ;
if (pc.getObject("condObs",0)) {
_projDeps.add(static_cast<RooArgSet&>(*pc.getObject("condObs",0))) ;
}
_dependents.add(observables) ;
_allDependents.add(_dependents) ;
_fitOptions = pc.getString("fitOpts") ;
_canAddFitResults = kTRUE ;
if (_extendedGen && _genProtoData && !_randProto) {
oocoutW(_fitModel,Generation) << "RooMCStudy::RooMCStudy: WARNING Using generator option 'e' (Poisson distribution of #events) together " << endl
<< " with a prototype dataset implies incomplete sampling or oversampling of proto data." << endl
<< " Use option \"r\" to randomize prototype dataset order and thus to randomize" << endl
<< " the set of over/undersampled prototype events for each generation cycle." << endl ;
}
_genParams = _genModel->getParameters(&_dependents) ;
if (!_binGenData) {
_genContext = _genModel->genContext(_dependents,_genProtoData,0,_verboseGen) ;
_genContext->attach(*_genParams) ;
} else {
_genContext = 0 ;
}
_genInitParams = (RooArgSet*) _genParams->snapshot(kFALSE) ;
_fitParams = _fitModel->getParameters(&_dependents) ;
_fitInitParams = (RooArgSet*) _fitParams->snapshot(kTRUE) ;
_nExpGen = _extendedGen ? _genModel->expectedEvents(&_dependents) : 0 ;
_nllVar = new RooRealVar("NLL","-log(Likelihood)",0) ;
_ngenVar = new RooRealVar("ngen","number of generated events",0) ;
RooArgSet tmp2(*_fitParams) ;
tmp2.add(*_nllVar) ;
tmp2.add(*_ngenVar) ;
tmp2.setAttribAll("StoreError",kTRUE) ;
tmp2.setAttribAll("StoreAsymError",kTRUE) ;
TString fpdName ;
if (_fitModel==_genModel) {
fpdName = Form("fitParData_%s",_fitModel->GetName()) ;
} else {
fpdName= Form("fitParData_%s_%s",_fitModel->GetName(),_genModel->GetName()) ;
}
_fitParData = new RooDataSet(fpdName.Data(),"Fit Parameters DataSet",tmp2) ;
tmp2.setAttribAll("StoreError",kFALSE) ;
tmp2.setAttribAll("StoreAsymError",kFALSE) ;
if (_perExptGenParams) {
_genParData = new RooDataSet("genParData","Generated Parameters dataset",*_genParams) ;
} else {
_genParData = 0 ;
}
if (_genProtoData) {
_allDependents.add(*_genProtoData->get(),kTRUE) ;
}
list<RooAbsMCStudyModule*>::iterator iter ;
for (iter=_modList.begin() ; iter!= _modList.end() ; ++iter) {
Bool_t ok = (*iter)->doInitializeInstance(*this) ;
if (!ok) {
oocoutE(_fitModel,Generation) << "RooMCStudy::ctor: removing study module " << (*iter)->GetName() << " from analysis chain because initialization failed" << endl ;
iter = _modList.erase(iter) ;
}
}
}
RooMCStudy::RooMCStudy(const RooAbsPdf& genModel, const RooAbsPdf& fitModel,
const RooArgSet& dependents, const char* genOptions,
const char* fitOptions, const RooDataSet* genProtoData,
const RooArgSet& projDeps) :
TNamed("mcstudy","mcstudy"),
_genModel((RooAbsPdf*)&genModel),
_genProtoData(genProtoData),
_projDeps(projDeps),
_constrPdf(0),
_constrGenContext(0),
_dependents(dependents),
_allDependents(dependents),
_fitModel((RooAbsPdf*)&fitModel),
_nllVar(0),
_ngenVar(0),
_genParData(0),
_fitOptions(fitOptions),
_canAddFitResults(kTRUE),
_perExptGenParams(0)
{
TString genOpt(genOptions) ;
genOpt.ToLower() ;
_verboseGen = genOpt.Contains("v") ;
_extendedGen = genOpt.Contains("e") ;
_binGenData = genOpt.Contains("b") ;
_randProto = genOpt.Contains("r") ;
if (_extendedGen && genProtoData && !_randProto) {
oocoutE(_fitModel,Generation) << "RooMCStudy::RooMCStudy: WARNING Using generator option 'e' (Poisson distribution of #events) together " << endl
<< " with a prototype dataset implies incomplete sampling or oversampling of proto data." << endl
<< " Use option \"r\" to randomize prototype dataset order and thus to randomize" << endl
<< " the set of over/undersampled prototype events for each generation cycle." << endl ;
}
if (!_binGenData) {
_genContext = genModel.genContext(dependents,genProtoData,0,_verboseGen) ;
} else {
_genContext = 0 ;
}
_genParams = _genModel->getParameters(&_dependents) ;
RooArgSet* tmp = genModel.getParameters(&dependents) ;
_genInitParams = (RooArgSet*) tmp->snapshot(kFALSE) ;
delete tmp ;
_fitParams = fitModel.getParameters(&dependents) ;
_fitInitParams = (RooArgSet*) _fitParams->snapshot(kTRUE) ;
_nExpGen = _extendedGen ? genModel.expectedEvents(&dependents) : 0 ;
_nllVar = new RooRealVar("NLL","-log(Likelihood)",0) ;
_ngenVar = new RooRealVar("ngen","number of generated events",0) ;
RooArgSet tmp2(*_fitParams) ;
tmp2.add(*_nllVar) ;
tmp2.add(*_ngenVar) ;
tmp2.setAttribAll("StoreError",kTRUE) ;
tmp2.setAttribAll("StoreAsymError",kTRUE) ;
_fitParData = new RooDataSet("fitParData","Fit Parameters DataSet",tmp2) ;
tmp2.setAttribAll("StoreError",kFALSE) ;
tmp2.setAttribAll("StoreAsymError",kFALSE) ;
if (genProtoData) {
_allDependents.add(*genProtoData->get(),kTRUE) ;
}
list<RooAbsMCStudyModule*>::iterator iter ;
for (iter=_modList.begin() ; iter!= _modList.end() ; ++iter) {
Bool_t ok = (*iter)->doInitializeInstance(*this) ;
if (!ok) {
oocoutE(_fitModel,Generation) << "RooMCStudy::ctor: removing study module " << (*iter)->GetName() << " from analysis chain because initialization failed" << endl ;
iter = _modList.erase(iter) ;
}
}
}
RooMCStudy::~RooMCStudy()
{
_genDataList.Delete() ;
_fitResList.Delete() ;
_fitOptList.Delete() ;
delete _ngenVar ;
delete _fitParData ;
delete _genParData ;
delete _fitInitParams ;
delete _fitParams ;
delete _genInitParams ;
delete _genParams ;
delete _genContext ;
delete _nllVar ;
delete _constrPdf ;
delete _constrGenContext ;
}
void RooMCStudy::addModule(RooAbsMCStudyModule& module)
{
module.doInitializeInstance(*this) ;
_modList.push_back(&module) ;
}
Bool_t RooMCStudy::run(Bool_t doGenerate, Bool_t DoFit, Int_t nSamples, Int_t nEvtPerSample, Bool_t keepGenData, const char* asciiFilePat)
{
RooFit::MsgLevel oldLevel(RooFit::FATAL) ;
if (_silence) {
oldLevel = RooMsgService::instance().globalKillBelow() ;
RooMsgService::instance().setGlobalKillBelow(RooFit::PROGRESS) ;
}
list<RooAbsMCStudyModule*>::iterator iter ;
for (iter=_modList.begin() ; iter!= _modList.end() ; ++iter) {
(*iter)->initializeRun(nSamples) ;
}
Int_t prescale = nSamples>100 ? Int_t(nSamples/100) : 1 ;
while(nSamples--) {
if (nSamples%prescale==0) {
oocoutP(_fitModel,Generation) << "RooMCStudy::run: " ;
if (doGenerate) ooccoutI(_fitModel,Generation) << "Generating " ;
if (doGenerate && DoFit) ooccoutI(_fitModel,Generation) << "and " ;
if (DoFit) ooccoutI(_fitModel,Generation) << "fitting " ;
ooccoutP(_fitModel,Generation) << "sample " << nSamples << endl ;
}
_genSample = 0;
Bool_t existingData = kFALSE ;
if (doGenerate) {
Int_t nEvt(nEvtPerSample) ;
*_genParams = *_genInitParams ;
if (_constrPdf) {
RooDataSet* tmp = _constrGenContext->generate(1) ;
*_genParams = *tmp->get() ;
delete tmp ;
}
if (_genParData) {
_genParData->add(*_genParams) ;
}
list<RooAbsMCStudyModule*>::iterator iter2 ;
for (iter2=_modList.begin() ; iter2!= _modList.end() ; ++iter2) {
(*iter2)->processBeforeGen(nSamples) ;
}
if (_binGenData) {
if (_extendedGen) {
_nExpGen = _genModel->expectedEvents(&_dependents) ;
nEvt = RooRandom::randomGenerator()->Poisson(nEvtPerSample==0?_nExpGen:nEvtPerSample) ;
}
_genSample = _genModel->generateBinned(_dependents,nEvt) ;
} else {
if (_extendedGen) {
_nExpGen = _genModel->expectedEvents(&_dependents) ;
nEvt = RooRandom::randomGenerator()->Poisson(nEvtPerSample==0?_nExpGen:nEvtPerSample) ;
}
if (_randProto && _genProtoData && _genProtoData->numEntries()!=nEvt) {
oocoutI(_fitModel,Generation) << "RooMCStudy: (Re)randomizing event order in prototype dataset (Nevt=" << nEvt << ")" << endl ;
Int_t* newOrder = _genModel->randomizeProtoOrder(_genProtoData->numEntries(),nEvt) ;
_genContext->setProtoDataOrder(newOrder) ;
delete[] newOrder ;
}
if (nEvt>0) {
_genSample = _genContext->generate(nEvt) ;
} else {
_genSample = new RooDataSet("emptySample","emptySample",_dependents) ;
}
}
} else if (asciiFilePat) {
char asciiFile[1024] ;
sprintf(asciiFile,asciiFilePat,nSamples) ;
RooArgList depList(_allDependents) ;
_genSample = RooDataSet::read(asciiFile,depList,"q") ;
} else {
_genSample = (RooDataSet*) _genDataList.At(nSamples) ;
existingData = kTRUE ;
if (!_genSample) {
oocoutW(_fitModel,Generation) << "RooMCStudy::run: WARNING: Sample #" << nSamples << " not loaded, skipping" << endl ;
continue ;
}
}
_ngenVar->setVal(_genSample->sumEntries()) ;
list<RooAbsMCStudyModule*>::iterator iter3 ;
for (iter3=_modList.begin() ; iter3!= _modList.end() ; ++iter3) {
(*iter3)->processBetweenGenAndFit(nSamples) ;
}
if (DoFit) fitSample(_genSample) ;
for (iter3=_modList.begin() ; iter3!= _modList.end() ; ++iter3) {
(*iter3)->processAfterFit(nSamples) ;
}
if (doGenerate && asciiFilePat && *asciiFilePat) {
char asciiFile[1024] ;
sprintf(asciiFile,asciiFilePat,nSamples) ;
RooDataSet* unbinnedData = dynamic_cast<RooDataSet*>(_genSample) ;
if (unbinnedData) {
unbinnedData->write(asciiFile) ;
} else {
coutE(InputArguments) << "RooMCStudy::run(" << GetName() << ") ERROR: ASCII writing of binned datasets is not supported" << endl ;
}
}
if (!existingData) {
if (keepGenData) {
_genDataList.Add(_genSample) ;
} else {
delete _genSample ;
}
}
}
for (iter=_modList.begin() ; iter!= _modList.end() ; ++iter) {
RooDataSet* auxData = (*iter)->finalizeRun() ;
if (auxData) {
_fitParData->merge(auxData) ;
}
}
_canAddFitResults = kFALSE ;
if (_genParData) {
const RooArgSet* genPars = _genParData->get() ;
TIterator* iter2 = genPars->createIterator() ;
RooAbsArg* arg ;
while((arg=(RooAbsArg*)iter2->Next())) {
_genParData->changeObservableName(arg->GetName(),Form("%s_gen",arg->GetName())) ;
}
delete iter2 ;
_fitParData->merge(_genParData) ;
}
if (DoFit) calcPulls() ;
if (_silence) {
RooMsgService::instance().setGlobalKillBelow(oldLevel) ;
}
return kFALSE ;
}
Bool_t RooMCStudy::generateAndFit(Int_t nSamples, Int_t nEvtPerSample, Bool_t keepGenData, const char* asciiFilePat)
{
_fitResList.Delete() ;
_genDataList.Delete() ;
_fitParData->reset() ;
return run(kTRUE,kTRUE,nSamples,nEvtPerSample,keepGenData,asciiFilePat) ;
}
Bool_t RooMCStudy::generate(Int_t nSamples, Int_t nEvtPerSample, Bool_t keepGenData, const char* asciiFilePat)
{
_genDataList.Delete() ;
return run(kTRUE,kFALSE,nSamples,nEvtPerSample,keepGenData,asciiFilePat) ;
}
Bool_t RooMCStudy::fit(Int_t nSamples, const char* asciiFilePat)
{
_fitResList.Delete() ;
_fitParData->reset() ;
return run(kFALSE,kTRUE,nSamples,0,kFALSE,asciiFilePat) ;
}
Bool_t RooMCStudy::fit(Int_t nSamples, TList& dataSetList)
{
_fitResList.Delete() ;
_genDataList.Delete() ;
_fitParData->reset() ;
TIterator* iter = dataSetList.MakeIterator() ;
RooAbsData* gset ;
while((gset=(RooAbsData*)iter->Next())) {
_genDataList.Add(gset) ;
}
delete iter ;
return run(kFALSE,kTRUE,nSamples,0,kTRUE,0) ;
}
void RooMCStudy::resetFitParams()
{
*_fitParams = *_fitInitParams ;
}
RooFitResult* RooMCStudy::doFit(RooAbsData* genSample)
{
TString fitOpt2(_fitOptions) ; fitOpt2.Append("r") ;
if (_silence) {
fitOpt2.Append("b") ;
}
RooAbsData* data ;
if (_binGenData) {
RooArgSet* depList = _fitModel->getObservables(genSample) ;
data = new RooDataHist(genSample->GetName(),genSample->GetTitle(),*depList,*genSample) ;
delete depList ;
} else {
data = genSample ;
}
RooFitResult* fr ;
if (_fitOptList.GetSize()==0) {
if (_projDeps.getSize()>0) {
fr = (RooFitResult*) _fitModel->fitTo(*data,RooFit::ConditionalObservables(_projDeps),RooFit::FitOptions(fitOpt2)) ;
} else {
fr = (RooFitResult*) _fitModel->fitTo(*data,RooFit::FitOptions(fitOpt2)) ;
}
} else {
RooCmdArg save = RooFit::Save() ;
RooCmdArg condo = RooFit::ConditionalObservables(_projDeps) ;
RooCmdArg plevel = RooFit::PrintLevel(-1) ;
RooLinkedList fitOptList(_fitOptList) ;
fitOptList.Add(&save) ;
if (_projDeps.getSize()>0) {
fitOptList.Add(&condo) ;
}
if (_silence) {
fitOptList.Add(&plevel) ;
}
fr = (RooFitResult*) _fitModel->fitTo(*data,fitOptList) ;
}
if (_binGenData) delete data ;
return fr ;
}
RooFitResult* RooMCStudy::refit(RooAbsData* genSample)
{
if (!genSample) {
genSample = _genSample ;
}
RooFitResult* fr(0) ;
if (genSample->sumEntries()>0) {
fr = doFit(genSample) ;
}
return fr ;
}
Bool_t RooMCStudy::fitSample(RooAbsData* genSample)
{
resetFitParams() ;
Bool_t ok ;
RooFitResult* fr(0) ;
if (genSample->sumEntries()>0) {
fr = doFit(genSample) ;
ok = (fr->status()==0) ;
} else {
ok = kFALSE ;
}
if (ok) {
_nllVar->setVal(fr->minNll()) ;
RooArgSet tmp(*_fitParams) ;
tmp.add(*_nllVar) ;
tmp.add(*_ngenVar) ;
_fitParData->add(tmp) ;
}
Bool_t userSaveRequest = kFALSE ;
if (_fitOptList.GetSize()>0) {
if (_fitOptList.FindObject("Save")) userSaveRequest = kTRUE ;
} else {
if (_fitOptions.Contains("r")) userSaveRequest = kTRUE ;
}
if (userSaveRequest) {
_fitResList.Add(fr) ;
} else {
delete fr ;
}
return !ok ;
}
Bool_t RooMCStudy::addFitResult(const RooFitResult& fr)
{
if (!_canAddFitResults) {
oocoutE(_fitModel,InputArguments) << "RooMCStudy::addFitResult: ERROR cannot add fit results in current state" << endl ;
return kTRUE ;
}
*_fitParams = RooArgSet(fr.floatParsFinal()) ;
Bool_t ok = (fr.status()==0) ;
if (ok) {
_nllVar->setVal(fr.minNll()) ;
RooArgSet tmp(*_fitParams) ;
tmp.add(*_nllVar) ;
tmp.add(*_ngenVar) ;
_fitParData->add(tmp) ;
}
if (_fitOptions.Contains("r")) {
_fitResList.Add((TObject*)&fr) ;
}
return kFALSE ;
}
void RooMCStudy::calcPulls()
{
TIterator* iter = _fitParams->createIterator() ;
RooRealVar* par ;
while((par=(RooRealVar*)iter->Next())) {
RooErrorVar* err = par->errorVar() ;
_fitParData->addColumn(*err) ;
delete err ;
TString name(par->GetName()), title(par->GetTitle()) ;
name.Append("pull") ;
title.Append(" Pull") ;
RooAbsReal* genParOrig = (RooAbsReal*) _fitParData->get()->find(Form("%s_gen",par->GetName())) ;
if (genParOrig && _perExptGenParams) {
RooPullVar pull(name,title,*par,*genParOrig) ;
_fitParData->addColumn(pull,kFALSE) ;
} else {
genParOrig = (RooAbsReal*)_genInitParams->find(par->GetName()) ;
if (genParOrig) {
RooAbsReal* genPar = (RooAbsReal*) genParOrig->Clone("truth") ;
RooPullVar pull(name,title,*par,*genPar) ;
_fitParData->addColumn(pull,kFALSE) ;
delete genPar ;
}
}
}
delete iter ;
}
const RooDataSet& RooMCStudy::fitParDataSet()
{
if (_canAddFitResults) {
calcPulls() ;
_canAddFitResults = kFALSE ;
}
return *_fitParData ;
}
const RooArgSet* RooMCStudy::fitParams(Int_t sampleNum) const
{
if (sampleNum<0 || sampleNum>=_fitParData->numEntries()) {
oocoutE(_fitModel,InputArguments) << "RooMCStudy::fitParams: ERROR, invalid sample number: " << sampleNum << endl ;
return 0 ;
}
return _fitParData->get(sampleNum) ;
}
const RooFitResult* RooMCStudy::fitResult(Int_t sampleNum) const
{
if (sampleNum<0 || sampleNum>=_fitResList.GetSize()) {
oocoutE(_fitModel,InputArguments) << "RooMCStudy::fitResult: ERROR, invalid sample number: " << sampleNum << endl ;
return 0 ;
}
const RooFitResult* fr = (RooFitResult*) _fitResList.At(sampleNum) ;
if (fr) {
return fr ;
} else {
oocoutE(_fitModel,InputArguments) << "RooMCStudy::fitResult: ERROR, no fit result saved for sample "
<< sampleNum << ", did you use the 'r; fit option?" << endl ;
}
return 0 ;
}
const RooDataSet* RooMCStudy::genData(Int_t sampleNum) const
{
if (_genDataList.GetSize()==0) {
oocoutE(_fitModel,InputArguments) << "RooMCStudy::genData() ERROR, generated data was not saved" << endl ;
return 0 ;
}
if (sampleNum<0 || sampleNum>=_genDataList.GetSize()) {
oocoutE(_fitModel,InputArguments) << "RooMCStudy::genData() ERROR, invalid sample number: " << sampleNum << endl ;
return 0 ;
}
return (RooDataSet*) _genDataList.At(sampleNum) ;
}
RooPlot* RooMCStudy::plotParamOn(RooPlot* frame, const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4,
const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
{
_fitParData->plotOn(frame,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8) ;
return frame ;
}
RooPlot* RooMCStudy::plotParam(const char* paramName, const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4,
const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
{
RooRealVar* param = static_cast<RooRealVar*>(_fitParData->get()->find(paramName)) ;
if (!param) {
oocoutE(_fitModel,InputArguments) << "RooMCStudy::plotParam: ERROR: no parameter defined with name " << paramName << endl ;
return 0 ;
}
return plotParam(*param,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8) ;
}
RooPlot* RooMCStudy::plotParam(const RooRealVar& param, const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4,
const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
{
RooLinkedList cmdList;
cmdList.Add(const_cast<RooCmdArg*>(&arg1)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg2)) ;
cmdList.Add(const_cast<RooCmdArg*>(&arg3)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg4)) ;
cmdList.Add(const_cast<RooCmdArg*>(&arg5)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg6)) ;
cmdList.Add(const_cast<RooCmdArg*>(&arg7)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg8)) ;
RooPlot* frame = makeFrameAndPlotCmd(param, cmdList) ;
if (frame) {
_fitParData->plotOn(frame, cmdList) ;
}
return frame ;
}
RooPlot* RooMCStudy::plotNLL(const RooCmdArg& arg1, const RooCmdArg& arg2,
const RooCmdArg& arg3, const RooCmdArg& arg4,
const RooCmdArg& arg5, const RooCmdArg& arg6,
const RooCmdArg& arg7, const RooCmdArg& arg8)
{
return plotParam(*_nllVar,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8) ;
}
RooPlot* RooMCStudy::plotError(const RooRealVar& param, const RooCmdArg& arg1, const RooCmdArg& arg2,
const RooCmdArg& arg3, const RooCmdArg& arg4,
const RooCmdArg& arg5, const RooCmdArg& arg6,
const RooCmdArg& arg7, const RooCmdArg& arg8)
{
if (_canAddFitResults) {
calcPulls() ;
_canAddFitResults=kFALSE ;
}
RooErrorVar* evar = param.errorVar() ;
RooRealVar* evar_rrv = static_cast<RooRealVar*>(evar->createFundamental()) ;
RooPlot* frame = plotParam(*evar_rrv,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8) ;
delete evar_rrv ;
delete evar ;
return frame ;
}
RooPlot* RooMCStudy::plotPull(const RooRealVar& param, const RooCmdArg& arg1, const RooCmdArg& arg2,
const RooCmdArg& arg3, const RooCmdArg& arg4,
const RooCmdArg& arg5, const RooCmdArg& arg6,
const RooCmdArg& arg7, const RooCmdArg& arg8)
{
RooLinkedList cmdList;
cmdList.Add(const_cast<RooCmdArg*>(&arg1)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg2)) ;
cmdList.Add(const_cast<RooCmdArg*>(&arg3)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg4)) ;
cmdList.Add(const_cast<RooCmdArg*>(&arg5)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg6)) ;
cmdList.Add(const_cast<RooCmdArg*>(&arg7)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg8)) ;
TString name(param.GetName()), title(param.GetTitle()) ;
name.Append("pull") ; title.Append(" Pull") ;
RooRealVar pvar(name,title,-100,100) ;
pvar.setBins(100) ;
RooPlot* frame = makeFrameAndPlotCmd(pvar, cmdList, kTRUE) ;
if (frame) {
RooCmdConfig pc(Form("RooMCStudy::plotPull(%s)",_genModel->GetName())) ;
pc.defineInt("fitGauss","FitGauss",0,0) ;
pc.allowUndefined() ;
pc.process(cmdList) ;
Bool_t fitGauss=pc.getInt("fitGauss") ;
pc.stripCmdList(cmdList,"FitGauss") ;
_fitParData->plotOn(frame,cmdList) ;
if (fitGauss) {
RooRealVar pullMean("pullMean","Mean of pull",0,-100,100) ;
RooRealVar pullSigma("pullSigma","Width of pull",1,0.1,5) ;
RooGenericPdf pullGauss("pullGauss","Gaussian of pull",
"exp(-0.5*(@0-@1)*(@0-@1)/(@2*@2))",
RooArgSet(pvar,pullMean,pullSigma)) ;
pullGauss.fitTo(*_fitParData,RooFit::Minos(0),RooFit::PrintLevel(-1)) ;
pullGauss.plotOn(frame) ;
pullGauss.paramOn(frame,_fitParData) ;
}
}
return frame ; ;
}
RooPlot* RooMCStudy::makeFrameAndPlotCmd(const RooRealVar& param, RooLinkedList& cmdList, Bool_t symRange) const
{
RooCmdConfig pc(Form("RooMCStudy::plotParam(%s)",_genModel->GetName())) ;
pc.defineInt("nbins","Bins",0,0) ;
pc.defineDouble("xlo","Range",0,0) ;
pc.defineDouble("xhi","Range",1,0) ;
pc.defineInt("dummy","FrameArgs",0,0) ;
pc.defineMutex("Bins","FrameArgs") ;
pc.defineMutex("Range","FrameArgs") ;
pc.allowUndefined() ;
pc.process(cmdList) ;
if (!pc.ok(kTRUE)) {
return 0 ;
}
Int_t nbins = pc.getInt("nbins") ;
Double_t xlo = pc.getDouble("xlo") ;
Double_t xhi = pc.getDouble("xhi") ;
RooPlot* frame ;
if (pc.hasProcessed("FrameArgs")) {
RooCmdArg* frameArg = static_cast<RooCmdArg*>(cmdList.FindObject("FrameArgs")) ;
frame = param.frame(frameArg->subArgs()) ;
} else {
RooCmdArg bins = RooFit::Bins(nbins) ;
RooCmdArg range = RooFit::Range(xlo,xhi) ;
RooCmdArg autor = symRange ? RooFit::AutoSymRange(*_fitParData,0.2) : RooFit::AutoRange(*_fitParData,0.2) ;
RooLinkedList frameCmdList ;
if (pc.hasProcessed("Bins")) frameCmdList.Add(&bins) ;
if (pc.hasProcessed("Range")) {
frameCmdList.Add(&range) ;
} else {
frameCmdList.Add(&autor) ;
}
frame = param.frame(frameCmdList) ;
}
pc.stripCmdList(cmdList,"FrameArgs,Bins,Range") ;
return frame ;
}
RooPlot* RooMCStudy::plotNLL(Double_t lo, Double_t hi, Int_t nBins)
{
RooPlot* frame = _nllVar->frame(lo,hi,nBins) ;
_fitParData->plotOn(frame) ;
return frame ;
}
RooPlot* RooMCStudy::plotError(const RooRealVar& param, Double_t lo, Double_t hi, Int_t nbins)
{
if (_canAddFitResults) {
calcPulls() ;
_canAddFitResults=kFALSE ;
}
RooErrorVar* evar = param.errorVar() ;
RooPlot* frame = evar->frame(lo,hi,nbins) ;
_fitParData->plotOn(frame) ;
delete evar ;
return frame ;
}
RooPlot* RooMCStudy::plotPull(const RooRealVar& param, Double_t lo, Double_t hi, Int_t nbins, Bool_t fitGauss)
{
if (_canAddFitResults) {
calcPulls() ;
_canAddFitResults=kFALSE ;
}
TString name(param.GetName()), title(param.GetTitle()) ;
name.Append("pull") ; title.Append(" Pull") ;
RooRealVar pvar(name,title,lo,hi) ;
pvar.setBins(nbins) ;
RooPlot* frame = pvar.frame() ;
_fitParData->plotOn(frame) ;
if (fitGauss) {
RooRealVar pullMean("pullMean","Mean of pull",0,lo,hi) ;
RooRealVar pullSigma("pullSigma","Width of pull",1,0,5) ;
RooGenericPdf pullGauss("pullGauss","Gaussian of pull",
"exp(-0.5*(@0-@1)*(@0-@1)/(@2*@2))",
RooArgSet(pvar,pullMean,pullSigma)) ;
pullGauss.fitTo(*_fitParData,"mh") ;
pullGauss.plotOn(frame) ;
pullGauss.paramOn(frame,_fitParData) ;
}
return frame ;
}