55 ClassImp(RooGenContext);
70 RooGenContext::RooGenContext(
const RooAbsPdf &model,
const RooArgSet &vars,
71 const RooDataSet *prototype,
const RooArgSet* auxProto,
72 Bool_t verbose,
const RooArgSet* forceDirect) :
73 RooAbsGenContext(model,vars,prototype,auxProto,verbose),
74 _cloneSet(0), _pdfClone(0), _acceptRejectFunc(0), _generator(0),
75 _maxVar(0), _uniIter(0), _updateFMaxPerEvent(0)
77 cxcoutI(Generation) <<
"RooGenContext::ctor() setting up event generator context for p.d.f. " << model.GetName()
78 <<
" for generation of observable(s) " << vars ;
79 if (prototype) ccxcoutI(Generation) <<
" with prototype data for " << *prototype->get() ;
80 if (auxProto && auxProto->getSize()>0) ccxcoutI(Generation) <<
" with auxiliary prototypes " << *auxProto ;
81 if (forceDirect && forceDirect->getSize()>0) ccxcoutI(Generation) <<
" with internal generation forced for observables " << *forceDirect ;
82 ccxcoutI(Generation) << endl ;
87 RooArgSet nodes(model,model.GetName());
88 _cloneSet= (RooArgSet*) nodes.snapshot(kTRUE);
90 coutE(Generation) <<
"RooGenContext::RooGenContext(" << GetName() <<
") Couldn't deep-clone PDF, abort," << endl ;
91 RooErrorHandler::softAbort() ;
95 _pdfClone = (RooAbsPdf*)_cloneSet->find(model.GetName());
96 _pdfClone->setOperMode(RooAbsArg::ADirty,kTRUE) ;
99 if (prototype&&_pdfClone->dependsOn(*prototype->get())) {
100 RooArgSet fullNormSet(vars) ;
101 fullNormSet.add(*prototype->get()) ;
102 _pdfClone->fixAddCoefNormalization(fullNormSet) ;
107 TIterator *iterator= vars.createIterator();
108 TIterator *servers= _pdfClone->serverIterator();
109 const RooAbsArg *tmp = 0;
110 const RooAbsArg *arg = 0;
111 while((_isValid && (tmp= (
const RooAbsArg*)iterator->Next()))) {
113 if(!tmp->isFundamental()) {
114 coutE(Generation) <<
"RooGenContext::ctor(): cannot generate values for derived \"" << tmp->GetName() <<
"\"" << endl;
119 arg= (
const RooAbsArg*)_cloneSet->find(tmp->GetName());
121 coutI(Generation) <<
"RooGenContext::ctor() WARNING model does not depend on \"" << tmp->GetName()
122 <<
"\" which will have uniform distribution" << endl;
123 _uniformVars.add(*tmp);
129 const RooAbsArg *direct= arg ;
130 if (forceDirect==0 || !forceDirect->find(direct->GetName())) {
131 if (!_pdfClone->isDirectGenSafe(*arg)) {
132 cxcoutD(Generation) <<
"RooGenContext::ctor() observable " << arg->GetName() <<
" has been determined to be unsafe for internal generation" << endl;
140 _directVars.add(*arg);
142 _otherVars.add(*arg);
149 coutE(Generation) <<
"RooGenContext::ctor() constructor failed with errors" << endl;
155 if (_directVars.getSize()>0) {
156 cxcoutD(Generation) <<
"RooGenContext::ctor() observables " << _directVars <<
" are safe for internal generator (if supported by p.d.f)" << endl ;
158 if (_otherVars.getSize()>0) {
159 cxcoutD(Generation) <<
"RooGenContext::ctor() observables " << _otherVars <<
" are NOT safe for internal generator (if supported by p.d.f)" << endl ;
164 Bool_t staticInitOK = !_pdfClone->dependsOn(_protoVars) ;
166 cxcoutD(Generation) <<
"RooGenContext::ctor() Model depends on supplied protodata observables, static initialization of internal generator will not be allowed" << endl ;
170 RooArgSet generatedVars;
171 if (_directVars.getSize()>0) {
172 _code= _pdfClone->getGenerator(_directVars,generatedVars,staticInitOK);
174 cxcoutD(Generation) <<
"RooGenContext::ctor() Model indicates that it can internally generate observables "
175 << generatedVars <<
" with configuration identifier " << _code << endl ;
181 _directVars.remove(generatedVars);
182 _otherVars.add(_directVars);
185 _directVars.removeAll();
186 _directVars.add(generatedVars);
188 cxcoutI(Generation) <<
"RooGenContext::ctor() Context will" ;
189 if (_directVars.getSize()>0) ccxcoutI(Generation) <<
" let model internally generate observables " << _directVars ;
190 if (_directVars.getSize()>0 && _otherVars.getSize()>0) ccxcoutI(Generation) <<
" and" ;
191 if (_otherVars.getSize()>0) ccxcoutI(Generation) <<
" generate variables " << _otherVars <<
" with accept/reject sampling" ;
192 if (_uniformVars.getSize()>0) ccxcoutI(Generation) <<
", non-observable variables " << _uniformVars <<
" will be generated with flat distribution" ;
193 ccxcoutI(Generation) << endl ;
196 RooArgSet *depList= _pdfClone->getObservables(_theEvent);
197 depList->remove(_otherVars);
199 TString nname(_pdfClone->GetName()) ;
200 nname.Append(
"_AccRej") ;
201 TString ntitle(_pdfClone->GetTitle()) ;
202 ntitle.Append(
" (Accept/Reject)") ;
205 RooArgSet* protoDeps = model.getObservables(_protoVars) ;
208 if (_protoVars.getSize()==0) {
212 if(depList->getSize()==0) {
216 Int_t maxFindCode = _pdfClone->getMaxVal(_otherVars) ;
217 if (maxFindCode != 0) {
218 coutI(Generation) <<
"RooGenContext::ctor() no prototype data provided, all observables are generated with numerically and "
219 <<
"model supports analytical maximum findin:, can provide analytical pdf maximum to numeric generator" << endl ;
220 Double_t maxVal = _pdfClone->maxVal(maxFindCode) / _pdfClone->getNorm(_theEvent) ;
221 _maxVar =
new RooRealVar(
"funcMax",
"function maximum",maxVal) ;
222 cxcoutD(Generation) <<
"RooGenContext::ctor() maximum value returned by RooAbsPdf::maxVal() is " << maxVal << endl ;
226 if (_otherVars.getSize()>0) {
227 _pdfClone->getVal(&vars) ;
228 _acceptRejectFunc=
new RooRealIntegral(nname,ntitle,*_pdfClone,*depList,&vars);
229 cxcoutI(Generation) <<
"RooGenContext::ctor() accept/reject sampling function is " << _acceptRejectFunc->GetName() << endl ;
231 _acceptRejectFunc = 0 ;
237 depList->remove(_protoVars,kTRUE,kTRUE) ;
238 _acceptRejectFunc= (RooRealIntegral*) _pdfClone->createIntegral(*depList,vars) ;
239 cxcoutI(Generation) <<
"RooGenContext::ctor() accept/reject sampling function is " << _acceptRejectFunc->GetName() << endl ;
242 RooArgSet allVars(_otherVars);
243 allVars.add(_directVars);
244 Int_t maxFindCode = _pdfClone->getMaxVal(allVars) ;
245 if (maxFindCode != 0) {
247 coutI(Generation) <<
"RooGenContext::ctor() prototype data provided, and "
248 <<
"model supports analytical maximum finding in the full phase space: "
249 <<
"can provide analytical pdf maximum to numeric generator" << endl ;
250 _maxVar =
new RooRealVar(
"funcMax",
"function maximum",_pdfClone->maxVal(maxFindCode)) ;
252 maxFindCode = _pdfClone->getMaxVal(_otherVars) ;
253 if (maxFindCode != 0) {
254 _updateFMaxPerEvent = maxFindCode ;
255 coutI(Generation) <<
"RooGenContext::ctor() prototype data provided, and "
256 <<
"model supports analytical maximum finding in the variables which are not"
257 <<
" internally generated. Can provide analytical pdf maximum to numeric "
258 <<
"generator" << endl;
259 cxcoutD(Generation) <<
"RooGenContext::ctor() maximum value must be reevaluated for each "
260 <<
"event with configuration code " << maxFindCode << endl ;
261 _maxVar =
new RooRealVar(
"funcMax",
"function maximum",1) ;
268 RooArgSet otherAndProto(_otherVars) ;
270 otherAndProto.add(*protoDeps) ;
272 if (_otherVars.getSize()>0) {
274 cxcoutD(Generation) <<
"RooGenContext::ctor() prototype data provided, observables are generated numericaly no "
275 <<
"analytical estimate of maximum function value provided by model, must determine maximum value through initial sampling space "
276 <<
"of accept/reject observables plus prototype observables: " << otherAndProto << endl ;
279 RooAbsNumGenerator* maxFinder = RooNumGenFactory::instance().createSampler(*_acceptRejectFunc,otherAndProto,RooArgSet(_protoVars),
280 *model.getGeneratorConfig(),_verbose) ;
282 Double_t max = maxFinder->getFuncMax() ;
283 _maxVar =
new RooRealVar(
"funcMax",
"function maximum",max) ;
286 coutE(Generation) <<
"RooGenContext::ctor(" << model.GetName()
287 <<
") ERROR: generating conditional p.d.f. which requires prior knowledge of function maximum, "
288 <<
"but chosen numeric generator (" << maxFinder->IsA()->GetName() <<
") does not support maximum finding" << endl ;
290 throw string(
"RooGenContext::ctor()") ;
294 cxcoutD(Generation) <<
"RooGenContext::ctor() maximum function value found through initial sampling is " << max << endl ;
300 if (_acceptRejectFunc && _otherVars.getSize()>0) {
301 _generator = RooNumGenFactory::instance().createSampler(*_acceptRejectFunc,_otherVars,RooArgSet(_protoVars),*model.getGeneratorConfig(),_verbose,_maxVar) ;
302 cxcoutD(Generation) <<
"RooGenContext::ctor() creating MC sampling generator " << _generator->IsA()->GetName() <<
" from function for observables " << _otherVars << endl ;
310 _otherVars.add(_uniformVars);
317 RooGenContext::~RooGenContext()
323 if (_generator)
delete _generator;
324 if (_acceptRejectFunc)
delete _acceptRejectFunc;
325 if (_maxVar)
delete _maxVar ;
326 if (_uniIter)
delete _uniIter ;
334 void RooGenContext::attach(
const RooArgSet& args)
336 _pdfClone->recursiveRedirectServers(args,kFALSE);
337 if (_acceptRejectFunc) {
338 _acceptRejectFunc->recursiveRedirectServers(args,kFALSE) ;
343 _generator->attachParameters(args) ;
353 void RooGenContext::initGenerator(
const RooArgSet &theEvent)
355 RooFIter iter = theEvent.fwdIterator() ;
357 while((arg=iter.next())) {
358 arg->setOperMode(RooAbsArg::ADirty) ;
364 _pdfClone->resetErrorCounters();
367 if (_directVars.getSize()>0) {
368 cxcoutD(Generation) <<
"RooGenContext::initGenerator() initializing internal generator of model with code " << _code << endl ;
369 _pdfClone->initGenerator(_code) ;
373 if (_uniformVars.getSize()>0) {
374 _uniIter = _uniformVars.createIterator() ;
383 void RooGenContext::generateEvent(RooArgSet &theEvent, Int_t remaining)
385 if(_otherVars.getSize() > 0) {
388 if (_updateFMaxPerEvent!=0) {
389 Double_t max = _pdfClone->maxVal(_updateFMaxPerEvent)/_pdfClone->getNorm(_otherVars) ;
390 cxcoutD(Generation) <<
"RooGenContext::initGenerator() reevaluation of maximum function value is required for each event, new value is " << max << endl ;
391 _maxVar->setVal(max) ;
395 Double_t resampleRatio(1) ;
396 const RooArgSet *subEvent= _generator->generateEvent(remaining,resampleRatio);
397 if (resampleRatio<1) {
398 coutI(Generation) <<
"RooGenContext::generateEvent INFO: accept/reject generator requests resampling of previously produced events by factor "
399 << resampleRatio <<
" due to increased maximum weight" << endl ;
400 resampleData(resampleRatio) ;
403 coutE(Generation) <<
"RooGenContext::generateEvent ERROR accept/reject generator failed" << endl ;
406 theEvent.assignValueOnly(*subEvent) ;
414 if(_directVars.getSize() > 0) {
415 _pdfClone->generateEvent(_code);
422 while((uniVar=(RooAbsArg*)_uniIter->Next())) {
423 RooAbsLValue* arglv =
dynamic_cast<RooAbsLValue*
>(uniVar) ;
425 coutE(Generation) <<
"RooGenContext::generateEvent(" << GetName() <<
") ERROR: uniform variable " << uniVar->GetName() <<
" is not an lvalue" << endl ;
426 RooErrorHandler::softAbort() ;
430 theEvent = _uniformVars ;
439 void RooGenContext::printMultiline(ostream &os, Int_t content, Bool_t verbose, TString indent)
const
441 RooAbsGenContext::printMultiline(os,content,verbose,indent);
442 os << indent <<
" --- RooGenContext --- " << endl ;
443 os << indent <<
"Using PDF ";
444 _pdfClone->printStream(os,kName|kArgs|kClassName,kSingleLine,indent);
447 os << indent <<
"Use PDF generator for " << _directVars << endl ;
448 os << indent <<
"Use MC sampling generator " << (_generator ? _generator->IsA()->GetName() :
"<none>") <<
" for " << _otherVars << endl ;
449 if (_protoVars.getSize()>0) {
450 os << indent <<
"Prototype observables are " << _protoVars << endl ;