44 ClassImp(RooConvGenContext);
59 RooConvGenContext::RooConvGenContext(
const RooAbsAnaConvPdf &model,
const RooArgSet &vars,
60 const RooDataSet *prototype,
const RooArgSet* auxProto, Bool_t verbose) :
61 RooAbsGenContext(model,vars,prototype,auxProto,verbose), _pdfVarsOwned(0), _modelVarsOwned(0)
63 cxcoutI(Generation) <<
"RooConvGenContext::ctor() setting up special generator context for analytical convolution p.d.f. " << model.GetName()
64 <<
" for generation of observable(s) " << vars << endl ;
67 _pdfCloneSet = (RooArgSet*) RooArgSet(model).snapshot(kTRUE) ;
69 coutE(Generation) <<
"RooConvGenContext::RooConvGenContext(" << GetName() <<
") Couldn't deep-clone PDF, abort," << endl ;
70 RooErrorHandler::softAbort() ;
73 RooAbsAnaConvPdf* pdfClone = (RooAbsAnaConvPdf*) _pdfCloneSet->find(model.GetName()) ;
74 RooTruthModel truthModel(
"truthModel",
"Truth resolution model",(RooRealVar&)*pdfClone->convVar()) ;
75 pdfClone->changeModel(truthModel) ;
76 ((RooRealVar*)pdfClone->convVar())->removeRange() ;
79 _pdfVars = (RooArgSet*) pdfClone->getObservables(&vars) ; ;
80 _pdfGen = pdfClone->genContext(*_pdfVars,prototype,auxProto,verbose) ;
83 _modelCloneSet = (RooArgSet*) RooArgSet(*model._convSet.at(0)).snapshot(kTRUE) ;
84 if (!_modelCloneSet) {
85 coutE(Generation) <<
"RooConvGenContext::RooConvGenContext(" << GetName() <<
") Couldn't deep-clone resolution model, abort," << endl ;
86 RooErrorHandler::softAbort() ;
88 RooResolutionModel* modelClone = (RooResolutionModel*)
89 _modelCloneSet->find(model._convSet.at(0)->GetName())->Clone(
"smearing") ;
90 _modelCloneSet->addOwned(*modelClone) ;
91 modelClone->changeBasis(0) ;
92 modelClone->convVar().removeRange() ;
95 _modelVars = (RooArgSet*) modelClone->getObservables(&vars) ;
97 _modelVars->add(modelClone->convVar()) ;
98 _convVarName = modelClone->convVar().GetName() ;
99 _modelGen = modelClone->genContext(*_modelVars,prototype,auxProto,verbose) ;
102 _pdfVars->add(*prototype->get()) ;
103 _modelVars->add(*prototype->get()) ;
108 _pdfVars->add(*auxProto) ;
109 _modelVars->add(*auxProto) ;
129 RooConvGenContext::RooConvGenContext(
const RooNumConvPdf &model,
const RooArgSet &vars,
130 const RooDataSet *prototype,
const RooArgSet* auxProto, Bool_t verbose) :
131 RooAbsGenContext(model,vars,prototype,auxProto,verbose)
133 cxcoutI(Generation) <<
"RooConvGenContext::ctor() setting up special generator context for numeric convolution p.d.f. " << model.GetName()
134 <<
" for generation of observable(s) " << vars << endl ;
137 _pdfVarsOwned = (RooArgSet*) model.conv().clonePdf().getObservables(&vars)->snapshot(kTRUE) ;
138 _pdfVars =
new RooArgSet(*_pdfVarsOwned) ;
139 _pdfGen = ((RooAbsPdf&)model.conv().clonePdf()).genContext(*_pdfVars,prototype,auxProto,verbose) ;
143 _modelVarsOwned = (RooArgSet*) model.conv().cloneModel().getObservables(&vars)->snapshot(kTRUE) ;
144 _modelVars =
new RooArgSet(*_modelVarsOwned) ;
145 _convVarName = model.conv().cloneVar().GetName() ;
146 _modelGen = ((RooAbsPdf&)model.conv().cloneModel()).genContext(*_modelVars,prototype,auxProto,verbose) ;
147 _modelCloneSet =
new RooArgSet ;
148 _modelCloneSet->add(model.conv().cloneModel()) ;
151 _pdfVars->add(*prototype->get()) ;
152 _modelVars->add(*prototype->get()) ;
169 RooConvGenContext::RooConvGenContext(
const RooFFTConvPdf &model,
const RooArgSet &vars,
170 const RooDataSet *prototype,
const RooArgSet* auxProto, Bool_t verbose) :
171 RooAbsGenContext(model,vars,prototype,auxProto,verbose)
173 cxcoutI(Generation) <<
"RooConvGenContext::ctor() setting up special generator context for fft convolution p.d.f. " << model.GetName()
174 <<
" for generation of observable(s) " << vars << endl ;
176 _convVarName = model._x.arg().GetName() ;
179 _pdfCloneSet = (RooArgSet*) RooArgSet(model._pdf1.arg()).snapshot(kTRUE) ;
180 RooAbsPdf* pdfClone = (RooAbsPdf*) _pdfCloneSet->find(model._pdf1.arg().GetName()) ;
181 RooRealVar* cvPdf = (RooRealVar*) _pdfCloneSet->find(model._x.arg().GetName()) ;
182 cvPdf->removeRange() ;
183 RooArgSet* tmp1 = pdfClone->getObservables(&vars) ;
184 _pdfVarsOwned = (RooArgSet*) tmp1->snapshot(kTRUE) ;
185 _pdfVars =
new RooArgSet(*_pdfVarsOwned) ;
186 _pdfGen = pdfClone->genContext(*_pdfVars,prototype,auxProto,verbose) ;
189 _modelCloneSet = (RooArgSet*) RooArgSet(model._pdf2.arg()).snapshot(kTRUE) ;
190 RooAbsPdf* modelClone = (RooAbsPdf*) _modelCloneSet->find(model._pdf2.arg().GetName()) ;
191 RooRealVar* cvModel = (RooRealVar*) _modelCloneSet->find(model._x.arg().GetName()) ;
192 cvModel->removeRange() ;
193 RooArgSet* tmp2 = modelClone->getObservables(&vars) ;
194 _modelVarsOwned = (RooArgSet*) tmp2->snapshot(kTRUE) ;
195 _modelVars =
new RooArgSet(*_modelVarsOwned) ;
196 _modelGen = modelClone->genContext(*_pdfVars,prototype,auxProto,verbose) ;
202 _pdfVars->add(*prototype->get()) ;
203 _modelVars->add(*prototype->get()) ;
212 RooConvGenContext::~RooConvGenContext()
217 delete _pdfCloneSet ;
218 delete _modelCloneSet ;
221 delete _pdfVarsOwned ;
222 delete _modelVarsOwned ;
231 void RooConvGenContext::attach(
const RooArgSet& args)
234 RooRealVar* cvModel = (RooRealVar*) _modelVars->find(_convVarName) ;
235 RooRealVar* cvPdf = (RooRealVar*) _pdfVars->find(_convVarName) ;
238 RooArgSet* pdfCommon = (RooArgSet*) args.selectCommon(*_pdfVars) ;
239 pdfCommon->remove(*cvPdf,kTRUE,kTRUE) ;
241 RooArgSet* modelCommon = (RooArgSet*) args.selectCommon(*_modelVars) ;
242 modelCommon->remove(*cvModel,kTRUE,kTRUE) ;
244 _pdfGen->attach(*pdfCommon) ;
245 _modelGen->attach(*modelCommon) ;
256 void RooConvGenContext::initGenerator(
const RooArgSet &theEvent)
259 _cvModel = (RooRealVar*) _modelVars->find(_convVarName) ;
260 _cvPdf = (RooRealVar*) _pdfVars->find(_convVarName) ;
261 _cvOut = (RooRealVar*) theEvent.find(_convVarName) ;
264 RooArgSet* pdfCommon = (RooArgSet*) theEvent.selectCommon(*_pdfVars) ;
265 pdfCommon->remove(*_cvPdf,kTRUE,kTRUE) ;
266 _pdfVars->replace(*pdfCommon) ;
269 RooArgSet* modelCommon = (RooArgSet*) theEvent.selectCommon(*_modelVars) ;
270 modelCommon->remove(*_cvModel,kTRUE,kTRUE) ;
271 _modelVars->replace(*modelCommon) ;
275 _pdfGen->initGenerator(*_pdfVars) ;
276 _modelGen->initGenerator(*_modelVars) ;
284 void RooConvGenContext::generateEvent(RooArgSet &theEvent, Int_t remaining)
289 _modelGen->generateEvent(*_modelVars,remaining) ;
290 _pdfGen->generateEvent(*_pdfVars,remaining) ;
293 Double_t convValSmeared = _cvPdf->getVal() + _cvModel->getVal() ;
294 if (_cvOut->isValidReal(convValSmeared)) {
296 theEvent = *_modelVars ;
297 theEvent = *_pdfVars ;
298 _cvOut->setVal(convValSmeared) ;
312 void RooConvGenContext::setProtoDataOrder(Int_t* lut)
314 RooAbsGenContext::setProtoDataOrder(lut) ;
315 _modelGen->setProtoDataOrder(lut) ;
316 _pdfGen->setProtoDataOrder(lut) ;
323 void RooConvGenContext::printMultiline(ostream &os, Int_t content, Bool_t verbose, TString indent)
const
325 RooAbsGenContext::printMultiline(os,content,verbose,indent) ;
326 os << indent <<
"--- RooConvGenContext ---" << endl ;
327 os << indent <<
"List of component generators" << endl ;
329 TString indent2(indent) ;
330 indent2.Append(
" ") ;
332 _modelGen->printMultiline(os,content,verbose,indent2);
333 _pdfGen->printMultiline(os,content,verbose,indent2);