Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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
32
33////////////////////////////////////////////////////////////////////////////////
34
35RooGaussian::RooGaussian(const char *name, const char *title,
37 RooAbsReal::Ref _sigma) :
38 RooAbsPdf(name,title),
39 x("x","Observable",this,_x),
40 mean("mean","Mean",this,_mean),
41 sigma("sigma","Width",this,_sigma)
42{
43 RooHelpers::checkRangeOfParameters(this, {&static_cast<RooAbsReal&>(_sigma)}, 0);
44}
45
46////////////////////////////////////////////////////////////////////////////////
47
49 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
50 sigma("sigma",this,other.sigma)
51{
52}
53
54////////////////////////////////////////////////////////////////////////////////
55
60
61
62////////////////////////////////////////////////////////////////////////////////
63/// Compute multiple values of Gaussian distribution.
65{
67 {ctx.at(x), ctx.at(mean), ctx.at(sigma)});
68}
69
70////////////////////////////////////////////////////////////////////////////////
71
72Int_t RooGaussian::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
73{
74 if (matchArgs(allVars,analVars,x)) return 1 ;
75 if (matchArgs(allVars,analVars,mean)) return 2 ;
76 return 0 ;
77}
78
79////////////////////////////////////////////////////////////////////////////////
80
81double RooGaussian::analyticalIntegral(Int_t code, const char* rangeName) const
82{
83 using namespace RooFit::Detail::MathFuncs;
84
85 auto& constant = code == 1 ? mean : x;
86 auto& integrand = code == 1 ? x : mean;
87
88 return gaussianIntegral(integrand.min(rangeName), integrand.max(rangeName), constant, sigma);
89}
90
91////////////////////////////////////////////////////////////////////////////////
92
94{
95 if (matchArgs(directVars,generateVars,x)) return 1 ;
96 if (matchArgs(directVars,generateVars,mean)) return 2 ;
97 return 0 ;
98}
99
100////////////////////////////////////////////////////////////////////////////////
101
103{
104 assert(code==1 || code==2) ;
105 double xgen ;
106 if(code==1){
107 while(true) {
109 if (xgen<x.max() && xgen>x.min()) {
110 x = xgen ;
111 break;
112 }
113 }
114 } else if(code==2){
115 while(true) {
117 if (xgen<mean.max() && xgen>mean.min()) {
118 mean = xgen ;
119 break;
120 }
121 }
122 } else {
123 std::cout << "error in RooGaussian generateEvent"<< std::endl;
124 }
125
126 return;
127}
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
RooRealProxy mean
Definition RooGaussian.h:56
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.
RooRealProxy sigma
Definition RooGaussian.h:57
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:55
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:47
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.
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:85
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.