Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
RooAbsGenContext.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 RooAbsGenContext.cxx
19 \class RooAbsGenContext
20 \ingroup Roofitcore
21 
22 RooAbsGenContext is the abstract base class for generator contexts of
23 RooAbsPdf objects. A generator context is an object that controls
24 the generation of events from a given p.d.f in one or more sessions.
25 This class defines the common interface for all such contexts and organizes
26 storage of common components, such as the observables definition, the
27 prototype data etc..
28 **/
29 
30 #include "RooFit.h"
31 
32 #include "TClass.h"
33 
34 #include "RooAbsGenContext.h"
35 #include "RooRandom.h"
36 #include "RooAbsPdf.h"
37 #include "RooDataSet.h"
38 #include "RooMsgService.h"
39 #include "RooGlobalFunc.h"
40 
41 #include "Riostream.h"
42 
43 
44 using namespace std;
45 
46 ClassImp(RooAbsGenContext);
47 ;
48 
49 
50 ////////////////////////////////////////////////////////////////////////////////
51 /// Constructor
52 
53 RooAbsGenContext::RooAbsGenContext(const RooAbsPdf& model, const RooArgSet &vars,
54  const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose) :
55  TNamed(model),
56  _prototype(prototype),
57  _theEvent(0),
58  _isValid(kTRUE),
59  _verbose(verbose),
60  _protoOrder(0),
61  _genData(0)
62 {
63  // Check PDF dependents
64  if (model.recursiveCheckObservables(&vars)) {
65  coutE(Generation) << "RooAbsGenContext::ctor: Error in PDF dependents" << endl ;
66  _isValid = kFALSE ;
67  return ;
68  }
69 
70  // Make a snapshot of the generated variables that we can overwrite.
71  _theEvent= (RooArgSet*)vars.snapshot(kFALSE);
72 
73  // Analyze the prototype dataset, if one is specified
74  _nextProtoIndex= 0;
75  if(0 != _prototype) {
76  TIterator *protoIterator= _prototype->get()->createIterator();
77  const RooAbsArg *proto = 0;
78  while((proto= (const RooAbsArg*)protoIterator->Next())) {
79  // is this variable being generated or taken from the prototype?
80  if(!_theEvent->contains(*proto)) {
81  _protoVars.add(*proto);
82  _theEvent->addClone(*proto);
83  }
84  }
85  delete protoIterator;
86  }
87 
88  // Add auxiliary protovars to _protoVars, if provided
89  if (auxProto) {
90  _protoVars.add(*auxProto) ;
91  _theEvent->addClone(*auxProto) ;
92  }
93 
94  // Remember the default number of events to generate when no prototype dataset is provided.
95  _extendMode = model.extendMode() ;
96  if (model.canBeExtended()) {
97  _expectedEvents= (Int_t)(model.expectedEvents(_theEvent) + 0.5);
98  } else {
99  _expectedEvents= 0 ;
100  }
101 
102  // Save normalization range
103  if (model.normRange()) {
104  _normRange = model.normRange() ;
105  }
106 }
107 
108 
109 
110 ////////////////////////////////////////////////////////////////////////////////
111 /// Destructor
112 
113 RooAbsGenContext::~RooAbsGenContext()
114 {
115  if(0 != _theEvent) delete _theEvent;
116  if (_protoOrder) delete[] _protoOrder ;
117 }
118 
119 
120 
121 ////////////////////////////////////////////////////////////////////////////////
122 /// Interface to attach given parameters to object in this context
123 
124 void RooAbsGenContext::attach(const RooArgSet& /*params*/)
125 {
126 }
127 
128 
129 
130 ////////////////////////////////////////////////////////////////////////////////
131 /// Create an empty dataset to hold the events that will be generated
132 
133 RooDataSet* RooAbsGenContext::createDataSet(const char* name, const char* title, const RooArgSet& obs)
134 {
135  RooDataSet* ret = new RooDataSet(name, title, obs);
136  ret->setDirtyProp(kFALSE) ;
137  return ret ;
138 }
139 
140 
141 ////////////////////////////////////////////////////////////////////////////////
142 /// Generate the specified number of events with nEvents>0 and
143 /// and return a dataset containing the generated events. With nEvents<=0,
144 /// generate the number of events in the prototype dataset, if available,
145 /// or else the expected number of events, if non-zero.
146 /// If extendedMode = true generate according to a Poisson(nEvents)
147 /// The returned dataset belongs to the caller. Return zero in case of an error.
148 /// Generation of individual events is delegated to a virtual generateEvent()
149 /// method. A virtual initGenerator() method is also called just before the
150 /// first call to generateEvent().
151 
152 RooDataSet *RooAbsGenContext::generate(Double_t nEvents, Bool_t skipInit, Bool_t extendedMode)
153 {
154  if(!isValid()) {
155  coutE(Generation) << ClassName() << "::" << GetName() << ": context is not valid" << endl;
156  return 0;
157  }
158 
159  // Calculate the expected number of events if necessary
160  if(nEvents <= 0) {
161  if(_prototype) {
162  nEvents= (Int_t)_prototype->numEntries();
163  }
164  else {
165  if (_extendMode == RooAbsPdf::CanNotBeExtended) {
166  coutE(Generation) << ClassName() << "::" << GetName()
167  << ":generate: PDF not extendable: cannot calculate expected number of events" << endl;
168  return 0;
169  }
170  nEvents= _expectedEvents;
171  }
172  if(nEvents <= 0) {
173  coutE(Generation) << ClassName() << "::" << GetName()
174  << ":generate: cannot calculate expected number of events" << endl;
175  return 0;
176  }
177  coutI(Generation) << ClassName() << "::" << GetName() << ":generate: will generate "
178  << nEvents << " events" << endl;
179 
180  }
181 
182  if (extendedMode) {
183  double nExpEvents = nEvents;
184  nEvents = RooRandom::randomGenerator()->Poisson(nEvents) ;
185  cxcoutI(Generation) << " Extended mode active, number of events generated (" << nEvents << ") is Poisson fluctuation on "
186  << GetName() << "::expectedEvents() = " << nExpEvents << endl ;
187  }
188 
189  // check that any prototype dataset still defines the variables we need
190  // (this is necessary since we never make a private clone, for efficiency)
191  if(_prototype) {
192  const RooArgSet *vars= _prototype->get();
193  TIterator *iterator= _protoVars.createIterator();
194  const RooAbsArg *arg = 0;
195  Bool_t ok(kTRUE);
196  while((arg= (const RooAbsArg*)iterator->Next())) {
197  if(vars->contains(*arg)) continue;
198  coutE(InputArguments) << ClassName() << "::" << GetName() << ":generate: prototype dataset is missing \""
199  << arg->GetName() << "\"" << endl;
200 
201  // WVE disable this for the moment
202  // ok= kFALSE;
203  }
204  delete iterator;
205  // coverity[DEADCODE]
206  if(!ok) return 0;
207  }
208 
209  if (_verbose) Print("v") ;
210 
211  // create a new dataset
212  TString name(GetName()),title(GetTitle());
213  name.Append("Data");
214  title.Prepend("Generated From ");
215 
216  // WVE need specialization here for simultaneous pdfs
217  _genData = createDataSet(name.Data(),title.Data(),*_theEvent) ;
218 
219  // Perform any subclass implementation-specific initialization
220  // Can be skipped if this is a rerun with an identical configuration
221  if (!skipInit) {
222  initGenerator(*_theEvent);
223  }
224 
225  // Loop over the events to generate
226  Int_t evt(0) ;
227  while(_genData->numEntries()<nEvents) {
228 
229  // first, load values from the prototype dataset, if one was provided
230  if(0 != _prototype) {
231  if(_nextProtoIndex >= _prototype->numEntries()) _nextProtoIndex= 0;
232 
233  Int_t actualProtoIdx = _protoOrder ? _protoOrder[_nextProtoIndex] : _nextProtoIndex ;
234 
235  const RooArgSet *subEvent= _prototype->get(actualProtoIdx);
236  _nextProtoIndex++;
237  if(0 != subEvent) {
238  *_theEvent= *subEvent;
239  }
240  else {
241  coutE(Generation) << ClassName() << "::" << GetName() << ":generate: cannot load event "
242  << actualProtoIdx << " from prototype dataset" << endl;
243  return 0;
244  }
245  }
246 
247  // delegate the generation of the rest of this event to our subclass implementation
248  generateEvent(*_theEvent, (Int_t)(nEvents - _genData->numEntries()));
249 
250 
251  // WVE add check that event is in normRange
252  if (_normRange.Length()>0 && !_theEvent->isInRange(_normRange.Data())) {
253  continue ;
254  }
255 
256  _genData->addFast(*_theEvent);
257  evt++ ;
258  }
259 
260  RooDataSet* output = _genData ;
261  _genData = 0 ;
262  output->setDirtyProp(kTRUE) ;
263 
264  return output;
265 }
266 
267 
268 
269 ////////////////////////////////////////////////////////////////////////////////
270 /// Interface function to initialize context for generation for given
271 /// set of observables
272 
273 void RooAbsGenContext::initGenerator(const RooArgSet&)
274 {
275 }
276 
277 
278 
279 ////////////////////////////////////////////////////////////////////////////////
280 /// Print name of context
281 
282 void RooAbsGenContext::printName(ostream& os) const
283 {
284  os << GetName() ;
285 }
286 
287 
288 
289 ////////////////////////////////////////////////////////////////////////////////
290 /// Print title of context
291 
292 void RooAbsGenContext::printTitle(ostream& os) const
293 {
294  os << GetTitle() ;
295 }
296 
297 
298 
299 ////////////////////////////////////////////////////////////////////////////////
300 /// Print class name of context
301 
302 void RooAbsGenContext::printClassName(ostream& os) const
303 {
304  os << IsA()->GetName() ;
305 }
306 
307 
308 
309 ////////////////////////////////////////////////////////////////////////////////
310 /// Print arguments of context, i.e. the observables being generated in this context
311 
312 void RooAbsGenContext::printArgs(ostream& os) const
313 {
314  os << "[ " ;
315  TIterator* iter = _theEvent->createIterator() ;
316  RooAbsArg* arg ;
317  Bool_t first(kTRUE) ;
318  while((arg=(RooAbsArg*)iter->Next())) {
319  if (first) {
320  first=kFALSE ;
321  } else {
322  os << "," ;
323  }
324  os << arg->GetName() ;
325  }
326  os << "]" ;
327  delete iter ;
328 }
329 
330 
331 
332 ////////////////////////////////////////////////////////////////////////////////
333 /// Interface for multi-line printing
334 
335 void RooAbsGenContext::printMultiline(ostream &/*os*/, Int_t /*contents*/, Bool_t /*verbose*/, TString /*indent*/) const
336 {
337 }
338 
339 
340 
341 
342 ////////////////////////////////////////////////////////////////////////////////
343 /// Set the traversal order of prototype data to that in the lookup tables
344 /// passed as argument. The LUT must be an array of integers with the same
345 /// size as the number of entries in the prototype dataset and must contain
346 /// integer values in the range [0,Nevt-1]
347 
348 void RooAbsGenContext::setProtoDataOrder(Int_t* lut)
349 {
350  // Delete any previous lookup table
351  if (_protoOrder) {
352  delete[] _protoOrder ;
353  _protoOrder = 0 ;
354  }
355 
356  // Copy new lookup table if provided and needed
357  if (lut && _prototype) {
358  Int_t n = _prototype->numEntries() ;
359  _protoOrder = new Int_t[n] ;
360  Int_t i ;
361  for (i=0 ; i<n ; i++) {
362  _protoOrder[i] = lut[i] ;
363  }
364  }
365 }
366 
367 
368 
369 
370 ////////////////////////////////////////////////////////////////////////////////
371 /// Rescale existing output buffer with given ratio
372 
373 void RooAbsGenContext::resampleData(Double_t& ratio)
374 {
375 
376  Int_t nOrig = _genData->numEntries() ;
377  Int_t nTarg = Int_t(nOrig*ratio+0.5) ;
378  RooDataSet* trimmedData = (RooDataSet*) _genData->reduce(RooFit::EventRange(0,nTarg)) ;
379 
380  cxcoutD(Generation) << "RooGenContext::resampleData*( existing production trimmed from " << nOrig << " to " << trimmedData->numEntries() << " events" << endl ;
381 
382  delete _genData ;
383  _genData = trimmedData ;
384 
385  if (_prototype) {
386  // Push back proto index by trimmed amount to force recycling of the
387  // proto entries that were trimmed away
388  _nextProtoIndex -= (nOrig-nTarg) ;
389  while (_nextProtoIndex<0) {
390  _nextProtoIndex += _prototype->numEntries() ;
391  }
392  }
393 
394 }
395 
396 
397 
398 
399 ////////////////////////////////////////////////////////////////////////////////
400 /// Define default contents when printing
401 
402 Int_t RooAbsGenContext::defaultPrintContents(Option_t* /*opt*/) const
403 {
404  return kName|kClassName|kValue ;
405 }
406 
407 
408 
409 ////////////////////////////////////////////////////////////////////////////////
410 /// Define default print style
411 
412 RooPrintable::StyleOption RooAbsGenContext::defaultPrintStyle(Option_t* opt) const
413 {
414  if (opt && TString(opt).Contains("v")) {
415  return kVerbose ;
416  }
417  return kStandard ;
418 }