Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
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 
22 RooConvGenContext is an efficient implementation of the generator context
23 specific for RooAbsAnaConvPdf objects. The physics model is generated
24 with a truth resolution model and the requested resolution model is generated
25 separately as a PDF. The convolution variable of the physics model is
26 subsequently 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 
42 using namespace std;
43 
44 ClassImp(RooConvGenContext);
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 
59 RooConvGenContext::RooConvGenContext(const RooAbsAnaConvPdf &model, const RooArgSet &vars,
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
67  _pdfCloneSet = (RooArgSet*) RooArgSet(model).snapshot(kTRUE) ;
68  if (!_pdfCloneSet) {
69  coutE(Generation) << "RooConvGenContext::RooConvGenContext(" << GetName() << ") Couldn't deep-clone PDF, abort," << endl ;
70  RooErrorHandler::softAbort() ;
71  }
72 
73  RooAbsAnaConvPdf* pdfClone = (RooAbsAnaConvPdf*) _pdfCloneSet->find(model.GetName()) ;
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 ;
86  RooErrorHandler::softAbort() ;
87  }
88  RooResolutionModel* modelClone = (RooResolutionModel*)
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 
129 RooConvGenContext::RooConvGenContext(const RooNumConvPdf &model, const RooArgSet &vars,
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
137  _pdfVarsOwned = (RooArgSet*) model.conv().clonePdf().getObservables(&vars)->snapshot(kTRUE) ;
138  _pdfVars = new RooArgSet(*_pdfVarsOwned) ;
139  _pdfGen = ((RooAbsPdf&)model.conv().clonePdf()).genContext(*_pdfVars,prototype,auxProto,verbose) ;
140  _pdfCloneSet = 0 ;
141 
142  // Create generator for resolution model as PDF
143  _modelVarsOwned = (RooArgSet*) model.conv().cloneModel().getObservables(&vars)->snapshot(kTRUE) ;
144  _modelVars = new RooArgSet(*_modelVarsOwned) ;
145  _convVarName = model.conv().cloneVar().GetName() ;
146  _modelGen = ((RooAbsPdf&)model.conv().cloneModel()).genContext(*_modelVars,prototype,auxProto,verbose) ;
147  _modelCloneSet = new RooArgSet ;
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 
169 RooConvGenContext::RooConvGenContext(const RooFFTConvPdf &model, const RooArgSet &vars,
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) ;
184  _pdfVarsOwned = (RooArgSet*) tmp1->snapshot(kTRUE) ;
185  _pdfVars = new RooArgSet(*_pdfVarsOwned) ;
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) ;
194  _modelVarsOwned = (RooArgSet*) tmp2->snapshot(kTRUE) ;
195  _modelVars = new RooArgSet(*_modelVarsOwned) ;
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 
212 RooConvGenContext::~RooConvGenContext()
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 
231 void RooConvGenContext::attach(const RooArgSet& args)
232 {
233  // Find convolution variable in input and output sets
234  RooRealVar* cvModel = (RooRealVar*) _modelVars->find(_convVarName) ;
235  RooRealVar* cvPdf = (RooRealVar*) _pdfVars->find(_convVarName) ;
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 
256 void RooConvGenContext::initGenerator(const RooArgSet &theEvent)
257 {
258  // Find convolution variable in input and output sets
259  _cvModel = (RooRealVar*) _modelVars->find(_convVarName) ;
260  _cvPdf = (RooRealVar*) _pdfVars->find(_convVarName) ;
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
275  _pdfGen->initGenerator(*_pdfVars) ;
276  _modelGen->initGenerator(*_modelVars) ;
277 }
278 
279 
280 
281 ////////////////////////////////////////////////////////////////////////////////
282 /// Generate a single event
283 
284 void RooConvGenContext::generateEvent(RooArgSet &theEvent, Int_t remaining)
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 
312 void RooConvGenContext::setProtoDataOrder(Int_t* lut)
313 {
314  RooAbsGenContext::setProtoDataOrder(lut) ;
315  _modelGen->setProtoDataOrder(lut) ;
316  _pdfGen->setProtoDataOrder(lut) ;
317 }
318 
319 
320 ////////////////////////////////////////////////////////////////////////////////
321 /// Print the details of this generator context
322 
323 void RooConvGenContext::printMultiline(ostream &os, Int_t content, Bool_t verbose, TString indent) const
324 {
325  RooAbsGenContext::printMultiline(os,content,verbose,indent) ;
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 }