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
22Efficient 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 std::ostream;
43
44
45////////////////////////////////////////////////////////////////////////////////
46/// Constructor for specialized generator context for analytical convolutions.
47///
48/// Builds a generator for the physics PDF convoluted with the truth model
49/// and a generator for the resolution model as PDF. Events are generated
50/// by sampling events from the p.d.f and smearings from the resolution model
51/// and adding these to obtain a distribution of events consistent with the
52/// convolution of these two. The advantage of this procedure is so that
53/// both p.d.f and resolution model can take advantage of any internal
54/// generators that may be defined.
55
57 const RooDataSet *prototype, const RooArgSet* auxProto, bool verbose) :
58 RooAbsGenContext(model,vars,prototype,auxProto,verbose)
59{
60 cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for analytical convolution p.d.f. " << model.GetName()
61 << " for generation of observable(s) " << vars << std::endl ;
62
63 // Clone PDF and change model to internal truth model
64 _pdfCloneSet = std::make_unique<RooArgSet>();
65 RooArgSet(model).snapshot(*_pdfCloneSet, true);
66 if (!_pdfCloneSet) {
67 coutE(Generation) << "RooConvGenContext::RooConvGenContext(" << GetName() << ") Couldn't deep-clone PDF, abort," << std::endl ;
69 }
70
71 RooAbsAnaConvPdf* pdfClone = static_cast<RooAbsAnaConvPdf*>(_pdfCloneSet->find(model.GetName())) ;
72 RooTruthModel truthModel("truthModel","Truth resolution model",*pdfClone->convVar()) ;
73 pdfClone->changeModel(truthModel) ;
74 auto convV = dynamic_cast<RooRealVar*>(pdfClone->convVar());
75 if (!convV) {
76 throw std::runtime_error("RooConvGenContext only works with convolution variables of type RooRealVar.");
77 }
78 convV->removeMin();
79 convV->removeMax();
80
81 // Create generator for physics X truth model
82 _pdfVars = std::unique_ptr<RooArgSet>{pdfClone->getObservables(&vars)};
83 _pdfGen.reset(pdfClone->genContext(*_pdfVars,prototype,auxProto,verbose));
84
85 // Clone resolution model and use as normal PDF
86 _modelCloneSet = std::make_unique<RooArgSet>();
87 RooArgSet(*model._convSet.at(0)).snapshot(*_modelCloneSet, true);
88 if (!_modelCloneSet) {
89 coutE(Generation) << "RooConvGenContext::RooConvGenContext(" << GetName() << ") Couldn't deep-clone resolution model, abort," << std::endl;
91 }
92 std::unique_ptr<RooResolutionModel> modelClone{static_cast<RooResolutionModel*>(_modelCloneSet->find(model._convSet.at(0)->GetName())->Clone("smearing"))};
93 modelClone->changeBasis(nullptr) ;
94 convV = dynamic_cast<RooRealVar*>(&modelClone->convVar());
95 if (!convV) {
96 throw std::runtime_error("RooConvGenContext only works with convolution variables of type RooRealVar.");
97 }
98 convV->removeMin();
99 convV->removeMax();
100
101 // Create generator for resolution model as PDF
102 _modelVars = std::unique_ptr<RooArgSet>{modelClone->getObservables(&vars)};
103
104 _modelVars->add(modelClone->convVar()) ;
105 _convVarName = modelClone->convVar().GetName() ;
106 _modelGen.reset(modelClone->genContext(*_modelVars,prototype,auxProto,verbose));
107
108 _modelCloneSet->addOwned(std::move(modelClone));
109
110 if (prototype) {
111 _pdfVars->add(*prototype->get()) ;
112 _modelVars->add(*prototype->get()) ;
113 }
114
115 // WVE ADD FOR DEBUGGING
116 if (auxProto) {
117 _pdfVars->add(*auxProto) ;
118 _modelVars->add(*auxProto) ;
119 }
120}
121
122
123
124////////////////////////////////////////////////////////////////////////////////
125/// Constructor for specialized generator context for numerical convolutions.
126///
127/// Builds a generator for the physics PDF convoluted with the truth model
128/// and a generator for the resolution model as PDF. Events are generated
129/// by sampling events from the p.d.f and smearings from the resolution model
130/// and adding these to obtain a distribution of events consistent with the
131/// convolution of these two. The advantage of this procedure is so that
132/// both p.d.f and resolution model can take advantage of any internal
133/// generators that may be defined.
134
136 const RooDataSet *prototype, const RooArgSet* auxProto, bool verbose) :
137 RooAbsGenContext(model,vars,prototype,auxProto,verbose)
138{
139 cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for numeric convolution p.d.f. " << model.GetName()
140 << " for generation of observable(s) " << vars << std::endl ;
141
142 // Create generator for physics X truth model
143 {
146 _pdfVarsOwned = std::make_unique<RooArgSet>();
147 clonedPdfObservables.snapshot(*_pdfVarsOwned, true);
148 }
149 _pdfVars = std::make_unique<RooArgSet>(*_pdfVarsOwned) ;
150 _pdfGen.reset(static_cast<RooAbsPdf&>(model.conv().clonePdf()).genContext(*_pdfVars,prototype,auxProto,verbose));
151
152 // Create generator for resolution model as PDF
153 {
156 _modelVarsOwned = std::make_unique<RooArgSet>();
158 }
159 _modelVars = std::make_unique<RooArgSet>(*_modelVarsOwned) ;
160 _convVarName = model.conv().cloneVar().GetName() ;
161 _modelGen.reset(static_cast<RooAbsPdf&>(model.conv().cloneModel()).genContext(*_modelVars,prototype,auxProto,verbose));
162 _modelCloneSet = std::make_unique<RooArgSet>();
163 _modelCloneSet->add(model.conv().cloneModel()) ;
164
165 if (prototype) {
166 _pdfVars->add(*prototype->get()) ;
167 _modelVars->add(*prototype->get()) ;
168 }
169}
170
171
172
173////////////////////////////////////////////////////////////////////////////////
174/// Constructor for specialized generator context for FFT numerical convolutions.
175///
176/// Builds a generator for the physics PDF convoluted with the truth model
177/// and a generator for the resolution model as PDF. Events are generated
178/// by sampling events from the p.d.f and smearings from the resolution model
179/// and adding these to obtain a distribution of events consistent with the
180/// convolution of these two. The advantage of this procedure is so that
181/// both p.d.f and resolution model can take advantage of any internal
182/// generators that may be defined.
183
185 const RooDataSet *prototype, const RooArgSet* auxProto, bool verbose) :
186 RooAbsGenContext(model,vars,prototype,auxProto,verbose)
187{
188 cxcoutI(Generation) << "RooConvGenContext::ctor() setting up special generator context for fft convolution p.d.f. " << model.GetName()
189 << " for generation of observable(s) " << vars << std::endl ;
190
191 _convVarName = model._x.arg().GetName() ;
192
193 // Create generator for physics model
194 _pdfCloneSet = std::make_unique<RooArgSet>();
195 RooArgSet(model._pdf1.arg()).snapshot(*_pdfCloneSet, true);
196 RooAbsPdf* pdfClone = static_cast<RooAbsPdf*>(_pdfCloneSet->find(model._pdf1.arg().GetName())) ;
197 RooRealVar* cvPdf = static_cast<RooRealVar*>(_pdfCloneSet->find(model._x.arg().GetName())) ;
198 cvPdf->removeMin() ;
199 cvPdf->removeMax() ;
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 = static_cast<RooAbsPdf*>(_modelCloneSet->find(model._pdf2.arg().GetName())) ;
211 RooRealVar* cvModel = static_cast<RooRealVar*>(_modelCloneSet->find(model._x.arg().GetName())) ;
212 cvModel->removeMin() ;
213 cvModel->removeMax() ;
215 modelClone->getObservables(&vars, tmp2) ;
216 _modelVarsOwned = std::make_unique<RooArgSet>();
217 tmp2.snapshot(*_modelVarsOwned, true);
218 _modelVars = std::make_unique<RooArgSet>(*_modelVarsOwned) ;
219 _modelGen.reset(modelClone->genContext(*_pdfVars,prototype,auxProto,verbose));
220
221 if (prototype) {
222 _pdfVars->add(*prototype->get()) ;
223 _modelVars->add(*prototype->get()) ;
224 }
225}
226
227
228////////////////////////////////////////////////////////////////////////////////
229/// Attach given set of arguments to internal clones of
230/// pdf and resolution model
231
233{
234 // Find convolution variable in input and output sets
235 auto* cvModel = static_cast<RooRealVar*>(_modelVars->find(_convVarName));
236 auto* cvPdf = static_cast<RooRealVar*>(_pdfVars->find(_convVarName));
237
238 // Replace all servers in _pdfVars and _modelVars with those in theEvent, except for the convolution variable
239 std::unique_ptr<RooArgSet> pdfCommon{args.selectCommon(*_pdfVars)};
240 pdfCommon->remove(*cvPdf,true,true) ;
241
242 std::unique_ptr<RooArgSet> modelCommon{args.selectCommon(*_modelVars)};
243 modelCommon->remove(*cvModel,true,true) ;
244
245 _pdfGen->attach(*pdfCommon) ;
246 _modelGen->attach(*modelCommon) ;
247}
248
249
250////////////////////////////////////////////////////////////////////////////////
251/// One-time initialization of generator context, attaches
252/// the context to the supplied event container
253
255{
256 // Find convolution variable in input and output sets
257 _cvModel = static_cast<RooRealVar*>(_modelVars->find(_convVarName)) ;
258 _cvPdf = static_cast<RooRealVar*>(_pdfVars->find(_convVarName)) ;
259 _cvOut = static_cast<RooRealVar*>(theEvent.find(_convVarName)) ;
260
261 // Replace all servers in _pdfVars and _modelVars with those in theEvent, except for the convolution variable
262 std::unique_ptr<RooArgSet> pdfCommon{theEvent.selectCommon(*_pdfVars)};
263 pdfCommon->remove(*_cvPdf,true,true) ;
264 _pdfVars->replace(*pdfCommon) ;
265
266 std::unique_ptr<RooArgSet> modelCommon{theEvent.selectCommon(*_modelVars)};
267 modelCommon->remove(*_cvModel,true,true) ;
268 _modelVars->replace(*modelCommon) ;
269
270 // Initialize component generators
271 _pdfGen->initGenerator(*_pdfVars) ;
272 _modelGen->initGenerator(*_modelVars) ;
273}
274
275
276
277////////////////////////////////////////////////////////////////////////////////
278/// Generate a single event
279
281{
282 while(true) {
283
284 // Generate pdf and model data
285 _modelGen->generateEvent(*_modelVars,remaining) ;
286 _pdfGen->generateEvent(*_pdfVars,remaining) ;
287
288 // Construct smeared convolution variable
289 double convValSmeared = _cvPdf->getVal() + _cvModel->getVal() ;
291 // Smeared value in acceptance range, transfer values to output set
292 theEvent.assign(*_modelVars) ;
293 theEvent.assign(*_pdfVars) ;
295 return ;
296 }
297 }
298}
299
300
301
302////////////////////////////////////////////////////////////////////////////////
303/// Set the traversal order for events in the prototype dataset
304/// The argument is a array of integers with a size identical
305/// to the number of events in the prototype dataset. Each element
306/// should contain an integer in the range 1-N.
307
309{
311 _modelGen->setProtoDataOrder(lut) ;
312 _pdfGen->setProtoDataOrder(lut) ;
313}
314
315
316////////////////////////////////////////////////////////////////////////////////
317/// Print the details of this generator context
318
319void RooConvGenContext::printMultiline(ostream &os, Int_t content, bool verbose, TString indent) const
320{
322 os << indent << "--- RooConvGenContext ---" << std::endl ;
323 os << indent << "List of component generators" << std::endl ;
324
326 indent2.Append(" ") ;
327
328 _modelGen->printMultiline(os,content,verbose,indent2);
329 _pdfGen->printMultiline(os,content,verbose,indent2);
330}
#define cxcoutI(a)
#define coutE(a)
static void indent(ostringstream &buf, int indent_level)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Base class for PDFs that represent a physics model that can be analytically convolved with a resoluti...
RooListProxy _convSet
Set of (resModel (x) basisFunc) convolution objects.
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.
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.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
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:107
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:159
RooArgSet * selectCommon(const RooAbsCollection &refColl) const
Use RooAbsCollection::selecCommon(), but return as RooArgSet.
Definition RooArgSet.h:154
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.
Container class to hold unbinned data.
Definition RooDataSet.h:32
static void softAbort()
Soft abort function that interrupts macro execution but doesn't kill ROOT.
PDF for the numerical (FFT) convolution of two PDFs.
RooRealProxy _pdf1
First input p.d.f.
RooRealProxy _x
Convolution observable.
RooRealProxy _pdf2
Second input p.d.f.
Numeric 1-dimensional convolution operator PDF.
RooNumConvolution & conv() const
RooAbsReal & cloneModel() const
RooRealVar & cloneVar() const
RooAbsReal & clonePdf() const
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...
const T & arg() const
Return reference to object held in proxy.
Implements a RooResolution model that corresponds to a delta function.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Basic string class.
Definition TString.h:138