Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooSimSplitGenContext.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 RooSimSplitGenContext.cxx
19\class RooSimSplitGenContext
20\ingroup Roofitcore
21
22Efficient implementation of the generator context
23specific for RooSimultaneous PDFs when generating more than one of the
24component pdfs.
25**/
26
28#include "RooSimultaneous.h"
29#include "RooRealProxy.h"
30#include "RooDataSet.h"
31#include "Roo1DTable.h"
32#include "RooCategory.h"
33#include "RooMsgService.h"
34#include "RooRandom.h"
35#include "RooGlobalFunc.h"
36
37#include <iostream>
38#include <string>
39
40
41
42////////////////////////////////////////////////////////////////////////////////
43/// Constructor of specialized generator context for RooSimultaneous p.d.f.s. This
44/// context creates a dedicated context for each component p.d.f.s and delegates
45/// generation of events to the appropriate component generator context
46
47RooSimSplitGenContext::RooSimSplitGenContext(const RooSimultaneous &model, const RooArgSet &vars, bool verbose, bool autoBinned, const char* binnedTag) :
48 RooAbsGenContext(model,vars,nullptr,nullptr,verbose), _pdf(&model)
49{
50 // Determine if we are requested to generate the index category
51 RooAbsCategoryLValue const& idxCat = model.indexCat();
52 RooArgSet pdfVars(vars) ;
53
55
57 allPdfVars.selectCommon(model.flattenedCatList(), catsAmongAllVars);
58
59 if(catsAmongAllVars.size() != model.flattenedCatList().size()) {
60 oocoutE(_pdf,Generation) << "RooSimSplitGenContext::ctor(" << GetName() << ") ERROR: This context must"
61 << " generate all components of the index category" << std::endl ;
62 _isValid = false ;
63 _numPdf = 0 ;
64 // coverity[UNINIT_CTOR]
65 return ;
66 }
67
68 // We must extended likelihood to determine the relative fractions of the components
69 _idxCatName = idxCat.GetName() ;
70 if (!model.canBeExtended()) {
71 oocoutE(_pdf,Generation) << "RooSimSplitGenContext::RooSimSplitGenContext(" << GetName() << "): All components of the simultaneous PDF "
72 << "must be extended PDFs. Otherwise, it is impossible to calculate the number of events to be generated per component." << std::endl ;
73 _isValid = false ;
74 _numPdf = 0 ;
75 // coverity[UNINIT_CTOR]
76 return ;
77 }
78
79 // Initialize fraction threshold array (used only in extended mode)
80 _numPdf = model._pdfProxyList.GetSize() ;
81 _fracThresh = new double[_numPdf+1] ;
82 _fracThresh[0] = 0 ;
83
84 // Generate index category and all registered PDFS
86 Int_t i(1) ;
88 auto pdf = static_cast<RooAbsPdf*>(proxy->absArg());
89
90 // Create generator context for this PDF
91 std::unique_ptr<RooArgSet> compVars{pdf->getObservables(pdfVars)};
92 RooAbsGenContext* cx = pdf->autoGenContext(*compVars,nullptr,nullptr,verbose,autoBinned,binnedTag) ;
93
94 const auto state = idxCat.lookupIndex(proxy->name());
95
96 cx->SetName(proxy->name()) ;
97 _gcList.push_back(cx) ;
98 _gcIndex.push_back(state);
99
100 // Fill fraction threshold array
101 _fracThresh[i] = _fracThresh[i-1] + pdf->expectedEvents(&allPdfVars) ;
102 i++ ;
103 }
104
105 for(i=0 ; i<_numPdf ; i++) {
107 }
108
109 // Clone the index category
110 if(RooArgSet(model.indexCat()).snapshot(_idxCatSet, true)) {
111 oocoutE(_pdf,Generation) << "RooSimSplitGenContext::RooSimSplitGenContext(" << GetName() << ") Couldn't deep-clone index category, abort," << std::endl ;
112 throw std::string("RooSimSplitGenContext::RooSimSplitGenContext() Couldn't deep-clone index category, abort") ;
113 }
114 _idxCat = static_cast<RooAbsCategoryLValue*>(_idxCatSet.find(model.indexCat().GetName()));
115}
116
117
118
119////////////////////////////////////////////////////////////////////////////////
120/// Destructor. Delete all owned subgenerator contexts
121
123{
124 delete[] _fracThresh ;
125 for (RooAbsGenContext *item : _gcList) {
126 delete item;
127 }
128}
129
130
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 (RooAbsGenContext *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 // Forward initGenerator call to all components
162 for (RooAbsGenContext *item : _gcList) {
163 item->initGenerator(theEvent) ;
164 }
165
166}
167
168
169
170////////////////////////////////////////////////////////////////////////////////
171
173{
174 if(!isValid()) {
175 coutE(Generation) << ClassName() << "::" << GetName() << ": context is not valid" << std::endl;
176 return nullptr;
177 }
178
179
180 // Calculate the expected number of events if necessary
181 if(nEvents <= 0) {
182 nEvents= _expectedEvents;
183 }
184 coutI(Generation) << ClassName() << "::" << GetName() << ":generate: will generate "
185 << nEvents << " events" << std::endl;
186
187 if (_verbose) Print("v") ;
188
189 // Perform any subclass implementation-specific initialization
190 // Can be skipped if this is a rerun with an identical configuration
191 if (!skipInit) {
193 }
194
195 // Generate lookup table from expected event counts
196 std::vector<double> nGen(_numPdf) ;
197 if (extendedMode ) {
198 Int_t i(0) ;
199 for(auto * proxy : static_range_cast<RooRealProxy*>(_pdf->_pdfProxyList)) {
200 RooAbsPdf* pdf=static_cast<RooAbsPdf*>(proxy->absArg()) ;
201 //nGen[i] = Int_t(pdf->expectedEvents(&_allVarsPdf)+0.5) ;
202 nGen[i] = pdf->expectedEvents(&_allVarsPdf) ;
203 i++ ;
204 }
205
206 } else {
207 Int_t i(1) ;
208 _fracThresh[0] = 0 ;
209 for(auto * proxy : static_range_cast<RooRealProxy*>(_pdf->_pdfProxyList)) {
210 RooAbsPdf* pdf=static_cast<RooAbsPdf*>(proxy->absArg()) ;
212 i++ ;
213 }
214 for(i=0 ; i<_numPdf ; i++) {
216 }
217
218 // Determine from that total number of events to be generated for each component
219 double nGenSoFar(0) ;
220 while (nGenSoFar<nEvents) {
221 double rand = RooRandom::uniform() ;
222 i=0 ;
223 for (i=0 ; i<_numPdf ; i++) {
224 if (rand>_fracThresh[i] && rand<_fracThresh[i+1]) {
225 nGen[i]++ ;
226 nGenSoFar++ ;
227 break ;
228 }
229 }
230 }
231 }
232
233
234
235 // Now loop over states
236 std::map<std::string,RooAbsData*> dataMap ;
237 Int_t icomp(0) ;
238 for(auto * proxy : static_range_cast<RooRealProxy*>(_pdf->_pdfProxyList)) {
239
240 // Calculate number of events to generate for this state
241 if (_gcList[icomp]) {
242 dataMap[proxy->GetName()] = _gcList[icomp]->generate(nGen[icomp],skipInit,extendedMode) ;
243 }
244
245 icomp++ ;
246 }
247
248 // Put all datasets together in a composite-store RooDataSet that links and owns the component datasets
250 return hmaster ;
251}
252
253
254
255////////////////////////////////////////////////////////////////////////////////
256/// Forward to components
257
259{
261 elem->setExpectedData(flag) ;
262 }
263}
264
265
266
267////////////////////////////////////////////////////////////////////////////////
268/// this method is empty because it is not used by this context
269
270RooDataSet* RooSimSplitGenContext::createDataSet(const char* , const char* , const RooArgSet& )
271{
272 return nullptr;
273}
274
275
276
277////////////////////////////////////////////////////////////////////////////////
278/// this method is empty because it is not used in this type of context
279
284
285
286
287
288////////////////////////////////////////////////////////////////////////////////
289/// this method is empty because proto datasets are not supported by this context
290
295
296
297////////////////////////////////////////////////////////////////////////////////
298/// Detailed printing interface
299
300void RooSimSplitGenContext::printMultiline(std::ostream &os, Int_t content, bool verbose, TString indent) const
301{
303 os << indent << "--- RooSimSplitGenContext ---" << std::endl ;
304 os << indent << "Using PDF ";
305 _pdf->printStream(os,kName|kArgs|kClassName,kSingleLine,indent);
306}
#define coutI(a)
#define oocoutE(o, 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.
bool recursiveRedirectServers(const RooAbsCollection &newServerList, bool mustReplaceAll=false, bool nameChange=false, bool recurseInNewSet=true)
Recursively replace all servers with the new servers in newSet.
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 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.
RooArgSet _theEvent
Pointer to observable event being generated.
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Interface for multi-line printing.
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
UInt_t _expectedEvents
Number of expected events from extended p.d.f.
bool _verbose
Verbose messaging?
bool _isValid
Is context in valid state?
bool isValid() const
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
virtual double expectedEvents(const RooArgSet *nset) const
Return expected number of events to be used in calculation of extended likelihood.
bool canBeExtended() const
If true, PDF can provide extended likelihood term.
Definition RooAbsPdf.h:218
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Object to represent discrete states.
Definition RooCategory.h:28
Container class to hold unbinned data.
Definition RooDataSet.h:34
static double uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition RooRandom.cxx:77
void setExpectedData(bool) override
Forward to components.
void generateEvent(RooArgSet &theEvent, Int_t remaining) override
this method is empty because it is not used in this type of context
std::vector< RooAbsGenContext * > _gcList
List of component generator contexts.
RooArgSet _allVarsPdf
All pdf variables.
void attach(const RooArgSet &params) override
Attach the index category clone to the given event buffer.
RooSimSplitGenContext(const RooSimultaneous &model, const RooArgSet &vars, bool _verbose=false, bool autoBinned=true, const char *binnedTag="")
Constructor of specialized generator context for RooSimultaneous p.d.f.s.
void initGenerator(const RooArgSet &theEvent) override
Perform one-time initialization of generator context.
double * _fracThresh
fraction thresholds
RooDataSet * createDataSet(const char *name, const char *title, const RooArgSet &obs) override
this method is empty because it is not used by this context
Int_t _numPdf
Number of generated PDFs.
void printMultiline(std::ostream &os, Int_t content, bool verbose=false, TString indent="") const override
Detailed printing interface.
RooArgSet _idxCatSet
Owner of index category components.
const RooSimultaneous * _pdf
Original PDF.
RooDataSet * generate(double nEvents=0, bool skipInit=false, bool extendedMode=false) override
Generate the specified number of events with nEvents>0 and and return a dataset containing the genera...
std::vector< int > _gcIndex
Index value corresponding to component.
void setProtoDataOrder(Int_t *lut) override
this method is empty because proto datasets are not supported by this context
~RooSimSplitGenContext() override
Destructor. Delete all owned subgenerator contexts.
RooAbsCategoryLValue * _idxCat
Clone of index category.
TString _idxCatName
Name of index category.
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.
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:49
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:226
Basic string class.
Definition TString.h:139
RooCmdArg OwnLinked()
RooCmdArg Index(RooCategory &icat)
RooCmdArg Link(const char *state, RooAbsData &data)