Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
RooUniform.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * @(#)root/roofit:$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 /** \class RooUniform
18  \ingroup Roofit
19 
20 Flat p.d.f. in N dimensions
21 **/
22 
23 #include "RooFit.h"
24 
25 #include "Riostream.h"
26 #include <math.h>
27 
28 #include "RooUniform.h"
29 #include "RooAbsReal.h"
30 #include "RooRealVar.h"
31 #include "RooRandom.h"
32 #include "RooMath.h"
33 #include "RooArgSet.h"
34 
35 using namespace std;
36 
37 ClassImp(RooUniform);
38 
39 ////////////////////////////////////////////////////////////////////////////////
40 
41 RooUniform::RooUniform(const char *name, const char *title, const RooArgSet& _x) :
42  RooAbsPdf(name,title),
43  x("x","Observables",this,kTRUE,kFALSE)
44 {
45  x.add(_x) ;
46 }
47 
48 ////////////////////////////////////////////////////////////////////////////////
49 
50 RooUniform::RooUniform(const RooUniform& other, const char* name) :
51  RooAbsPdf(other,name), x("x",this,other.x)
52 {
53 }
54 
55 ////////////////////////////////////////////////////////////////////////////////
56 
57 Double_t RooUniform::evaluate() const
58 {
59  return 1 ;
60 }
61 
62 ////////////////////////////////////////////////////////////////////////////////
63 /// Advertise analytical integral
64 
65 Int_t RooUniform::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
66 {
67  Int_t nx = x.getSize() ;
68  if (nx>31) {
69  // Warn that analytical integration is only provided for the first 31 observables
70  coutW(Integration) << "RooUniform::getAnalyticalIntegral(" << GetName() << ") WARNING: p.d.f. has " << x.getSize()
71  << " observables, analytical integration is only implemented for the first 31 observables" << endl ;
72  nx=31 ;
73  }
74 
75  Int_t code(0) ;
76  for (int i=0 ; i<x.getSize() ; i++) {
77  if (allVars.find(x.at(i)->GetName())) {
78  code |= (1<<i) ;
79  analVars.add(*allVars.find(x.at(i)->GetName())) ;
80  }
81  }
82  return code ;
83 }
84 
85 ////////////////////////////////////////////////////////////////////////////////
86 /// Implement analytical integral
87 
88 Double_t RooUniform::analyticalIntegral(Int_t code, const char* rangeName) const
89 {
90  Double_t ret(1) ;
91  for (int i=0 ; i<32 ; i++) {
92  if (code&(1<<i)) {
93  RooAbsRealLValue* var = (RooAbsRealLValue*)x.at(i) ;
94  ret *= (var->getMax(rangeName) - var->getMin(rangeName)) ;
95  }
96  }
97  return ret ;
98 }
99 
100 ////////////////////////////////////////////////////////////////////////////////
101 /// Advertise internal generator
102 
103 Int_t RooUniform::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
104 {
105  Int_t nx = x.getSize() ;
106  if (nx>31) {
107  // Warn that analytical integration is only provided for the first 31 observables
108  coutW(Integration) << "RooUniform::getGenerator(" << GetName() << ") WARNING: p.d.f. has " << x.getSize()
109  << " observables, internal integrator is only implemented for the first 31 observables" << endl ;
110  nx=31 ;
111  }
112 
113  Int_t code(0) ;
114  for (int i=0 ; i<x.getSize() ; i++) {
115  if (directVars.find(x.at(i)->GetName())) {
116  code |= (1<<i) ;
117  generateVars.add(*directVars.find(x.at(i)->GetName())) ;
118  }
119  }
120  return code ;
121  return 0 ;
122 }
123 
124 ////////////////////////////////////////////////////////////////////////////////
125 /// Implement internal generator
126 
127 void RooUniform::generateEvent(Int_t code)
128 {
129  // Fast-track handling of one-observable case
130  if (code==1) {
131  ((RooAbsRealLValue*)x.at(0))->randomize() ;
132  return ;
133  }
134 
135  for (int i=0 ; i<32 ; i++) {
136  if (code&(1<<i)) {
137  RooAbsRealLValue* var = (RooAbsRealLValue*)x.at(i) ;
138  var->randomize() ;
139  }
140  }
141 }