// RooProdGenContext is an efficient implementation of the generator context
// specific for RooProdPdf PDFs. The sim-context owns a list of
// component generator contexts that are used to generate the dependents
// for each component PDF sequentially.
// END_HTML
#include "RooFit.h"
#include "Riostream.h"
#include "RooMsgService.h"
#include "RooProdGenContext.h"
#include "RooProdGenContext.h"
#include "RooProdPdf.h"
#include "RooDataSet.h"
#include "RooRealVar.h"
#include "RooGlobalFunc.h"
ClassImp(RooProdGenContext)
;
RooProdGenContext::RooProdGenContext(const RooProdPdf &model, const RooArgSet &vars,
const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
RooAbsGenContext(model,vars,prototype,auxProto,verbose), _pdf(&model)
{
cxcoutI(Generation) << "RooProdGenContext::ctor() setting up event special generator context for product p.d.f. " << model.GetName()
<< " for generation of observable(s) " << vars ;
if (prototype) ccxcoutI(Generation) << " with prototype data for " << *prototype->get() ;
if (auxProto && auxProto->getSize()>0) ccxcoutI(Generation) << " with auxiliary prototypes " << *auxProto ;
ccxcoutI(Generation) << endl ;
RooArgSet deps(vars) ;
if (prototype) {
RooArgSet* protoDeps = model.getObservables(*prototype->get()) ;
deps.remove(*protoDeps,kTRUE,kTRUE) ;
delete protoDeps ;
}
RooLinkedList termList,depsList,impDepList,crossDepList,intList ;
model.factorizeProduct(deps,RooArgSet(),termList,depsList,impDepList,crossDepList,intList) ;
TIterator* termIter = termList.MakeIterator() ;
TIterator* normIter = depsList.MakeIterator() ;
TIterator* impIter = impDepList.MakeIterator() ;
if (dologD(Generation)) {
cxcoutD(Generation) << "RooProdGenContext::ctor() factorizing product expression in irriducible terms " ;
while(RooArgSet* t=(RooArgSet*)termIter->Next()) {
ccxcoutD(Generation) << *t ;
}
ccxcoutD(Generation) << endl ;
}
RooArgSet genDeps ;
Bool_t anyAction = kTRUE ;
Bool_t go=kTRUE ;
while(go) {
RooAbsPdf* pdf ;
RooArgSet* term ;
RooArgSet* impDeps ;
RooArgSet* termDeps ;
termIter->Reset() ;
impIter->Reset() ;
normIter->Reset() ;
Bool_t anyPrevAction=anyAction ;
anyAction=kFALSE ;
if (termList.GetSize()==0) {
break ;
}
while((term=(RooArgSet*)termIter->Next())) {
impDeps = (RooArgSet*)impIter->Next() ;
termDeps = (RooArgSet*)normIter->Next() ;
cxcoutD(Generation) << "RooProdGenContext::ctor() analyzing product term " << *term << " with observable(s) " << *termDeps ;
if (impDeps->getSize()>0) {
ccxcoutD(Generation) << " which has dependence of external observable(s) " << *impDeps << " that to be generated first by other terms" ;
}
ccxcoutD(Generation) << endl ;
RooArgSet neededDeps(*impDeps) ;
neededDeps.remove(genDeps,kTRUE,kTRUE) ;
if (neededDeps.getSize()>0) {
if (!anyPrevAction) {
cxcoutD(Generation) << "RooProdGenContext::ctor() no convergence in single term analysis loop, terminating loop and process remainder of terms as single unit " << endl ;
go=kFALSE ;
break ;
}
cxcoutD(Generation) << "RooProdGenContext::ctor() skipping this term for now because it needs imported dependents that are not generated yet" << endl ;
continue ;
}
if (termDeps->getSize()==0) {
cxcoutD(Generation) << "RooProdGenContext::ctor() term has no observables requested to be generated, removing it" << endl ;
termList.Remove(term) ;
depsList.Remove(termDeps) ;
impDepList.Remove(impDeps) ;
delete term ;
delete termDeps ;
delete impDeps ;
anyAction=kTRUE ;
continue ;
}
TIterator* pdfIter = term->createIterator() ;
if (term->getSize()==1) {
pdf = (RooAbsPdf*) pdfIter->Next() ;
RooArgSet* pdfDep = pdf->getObservables(termDeps) ;
if (pdfDep->getSize()>0) {
coutI(Generation) << "RooProdGenContext::ctor() creating subcontext for generation of observables " << *pdfDep << " from model " << pdf->GetName() << endl ;
RooArgSet* auxProto2 = impDeps ? pdf->getObservables(impDeps) : 0 ;
RooAbsGenContext* cx = pdf->genContext(*pdfDep,prototype,auxProto2,verbose) ;
delete auxProto2 ;
_gcList.Add(cx) ;
}
genDeps.add(*pdfDep) ;
delete pdfDep ;
} else {
if (termDeps->getSize()>0) {
const char* name = model.makeRGPPName("PRODGEN_",*term,RooArgSet(),RooArgSet(),0) ;
RooLinkedList cmdList ;
RooLinkedList pdfSetList ;
pdfIter->Reset() ;
RooArgSet fullPdfSet ;
while((pdf=(RooAbsPdf*)pdfIter->Next())) {
RooArgSet* pdfnset = model.findPdfNSet(*pdf) ;
RooArgSet* pdfSet = new RooArgSet(*pdf) ;
pdfSetList.Add(pdfSet) ;
if (pdfnset && pdfnset->getSize()>0) {
cmdList.Add(RooFit::Conditional(*pdfSet,*pdfnset).Clone()) ;
} else {
fullPdfSet.add(*pdfSet) ;
}
}
RooProdPdf* multiPdf = new RooProdPdf(name,name,fullPdfSet,cmdList) ;
cmdList.Delete() ;
pdfSetList.Delete() ;
multiPdf->useDefaultGen(kTRUE) ;
_ownedMultiProds.addOwned(*multiPdf) ;
coutI(Generation) << "RooProdGenContext()::ctor creating subcontext for generation of observables " << *termDeps
<< "for irriducuble composite term using sub-product object " << multiPdf->GetName() ;
RooAbsGenContext* cx = multiPdf->genContext(*termDeps,prototype,auxProto,verbose) ;
_gcList.Add(cx) ;
genDeps.add(*termDeps) ;
}
}
delete pdfIter ;
termList.Remove(term) ;
depsList.Remove(termDeps) ;
impDepList.Remove(impDeps) ;
delete term ;
delete termDeps ;
delete impDeps ;
anyAction=kTRUE ;
}
}
if (termList.GetSize()>0) {
cxcoutD(Generation) << "RooProdGenContext::ctor() there are left-over terms that need to be generated separately" << endl ;
RooAbsPdf* pdf ;
RooArgSet* term ;
termIter->Reset() ;
normIter->Reset() ;
RooArgSet trailerTerm ;
RooArgSet trailerTermDeps ;
while((term=(RooArgSet*)termIter->Next())) {
RooArgSet* termDeps = (RooArgSet*)normIter->Next() ;
trailerTerm.add(*term) ;
trailerTermDeps.add(*termDeps) ;
}
const char* name = model.makeRGPPName("PRODGEN_",trailerTerm,RooArgSet(),RooArgSet(),0) ;
RooLinkedList cmdList ;
RooLinkedList pdfSetList ;
RooArgSet fullPdfSet ;
TIterator* pdfIter = trailerTerm.createIterator() ;
while((pdf=(RooAbsPdf*)pdfIter->Next())) {
RooArgSet* pdfnset = model.findPdfNSet(*pdf) ;
RooArgSet* pdfSet = new RooArgSet(*pdf) ;
pdfSetList.Add(pdfSet) ;
if (pdfnset && pdfnset->getSize()>0) {
cmdList.Add(RooFit::Conditional(*pdfSet,*pdfnset).Clone()) ;
} else {
fullPdfSet.add(*pdfSet) ;
}
}
RooProdPdf* multiPdf = new RooProdPdf(name,name,fullPdfSet,cmdList) ;
cmdList.Delete() ;
pdfSetList.Delete() ;
multiPdf->useDefaultGen(kTRUE) ;
_ownedMultiProds.addOwned(*multiPdf) ;
cxcoutD(Generation) << "RooProdGenContext(" << model.GetName() << "): creating context for irreducible composite trailer term "
<< multiPdf->GetName() << " that generates observables " << trailerTermDeps << endl ;
RooAbsGenContext* cx = multiPdf->genContext(trailerTermDeps,prototype,auxProto,verbose) ;
_gcList.Add(cx) ;
}
delete termIter ;
delete impIter ;
delete normIter ;
_gcIter = _gcList.MakeIterator() ;
termList.Delete() ;
depsList.Delete() ;
impDepList.Delete() ;
crossDepList.Delete() ;
intList.Delete() ;
}
RooProdGenContext::~RooProdGenContext()
{
delete _gcIter ;
_gcList.Delete() ;
}
void RooProdGenContext::attach(const RooArgSet& args)
{
RooAbsGenContext* gc ;
_gcIter->Reset() ;
while((gc=(RooAbsGenContext*)_gcIter->Next())){
gc->attach(args) ;
}
}
void RooProdGenContext::initGenerator(const RooArgSet &theEvent)
{
RooAbsGenContext* gc ;
_gcIter->Reset() ;
while((gc=(RooAbsGenContext*)_gcIter->Next())){
gc->initGenerator(theEvent) ;
}
}
void RooProdGenContext::generateEvent(RooArgSet &theEvent, Int_t remaining)
{
TList compData ;
RooAbsGenContext* gc ;
_gcIter->Reset() ;
while((gc=(RooAbsGenContext*)_gcIter->Next())) {
gc->generateEvent(theEvent,remaining) ;
}
}
void RooProdGenContext::setProtoDataOrder(Int_t* lut)
{
RooAbsGenContext::setProtoDataOrder(lut) ;
_gcIter->Reset() ;
RooAbsGenContext* gc ;
while((gc=(RooAbsGenContext*)_gcIter->Next())) {
gc->setProtoDataOrder(lut) ;
}
}
void RooProdGenContext::printMultiline(ostream &os, Int_t content, Bool_t verbose, TString indent) const
{
RooAbsGenContext::printMultiline(os,content,verbose,indent) ;
os << indent << "--- RooProdGenContext ---" << endl ;
os << indent << "Using PDF ";
_pdf->printStream(os,kName|kArgs|kClassName,kSingleLine,indent);
os << indent << "List of component generators" << endl ;
TString indent2(indent) ;
indent2.Append(" ") ;
RooAbsGenContext* gc ;
_gcIter->Reset() ;
while((gc=(RooAbsGenContext*)_gcIter->Next())) {
gc->printMultiline(os,content,verbose,indent2) ;
}
}
Last change: Wed Jun 25 08:33:49 2008
Last generated: 2008-06-25 08:33
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.