Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooSimGenContext.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 RooSimGenContext.cxx
19\class RooSimGenContext
20\ingroup Roofitcore
21
22Efficient implementation of the generator context
23specific for RooSimultaneous PDFs when generating more than one of the
24component pdfs.
25It runs in two modes:
26- Proto data with category entries are given: An event from the same category as
27in the proto data is created.
28- No proto data: A category is chosen randomly.
29\note This requires that the PDFs are extended, to determine the relative probabilities
30that an event originates from a certain category.
31**/
32
33#include "RooSimGenContext.h"
34#include "RooSimultaneous.h"
35#include "RooRealProxy.h"
36#include "RooDataSet.h"
37#include "Roo1DTable.h"
38#include "RooCategory.h"
39#include "RooMsgService.h"
40#include "RooRandom.h"
41#include "RooGlobalFunc.h"
42
43#include <iostream>
44#include <string>
45
47
48////////////////////////////////////////////////////////////////////////////////
49/// Constructor of specialized generator context for RooSimultaneous p.d.f.s. This
50/// context creates a dedicated context for each component p.d.f.s and delegates
51/// generation of events to the appropriate component generator context
52
54 const RooDataSet *prototype, const RooArgSet* auxProto, bool verbose) :
55 RooAbsGenContext(model,vars,prototype,auxProto,verbose), _pdf(&model)
56{
57 // Determine if we are requested to generate the index category
58 RooAbsCategoryLValue const& idxCat = model.indexCat();
59 RooArgSet pdfVars(vars) ;
60
61 RooArgSet allPdfVars(pdfVars) ;
62 if (prototype) allPdfVars.add(*prototype->get(),true) ;
63
64 RooArgSet catsAmongAllVars;
65 allPdfVars.selectCommon(model.flattenedCatList(), catsAmongAllVars);
66
67 if(catsAmongAllVars.size() != model.flattenedCatList().size()) {
68 oocoutE(_pdf,Generation) << "RooSimGenContext::ctor(" << GetName() << ") ERROR: This context must"
69 << " generate all components of the index category" << std::endl ;
70 _isValid = false ;
71 _numPdf = 0 ;
72 _haveIdxProto = false ;
73 return ;
74 }
75
76 // We must either have the prototype or extended likelihood to determined
77 // the relative fractions of the components
78 _haveIdxProto = prototype ? true : false ;
79 _idxCatName = idxCat.GetName() ;
80 if (!_haveIdxProto && !model.canBeExtended()) {
81 oocoutE(_pdf,Generation) << "RooSimGenContext::ctor(" << GetName() << ") ERROR: Need either extended mode"
82 << " or prototype data to calculate number of events per category" << std::endl ;
83 _isValid = false ;
84 _numPdf = 0 ;
85 return ;
86 }
87
88 // Initialize fraction threshold array (used only in extended mode)
89 _numPdf = model._pdfProxyList.GetSize() ;
90 _fracThresh.resize(_numPdf+1);
91 _fracThresh[0] = 0 ;
92
93 // Generate index category and all registered PDFS
94 _allVarsPdf.add(allPdfVars) ;
95 Int_t i(1) ;
96 for(auto * proxy : static_range_cast<RooRealProxy*>(model._pdfProxyList)) {
97 auto* pdf = static_cast<RooAbsPdf*>(proxy->absArg());
98
99 // Name the context after the associated state and add to list
100 _gcList.emplace_back(pdf->genContext(pdfVars,prototype,auxProto,verbose)) ;
101
102 // Create generator context for this PDF
103 _gcList.back()->SetName(proxy->name()) ;
104
105 _gcIndex.push_back(idxCat.lookupIndex(proxy->name()));
106
107 // Fill fraction threshold array
108 _fracThresh[i] = _fracThresh[i-1] + (_haveIdxProto?0:pdf->expectedEvents(&allPdfVars)) ;
109 i++ ;
110 }
111
112 // Normalize fraction threshold array
113 if (!_haveIdxProto) {
114 for(i=0 ; i<_numPdf ; i++)
116 }
117
118
119 // Clone the index category
120 _idxCatSet = std::make_unique<RooArgSet>();
121 RooArgSet(model.indexCat()).snapshot(*_idxCatSet, true);
122 if (!_idxCatSet) {
123 oocoutE(_pdf,Generation) << "RooSimGenContext::RooSimGenContext(" << GetName() << ") Couldn't deep-clone index category, abort," << std::endl ;
124 throw std::string("RooSimGenContext::RooSimGenContext() Couldn't deep-clone index category, abort") ;
125 }
126
127 _idxCat = static_cast<RooAbsCategoryLValue*>(_idxCatSet->find(model.indexCat().GetName()));
128}
129
131
132////////////////////////////////////////////////////////////////////////////////
133/// Attach the index category clone to the given event buffer
134
136{
137 if (_idxCat->isDerived()) {
139 }
140
141 // Forward initGenerator call to all components
142 for(auto& item : _gcList) {
143 item->attach(args) ;
144 }
145
146}
147
148
149////////////////////////////////////////////////////////////////////////////////
150/// Perform one-time initialization of generator context
151
153{
154 // Attach the index category clone to the event
155 if (_idxCat->isDerived()) {
157 } else {
158 _idxCat = static_cast<RooAbsCategoryLValue*>(theEvent.find(_idxCat->GetName())) ;
159 }
160
161 // Update fractions reflecting possible new parameter values
163
164 // Forward initGenerator call to all components
165 for(auto& item : _gcList) {
166 item->initGenerator(theEvent) ;
167 }
168
169}
170
171
172////////////////////////////////////////////////////////////////////////////////
173/// Create an empty dataset to hold the events that will be generated
174
175RooDataSet* RooSimGenContext::createDataSet(const char* name, const char* title, const RooArgSet& obs)
176{
177
178 // If the observables do not contain the index, make a plain dataset
179 if (!obs.contains(*_idxCat)) {
180 return new RooDataSet(name,title,obs) ;
181 }
182
183 if (!_protoData) {
184 std::map<std::string,RooAbsData*> dmap ;
185 for (const auto& nameIdx : *_idxCat) {
186 RooAbsPdf* slicePdf = _pdf->getPdf(nameIdx.first.c_str());
187 std::unique_ptr<RooArgSet> sliceObs{slicePdf->getObservables(obs)};
188 std::string sliceName = name + ("_slice_" + nameIdx.first);
189 std::string sliceTitle = title + (" (index slice " + nameIdx.first + ")");
190 dmap[nameIdx.first] = new RooDataSet(sliceName,sliceTitle,*sliceObs);
191 }
192 using namespace RooFit;
193 _protoData = std::make_unique<RooDataSet>(name, title, obs, Index(static_cast<RooCategory&>(*_idxCat)), Link(dmap), OwnLinked()) ;
194 }
195
196 return new RooDataSet(*_protoData,name) ;
197}
198
199
200
201
202
203////////////////////////////////////////////////////////////////////////////////
204/// Generate event appropriate for current index state.
205/// The index state is taken either from the prototype
206/// or is generated from the fraction threshold table.
207
209{
210 if (_haveIdxProto) {
211
212 // Lookup pdf from selected prototype index state
213 Int_t gidx(0);
214 Int_t cidx = _idxCat->getCurrentIndex();
215 for (Int_t i=0 ; i<(Int_t)_gcIndex.size() ; i++) {
216 if (_gcIndex[i]==cidx) { gidx = i ; break ; }
217 }
218 RooAbsGenContext* cx = _gcList[gidx].get();
219 if (cx) {
220 cx->generateEvent(theEvent,remaining) ;
221 } else {
222 oocoutW(_pdf,Generation) << "RooSimGenContext::generateEvent: WARNING, no PDF to generate event of type " << cidx << std::endl ;
223 }
224
225
226 } else {
227
228 // Throw a random number and select PDF from fraction threshold table
229 double rand = RooRandom::uniform() ;
230 Int_t i=0 ;
231 for (i=0 ; i<_numPdf ; i++) {
232 if (rand>_fracThresh[i] && rand<_fracThresh[i+1]) {
233 RooAbsGenContext* gen=_gcList[i].get();
234 gen->generateEvent(theEvent,remaining) ;
235 //Write through to sub-categories because they might be written to dataset:
237 return ;
238 }
239 }
240
241 }
242}
243
244
245
246////////////////////////////////////////////////////////////////////////////////
247/// No action needed if we have a proto index
248
250{
251 if (_haveIdxProto) return ;
252
253 // Generate index category and all registered PDFS
254 Int_t i(1) ;
255 for(auto * proxy : static_range_cast<RooRealProxy*>(_pdf->_pdfProxyList)) {
256 auto* pdf = static_cast<RooAbsPdf*>(proxy->absArg());
257
258 // Fill fraction threshold array
259 _fracThresh[i] = _fracThresh[i-1] + (_haveIdxProto?0:pdf->expectedEvents(&_allVarsPdf)) ;
260 i++ ;
261 }
262
263 // Normalize fraction threshold array
264 if (!_haveIdxProto) {
265 for(i=0 ; i<_numPdf ; i++)
267 }
268
269}
270
271
272
273////////////////////////////////////////////////////////////////////////////////
274/// Set the traversal order of the prototype data to that in the
275/// given lookup table. This information is passed to all
276/// component generator contexts
277
279{
281
282 for (auto& item : _gcList) {
283 item->setProtoDataOrder(lut) ;
284 }
285}
286
287
288////////////////////////////////////////////////////////////////////////////////
289/// Detailed printing interface
290
291void RooSimGenContext::printMultiline(std::ostream &os, Int_t content, bool verbose, TString indent) const
292{
293 RooAbsGenContext::printMultiline(os,content,verbose,indent) ;
294 os << indent << "--- RooSimGenContext ---" << std::endl ;
295 os << indent << "Using PDF ";
297 os << indent << "List of component generators" << std::endl ;
298
299 TString indent2(indent) ;
300 indent2.Append(" ") ;
301
302 for (auto& item : _gcList) {
303 item->printMultiline(os,content,verbose,indent2);
304 }
305}
#define oocoutW(o, a)
#define oocoutE(o, a)
int Int_t
Definition RtypesCore.h:45
#define ClassImp(name)
Definition Rtypes.h:377
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition TGX11.cxx:110
bool recursiveRedirectServers(const RooAbsCollection &newServerList, bool mustReplaceAll=false, bool nameChange=false, bool recurseInNewSet=true)
Recursively replace all servers with the new servers in newSet.
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.
void SetName(const char *name) override
Set the name of the TNamed.
virtual bool isDerived() const
Does value or shape of this arg depend on any other arg?
Definition RooAbsArg.h:97
Abstract base class for objects that represent a discrete value that can be set from the outside,...
virtual bool setIndex(value_type index, bool printError=true)=0
Change category state by specifying the index code of the desired state.
virtual value_type getCurrentIndex() const
Return index number of current state.
value_type lookupIndex(const std::string &stateName) const
Find the index number corresponding to the state name.
bool contains(const RooAbsArg &var) const
Check if collection contains an argument with the same name as var.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
RooAbsArg * find(const char *name) const
Find object with given name in list.
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 generateEvent(RooArgSet &theEvent, Int_t remaining)=0
bool _isValid
Is context in valid state?
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:40
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:219
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:191
RooArgSet * selectCommon(const RooAbsCollection &refColl) const
Use RooAbsCollection::selecCommon(), but return as RooArgSet.
Definition RooArgSet.h:186
Object to represent discrete states.
Definition RooCategory.h:28
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'.
virtual void printStream(std::ostream &os, Int_t contents, StyleOption style, TString indent="") const
Print description of object on ostream, printing contents set by contents integer,...
static double uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition RooRandom.cxx:78
Efficient implementation of the generator context specific for RooSimultaneous PDFs when generating m...
RooAbsCategoryLValue * _idxCat
Clone of index category.
void initGenerator(const RooArgSet &theEvent) override
Perform one-time initialization of generator context.
std::unique_ptr< RooDataSet > _protoData
! Prototype dataset
void printMultiline(std::ostream &os, Int_t content, bool verbose=false, TString indent="") const override
Detailed printing interface.
std::unique_ptr< RooArgSet > _idxCatSet
Owner of index category components.
RooArgSet _allVarsPdf
All pdf variables.
void setProtoDataOrder(Int_t *lut) override
Set the traversal order of the prototype data to that in the given lookup table.
void updateFractions()
No action needed if we have a proto index.
std::vector< std::unique_ptr< RooAbsGenContext > > _gcList
List of component generator contexts.
const RooSimultaneous * _pdf
Original PDF.
RooSimGenContext(const RooSimultaneous &model, const RooArgSet &vars, const RooDataSet *prototype=nullptr, const RooArgSet *auxProto=nullptr, bool _verbose=false)
Constructor of specialized generator context for RooSimultaneous p.d.f.s.
RooDataSet * createDataSet(const char *name, const char *title, const RooArgSet &obs) override
Create an empty dataset to hold the events that will be generated.
TString _idxCatName
Name of index category.
void generateEvent(RooArgSet &theEvent, Int_t remaining) override
Generate event appropriate for current index state.
Int_t _numPdf
Number of generated PDFs.
std::vector< double > _fracThresh
[_numPdf] Fraction threshold array
std::vector< int > _gcIndex
Index value corresponding to component.
void attach(const RooArgSet &params) override
Attach the index category clone to the given event buffer.
~RooSimGenContext() override
bool _haveIdxProto
Flag set if generation of index is requested.
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
TList _pdfProxyList
List of PDF proxies (named after applicable category state)
RooArgSet const & flattenedCatList() const
Internal utility function to get a list of all category components for this RooSimultaneous.
RooAbsPdf * getPdf(RooStringView catName) const
Return the p.d.f associated with the given index category name.
const RooAbsCategoryLValue & indexCat() const
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
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:572
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26