Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 "RooGaussian.h"
24#include "RooBatchCompute.h"
25#include "RooHelpers.h"
26#include "RooRandom.h"
27
29
30#include <vector>
31
33
34////////////////////////////////////////////////////////////////////////////////
35
36RooGaussian::RooGaussian(const char *name, const char *title,
38 RooAbsReal::Ref _sigma) :
39 RooAbsPdf(name,title),
40 x("x","Observable",this,_x),
41 mean("mean","Mean",this,_mean),
42 sigma("sigma","Width",this,_sigma)
43{
44 RooHelpers::checkRangeOfParameters(this, {&static_cast<RooAbsReal&>(_sigma)}, 0);
45}
46
47////////////////////////////////////////////////////////////////////////////////
48
49RooGaussian::RooGaussian(const RooGaussian& other, const char* name) :
50 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
51 sigma("sigma",this,other.sigma)
52{
53}
54
55////////////////////////////////////////////////////////////////////////////////
56
58{
60}
61
62
63////////////////////////////////////////////////////////////////////////////////
64/// Compute multiple values of Gaussian distribution.
66{
68 {ctx.at(x), ctx.at(mean), ctx.at(sigma)});
69}
70
71////////////////////////////////////////////////////////////////////////////////
72
73Int_t RooGaussian::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
74{
75 if (matchArgs(allVars,analVars,x)) return 1 ;
76 if (matchArgs(allVars,analVars,mean)) return 2 ;
77 return 0 ;
78}
79
80////////////////////////////////////////////////////////////////////////////////
81
82double RooGaussian::analyticalIntegral(Int_t code, const char* rangeName) const
83{
84 using namespace RooFit::Detail::MathFuncs;
85
86 auto& constant = code == 1 ? mean : x;
87 auto& integrand = code == 1 ? x : mean;
88
89 return gaussianIntegral(integrand.min(rangeName), integrand.max(rangeName), constant, sigma);
90}
91
92////////////////////////////////////////////////////////////////////////////////
93
94Int_t RooGaussian::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool /*staticInitOK*/) const
95{
96 if (matchArgs(directVars,generateVars,x)) return 1 ;
97 if (matchArgs(directVars,generateVars,mean)) return 2 ;
98 return 0 ;
99}
100
101////////////////////////////////////////////////////////////////////////////////
102
104{
105 assert(code==1 || code==2) ;
106 double xgen ;
107 if(code==1){
108 while(true) {
110 if (xgen<x.max() && xgen>x.min()) {
111 x = xgen ;
112 break;
113 }
114 }
115 } else if(code==2){
116 while(true) {
118 if (xgen<mean.max() && xgen>mean.min()) {
119 mean = xgen ;
120 break;
121 }
122 }
123 } else {
124 std::cout << "error in RooGaussian generateEvent"<< std::endl;
125 }
126
127 return;
128}
129
130////////////////////////////////////////////////////////////////////////////////
131
133{
134 // Build a call to the stateless gaussian defined later.
135 ctx.addResult(this, ctx.buildCall("RooFit::Detail::MathFuncs::gaussian", x, mean, sigma));
136}
137
138////////////////////////////////////////////////////////////////////////////////
139
140std::string RooGaussian::buildCallToAnalyticIntegral(Int_t code, const char *rangeName,
142{
143 auto& constant = code == 1 ? mean : x;
144 auto& integrand = code == 1 ? x : mean;
145
146 return ctx.buildCall("RooFit::Detail::MathFuncs::gaussianIntegral",
147 integrand.min(rangeName), integrand.max(rangeName), constant, sigma);
148}
#define ClassImp(name)
Definition Rtypes.h:382
char name[80]
Definition TGX11.cxx:110
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
A RooAbsReal::Ref can be constructed from a RooAbsReal& or a double that will be implicitly converted...
Definition RooAbsReal.h:68
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
bool 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:24
A class to maintain the context for squashing of RooFit models into code.
std::string buildCall(std::string const &funcname, Args_t const &...args)
Build the code to call the function with name funcname, passing some arguments.
void addResult(RooAbsArg const *key, std::string const &value)
A function to save an expression that includes/depends on the result of the input node.
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
std::string buildCallToAnalyticIntegral(Int_t code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const override
This function defines the analytical integral translation for the class.
RooRealProxy mean
Definition RooGaussian.h:60
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
void generateEvent(Int_t code) override
Interface for generation of an event using the algorithm corresponding to the specified code.
void translate(RooFit::Detail::CodeSquashContext &ctx) const override
This function defines a translation for each RooAbsReal based object that can be used to express the ...
RooRealProxy sigma
Definition RooGaussian.h:61
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const override
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
void doEval(RooFit::EvalContext &) const override
Compute multiple values of Gaussian distribution.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooRealProxy x
Definition RooGaussian.h:59
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:48
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:275
const Double_t sigma
Double_t x[n]
Definition legend1.C:17
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})
double gaussian(double x, double mean, double sigma)
Function to evaluate an un-normalized RooGaussian.
Definition MathFuncs.h:86
void checkRangeOfParameters(const RooAbsReal *callingClass, std::initializer_list< const RooAbsReal * > pars, double min=-std::numeric_limits< double >::max(), double max=std::numeric_limits< double >::max(), bool limitsInAllowedRange=false, std::string const &extraMessage="")
Check if the parameters have a range, and warn if the range extends below / above the set limits.