00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024 #include "RooFit.h"
00025
00026 #include "Riostream.h"
00027 #include <math.h>
00028
00029 #include "RooUniform.h"
00030 #include "RooAbsReal.h"
00031 #include "RooRealVar.h"
00032 #include "RooRandom.h"
00033 #include "RooMath.h"
00034 #include "RooArgSet.h"
00035
00036 ClassImp(RooUniform)
00037
00038
00039
00040 RooUniform::RooUniform(const char *name, const char *title, const RooArgSet& _x) :
00041 RooAbsPdf(name,title),
00042 x("x","Observables",this,kTRUE,kFALSE)
00043 {
00044 x.add(_x) ;
00045 }
00046
00047
00048
00049
00050 RooUniform::RooUniform(const RooUniform& other, const char* name) :
00051 RooAbsPdf(other,name), x("x",this,other.x)
00052 {
00053 }
00054
00055
00056
00057
00058 Double_t RooUniform::evaluate() const
00059 {
00060 return 1 ;
00061 }
00062
00063
00064
00065
00066 Int_t RooUniform::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* ) const
00067 {
00068
00069
00070 Int_t nx = x.getSize() ;
00071 if (nx>31) {
00072
00073 coutW(Integration) << "RooUniform::getAnalyticalIntegral(" << GetName() << ") WARNING: p.d.f. has " << x.getSize()
00074 << " observables, analytical integration is only implemented for the first 31 observables" << endl ;
00075 nx=31 ;
00076 }
00077
00078 Int_t code(0) ;
00079 for (int i=0 ; i<x.getSize() ; i++) {
00080 if (allVars.find(x.at(i)->GetName())) {
00081 code |= (1<<i) ;
00082 analVars.add(*allVars.find(x.at(i)->GetName())) ;
00083 }
00084 }
00085 return code ;
00086 }
00087
00088
00089
00090
00091 Double_t RooUniform::analyticalIntegral(Int_t code, const char* rangeName) const
00092 {
00093
00094 Double_t ret(1) ;
00095 for (int i=0 ; i<32 ; i++) {
00096 if (code&(1<<i)) {
00097 RooAbsRealLValue* var = (RooAbsRealLValue*)x.at(i) ;
00098 ret *= (var->getMax(rangeName) - var->getMin(rangeName)) ;
00099 }
00100 }
00101 return ret ;
00102 }
00103
00104
00105
00106
00107
00108 Int_t RooUniform::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t ) const
00109 {
00110
00111
00112 Int_t nx = x.getSize() ;
00113 if (nx>31) {
00114
00115 coutW(Integration) << "RooUniform::getGenerator(" << GetName() << ") WARNING: p.d.f. has " << x.getSize()
00116 << " observables, internal integrator is only implemented for the first 31 observables" << endl ;
00117 nx=31 ;
00118 }
00119
00120 Int_t code(0) ;
00121 for (int i=0 ; i<x.getSize() ; i++) {
00122 if (directVars.find(x.at(i)->GetName())) {
00123 code |= (1<<i) ;
00124 generateVars.add(*directVars.find(x.at(i)->GetName())) ;
00125 }
00126 }
00127 return code ;
00128 return 0 ;
00129 }
00130
00131
00132
00133
00134 void RooUniform::generateEvent(Int_t code)
00135 {
00136
00137
00138 for (int i=0 ; i<32 ; i++) {
00139 if (code&(1<<i)) {
00140 RooAbsRealLValue* var = (RooAbsRealLValue*)x.at(i) ;
00141 var->randomize() ;
00142 }
00143 }
00144 }
00145
00146