34 using namespace BatchHelpers;
37 ClassImp(RooGaussian);
41 RooGaussian::RooGaussian(
const char *name,
const char *title,
42 RooAbsReal& _x, RooAbsReal& _mean,
44 RooAbsPdf(name,title),
45 x(
"x",
"Observable",this,_x),
46 mean(
"mean",
"Mean",this,_mean),
47 sigma(
"sigma",
"Width",this,_sigma)
53 RooGaussian::RooGaussian(
const RooGaussian& other,
const char* name) :
54 RooAbsPdf(other,name), x(
"x",this,other.x), mean(
"mean",this,other.mean),
55 sigma(
"sigma",this,other.sigma)
61 Double_t RooGaussian::evaluate()
const
63 const double arg = x - mean;
64 const double sig = sigma;
65 return exp(-0.5*arg*arg/(sig*sig));
75 template<
class Tx,
class TMean,
class TSig>
76 void compute(RooSpan<double> output, Tx x, TMean mean, TSig sigma) {
77 const int n = output.size();
79 for (
int i = 0; i < n; ++i) {
80 const double arg = x[i] - mean[i];
81 const double halfBySigmaSq = -0.5 / (sigma[i] * sigma[i]);
83 output[i] = _rf_fast_exp(arg*arg * halfBySigmaSq);
99 RooSpan<double> RooGaussian::evaluateBatch(std::size_t begin, std::size_t batchSize)
const {
100 auto xData = x.getValBatch(begin, batchSize);
101 auto meanData = mean.getValBatch(begin, batchSize);
102 auto sigmaData = sigma.getValBatch(begin, batchSize);
105 const bool batchX = !xData.empty();
106 const bool batchMean = !meanData.empty();
107 const bool batchSigma = !sigmaData.empty();
109 if (!(batchX || batchMean || batchSigma)) {
113 auto output = _batchData.makeWritableBatchUnInit(begin, batchSize);
115 if (batchX && !batchMean && !batchSigma) {
116 compute(output, xData, BracketAdapter<double>(mean), BracketAdapter<double>(sigma));
118 else if (batchX && batchMean && !batchSigma) {
119 compute(output, xData, meanData, BracketAdapter<double>(sigma));
121 else if (batchX && !batchMean && batchSigma) {
122 compute(output, xData, BracketAdapter<double>(mean), sigmaData);
124 else if (batchX && batchMean && batchSigma) {
125 compute(output, xData, meanData, sigmaData);
127 else if (!batchX && batchMean && !batchSigma) {
128 compute(output, BracketAdapter<double>(x), meanData, BracketAdapter<double>(sigma));
130 else if (!batchX && !batchMean && batchSigma) {
131 compute(output, BracketAdapter<double>(x), BracketAdapter<double>(mean), sigmaData);
133 else if (!batchX && batchMean && batchSigma) {
134 compute(output, BracketAdapter<double>(x), meanData, sigmaData);
142 Int_t RooGaussian::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars,
const char* )
const
144 if (matchArgs(allVars,analVars,x))
return 1 ;
145 if (matchArgs(allVars,analVars,mean))
return 2 ;
151 Double_t RooGaussian::analyticalIntegral(Int_t code,
const char* rangeName)
const
153 assert(code==1 || code==2);
158 const double resultScale = sqrt(TMath::TwoPi()) * sigma;
161 const double xscale = TMath::Sqrt2() * sigma;
165 max = (x.max(rangeName)-mean)/xscale;
166 min = (x.min(rangeName)-mean)/xscale;
168 max = (mean.max(rangeName)-x)/xscale;
169 min = (mean.min(rangeName)-x)/xscale;
177 const double ecmin = std::erfc(std::abs(min));
178 const double ecmax = std::erfc(std::abs(max));
181 return resultScale * 0.5 * (
182 min*max < 0.0 ? 2.0 - (ecmin + ecmax)
183 : max <= 0. ? ecmax - ecmin : ecmin - ecmax
189 Int_t RooGaussian::getGenerator(
const RooArgSet& directVars, RooArgSet &generateVars, Bool_t )
const
191 if (matchArgs(directVars,generateVars,x))
return 1 ;
192 if (matchArgs(directVars,generateVars,mean))
return 2 ;
198 void RooGaussian::generateEvent(Int_t code)
200 assert(code==1 || code==2) ;
204 xgen = RooRandom::randomGenerator()->Gaus(mean,sigma);
205 if (xgen<x.max() && xgen>x.min()) {
212 xgen = RooRandom::randomGenerator()->Gaus(x,sigma);
213 if (xgen<mean.max() && xgen>mean.min()) {
219 cout <<
"error in RooGaussian generateEvent"<< endl;