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