ROOT  6.06/09
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 //
19 // BEGIN_HTML
20 // RooConvGenContext is an efficient implementation of the generator context
21 // specific for RooAbsAnaConvPdf objects. The physics model is generated
22 // with a truth resolution model and the requested resolution model is generated
23 // separately as a PDF. The convolution variable of the physics model is
24 // subsequently explicitly smeared with the resolution model distribution.
25 // END_HTML
26 //
27 
28 #include "RooMsgService.h"
29 #include "RooErrorHandler.h"
30 #include "RooConvGenContext.h"
31 #include "RooAbsAnaConvPdf.h"
32 #include "RooNumConvPdf.h"
33 #include "RooFFTConvPdf.h"
34 #include "RooProdPdf.h"
35 #include "RooDataSet.h"
36 #include "RooArgSet.h"
37 #include "RooTruthModel.h"
38 #include "Riostream.h"
39 
40 
41 using namespace std;
42 
44 ;
45 
46 
47 ////////////////////////////////////////////////////////////////////////////////
48 /// Constructor for specialized generator context for analytical convolutions.
49 ///
50 /// Builds a generator for the physics PDF convoluted with the truth model
51 /// and a generator for the resolution model as PDF. Events are generated
52 /// by sampling events from the p.d.f and smearings from the resolution model
53 /// and adding these to obtain a distribution of events consistent with the
54 /// convolution of these two. The advantage of this procedure is so that
55 /// both p.d.f and resolution model can take advantage of any internal
56 /// generators that may be defined.
57 
59  const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
60  RooAbsGenContext(model,vars,prototype,auxProto,verbose), _pdfVarsOwned(0), _modelVarsOwned(0)
61 {
62  cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for analytical convolution p.d.f. " << model.GetName()
63  << " for generation of observable(s) " << vars << endl ;
64 
65  // Clone PDF and change model to internal truth model
67  if (!_pdfCloneSet) {
68  coutE(Generation) << "RooConvGenContext::RooConvGenContext(" << GetName() << ") Couldn't deep-clone PDF, abort," << endl ;
70  }
71 
72  RooAbsAnaConvPdf* pdfClone = (RooAbsAnaConvPdf*) _pdfCloneSet->find(model.GetName()) ;
73  RooTruthModel truthModel("truthModel","Truth resolution model",(RooRealVar&)*pdfClone->convVar()) ;
74  pdfClone->changeModel(truthModel) ;
75  ((RooRealVar*)pdfClone->convVar())->removeRange() ;
76 
77  // Create generator for physics X truth model
78  _pdfVars = (RooArgSet*) pdfClone->getObservables(&vars) ; ;
79  _pdfGen = pdfClone->genContext(*_pdfVars,prototype,auxProto,verbose) ;
80 
81  // Clone resolution model and use as normal PDF
82  _modelCloneSet = (RooArgSet*) RooArgSet(*model._convSet.at(0)).snapshot(kTRUE) ;
83  if (!_modelCloneSet) {
84  coutE(Generation) << "RooConvGenContext::RooConvGenContext(" << GetName() << ") Couldn't deep-clone resolution model, abort," << endl ;
86  }
87  RooResolutionModel* modelClone = (RooResolutionModel*)
88  _modelCloneSet->find(model._convSet.at(0)->GetName())->Clone("smearing") ;
89  _modelCloneSet->addOwned(*modelClone) ;
90  modelClone->changeBasis(0) ;
91  modelClone->convVar().removeRange() ;
92 
93  // Create generator for resolution model as PDF
94  _modelVars = (RooArgSet*) modelClone->getObservables(&vars) ;
95 
96  _modelVars->add(modelClone->convVar()) ;
97  _convVarName = modelClone->convVar().GetName() ;
98  _modelGen = modelClone->genContext(*_modelVars,prototype,auxProto,verbose) ;
99 
100  if (prototype) {
101  _pdfVars->add(*prototype->get()) ;
102  _modelVars->add(*prototype->get()) ;
103  }
104 
105  // WVE ADD FOR DEBUGGING
106  if (auxProto) {
107  _pdfVars->add(*auxProto) ;
108  _modelVars->add(*auxProto) ;
109  }
110 
111 // cout << "RooConvGenContext::ctor(" << this << "," << GetName() << ") _pdfVars = " << _pdfVars << " " ; _pdfVars->Print("1") ;
112 // cout << "RooConvGenContext::ctor(" << this << "," << GetName() << ") _modelVars = " << _modelVars << " " ; _modelVars->Print("1") ;
113 }
114 
115 
116 
117 ////////////////////////////////////////////////////////////////////////////////
118 /// Constructor for specialized generator context for numerical convolutions.
119 ///
120 /// Builds a generator for the physics PDF convoluted with the truth model
121 /// and a generator for the resolution model as PDF. Events are generated
122 /// by sampling events from the p.d.f and smearings from the resolution model
123 /// and adding these to obtain a distribution of events consistent with the
124 /// convolution of these two. The advantage of this procedure is so that
125 /// both p.d.f and resolution model can take advantage of any internal
126 /// generators that may be defined.
127 
129  const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
130  RooAbsGenContext(model,vars,prototype,auxProto,verbose)
131 {
132  cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for numeric convolution p.d.f. " << model.GetName()
133  << " for generation of observable(s) " << vars << endl ;
134 
135  // Create generator for physics X truth model
138  _pdfGen = ((RooAbsPdf&)model.conv().clonePdf()).genContext(*_pdfVars,prototype,auxProto,verbose) ;
139  _pdfCloneSet = 0 ;
140 
141  // Create generator for resolution model as PDF
144  _convVarName = model.conv().cloneVar().GetName() ;
145  _modelGen = ((RooAbsPdf&)model.conv().cloneModel()).genContext(*_modelVars,prototype,auxProto,verbose) ;
146  _modelCloneSet = new RooArgSet ;
147  _modelCloneSet->add(model.conv().cloneModel()) ;
148 
149  if (prototype) {
150  _pdfVars->add(*prototype->get()) ;
151  _modelVars->add(*prototype->get()) ;
152  }
153 }
154 
155 
156 
157 ////////////////////////////////////////////////////////////////////////////////
158 /// Constructor for specialized generator context for FFT numerical convolutions.
159 ///
160 /// Builds a generator for the physics PDF convoluted with the truth model
161 /// and a generator for the resolution model as PDF. Events are generated
162 /// by sampling events from the p.d.f and smearings from the resolution model
163 /// and adding these to obtain a distribution of events consistent with the
164 /// convolution of these two. The advantage of this procedure is so that
165 /// both p.d.f and resolution model can take advantage of any internal
166 /// generators that may be defined.
167 
169  const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
170  RooAbsGenContext(model,vars,prototype,auxProto,verbose)
171 {
172  cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for fft convolution p.d.f. " << model.GetName()
173  << " for generation of observable(s) " << vars << endl ;
174 
175  _convVarName = model._x.arg().GetName() ;
176 
177  // Create generator for physics model
178  _pdfCloneSet = (RooArgSet*) RooArgSet(model._pdf1.arg()).snapshot(kTRUE) ;
179  RooAbsPdf* pdfClone = (RooAbsPdf*) _pdfCloneSet->find(model._pdf1.arg().GetName()) ;
180  RooRealVar* cvPdf = (RooRealVar*) _pdfCloneSet->find(model._x.arg().GetName()) ;
181  cvPdf->removeRange() ;
182  RooArgSet* tmp1 = pdfClone->getObservables(&vars) ;
183  _pdfVarsOwned = (RooArgSet*) tmp1->snapshot(kTRUE) ;
185  _pdfGen = pdfClone->genContext(*_pdfVars,prototype,auxProto,verbose) ;
186 
187  // Create generator for resolution model
188  _modelCloneSet = (RooArgSet*) RooArgSet(model._pdf2.arg()).snapshot(kTRUE) ;
189  RooAbsPdf* modelClone = (RooAbsPdf*) _modelCloneSet->find(model._pdf2.arg().GetName()) ;
190  RooRealVar* cvModel = (RooRealVar*) _modelCloneSet->find(model._x.arg().GetName()) ;
191  cvModel->removeRange() ;
192  RooArgSet* tmp2 = modelClone->getObservables(&vars) ;
195  _modelGen = modelClone->genContext(*_pdfVars,prototype,auxProto,verbose) ;
196 
197  delete tmp1 ;
198  delete tmp2 ;
199 
200  if (prototype) {
201  _pdfVars->add(*prototype->get()) ;
202  _modelVars->add(*prototype->get()) ;
203  }
204 }
205 
206 
207 
208 ////////////////////////////////////////////////////////////////////////////////
209 /// Destructor
210 
212 {
213  // Destructor. Delete all owned subgenerator contexts
214  delete _pdfGen ;
215  delete _modelGen ;
216  delete _pdfCloneSet ;
217  delete _modelCloneSet ;
218  delete _modelVars ;
219  delete _pdfVars ;
220  delete _pdfVarsOwned ;
221  delete _modelVarsOwned ;
222 }
223 
224 
225 
226 ////////////////////////////////////////////////////////////////////////////////
227 /// Attach given set of arguments to internal clones of
228 /// pdf and resolution model
229 
231 {
232  // Find convolution variable in input and output sets
235 
236  // Replace all servers in _pdfVars and _modelVars with those in theEvent, except for the convolution variable
237  RooArgSet* pdfCommon = (RooArgSet*) args.selectCommon(*_pdfVars) ;
238  pdfCommon->remove(*cvPdf,kTRUE,kTRUE) ;
239 
240  RooArgSet* modelCommon = (RooArgSet*) args.selectCommon(*_modelVars) ;
241  modelCommon->remove(*cvModel,kTRUE,kTRUE) ;
242 
243  _pdfGen->attach(*pdfCommon) ;
244  _modelGen->attach(*modelCommon) ;
245 
246  delete pdfCommon ;
247  delete modelCommon ;
248 }
249 
250 
251 ////////////////////////////////////////////////////////////////////////////////
252 /// One-time initialization of generator context, attaches
253 /// the context to the supplied event container
254 
256 {
257  // Find convolution variable in input and output sets
260  _cvOut = (RooRealVar*) theEvent.find(_convVarName) ;
261 
262  // Replace all servers in _pdfVars and _modelVars with those in theEvent, except for the convolution variable
263  RooArgSet* pdfCommon = (RooArgSet*) theEvent.selectCommon(*_pdfVars) ;
264  pdfCommon->remove(*_cvPdf,kTRUE,kTRUE) ;
265  _pdfVars->replace(*pdfCommon) ;
266  delete pdfCommon ;
267 
268  RooArgSet* modelCommon = (RooArgSet*) theEvent.selectCommon(*_modelVars) ;
269  modelCommon->remove(*_cvModel,kTRUE,kTRUE) ;
270  _modelVars->replace(*modelCommon) ;
271  delete modelCommon ;
272 
273  // Initialize component generators
276 }
277 
278 
279 
280 ////////////////////////////////////////////////////////////////////////////////
281 /// Generate a single event
282 
284 {
285  while(1) {
286 
287  // Generate pdf and model data
288  _modelGen->generateEvent(*_modelVars,remaining) ;
289  _pdfGen->generateEvent(*_pdfVars,remaining) ;
290 
291  // Construct smeared convolution variable
292  Double_t convValSmeared = _cvPdf->getVal() + _cvModel->getVal() ;
293  if (_cvOut->isValidReal(convValSmeared)) {
294  // Smeared value in acceptance range, transfer values to output set
295  theEvent = *_modelVars ;
296  theEvent = *_pdfVars ;
297  _cvOut->setVal(convValSmeared) ;
298  return ;
299  }
300  }
301 }
302 
303 
304 
305 ////////////////////////////////////////////////////////////////////////////////
306 /// Set the traversal order for events in the prototype dataset
307 /// The argument is a array of integers with a size identical
308 /// to the number of events in the prototype dataset. Each element
309 /// should contain an integer in the range 1-N.
310 
312 {
315  _pdfGen->setProtoDataOrder(lut) ;
316 }
317 
318 
319 ////////////////////////////////////////////////////////////////////////////////
320 /// Print the details of this generator context
321 
323 {
324  RooAbsGenContext::printMultiline(os,content,verbose,indent) ;
325  os << indent << "--- RooConvGenContext ---" << endl ;
326  os << indent << "List of component generators" << endl ;
327 
328  TString indent2(indent) ;
329  indent2.Append(" ") ;
330 
331  _modelGen->printMultiline(os,content,verbose,indent2);
332  _pdfGen->printMultiline(os,content,verbose,indent2);
333 }
RooRealProxy _pdf1
Definition: RooFFTConvPdf.h:69
#define coutE(a)
Definition: RooMsgService.h:35
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents: An owning collection is returned containing clones o...
RooAbsReal & cloneModel() const
static void softAbort()
#define cxcoutI(a)
Definition: RooMsgService.h:84
RooArgSet * _modelCloneSet
RooAbsCollection * selectCommon(const RooAbsCollection &refColl) const
Create a subset of the current collection, consisting only of those elements that are contained as we...
RooAbsReal & clonePdf() const
virtual Bool_t replace(const RooAbsArg &var1, const RooAbsArg &var2)
Replace var1 with var2 and return kTRUE for success.
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Definition: RooAbsArg.h:194
RooRealProxy _pdf2
Definition: RooFFTConvPdf.h:70
const RooRealVar * convVar() const
Return a pointer to the convolution variable instance used in the resolution model.
Basic string class.
Definition: TString.h:137
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
virtual Bool_t isValidReal(Double_t value, Bool_t printError=kFALSE) const
Check if given value is valid.
STL namespace.
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
virtual void setProtoDataOrder(Int_t *lut)
Set the traversal order for events in the prototype dataset The argument is a array of integers with ...
virtual void initGenerator(const RooArgSet &theEvent)
One-time initialization of generator context, attaches the context to the supplied event container...
virtual ~RooConvGenContext()
Destructor.
virtual void attach(const RooArgSet &params)
Attach given set of arguments to internal clones of pdf and resolution model.
RooArgSet * _pdfCloneSet
void removeRange(const char *name=0)
Definition: RooRealVar.h:90
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
Add element to an owning set.
Definition: RooArgSet.cxx:461
RooListProxy _convSet
virtual void generateEvent(RooArgSet &theEvent, Int_t remaining)
Generate a single event.
RooRealVar & cloneVar() const
TString & Append(const char *cs)
Definition: TString.h:492
RooArgSet * _modelVars
virtual Bool_t changeModel(const RooResolutionModel &newModel)
Change the current resolution model to newModel.
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition: TNamed.cxx:63
RooAbsGenContext * _pdfGen
RooRealProxy _x
Definition: RooFFTConvPdf.h:67
virtual void attach(const RooArgSet &params)
Interface to attach given parameters to object in this context.
virtual void changeBasis(RooFormulaVar *basis)
Change the basis function we convolute with.
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Definition: RooRealVar.cxx:202
RooAbsArg * find(const char *name) const
Find object with given name in list.
return
Definition: TBase64.cxx:62
RooAbsGenContext * _modelGen
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.
RooRealVar & convVar() const
Return the convolution variable of the resolution model.
bool verbose
virtual void generateEvent(RooArgSet &theEvent, Int_t remaining)=0
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
static void indent(ostringstream &buf, int indent_level)
RooArgSet * _modelVarsOwned
RooArgSet * _pdfVarsOwned
RooNumConvolution & conv() const
Definition: RooNumConvPdf.h:63
RooAbsArg * at(Int_t idx) const
Definition: RooArgList.h:84
double Double_t
Definition: RtypesCore.h:55
virtual void printMultiline(std::ostream &os, Int_t content, Bool_t verbose=kFALSE, TString indent="") const
Print the details of this generator context.
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Interface for multi-line printing.
virtual void setProtoDataOrder(Int_t *lut)
Set the traversal order of prototype data to that in the lookup tables passed as argument.
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
ClassImp(RooConvGenContext)
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
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.
RooRealVar * _cvModel
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual void initGenerator(const RooArgSet &theEvent)
Interface function to initialize context for generation for given set of observables.
virtual const RooArgSet * get(Int_t index) const
Return RooArgSet with coordinates of event 'index'.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add element to non-owning set.
Definition: RooArgSet.cxx:448
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:1638