Logo ROOT   6.18/05
Reference Guide
RooGaussian.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/** \class RooGaussian
18 \ingroup Roofit
19
20Plain Gaussian p.d.f
21**/
22
23#include "RooFit.h"
24
25#include "Riostream.h"
26#include "Riostream.h"
27#include <math.h>
28
29#include "RooGaussian.h"
30#include "RooAbsReal.h"
31#include "RooRealVar.h"
32#include "RooRandom.h"
33#include "RooMath.h"
34
35using namespace std;
36
38
39////////////////////////////////////////////////////////////////////////////////
40
41RooGaussian::RooGaussian(const char *name, const char *title,
42 RooAbsReal& _x, RooAbsReal& _mean,
43 RooAbsReal& _sigma) :
44 RooAbsPdf(name,title),
45 x("x","Observable",this,_x),
46 mean("mean","Mean",this,_mean),
47 sigma("sigma","Width",this,_sigma)
48{
49}
50
51////////////////////////////////////////////////////////////////////////////////
52
53RooGaussian::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)
56{
57}
58
59////////////////////////////////////////////////////////////////////////////////
60
62{
63 const double arg = x - mean;
64 const double sig = sigma;
65 return exp(-0.5*arg*arg/(sig*sig)) ;
66}
67
68////////////////////////////////////////////////////////////////////////////////
69/// calculate and return the negative log-likelihood of the Poisson
70
72{
73 return RooAbsPdf::getLogVal(set) ;
74// Double_t prob = getVal(set) ;
75// return log(prob) ;
76
77 Double_t arg= x - mean;
78 Double_t sig = sigma ;
79
80 //static const Double_t rootPiBy2 = sqrt(atan2(0.0,-1.0)/2.0);
81 //Double_t extra = -0.5*arg*arg/(sig*sig) - log(2*rootPiBy2*sig) ;
82 Double_t extra = -0.5*arg*arg/(sig*sig) - log(analyticalIntegral(1,0)) ;
83
84 return extra ;
85
86}
87
88////////////////////////////////////////////////////////////////////////////////
89
90Int_t RooGaussian::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
91{
92 if (matchArgs(allVars,analVars,x)) return 1 ;
93 if (matchArgs(allVars,analVars,mean)) return 2 ;
94 return 0 ;
95}
96
97////////////////////////////////////////////////////////////////////////////////
98
99Double_t RooGaussian::analyticalIntegral(Int_t code, const char* rangeName) const
100{
101 assert(code==1 || code==2);
102
103 //The normalisation constant 1./sqrt(2*pi*sigma^2) is left out in evaluate().
104 //Therefore, the integral is scaled up by that amount to make RooFit normalise
105 //correctly.
106 const double resultScale = sqrt(TMath::TwoPi()) * sigma;
107
108 //Here everything is scaled and shifted into a standard normal distribution:
109 const double xscale = TMath::Sqrt2() * sigma;
110 double max = 0.;
111 double min = 0.;
112 if (code == 1){
113 max = (x.max(rangeName)-mean)/xscale;
114 min = (x.min(rangeName)-mean)/xscale;
115 } else { //No == 2 test because of assert
116 max = (mean.max(rangeName)-x)/xscale;
117 min = (mean.min(rangeName)-x)/xscale;
118 }
119
120
121 //Here we go for maximum precision: We compute all integrals in the UPPER
122 //tail of the Gaussian, because erfc has the highest precision there.
123 //Therefore, the different cases for range limits in the negative hemisphere are mapped onto
124 //the equivalent points in the upper hemisphere using erfc(-x) = 2. - erfc(x)
125 const double ecmin = std::erfc(std::abs(min));
126 const double ecmax = std::erfc(std::abs(max));
127
128
129 return resultScale * 0.5 * (
130 min*max < 0.0 ? 2.0 - (ecmin + ecmax)
131 : max <= 0. ? ecmax - ecmin : ecmin - ecmax
132 );
133}
134
135////////////////////////////////////////////////////////////////////////////////
136
137Int_t RooGaussian::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
138{
139 if (matchArgs(directVars,generateVars,x)) return 1 ;
140 if (matchArgs(directVars,generateVars,mean)) return 2 ;
141 return 0 ;
142}
143
144////////////////////////////////////////////////////////////////////////////////
145
147{
148 assert(code==1 || code==2) ;
149 Double_t xgen ;
150 if(code==1){
151 while(1) {
153 if (xgen<x.max() && xgen>x.min()) {
154 x = xgen ;
155 break;
156 }
157 }
158 } else if(code==2){
159 while(1) {
161 if (xgen<mean.max() && xgen>mean.min()) {
162 mean = xgen ;
163 break;
164 }
165 }
166 } else {
167 cout << "error in RooGaussian generateEvent"<< endl;
168 }
169
170 return;
171}
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
#define ClassImp(name)
Definition: Rtypes.h:365
char name[80]
Definition: TGX11.cxx:109
double sqrt(double)
double exp(double)
double log(double)
virtual Double_t getLogVal(const RooArgSet *set=0) const
Return the log of the current value with given normalization An error message is printed if the argum...
Definition: RooAbsPdf.cxx:621
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
Plain Gaussian p.d.f.
Definition: RooGaussian.h:25
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
Definition: RooGaussian.cxx:90
RooRealProxy mean
Definition: RooGaussian.h:45
Double_t getLogVal(const RooArgSet *set) const
calculate and return the negative log-likelihood of the Poisson
Definition: RooGaussian.cxx:71
void generateEvent(Int_t code)
Interface for generation of an event using the algorithm corresponding to the specified code.
RooRealProxy sigma
Definition: RooGaussian.h:46
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Definition: RooGaussian.cxx:99
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition: RooGaussian.cxx:61
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
RooRealProxy x
Definition: RooGaussian.h:44
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition: RooRandom.cxx:54
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:263
double erfc(double x)
Complementary error function.
const Double_t sigma
Double_t x[n]
Definition: legend1.C:17
constexpr Double_t Sqrt2()
Definition: TMath.h:89
constexpr Double_t TwoPi()
Definition: TMath.h:45