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 "Riostream.h"
00028 #include <math.h>
00029
00030 #include "RooGaussian.h"
00031 #include "RooAbsReal.h"
00032 #include "RooRealVar.h"
00033 #include "RooRandom.h"
00034 #include "RooMath.h"
00035
00036 ClassImp(RooGaussian)
00037
00038
00039
00040 RooGaussian::RooGaussian(const char *name, const char *title,
00041 RooAbsReal& _x, RooAbsReal& _mean,
00042 RooAbsReal& _sigma) :
00043 RooAbsPdf(name,title),
00044 x("x","Observable",this,_x),
00045 mean("mean","Mean",this,_mean),
00046 sigma("sigma","Width",this,_sigma)
00047 {
00048 }
00049
00050
00051
00052
00053 RooGaussian::RooGaussian(const RooGaussian& other, const char* name) :
00054 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
00055 sigma("sigma",this,other.sigma)
00056 {
00057 }
00058
00059
00060
00061
00062 Double_t RooGaussian::evaluate() const
00063 {
00064 Double_t arg= x - mean;
00065 Double_t ret =exp(-0.5*arg*arg/(sigma*sigma)) ;
00066
00067 return ret ;
00068 }
00069
00070
00071
00072
00073 Int_t RooGaussian::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* ) const
00074 {
00075 if (matchArgs(allVars,analVars,x)) return 1 ;
00076 if (matchArgs(allVars,analVars,mean)) return 2 ;
00077 return 0 ;
00078 }
00079
00080
00081
00082
00083 Double_t RooGaussian::analyticalIntegral(Int_t code, const char* rangeName) const
00084 {
00085 assert(code==1 || code==2) ;
00086
00087 static const Double_t root2 = sqrt(2.) ;
00088 static const Double_t rootPiBy2 = sqrt(atan2(0.0,-1.0)/2.0);
00089 Double_t xscale = root2*sigma;
00090 Double_t ret = 0;
00091 if(code==1){
00092 ret = rootPiBy2*sigma*(RooMath::erf((x.max(rangeName)-mean)/xscale)-RooMath::erf((x.min(rangeName)-mean)/xscale));
00093
00094 } else if(code==2) {
00095 ret = rootPiBy2*sigma*(RooMath::erf((mean.max(rangeName)-mean)/xscale)-RooMath::erf((mean.min(rangeName)-mean)/xscale));
00096 } else{
00097 cout << "Error in RooGaussian::analyticalIntegral" << endl;
00098 }
00099 return ret ;
00100
00101 }
00102
00103
00104
00105
00106
00107 Int_t RooGaussian::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t ) const
00108 {
00109 if (matchArgs(directVars,generateVars,x)) return 1 ;
00110 if (matchArgs(directVars,generateVars,mean)) return 2 ;
00111 return 0 ;
00112 }
00113
00114
00115
00116
00117 void RooGaussian::generateEvent(Int_t code)
00118 {
00119 assert(code==1 || code==2) ;
00120 Double_t xgen ;
00121 if(code==1){
00122 while(1) {
00123 xgen = RooRandom::randomGenerator()->Gaus(mean,sigma);
00124 if (xgen<x.max() && xgen>x.min()) {
00125 x = xgen ;
00126 break;
00127 }
00128 }
00129 } else if(code==2){
00130 while(1) {
00131 xgen = RooRandom::randomGenerator()->Gaus(x,sigma);
00132 if (xgen<mean.max() && xgen>mean.min()) {
00133 mean = xgen ;
00134 break;
00135 }
00136 }
00137 } else {
00138 cout << "error in RooGaussian generateEvent"<< endl;
00139 }
00140
00141 return;
00142 }
00143
00144