Logo ROOT  
Reference Guide
RooConvGenContext.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17/**
18\file RooConvGenContext.cxx
19\class RooConvGenContext
20\ingroup Roofitcore
21
22RooConvGenContext is an efficient implementation of the generator context
23specific for RooAbsAnaConvPdf objects. The physics model is generated
24with a truth resolution model and the requested resolution model is generated
25separately as a PDF. The convolution variable of the physics model is
26subsequently explicitly smeared with the resolution model distribution.
27**/
28
29#include "RooMsgService.h"
30#include "RooErrorHandler.h"
31#include "RooConvGenContext.h"
32#include "RooAbsAnaConvPdf.h"
33#include "RooNumConvPdf.h"
34#include "RooFFTConvPdf.h"
35#include "RooProdPdf.h"
36#include "RooDataSet.h"
37#include "RooArgSet.h"
38#include "RooTruthModel.h"
39#include "Riostream.h"
40
41
42using namespace std;
43
45;
46
47
48////////////////////////////////////////////////////////////////////////////////
49/// Constructor for specialized generator context for analytical convolutions.
50///
51/// Builds a generator for the physics PDF convoluted with the truth model
52/// and a generator for the resolution model as PDF. Events are generated
53/// by sampling events from the p.d.f and smearings from the resolution model
54/// and adding these to obtain a distribution of events consistent with the
55/// convolution of these two. The advantage of this procedure is so that
56/// both p.d.f and resolution model can take advantage of any internal
57/// generators that may be defined.
58
60 const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
61 RooAbsGenContext(model,vars,prototype,auxProto,verbose), _pdfVarsOwned(0), _modelVarsOwned(0)
62{
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 ;
65
66 // Clone PDF and change model to internal truth model
68 if (!_pdfCloneSet) {
69 coutE(Generation) << "RooConvGenContext::RooConvGenContext(" << GetName() << ") Couldn't deep-clone PDF, abort," << endl ;
71 }
72
74 RooTruthModel truthModel("truthModel","Truth resolution model",(RooRealVar&)*pdfClone->convVar()) ;
75 pdfClone->changeModel(truthModel) ;
76 ((RooRealVar*)pdfClone->convVar())->removeRange() ;
77
78 // Create generator for physics X truth model
79 _pdfVars = (RooArgSet*) pdfClone->getObservables(&vars) ; ;
80 _pdfGen = pdfClone->genContext(*_pdfVars,prototype,auxProto,verbose) ;
81
82 // Clone resolution model and use as normal PDF
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 ;
87 }
89 _modelCloneSet->find(model._convSet.at(0)->GetName())->Clone("smearing") ;
90 _modelCloneSet->addOwned(*modelClone) ;
91 modelClone->changeBasis(0) ;
92 modelClone->convVar().removeRange() ;
93
94 // Create generator for resolution model as PDF
95 _modelVars = (RooArgSet*) modelClone->getObservables(&vars) ;
96
97 _modelVars->add(modelClone->convVar()) ;
98 _convVarName = modelClone->convVar().GetName() ;
99 _modelGen = modelClone->genContext(*_modelVars,prototype,auxProto,verbose) ;
100
101 if (prototype) {
102 _pdfVars->add(*prototype->get()) ;
103 _modelVars->add(*prototype->get()) ;
104 }
105
106 // WVE ADD FOR DEBUGGING
107 if (auxProto) {
108 _pdfVars->add(*auxProto) ;
109 _modelVars->add(*auxProto) ;
110 }
111
112// cout << "RooConvGenContext::ctor(" << this << "," << GetName() << ") _pdfVars = " << _pdfVars << " " ; _pdfVars->Print("1") ;
113// cout << "RooConvGenContext::ctor(" << this << "," << GetName() << ") _modelVars = " << _modelVars << " " ; _modelVars->Print("1") ;
114}
115
116
117
118////////////////////////////////////////////////////////////////////////////////
119/// Constructor for specialized generator context for numerical convolutions.
120///
121/// Builds a generator for the physics PDF convoluted with the truth model
122/// and a generator for the resolution model as PDF. Events are generated
123/// by sampling events from the p.d.f and smearings from the resolution model
124/// and adding these to obtain a distribution of events consistent with the
125/// convolution of these two. The advantage of this procedure is so that
126/// both p.d.f and resolution model can take advantage of any internal
127/// generators that may be defined.
128
130 const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
131 RooAbsGenContext(model,vars,prototype,auxProto,verbose)
132{
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 ;
135
136 // Create generator for physics X truth model
139 _pdfGen = ((RooAbsPdf&)model.conv().clonePdf()).genContext(*_pdfVars,prototype,auxProto,verbose) ;
140 _pdfCloneSet = 0 ;
141
142 // Create generator for resolution model as PDF
145 _convVarName = model.conv().cloneVar().GetName() ;
146 _modelGen = ((RooAbsPdf&)model.conv().cloneModel()).genContext(*_modelVars,prototype,auxProto,verbose) ;
148 _modelCloneSet->add(model.conv().cloneModel()) ;
149
150 if (prototype) {
151 _pdfVars->add(*prototype->get()) ;
152 _modelVars->add(*prototype->get()) ;
153 }
154}
155
156
157
158////////////////////////////////////////////////////////////////////////////////
159/// Constructor for specialized generator context for FFT numerical convolutions.
160///
161/// Builds a generator for the physics PDF convoluted with the truth model
162/// and a generator for the resolution model as PDF. Events are generated
163/// by sampling events from the p.d.f and smearings from the resolution model
164/// and adding these to obtain a distribution of events consistent with the
165/// convolution of these two. The advantage of this procedure is so that
166/// both p.d.f and resolution model can take advantage of any internal
167/// generators that may be defined.
168
170 const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
171 RooAbsGenContext(model,vars,prototype,auxProto,verbose)
172{
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 ;
175
176 _convVarName = model._x.arg().GetName() ;
177
178 // Create generator for physics model
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) ;
186 _pdfGen = pdfClone->genContext(*_pdfVars,prototype,auxProto,verbose) ;
187
188 // Create generator for resolution model
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) ;
196 _modelGen = modelClone->genContext(*_pdfVars,prototype,auxProto,verbose) ;
197
198 delete tmp1 ;
199 delete tmp2 ;
200
201 if (prototype) {
202 _pdfVars->add(*prototype->get()) ;
203 _modelVars->add(*prototype->get()) ;
204 }
205}
206
207
208
209////////////////////////////////////////////////////////////////////////////////
210/// Destructor
211
213{
214 // Destructor. Delete all owned subgenerator contexts
215 delete _pdfGen ;
216 delete _modelGen ;
217 delete _pdfCloneSet ;
218 delete _modelCloneSet ;
219 delete _modelVars ;
220 delete _pdfVars ;
221 delete _pdfVarsOwned ;
222 delete _modelVarsOwned ;
223}
224
225
226
227////////////////////////////////////////////////////////////////////////////////
228/// Attach given set of arguments to internal clones of
229/// pdf and resolution model
230
232{
233 // Find convolution variable in input and output sets
236
237 // Replace all servers in _pdfVars and _modelVars with those in theEvent, except for the convolution variable
238 RooArgSet* pdfCommon = (RooArgSet*) args.selectCommon(*_pdfVars) ;
239 pdfCommon->remove(*cvPdf,kTRUE,kTRUE) ;
240
241 RooArgSet* modelCommon = (RooArgSet*) args.selectCommon(*_modelVars) ;
242 modelCommon->remove(*cvModel,kTRUE,kTRUE) ;
243
244 _pdfGen->attach(*pdfCommon) ;
245 _modelGen->attach(*modelCommon) ;
246
247 delete pdfCommon ;
248 delete modelCommon ;
249}
250
251
252////////////////////////////////////////////////////////////////////////////////
253/// One-time initialization of generator context, attaches
254/// the context to the supplied event container
255
257{
258 // Find convolution variable in input and output sets
261 _cvOut = (RooRealVar*) theEvent.find(_convVarName) ;
262
263 // Replace all servers in _pdfVars and _modelVars with those in theEvent, except for the convolution variable
264 RooArgSet* pdfCommon = (RooArgSet*) theEvent.selectCommon(*_pdfVars) ;
265 pdfCommon->remove(*_cvPdf,kTRUE,kTRUE) ;
266 _pdfVars->replace(*pdfCommon) ;
267 delete pdfCommon ;
268
269 RooArgSet* modelCommon = (RooArgSet*) theEvent.selectCommon(*_modelVars) ;
270 modelCommon->remove(*_cvModel,kTRUE,kTRUE) ;
271 _modelVars->replace(*modelCommon) ;
272 delete modelCommon ;
273
274 // Initialize component generators
277}
278
279
280
281////////////////////////////////////////////////////////////////////////////////
282/// Generate a single event
283
285{
286 while(1) {
287
288 // Generate pdf and model data
289 _modelGen->generateEvent(*_modelVars,remaining) ;
290 _pdfGen->generateEvent(*_pdfVars,remaining) ;
291
292 // Construct smeared convolution variable
293 Double_t convValSmeared = _cvPdf->getVal() + _cvModel->getVal() ;
294 if (_cvOut->isValidReal(convValSmeared)) {
295 // Smeared value in acceptance range, transfer values to output set
296 theEvent = *_modelVars ;
297 theEvent = *_pdfVars ;
298 _cvOut->setVal(convValSmeared) ;
299 return ;
300 }
301 }
302}
303
304
305
306////////////////////////////////////////////////////////////////////////////////
307/// Set the traversal order for events in the prototype dataset
308/// The argument is a array of integers with a size identical
309/// to the number of events in the prototype dataset. Each element
310/// should contain an integer in the range 1-N.
311
313{
317}
318
319
320////////////////////////////////////////////////////////////////////////////////
321/// Print the details of this generator context
322
324{
326 os << indent << "--- RooConvGenContext ---" << endl ;
327 os << indent << "List of component generators" << endl ;
328
329 TString indent2(indent) ;
330 indent2.Append(" ") ;
331
332 _modelGen->printMultiline(os,content,verbose,indent2);
333 _pdfGen->printMultiline(os,content,verbose,indent2);
334}
#define cxcoutI(a)
Definition: RooMsgService.h:86
#define coutE(a)
Definition: RooMsgService.h:34
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:365
static void indent(ostringstream &buf, int indent_level)
RooAbsAnaConvPdf is the base class for PDFs that represent a physics model that can be analytically c...
virtual Bool_t changeModel(const RooResolutionModel &newModel)
Change the current resolution model to newModel.
const RooRealVar * convVar() const
Return a pointer to the convolution variable instance used in the resolution model.
virtual RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=0, const RooArgSet *auxProto=0, Bool_t verbose=kFALSE) const
Create a generator context for this p.d.f.
RooListProxy _convSet
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Return the observables of this pdf given a set of observables.
Definition: RooAbsArg.h:240
RooAbsCollection * selectCommon(const RooAbsCollection &refColl) const
Create a subset of the current collection, consisting only of those elements that are contained as we...
virtual Bool_t replace(const RooAbsArg &var1, const RooAbsArg &var2)
Replace var1 with var2 and return kTRUE for success.
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooAbsGenContext is the abstract base class for generator contexts of RooAbsPdf objects.
virtual void attach(const RooArgSet &params)
Interface to attach given parameters to object in this context.
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Interface for multi-line printing.
virtual void initGenerator(const RooArgSet &theEvent)
Interface function to initialize context for generation for given set of observables.
virtual void generateEvent(RooArgSet &theEvent, Int_t remaining)=0
virtual void setProtoDataOrder(Int_t *lut)
Set the traversal order of prototype data to that in the lookup tables passed as argument.
virtual RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=0, const RooArgSet *auxProto=0, Bool_t verbose=kFALSE) const
Interface function to create a generator context from a p.d.f.
Definition: RooAbsPdf.cxx:2033
virtual Bool_t isValidReal(Double_t value, Bool_t printError=kFALSE) const
Check if given value is valid.
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:87
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition: RooArgList.h:74
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition: RooArgSet.h:134
virtual Bool_t addOwned(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling addOwned() for each element in the source...
Definition: RooArgSet.h:92
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:88
RooConvGenContext is an efficient implementation of the generator context specific for RooAbsAnaConvP...
virtual ~RooConvGenContext()
Destructor.
RooArgSet * _modelCloneSet
RooArgSet * _modelVarsOwned
RooConvGenContext(const RooFFTConvPdf &model, const RooArgSet &vars, const RooDataSet *prototype=0, const RooArgSet *auxProto=0, Bool_t _verbose=kFALSE)
Constructor for specialized generator context for FFT numerical convolutions.
RooAbsGenContext * _pdfGen
RooAbsGenContext * _modelGen
RooArgSet * _pdfCloneSet
virtual void printMultiline(std::ostream &os, Int_t content, Bool_t verbose=kFALSE, TString indent="") const
Print the details of this generator context.
RooArgSet * _pdfVarsOwned
RooRealVar * _cvModel
virtual void generateEvent(RooArgSet &theEvent, Int_t remaining)
Generate a single event.
virtual void initGenerator(const RooArgSet &theEvent)
One-time initialization of generator context, attaches the context to the supplied event container.
virtual void attach(const RooArgSet &params)
Attach given set of arguments to internal clones of pdf and resolution model.
RooArgSet * _modelVars
virtual void setProtoDataOrder(Int_t *lut)
Set the traversal order for events in the prototype dataset The argument is a array of integers with ...
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
virtual const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
static void softAbort()
PDF for the numerical (FFT) convolution of two PDFs.
Definition: RooFFTConvPdf.h:26
RooRealProxy _pdf1
Definition: RooFFTConvPdf.h:70
RooRealProxy _x
Definition: RooFFTConvPdf.h:68
RooRealProxy _pdf2
Definition: RooFFTConvPdf.h:71
Numeric 1-dimensional convolution operator PDF.
Definition: RooNumConvPdf.h:26
RooNumConvolution & conv() const
Definition: RooNumConvPdf.h:63
RooAbsReal & cloneModel() const
RooRealVar & cloneVar() const
RooAbsReal & clonePdf() const
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
void removeRange(const char *name=0)
Remove range limits for binning with given name. Empty name means default range.
Definition: RooRealVar.cxx:421
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Definition: RooRealVar.cxx:278
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...
virtual void changeBasis(RooFormulaVar *basis)
Change the basis function we convolute with.
RooRealVar & convVar() const
Return the convolution variable of the resolution model.
const T & arg() const
Return reference to object held in proxy.
RooTruthModel is an implementation of RooResolution model that provides a delta-function resolution m...
Definition: RooTruthModel.h:21
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition: TNamed.cxx:74
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Basic string class.
Definition: TString.h:131
TString & Append(const char *cs)
Definition: TString.h:559
@ Generation
Definition: RooGlobalFunc.h:67