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
30
31#include <vector>
32
34
35////////////////////////////////////////////////////////////////////////////////
36
37RooGaussian::RooGaussian(const char *name, const char *title,
39 RooAbsReal::Ref _sigma) :
40 RooAbsPdf(name,title),
41 x("x","Observable",this,_x),
42 mean("mean","Mean",this,_mean),
43 sigma("sigma","Width",this,_sigma)
44{
45 RooHelpers::checkRangeOfParameters(this, {&static_cast<RooAbsReal&>(_sigma)}, 0);
46}
47
48////////////////////////////////////////////////////////////////////////////////
49
50RooGaussian::RooGaussian(const RooGaussian& other, const char* name) :
51 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
52 sigma("sigma",this,other.sigma)
53{
54}
55
56////////////////////////////////////////////////////////////////////////////////
57
59{
61}
62
63
64////////////////////////////////////////////////////////////////////////////////
65/// Compute multiple values of Gaussian distribution.
66void RooGaussian::computeBatch(cudaStream_t* stream, double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
67{
69 dispatch->compute(stream, RooBatchCompute::Gaussian, output, nEvents,
70 {dataMap.at(x), dataMap.at(mean), dataMap.at(sigma)});
71}
72
73////////////////////////////////////////////////////////////////////////////////
74
75Int_t RooGaussian::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
76{
77 if (matchArgs(allVars,analVars,x)) return 1 ;
78 if (matchArgs(allVars,analVars,mean)) return 2 ;
79 return 0 ;
80}
81
82////////////////////////////////////////////////////////////////////////////////
83
84double RooGaussian::analyticalIntegral(Int_t code, const char* rangeName) const
85{
87
88 auto& constant = code == 1 ? mean : x;
89 auto& integrand = code == 1 ? x : mean;
90
91 return gaussianIntegral(integrand.min(rangeName), integrand.max(rangeName), constant, sigma);
92}
93
94////////////////////////////////////////////////////////////////////////////////
95
96Int_t RooGaussian::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool /*staticInitOK*/) const
97{
98 if (matchArgs(directVars,generateVars,x)) return 1 ;
99 if (matchArgs(directVars,generateVars,mean)) return 2 ;
100 return 0 ;
101}
102
103////////////////////////////////////////////////////////////////////////////////
104
106{
107 assert(code==1 || code==2) ;
108 double xgen ;
109 if(code==1){
110 while(1) {
112 if (xgen<x.max() && xgen>x.min()) {
113 x = xgen ;
114 break;
115 }
116 }
117 } else if(code==2){
118 while(1) {
120 if (xgen<mean.max() && xgen>mean.min()) {
121 mean = xgen ;
122 break;
123 }
124 }
125 } else {
126 std::cout << "error in RooGaussian generateEvent"<< std::endl;
127 }
128
129 return;
130}
131
132////////////////////////////////////////////////////////////////////////////////
133
135{
136 // Build a call to the stateless gaussian defined later.
137 ctx.addResult(this, ctx.buildCall("RooFit::Detail::EvaluateFuncs::gaussianEvaluate", x, mean, sigma));
138}
139
140////////////////////////////////////////////////////////////////////////////////
141
142std::string RooGaussian::buildCallToAnalyticIntegral(Int_t code, const char *rangeName,
144{
145 auto& constant = code == 1 ? mean : x;
146 auto& integrand = code == 1 ? x : mean;
147
148 return ctx.buildCall("RooFit::Detail::AnalyticalIntegrals::gaussianIntegral",
149 integrand.min(rangeName), integrand.max(rangeName), constant, sigma);
150}
#define ClassImp(name)
Definition Rtypes.h:375
char name[80]
Definition TGX11.cxx:110
A RooAbsReal::Ref can be constructed from a RooAbsReal& or a double that will be implicitly converted...
Definition RooAbsReal.h:70
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
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:55
virtual void compute(cudaStream_t *, Computer, RestrictArr, size_t, const VarVector &, ArgVector &)=0
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.
auto & at(RooAbsArg const *arg, RooAbsArg const *=nullptr)
Definition DataMap.h:88
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:61
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 computeBatch(cudaStream_t *, double *output, size_t size, RooFit::Detail::DataMap const &) const override
Compute multiple values of Gaussian distribution.
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:62
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...
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:60
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:51
double max(const char *rname=nullptr) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
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:274
const Double_t sigma
Double_t x[n]
Definition legend1.C:17
R__EXTERN RooBatchComputeInterface * dispatchCUDA
R__EXTERN RooBatchComputeInterface * dispatchCPU
This dispatch pointer points to an implementation of the compute library, provided one has been loade...
double gaussianEvaluate(double x, double mean, double sigma)
Function to evaluate an un-normalized RooGaussian.
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.
static void output()