#include "RooFit.h"
#include "RooMsgService.h"
#include "RooConvGenContext.h"
#include "RooAbsAnaConvPdf.h"
#include "RooNumConvPdf.h"
#include "RooFFTConvPdf.h"
#include "RooProdPdf.h"
#include "RooDataSet.h"
#include "RooArgSet.h"
#include "RooTruthModel.h"
#include "Riostream.h"
ClassImp(RooConvGenContext)
;
RooConvGenContext::RooConvGenContext(const RooAbsAnaConvPdf &model, const RooArgSet &vars,
const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
RooAbsGenContext(model,vars,prototype,auxProto,verbose), _pdfVarsOwned(0), _modelVarsOwned(0)
{
cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for analytical convolution p.d.f. " << model.GetName()
<< " for generation of observable(s) " << vars ;
_pdfCloneSet = (RooArgSet*) RooArgSet(model).snapshot(kTRUE) ;
if (!_pdfCloneSet) {
coutE(Generation) << "RooConvGenContext::RooConvGenContext(" << GetName() << ") Couldn't deep-clone PDF, abort," << endl ;
RooErrorHandler::softAbort() ;
}
RooAbsAnaConvPdf* pdfClone = (RooAbsAnaConvPdf*) _pdfCloneSet->find(model.GetName()) ;
RooTruthModel truthModel("truthModel","Truth resolution model",(RooRealVar&)*pdfClone->convVar()) ;
pdfClone->changeModel(truthModel) ;
((RooRealVar*)pdfClone->convVar())->removeRange() ;
_pdfVars = (RooArgSet*) pdfClone->getObservables(&vars) ; ;
_pdfGen = pdfClone->genContext(*_pdfVars,prototype,auxProto,verbose) ;
_modelCloneSet = (RooArgSet*) RooArgSet(*model._convSet.at(0)).snapshot(kTRUE) ;
if (!_modelCloneSet) {
coutE(Generation) << "RooConvGenContext::RooConvGenContext(" << GetName() << ") Couldn't deep-clone resolution model, abort," << endl ;
RooErrorHandler::softAbort() ;
}
RooResolutionModel* modelClone = (RooResolutionModel*)
_modelCloneSet->find(model._convSet.at(0)->GetName())->Clone("smearing") ;
_modelCloneSet->addOwned(*modelClone) ;
modelClone->changeBasis(0) ;
modelClone->convVar().removeRange() ;
_modelVars = (RooArgSet*) modelClone->getObservables(&vars) ;
_modelVars->add(modelClone->convVar()) ;
_convVarName = modelClone->convVar().GetName() ;
_modelGen = modelClone->genContext(*_modelVars,prototype,auxProto,verbose) ;
if (prototype) {
_pdfVars->add(*prototype->get()) ;
_modelVars->add(*prototype->get()) ;
}
if (auxProto) {
_pdfVars->add(*auxProto) ;
_modelVars->add(*auxProto) ;
}
}
RooConvGenContext::RooConvGenContext(const RooNumConvPdf &model, const RooArgSet &vars,
const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
RooAbsGenContext(model,vars,prototype,auxProto,verbose)
{
cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for numeric convolution p.d.f. " << model.GetName()
<< " for generation of observable(s) " << vars << endl ;
_pdfVarsOwned = (RooArgSet*) model.conv().clonePdf().getObservables(&vars)->snapshot(kTRUE) ;
_pdfVars = new RooArgSet(*_pdfVarsOwned) ;
_pdfGen = ((RooAbsPdf&)model.conv().clonePdf()).genContext(*_pdfVars,prototype,auxProto,verbose) ;
_pdfCloneSet = 0 ;
_modelVarsOwned = (RooArgSet*) model.conv().cloneModel().getObservables(&vars)->snapshot(kTRUE) ;
_modelVars = new RooArgSet(*_modelVarsOwned) ;
_convVarName = model.conv().cloneVar().GetName() ;
_modelGen = ((RooAbsPdf&)model.conv().cloneModel()).genContext(*_modelVars,prototype,auxProto,verbose) ;
_modelCloneSet = 0 ;
if (prototype) {
_pdfVars->add(*prototype->get()) ;
_modelVars->add(*prototype->get()) ;
}
}
RooConvGenContext::RooConvGenContext(const RooFFTConvPdf &model, const RooArgSet &vars,
const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
RooAbsGenContext(model,vars,prototype,auxProto,verbose)
{
cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for fft convolution p.d.f. " << model.GetName()
<< " for generation of observable(s) " << vars << endl ;
_pdfVarsOwned = (RooArgSet*) model._pdf1.arg().getObservables(&vars)->snapshot(kTRUE) ;
_pdfVars = new RooArgSet(*_pdfVarsOwned) ;
_pdfGen = ((RooAbsPdf&)model._pdf1.arg()).genContext(*_pdfVars,prototype,auxProto,verbose) ;
_pdfCloneSet = 0 ;
_modelVarsOwned = (RooArgSet*) model._pdf2.arg().getObservables(&vars)->snapshot(kTRUE) ;
_modelVars = new RooArgSet(*_modelVarsOwned) ;
_convVarName = model._x.arg().GetName() ;
_modelGen = ((RooAbsPdf&)model._pdf2.arg()).genContext(*_modelVars,prototype,auxProto,verbose) ;
_modelCloneSet = 0 ;
if (prototype) {
_pdfVars->add(*prototype->get()) ;
_modelVars->add(*prototype->get()) ;
}
}
RooConvGenContext::~RooConvGenContext()
{
delete _pdfGen ;
delete _modelGen ;
delete _pdfCloneSet ;
delete _modelCloneSet ;
delete _modelVars ;
delete _pdfVars ;
delete _pdfVarsOwned ;
delete _modelVarsOwned ;
}
void RooConvGenContext::printToStream(ostream &os, PrintOption opt, TString indent) const
{
RooAbsGenContext::printToStream(os,opt,indent);
}
void RooConvGenContext::initGenerator(const RooArgSet &theEvent)
{
_cvModel = (RooRealVar*) _modelVars->find(_convVarName) ;
_cvPdf = (RooRealVar*) _pdfVars->find(_convVarName) ;
_cvOut = (RooRealVar*) theEvent.find(_convVarName) ;
RooArgSet* pdfCommon = (RooArgSet*) theEvent.selectCommon(*_pdfVars) ;
pdfCommon->remove(*_cvPdf,kTRUE,kTRUE) ;
_pdfVars->replace(*pdfCommon) ;
delete pdfCommon ;
RooArgSet* modelCommon = (RooArgSet*) theEvent.selectCommon(*_modelVars) ;
modelCommon->remove(*_cvModel,kTRUE,kTRUE) ;
_modelVars->replace(*modelCommon) ;
delete modelCommon ;
_pdfGen->initGenerator(*_pdfVars) ;
_modelGen->initGenerator(*_modelVars) ;
}
void RooConvGenContext::generateEvent(RooArgSet &theEvent, Int_t remaining)
{
while(1) {
_modelGen->generateEvent(*_modelVars,remaining) ;
_pdfGen->generateEvent(*_pdfVars,remaining) ;
Double_t convValSmeared = _cvPdf->getVal() + _cvModel->getVal() ;
if (_cvOut->isValidReal(convValSmeared)) {
theEvent = *_modelVars ;
theEvent = *_pdfVars ;
_cvOut->setVal(convValSmeared) ;
return ;
}
}
}
void RooConvGenContext::setProtoDataOrder(Int_t* lut)
{
RooAbsGenContext::setProtoDataOrder(lut) ;
_modelGen->setProtoDataOrder(lut) ;
_pdfGen->setProtoDataOrder(lut) ;
}
Last update: Thu Jan 17 08:44:20 2008
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.