ROOT  6.06/09
Reference Guide
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 //////////////////////////////////////////////////////////////////////////////
18 //
19 // BEGIN_HTML
20 // Flat p.d.f. in N dimensions
21 // END_HTML
22 //
23 
24 #include "RooFit.h"
25 
26 #include "Riostream.h"
27 #include <math.h>
28 
29 #include "RooUniform.h"
30 #include "RooAbsReal.h"
31 #include "RooRealVar.h"
32 #include "RooRandom.h"
33 #include "RooMath.h"
34 #include "RooArgSet.h"
35 
36 using namespace std;
37 
39 
40 
41 ////////////////////////////////////////////////////////////////////////////////
42 
43 RooUniform::RooUniform(const char *name, const char *title, const RooArgSet& _x) :
44  RooAbsPdf(name,title),
45  x("x","Observables",this,kTRUE,kFALSE)
46 {
47  x.add(_x) ;
48 }
49 
50 
51 
52 ////////////////////////////////////////////////////////////////////////////////
53 
54 RooUniform::RooUniform(const RooUniform& other, const char* name) :
55  RooAbsPdf(other,name), x("x",this,other.x)
56 {
57 }
58 
59 
60 
61 ////////////////////////////////////////////////////////////////////////////////
62 
64 {
65  return 1 ;
66 }
67 
68 
69 
70 ////////////////////////////////////////////////////////////////////////////////
71 /// Advertise analytical integral
72 
73 Int_t RooUniform::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
74 {
75  Int_t nx = x.getSize() ;
76  if (nx>31) {
77  // Warn that analytical integration is only provided for the first 31 observables
78  coutW(Integration) << "RooUniform::getAnalyticalIntegral(" << GetName() << ") WARNING: p.d.f. has " << x.getSize()
79  << " observables, analytical integration is only implemented for the first 31 observables" << endl ;
80  nx=31 ;
81  }
82 
83  Int_t code(0) ;
84  for (int i=0 ; i<x.getSize() ; i++) {
85  if (allVars.find(x.at(i)->GetName())) {
86  code |= (1<<i) ;
87  analVars.add(*allVars.find(x.at(i)->GetName())) ;
88  }
89  }
90  return code ;
91 }
92 
93 
94 
95 ////////////////////////////////////////////////////////////////////////////////
96 /// Implement analytical integral
97 
98 Double_t RooUniform::analyticalIntegral(Int_t code, const char* rangeName) const
99 {
100  Double_t ret(1) ;
101  for (int i=0 ; i<32 ; i++) {
102  if (code&(1<<i)) {
103  RooAbsRealLValue* var = (RooAbsRealLValue*)x.at(i) ;
104  ret *= (var->getMax(rangeName) - var->getMin(rangeName)) ;
105  }
106  }
107  return ret ;
108 }
109 
110 
111 
112 
113 ////////////////////////////////////////////////////////////////////////////////
114 /// Advertise internal generator
115 
116 Int_t RooUniform::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
117 {
118  Int_t nx = x.getSize() ;
119  if (nx>31) {
120  // Warn that analytical integration is only provided for the first 31 observables
121  coutW(Integration) << "RooUniform::getGenerator(" << GetName() << ") WARNING: p.d.f. has " << x.getSize()
122  << " observables, internal integrator is only implemented for the first 31 observables" << endl ;
123  nx=31 ;
124  }
125 
126  Int_t code(0) ;
127  for (int i=0 ; i<x.getSize() ; i++) {
128  if (directVars.find(x.at(i)->GetName())) {
129  code |= (1<<i) ;
130  generateVars.add(*directVars.find(x.at(i)->GetName())) ;
131  }
132  }
133  return code ;
134  return 0 ;
135 }
136 
137 
138 
139 ////////////////////////////////////////////////////////////////////////////////
140 /// Implement internal generator
141 
143 {
144  // Fast-track handling of one-observable case
145  if (code==1) {
146  ((RooAbsRealLValue*)x.at(0))->randomize() ;
147  return ;
148  }
149 
150  for (int i=0 ; i<32 ; i++) {
151  if (code&(1<<i)) {
152  RooAbsRealLValue* var = (RooAbsRealLValue*)x.at(i) ;
153  var->randomize() ;
154  }
155  }
156 }
157 
158 
const int nx
Definition: kalman.C:16
void generateEvent(Int_t code)
Implement internal generator.
Definition: RooUniform.cxx:142
virtual void randomize(const char *rangeName=0)
Set a new value sampled from a uniform distribution over the fit range.
virtual Double_t getMin(const char *name=0) const
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
STL namespace.
#define coutW(a)
Definition: RooMsgService.h:34
Double_t x[n]
Definition: legend1.C:17
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implement analytical integral.
Definition: RooUniform.cxx:98
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Advertise analytical integral.
Definition: RooUniform.cxx:73
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Advertise internal generator.
Definition: RooUniform.cxx:116
RooListProxy x
Definition: RooUniform.h:40
RooAbsArg * find(const char *name) const
Find object with given name in list.
return
Definition: TBase64.cxx:62
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
RooAbsArg * at(Int_t idx) const
Definition: RooArgList.h:84
double Double_t
Definition: RtypesCore.h:55
#define name(a, b)
Definition: linkTestLib0.cpp:5
virtual Double_t getMax(const char *name=0) const
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Int_t getSize() const
Double_t evaluate() const
Definition: RooUniform.cxx:63
const Bool_t kTRUE
Definition: Rtypes.h:91
ClassImp(RooUniform) RooUniform
Definition: RooUniform.cxx:38
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add element to non-owning set.
Definition: RooArgSet.cxx:448