Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
RooBinnedGenContext.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 RooBinnedGenContext.cxx
19 \class RooBinnedGenContext
20 \ingroup Roofitcore
21 
22 RooBinnedGenContext is an efficient implementation of the
23 generator context specific for binned pdfs.
24 **/
25 
26 
27 #include "RooFit.h"
28 
29 #include "Riostream.h"
30 
31 
32 #include "RooMsgService.h"
33 #include "RooBinnedGenContext.h"
34 #include "RooAbsPdf.h"
35 #include "RooRealVar.h"
36 #include "RooDataHist.h"
37 #include "RooDataSet.h"
38 #include "RooRandom.h"
39 
40 using namespace std;
41 
42 ClassImp(RooBinnedGenContext);
43 ;
44 
45 
46 ////////////////////////////////////////////////////////////////////////////////
47 /// Constructor
48 
49 RooBinnedGenContext::RooBinnedGenContext(const RooAbsPdf &model, const RooArgSet &vars,
50  const RooDataSet *prototype, const RooArgSet* auxProto,
51  Bool_t verbose) :
52  RooAbsGenContext(model,vars,prototype,auxProto,verbose)
53 {
54  cxcoutI(Generation) << "RooBinnedGenContext::ctor() setting up event special generator context for sum p.d.f. " << model.GetName()
55  << " for generation of observable(s) " << vars ;
56  if (prototype) ccxcoutI(Generation) << " with prototype data for " << *prototype->get() ;
57  if (auxProto && auxProto->getSize()>0) ccxcoutI(Generation) << " with auxiliary prototypes " << *auxProto ;
58  ccxcoutI(Generation) << endl ;
59 
60  // Constructor. Build an array of generator contexts for each product component PDF
61  _pdfSet = (RooArgSet*) RooArgSet(model).snapshot(kTRUE) ;
62  _pdf = (RooAbsPdf*) _pdfSet->find(model.GetName()) ;
63  _pdf->setOperMode(RooAbsArg::ADirty,kTRUE) ;
64 
65  // Fix normalization set of this RooAddPdf
66  if (prototype)
67  {
68  RooArgSet coefNSet(vars) ;
69  coefNSet.add(*prototype->get()) ;
70  _pdf->fixAddCoefNormalization(coefNSet) ;
71  }
72 
73  _pdf->recursiveRedirectServers(*_theEvent) ;
74  _vars = _pdf->getObservables(vars) ;
75 
76  // If pdf has boundary definitions, follow those for the binning
77  RooFIter viter = _vars->fwdIterator() ;
78  RooAbsArg* var ;
79  while((var=viter.next())) {
80  RooRealVar* rvar = dynamic_cast<RooRealVar*>(var) ;
81  if (rvar) {
82  list<Double_t>* binb = model.binBoundaries(*rvar,rvar->getMin(),rvar->getMax()) ;
83  delete binb ;
84  }
85  }
86 
87 
88  // Create empty RooDataHist
89  _hist = new RooDataHist("genData","genData",*_vars) ;
90 
91  _expectedData = kFALSE ;
92 }
93 
94 
95 ////////////////////////////////////////////////////////////////////////////////
96 /// Destructor. Delete all owned subgenerator contexts
97 
98 RooBinnedGenContext::~RooBinnedGenContext()
99 {
100  delete _vars ;
101  delete _pdfSet ;
102  delete _hist ;
103 }
104 
105 
106 
107 ////////////////////////////////////////////////////////////////////////////////
108 /// Attach given set of variables to internal p.d.f. clone
109 
110 void RooBinnedGenContext::attach(const RooArgSet& args)
111 {
112  _pdf->recursiveRedirectServers(args) ;
113 }
114 
115 
116 
117 ////////////////////////////////////////////////////////////////////////////////
118 /// One-time initialization of generator contex. Attach theEvent
119 /// to internal p.d.f clone and forward initialization call to
120 /// the component generators
121 
122 void RooBinnedGenContext::initGenerator(const RooArgSet &theEvent)
123 {
124  _pdf->recursiveRedirectServers(theEvent) ;
125 
126 }
127 
128 
129 ////////////////////////////////////////////////////////////////////////////////
130 
131 void RooBinnedGenContext::setExpectedData(Bool_t flag)
132 {
133  _expectedData = flag ;
134 }
135 
136 
137 ////////////////////////////////////////////////////////////////////////////////
138 
139 RooDataSet *RooBinnedGenContext::generate(Double_t nEvt, Bool_t /*skipInit*/, Bool_t extended)
140 {
141  // Scale to number of events and introduce Poisson fluctuations
142  _hist->reset() ;
143 
144  Double_t nEvents = nEvt ;
145 
146  if (nEvents<=0) {
147  if (!_pdf->canBeExtended()) {
148  coutE(InputArguments) << "RooAbsPdf::generateBinned(" << GetName()
149  << ") ERROR: No event count provided and p.d.f does not provide expected number of events" << endl ;
150  return 0 ;
151  } else {
152  // Don't round in expectedData mode
153  if (_expectedData || extended) {
154  nEvents = _pdf->expectedEvents(_vars) ;
155  } else {
156  nEvents = Int_t(_pdf->expectedEvents(_vars)+0.5) ;
157  }
158  }
159  }
160 
161  // Sample p.d.f. distribution
162  _pdf->fillDataHist(_hist,_vars,1,kTRUE) ;
163 
164  // Output container
165  RooRealVar weight("weight","weight",0,1e9) ;
166  RooArgSet tmp(*_vars) ;
167  tmp.add(weight) ;
168  RooDataSet* wudata = new RooDataSet("wu","wu",tmp,RooFit::WeightVar("weight")) ;
169 
170  vector<int> histOut(_hist->numEntries()) ;
171  Double_t histMax(-1) ;
172  Int_t histOutSum(0) ;
173  for (int i=0 ; i<_hist->numEntries() ; i++) {
174  _hist->get(i) ;
175  if (_expectedData) {
176 
177  // Expected data, multiply p.d.f by nEvents
178  Double_t w=_hist->weight()*nEvents ;
179  wudata->add(*_hist->get(),w) ;
180 
181  } else if (extended) {
182 
183  // Extended mode, set contents to Poisson(pdf*nEvents)
184  Double_t w = RooRandom::randomGenerator()->Poisson(_hist->weight()*nEvents) ;
185  wudata->add(*_hist->get(),w) ;
186 
187  } else {
188 
189  // Regular mode, fill array of weights with Poisson(pdf*nEvents), but to not fill
190  // histogram yet.
191  if (_hist->weight()>histMax) {
192  histMax = _hist->weight() ;
193  }
194  histOut[i] = RooRandom::randomGenerator()->Poisson(_hist->weight()*nEvents) ;
195  histOutSum += histOut[i] ;
196  }
197  }
198 
199 
200  if (!_expectedData && !extended) {
201 
202  // Second pass for regular mode - Trim/Extend dataset to exact number of entries
203 
204  // Calculate difference between what is generated so far and what is requested
205  Int_t nEvtExtra = abs(Int_t(nEvents)-histOutSum) ;
206  Int_t wgt = (histOutSum>nEvents) ? -1 : 1 ;
207 
208  // Perform simple binned accept/reject procedure to get to exact event count
209  while(nEvtExtra>0) {
210 
211  Int_t ibinRand = RooRandom::randomGenerator()->Integer(_hist->numEntries()) ;
212  _hist->get(ibinRand) ;
213  Double_t ranY = RooRandom::randomGenerator()->Uniform(histMax) ;
214 
215  if (ranY<_hist->weight()) {
216  if (wgt==1) {
217  histOut[ibinRand]++ ;
218  } else {
219  // If weight is negative, prior bin content must be at least 1
220  if (histOut[ibinRand]>0) {
221  histOut[ibinRand]-- ;
222  } else {
223  continue ;
224  }
225  }
226  nEvtExtra-- ;
227  }
228  }
229 
230  // Transfer working array to histogram
231  for (int i=0 ; i<_hist->numEntries() ; i++) {
232  _hist->get(i) ;
233  wudata->add(*_hist->get(),histOut[i]) ;
234  }
235 
236  }
237 
238  return wudata ;
239 
240 }
241 
242 
243 
244 ////////////////////////////////////////////////////////////////////////////////
245 /// this method is not implemented for this context
246 
247 void RooBinnedGenContext::generateEvent(RooArgSet&, Int_t)
248 {
249  assert(0) ;
250 }
251 
252 
253 
254 ////////////////////////////////////////////////////////////////////////////////
255 /// Print the details of the context
256 
257 void RooBinnedGenContext::printMultiline(ostream &os, Int_t content, Bool_t verbose, TString indent) const
258 {
259  RooAbsGenContext::printMultiline(os,content,verbose,indent) ;
260  os << indent << "--- RooBinnedGenContext ---" << endl ;
261  os << indent << "Using PDF ";
262  _pdf->printStream(os,kName|kArgs|kClassName,kSingleLine,indent);
263 }