00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 #include "RooFit.h"
00029
00030 #include "RooMsgService.h"
00031 #include "RooConvGenContext.h"
00032 #include "RooAbsAnaConvPdf.h"
00033 #include "RooNumConvPdf.h"
00034 #include "RooFFTConvPdf.h"
00035 #include "RooProdPdf.h"
00036 #include "RooDataSet.h"
00037 #include "RooArgSet.h"
00038 #include "RooTruthModel.h"
00039 #include "Riostream.h"
00040
00041
00042 ClassImp(RooConvGenContext)
00043 ;
00044
00045
00046
00047 RooConvGenContext::RooConvGenContext(const RooAbsAnaConvPdf &model, const RooArgSet &vars,
00048 const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
00049 RooAbsGenContext(model,vars,prototype,auxProto,verbose), _pdfVarsOwned(0), _modelVarsOwned(0)
00050 {
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061 cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for analytical convolution p.d.f. " << model.GetName()
00062 << " for generation of observable(s) " << vars << endl ;
00063
00064
00065 _pdfCloneSet = (RooArgSet*) RooArgSet(model).snapshot(kTRUE) ;
00066 if (!_pdfCloneSet) {
00067 coutE(Generation) << "RooConvGenContext::RooConvGenContext(" << GetName() << ") Couldn't deep-clone PDF, abort," << endl ;
00068 RooErrorHandler::softAbort() ;
00069 }
00070
00071 RooAbsAnaConvPdf* pdfClone = (RooAbsAnaConvPdf*) _pdfCloneSet->find(model.GetName()) ;
00072 RooTruthModel truthModel("truthModel","Truth resolution model",(RooRealVar&)*pdfClone->convVar()) ;
00073 pdfClone->changeModel(truthModel) ;
00074 ((RooRealVar*)pdfClone->convVar())->removeRange() ;
00075
00076
00077 _pdfVars = (RooArgSet*) pdfClone->getObservables(&vars) ; ;
00078 _pdfGen = pdfClone->genContext(*_pdfVars,prototype,auxProto,verbose) ;
00079
00080
00081 _modelCloneSet = (RooArgSet*) RooArgSet(*model._convSet.at(0)).snapshot(kTRUE) ;
00082 if (!_modelCloneSet) {
00083 coutE(Generation) << "RooConvGenContext::RooConvGenContext(" << GetName() << ") Couldn't deep-clone resolution model, abort," << endl ;
00084 RooErrorHandler::softAbort() ;
00085 }
00086 RooResolutionModel* modelClone = (RooResolutionModel*)
00087 _modelCloneSet->find(model._convSet.at(0)->GetName())->Clone("smearing") ;
00088 _modelCloneSet->addOwned(*modelClone) ;
00089 modelClone->changeBasis(0) ;
00090 modelClone->convVar().removeRange() ;
00091
00092
00093 _modelVars = (RooArgSet*) modelClone->getObservables(&vars) ;
00094
00095 _modelVars->add(modelClone->convVar()) ;
00096 _convVarName = modelClone->convVar().GetName() ;
00097 _modelGen = modelClone->genContext(*_modelVars,prototype,auxProto,verbose) ;
00098
00099 if (prototype) {
00100 _pdfVars->add(*prototype->get()) ;
00101 _modelVars->add(*prototype->get()) ;
00102 }
00103
00104
00105 if (auxProto) {
00106 _pdfVars->add(*auxProto) ;
00107 _modelVars->add(*auxProto) ;
00108 }
00109
00110
00111
00112 }
00113
00114
00115
00116
00117 RooConvGenContext::RooConvGenContext(const RooNumConvPdf &model, const RooArgSet &vars,
00118 const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
00119 RooAbsGenContext(model,vars,prototype,auxProto,verbose)
00120 {
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131 cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for numeric convolution p.d.f. " << model.GetName()
00132 << " for generation of observable(s) " << vars << endl ;
00133
00134
00135 _pdfVarsOwned = (RooArgSet*) model.conv().clonePdf().getObservables(&vars)->snapshot(kTRUE) ;
00136 _pdfVars = new RooArgSet(*_pdfVarsOwned) ;
00137 _pdfGen = ((RooAbsPdf&)model.conv().clonePdf()).genContext(*_pdfVars,prototype,auxProto,verbose) ;
00138 _pdfCloneSet = 0 ;
00139
00140
00141 _modelVarsOwned = (RooArgSet*) model.conv().cloneModel().getObservables(&vars)->snapshot(kTRUE) ;
00142 _modelVars = new RooArgSet(*_modelVarsOwned) ;
00143 _convVarName = model.conv().cloneVar().GetName() ;
00144 _modelGen = ((RooAbsPdf&)model.conv().cloneModel()).genContext(*_modelVars,prototype,auxProto,verbose) ;
00145 _modelCloneSet = 0 ;
00146
00147 if (prototype) {
00148 _pdfVars->add(*prototype->get()) ;
00149 _modelVars->add(*prototype->get()) ;
00150 }
00151 }
00152
00153
00154
00155
00156 RooConvGenContext::RooConvGenContext(const RooFFTConvPdf &model, const RooArgSet &vars,
00157 const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
00158 RooAbsGenContext(model,vars,prototype,auxProto,verbose)
00159 {
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170 cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for fft convolution p.d.f. " << model.GetName()
00171 << " for generation of observable(s) " << vars << endl ;
00172
00173 _convVarName = model._x.arg().GetName() ;
00174
00175
00176 _pdfCloneSet = (RooArgSet*) RooArgSet(model._pdf1.arg()).snapshot(kTRUE) ;
00177 RooAbsPdf* pdfClone = (RooAbsPdf*) _pdfCloneSet->find(model._pdf1.arg().GetName()) ;
00178 RooRealVar* cvPdf = (RooRealVar*) _pdfCloneSet->find(model._x.arg().GetName()) ;
00179 cvPdf->removeRange() ;
00180 RooArgSet* tmp1 = pdfClone->getObservables(&vars) ;
00181 _pdfVarsOwned = (RooArgSet*) tmp1->snapshot(kTRUE) ;
00182 _pdfVars = new RooArgSet(*_pdfVarsOwned) ;
00183 _pdfGen = pdfClone->genContext(*_pdfVars,prototype,auxProto,verbose) ;
00184
00185
00186 _modelCloneSet = (RooArgSet*) RooArgSet(model._pdf2.arg()).snapshot(kTRUE) ;
00187 RooAbsPdf* modelClone = (RooAbsPdf*) _modelCloneSet->find(model._pdf2.arg().GetName()) ;
00188 RooRealVar* cvModel = (RooRealVar*) _modelCloneSet->find(model._x.arg().GetName()) ;
00189 cvModel->removeRange() ;
00190 RooArgSet* tmp2 = modelClone->getObservables(&vars) ;
00191 _modelVarsOwned = (RooArgSet*) tmp2->snapshot(kTRUE) ;
00192 _modelVars = new RooArgSet(*_modelVarsOwned) ;
00193 _modelGen = modelClone->genContext(*_pdfVars,prototype,auxProto,verbose) ;
00194
00195 delete tmp1 ;
00196 delete tmp2 ;
00197
00198 if (prototype) {
00199 _pdfVars->add(*prototype->get()) ;
00200 _modelVars->add(*prototype->get()) ;
00201 }
00202 }
00203
00204
00205
00206
00207 RooConvGenContext::~RooConvGenContext()
00208 {
00209
00210
00211
00212 delete _pdfGen ;
00213 delete _modelGen ;
00214 delete _pdfCloneSet ;
00215 delete _modelCloneSet ;
00216 delete _modelVars ;
00217 delete _pdfVars ;
00218 delete _pdfVarsOwned ;
00219 delete _modelVarsOwned ;
00220 }
00221
00222
00223
00224
00225 void RooConvGenContext::attach(const RooArgSet& args)
00226 {
00227
00228
00229
00230
00231 RooRealVar* cvModel = (RooRealVar*) _modelVars->find(_convVarName) ;
00232 RooRealVar* cvPdf = (RooRealVar*) _pdfVars->find(_convVarName) ;
00233
00234
00235 RooArgSet* pdfCommon = (RooArgSet*) args.selectCommon(*_pdfVars) ;
00236 pdfCommon->remove(*cvPdf,kTRUE,kTRUE) ;
00237
00238 RooArgSet* modelCommon = (RooArgSet*) args.selectCommon(*_modelVars) ;
00239 modelCommon->remove(*cvModel,kTRUE,kTRUE) ;
00240
00241 _pdfGen->attach(*pdfCommon) ;
00242 _modelGen->attach(*modelCommon) ;
00243
00244 delete pdfCommon ;
00245 delete modelCommon ;
00246 }
00247
00248
00249
00250 void RooConvGenContext::initGenerator(const RooArgSet &theEvent)
00251 {
00252
00253
00254
00255
00256 _cvModel = (RooRealVar*) _modelVars->find(_convVarName) ;
00257 _cvPdf = (RooRealVar*) _pdfVars->find(_convVarName) ;
00258 _cvOut = (RooRealVar*) theEvent.find(_convVarName) ;
00259
00260
00261 RooArgSet* pdfCommon = (RooArgSet*) theEvent.selectCommon(*_pdfVars) ;
00262 pdfCommon->remove(*_cvPdf,kTRUE,kTRUE) ;
00263 _pdfVars->replace(*pdfCommon) ;
00264 delete pdfCommon ;
00265
00266 RooArgSet* modelCommon = (RooArgSet*) theEvent.selectCommon(*_modelVars) ;
00267 modelCommon->remove(*_cvModel,kTRUE,kTRUE) ;
00268 _modelVars->replace(*modelCommon) ;
00269 delete modelCommon ;
00270
00271
00272 _pdfGen->initGenerator(*_pdfVars) ;
00273 _modelGen->initGenerator(*_modelVars) ;
00274 }
00275
00276
00277
00278
00279 void RooConvGenContext::generateEvent(RooArgSet &theEvent, Int_t remaining)
00280 {
00281
00282
00283 while(1) {
00284
00285
00286 _modelGen->generateEvent(*_modelVars,remaining) ;
00287 _pdfGen->generateEvent(*_pdfVars,remaining) ;
00288
00289
00290 Double_t convValSmeared = _cvPdf->getVal() + _cvModel->getVal() ;
00291 if (_cvOut->isValidReal(convValSmeared)) {
00292
00293 theEvent = *_modelVars ;
00294 theEvent = *_pdfVars ;
00295 _cvOut->setVal(convValSmeared) ;
00296 return ;
00297 }
00298 }
00299 }
00300
00301
00302
00303
00304 void RooConvGenContext::setProtoDataOrder(Int_t* lut)
00305 {
00306
00307
00308
00309
00310
00311 RooAbsGenContext::setProtoDataOrder(lut) ;
00312 _modelGen->setProtoDataOrder(lut) ;
00313 _pdfGen->setProtoDataOrder(lut) ;
00314 }
00315
00316
00317
00318 void RooConvGenContext::printMultiline(ostream &os, Int_t content, Bool_t verbose, TString indent) const
00319 {
00320
00321
00322 RooAbsGenContext::printMultiline(os,content,verbose,indent) ;
00323 os << indent << "--- RooConvGenContext ---" << endl ;
00324 os << indent << "List of component generators" << endl ;
00325
00326 TString indent2(indent) ;
00327 indent2.Append(" ") ;
00328
00329 _modelGen->printMultiline(os,content,verbose,indent2);
00330 _pdfGen->printMultiline(os,content,verbose,indent2);
00331 }