Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 verbose) :
61 RooAbsGenContext(model,vars,prototype,auxProto,verbose)
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
67 _pdfCloneSet = std::make_unique<RooArgSet>();
68 RooArgSet(model).snapshot(*_pdfCloneSet, true);
69 if (!_pdfCloneSet) {
70 coutE(Generation) << "RooConvGenContext::RooConvGenContext(" << GetName() << ") Couldn't deep-clone PDF, abort," << endl ;
72 }
73
74 RooAbsAnaConvPdf* pdfClone = (RooAbsAnaConvPdf*) _pdfCloneSet->find(model.GetName()) ;
75 RooTruthModel truthModel("truthModel","Truth resolution model",*pdfClone->convVar()) ;
76 pdfClone->changeModel(truthModel) ;
77 auto convV = dynamic_cast<RooRealVar*>(pdfClone->convVar());
78 if (!convV) {
79 throw std::runtime_error("RooConvGenContext only works with convolution variables of type RooRealVar.");
80 }
81 convV->removeRange();
82
83 // Create generator for physics X truth model
84 _pdfVars = std::unique_ptr<RooArgSet>{pdfClone->getObservables(&vars)};
85 _pdfGen.reset(pdfClone->genContext(*_pdfVars,prototype,auxProto,verbose));
86
87 // Clone resolution model and use as normal PDF
88 _modelCloneSet = std::make_unique<RooArgSet>();
89 RooArgSet(*model._convSet.at(0)).snapshot(*_modelCloneSet, true);
90 if (!_modelCloneSet) {
91 coutE(Generation) << "RooConvGenContext::RooConvGenContext(" << GetName() << ") Couldn't deep-clone resolution model, abort," << std::endl;
93 }
94 std::unique_ptr<RooResolutionModel> modelClone{static_cast<RooResolutionModel*>(_modelCloneSet->find(model._convSet.at(0)->GetName())->Clone("smearing"))};
95 modelClone->changeBasis(nullptr) ;
96 convV = dynamic_cast<RooRealVar*>(&modelClone->convVar());
97 if (!convV) {
98 throw std::runtime_error("RooConvGenContext only works with convolution variables of type RooRealVar.");
99 }
100 convV->removeRange();
101
102 // Create generator for resolution model as PDF
103 _modelVars = std::unique_ptr<RooArgSet>{modelClone->getObservables(&vars)};
104
105 _modelVars->add(modelClone->convVar()) ;
106 _convVarName = modelClone->convVar().GetName() ;
107 _modelGen.reset(modelClone->genContext(*_modelVars,prototype,auxProto,verbose));
108
109 _modelCloneSet->addOwned(std::move(modelClone));
110
111 if (prototype) {
112 _pdfVars->add(*prototype->get()) ;
113 _modelVars->add(*prototype->get()) ;
114 }
115
116 // WVE ADD FOR DEBUGGING
117 if (auxProto) {
118 _pdfVars->add(*auxProto) ;
119 _modelVars->add(*auxProto) ;
120 }
121}
122
123
124
125////////////////////////////////////////////////////////////////////////////////
126/// Constructor for specialized generator context for numerical convolutions.
127///
128/// Builds a generator for the physics PDF convoluted with the truth model
129/// and a generator for the resolution model as PDF. Events are generated
130/// by sampling events from the p.d.f and smearings from the resolution model
131/// and adding these to obtain a distribution of events consistent with the
132/// convolution of these two. The advantage of this procedure is so that
133/// both p.d.f and resolution model can take advantage of any internal
134/// generators that may be defined.
135
137 const RooDataSet *prototype, const RooArgSet* auxProto, bool verbose) :
138 RooAbsGenContext(model,vars,prototype,auxProto,verbose)
139{
140 cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for numeric convolution p.d.f. " << model.GetName()
141 << " for generation of observable(s) " << vars << endl ;
142
143 // Create generator for physics X truth model
144 {
145 RooArgSet clonedPdfObservables;
146 model.conv().clonePdf().getObservables(&vars, clonedPdfObservables);
147 _pdfVarsOwned = std::make_unique<RooArgSet>();
148 clonedPdfObservables.snapshot(*_pdfVarsOwned, true);
149 }
150 _pdfVars = std::make_unique<RooArgSet>(*_pdfVarsOwned) ;
151 _pdfGen.reset(static_cast<RooAbsPdf&>(model.conv().clonePdf()).genContext(*_pdfVars,prototype,auxProto,verbose));
152
153 // Create generator for resolution model as PDF
154 {
155 RooArgSet clonedModelObservables;
156 model.conv().cloneModel().getObservables(&vars, clonedModelObservables);
157 _modelVarsOwned = std::make_unique<RooArgSet>();
158 clonedModelObservables.snapshot(*_modelVarsOwned, true);
159 }
160 _modelVars = std::make_unique<RooArgSet>(*_modelVarsOwned) ;
161 _convVarName = model.conv().cloneVar().GetName() ;
162 _modelGen.reset(static_cast<RooAbsPdf&>(model.conv().cloneModel()).genContext(*_modelVars,prototype,auxProto,verbose));
163 _modelCloneSet = std::make_unique<RooArgSet>();
164 _modelCloneSet->add(model.conv().cloneModel()) ;
165
166 if (prototype) {
167 _pdfVars->add(*prototype->get()) ;
168 _modelVars->add(*prototype->get()) ;
169 }
170}
171
172
173
174////////////////////////////////////////////////////////////////////////////////
175/// Constructor for specialized generator context for FFT numerical convolutions.
176///
177/// Builds a generator for the physics PDF convoluted with the truth model
178/// and a generator for the resolution model as PDF. Events are generated
179/// by sampling events from the p.d.f and smearings from the resolution model
180/// and adding these to obtain a distribution of events consistent with the
181/// convolution of these two. The advantage of this procedure is so that
182/// both p.d.f and resolution model can take advantage of any internal
183/// generators that may be defined.
184
186 const RooDataSet *prototype, const RooArgSet* auxProto, bool verbose) :
187 RooAbsGenContext(model,vars,prototype,auxProto,verbose)
188{
189 cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for fft convolution p.d.f. " << model.GetName()
190 << " for generation of observable(s) " << vars << endl ;
191
192 _convVarName = model._x.arg().GetName() ;
193
194 // Create generator for physics model
195 _pdfCloneSet = std::make_unique<RooArgSet>();
196 RooArgSet(model._pdf1.arg()).snapshot(*_pdfCloneSet, true);
197 RooAbsPdf* pdfClone = (RooAbsPdf*) _pdfCloneSet->find(model._pdf1.arg().GetName()) ;
198 RooRealVar* cvPdf = (RooRealVar*) _pdfCloneSet->find(model._x.arg().GetName()) ;
199 cvPdf->removeRange() ;
200 RooArgSet tmp1;
201 pdfClone->getObservables(&vars, tmp1) ;
202 _pdfVarsOwned = std::make_unique<RooArgSet>();
203 tmp1.snapshot(*_pdfVarsOwned, true);
204 _pdfVars = std::make_unique<RooArgSet>(*_pdfVarsOwned) ;
205 _pdfGen.reset(pdfClone->genContext(*_pdfVars,prototype,auxProto,verbose));
206
207 // Create generator for resolution model
208 _modelCloneSet = std::make_unique<RooArgSet>();
209 RooArgSet(model._pdf2.arg()).snapshot(*_modelCloneSet, true);
210 RooAbsPdf* modelClone = (RooAbsPdf*) _modelCloneSet->find(model._pdf2.arg().GetName()) ;
211 RooRealVar* cvModel = (RooRealVar*) _modelCloneSet->find(model._x.arg().GetName()) ;
212 cvModel->removeRange() ;
213 RooArgSet tmp2;
214 modelClone->getObservables(&vars, tmp2) ;
215 _modelVarsOwned = std::make_unique<RooArgSet>();
216 tmp2.snapshot(*_modelVarsOwned, true);
217 _modelVars = std::make_unique<RooArgSet>(*_modelVarsOwned) ;
218 _modelGen.reset(modelClone->genContext(*_pdfVars,prototype,auxProto,verbose));
219
220 if (prototype) {
221 _pdfVars->add(*prototype->get()) ;
222 _modelVars->add(*prototype->get()) ;
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
234 auto* cvModel = static_cast<RooRealVar*>(_modelVars->find(_convVarName));
235 auto* cvPdf = static_cast<RooRealVar*>(_pdfVars->find(_convVarName));
236
237 // Replace all servers in _pdfVars and _modelVars with those in theEvent, except for the convolution variable
238 std::unique_ptr<RooArgSet> pdfCommon{static_cast<RooArgSet*>(args.selectCommon(*_pdfVars))};
239 pdfCommon->remove(*cvPdf,true,true) ;
240
241 std::unique_ptr<RooArgSet> modelCommon{static_cast<RooArgSet*>(args.selectCommon(*_modelVars))};
242 modelCommon->remove(*cvModel,true,true) ;
243
244 _pdfGen->attach(*pdfCommon) ;
245 _modelGen->attach(*modelCommon) ;
246}
247
248
249////////////////////////////////////////////////////////////////////////////////
250/// One-time initialization of generator context, attaches
251/// the context to the supplied event container
252
254{
255 // Find convolution variable in input and output sets
258 _cvOut = (RooRealVar*) theEvent.find(_convVarName) ;
259
260 // Replace all servers in _pdfVars and _modelVars with those in theEvent, except for the convolution variable
261 std::unique_ptr<RooArgSet> pdfCommon{static_cast<RooArgSet*>(theEvent.selectCommon(*_pdfVars))};
262 pdfCommon->remove(*_cvPdf,true,true) ;
263 _pdfVars->replace(*pdfCommon) ;
264
265 std::unique_ptr<RooArgSet> modelCommon{static_cast<RooArgSet*>(theEvent.selectCommon(*_modelVars))};
266 modelCommon->remove(*_cvModel,true,true) ;
267 _modelVars->replace(*modelCommon) ;
268
269 // Initialize component generators
270 _pdfGen->initGenerator(*_pdfVars) ;
271 _modelGen->initGenerator(*_modelVars) ;
272}
273
274
275
276////////////////////////////////////////////////////////////////////////////////
277/// Generate a single event
278
280{
281 while(true) {
282
283 // Generate pdf and model data
284 _modelGen->generateEvent(*_modelVars,remaining) ;
285 _pdfGen->generateEvent(*_pdfVars,remaining) ;
286
287 // Construct smeared convolution variable
288 double convValSmeared = _cvPdf->getVal() + _cvModel->getVal() ;
289 if (_cvOut->isValidReal(convValSmeared)) {
290 // Smeared value in acceptance range, transfer values to output set
291 theEvent.assign(*_modelVars) ;
292 theEvent.assign(*_pdfVars) ;
293 _cvOut->setVal(convValSmeared) ;
294 return ;
295 }
296 }
297}
298
299
300
301////////////////////////////////////////////////////////////////////////////////
302/// Set the traversal order for events in the prototype dataset
303/// The argument is a array of integers with a size identical
304/// to the number of events in the prototype dataset. Each element
305/// should contain an integer in the range 1-N.
306
308{
310 _modelGen->setProtoDataOrder(lut) ;
311 _pdfGen->setProtoDataOrder(lut) ;
312}
313
314
315////////////////////////////////////////////////////////////////////////////////
316/// Print the details of this generator context
317
318void RooConvGenContext::printMultiline(ostream &os, Int_t content, bool verbose, TString indent) const
319{
320 RooAbsGenContext::printMultiline(os,content,verbose,indent) ;
321 os << indent << "--- RooConvGenContext ---" << endl ;
322 os << indent << "List of component generators" << endl ;
323
324 TString indent2(indent) ;
325 indent2.Append(" ") ;
326
327 _modelGen->printMultiline(os,content,verbose,indent2);
328 _pdfGen->printMultiline(os,content,verbose,indent2);
329}
#define cxcoutI(a)
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:377
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 changeModel(const RooResolutionModel &newModel)
Change the current resolution model to newModel.
RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false) const override
Create a generator context for this p.d.f.
RooAbsRealLValue * convVar()
Retrieve the convolution variable.
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
void assign(const RooAbsCollection &other) const
Sets the value, cache and constant attribute of any argument in our set that also appears in the othe...
bool selectCommon(const RooAbsCollection &refColl, RooAbsCollection &outColl) const
Create a subset of the current collection, consisting only of those elements that are contained as we...
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.
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
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 RooAbsGenContext * genContext(const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool verbose=false) const
Interface function to create a generator context from a p.d.f.
bool isValidReal(double value, bool printError=false) const override
Check if given value is valid.
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:178
RooConvGenContext is an efficient implementation of the generator context specific for RooAbsAnaConvP...
std::unique_ptr< RooAbsGenContext > _modelGen
Resolution model generator context.
RooConvGenContext(const RooFFTConvPdf &model, const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool _verbose=false)
Constructor for specialized generator context for FFT numerical convolutions.
std::unique_ptr< RooArgSet > _pdfCloneSet
Owner of PDF clone.
void setProtoDataOrder(Int_t *lut) override
Set the traversal order for events in the prototype dataset The argument is a array of integers with ...
std::unique_ptr< RooArgSet > _pdfVars
Holder of PDF x truth event.
std::unique_ptr< RooArgSet > _modelVars
Holder of resModel event.
RooRealVar * _cvPdf
Convolution variable in PDFxTruth event.
std::unique_ptr< RooArgSet > _modelVarsOwned
Owning version of modelVars ;.
void attach(const RooArgSet &params) override
Attach given set of arguments to internal clones of pdf and resolution model.
std::unique_ptr< RooArgSet > _modelCloneSet
Owner of resModel clone.
RooRealVar * _cvOut
Convolution variable in output event.
RooRealVar * _cvModel
Convolution variable in resModel event.
void generateEvent(RooArgSet &theEvent, Int_t remaining) override
Generate a single event.
void initGenerator(const RooArgSet &theEvent) override
One-time initialization of generator context, attaches the context to the supplied event container.
std::unique_ptr< RooArgSet > _pdfVarsOwned
Owning version of pdfVars ;.
TString _convVarName
Name of convolution variable.
void printMultiline(std::ostream &os, Int_t content, bool verbose=false, TString indent="") const override
Print the details of this generator context.
std::unique_ptr< RooAbsGenContext > _pdfGen
Physics model generator context.
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:57
const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
static void softAbort()
Soft abort function that interrupts macro execution but doesn't kill ROOT.
PDF for the numerical (FFT) convolution of two PDFs.
Numeric 1-dimensional convolution operator PDF.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
void removeRange(const char *name=nullptr)
Remove range limits for binning with given name. Empty name means default range.
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.
RooTruthModel is an implementation of RooResolution model that provides a delta-function resolution m...
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Basic string class.
Definition TString.h:139
TString & Append(const char *cs)
Definition TString.h:576