Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooGaussian.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * File: $Id: RooGaussian.h,v 1.16 2007/07/12 20:30:49 wouter Exp $
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#ifndef ROO_GAUSSIAN
17#define ROO_GAUSSIAN
18
19#include "RooAbsPdf.h"
20#include "RooRealProxy.h"
21
22class RooAbsReal;
23
24class RooGaussian : public RooAbsPdf {
25public:
27 // Original constructor without RooAbsReal::Ref for backwards compatibility.
28 inline RooGaussian(const char *name, const char *title,
29 RooAbsReal& _x, RooAbsReal& _mean, RooAbsReal& _sigma)
30 : RooGaussian{name, title, RooAbsReal::Ref{_x}, RooAbsReal::Ref{_mean}, RooAbsReal::Ref{_sigma}} {}
31 RooGaussian(const char *name, const char *title,
33 RooGaussian(const RooGaussian& other, const char* name=nullptr);
34 TObject* clone(const char* newname) const override {
35 return new RooGaussian(*this,newname);
36 }
37
38 Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=nullptr) const override;
39 double analyticalIntegral(Int_t code, const char* rangeName=nullptr) const override;
40
41 Int_t getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool staticInitOK=true) const override;
42 void generateEvent(Int_t code) override;
43
44 /// Get the x variable.
45 RooAbsReal const& getX() const { return x.arg(); }
46
47 /// Get the mean parameter.
48 RooAbsReal const& getMean() const { return mean.arg(); }
49
50 /// Get the sigma parameter.
51 RooAbsReal const& getSigma() const { return sigma.arg(); }
52
53 void translate(RooFit::Detail::CodeSquashContext &ctx) const override;
54 std::string
55 buildCallToAnalyticIntegral(Int_t code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const override;
56
57protected:
58
62
63 double evaluate() const override;
64 void doEval(RooFit::EvalContext &) const override;
65 inline bool canComputeBatchWithCuda() const override { return true; }
66
67private:
68
69 ClassDefOverride(RooGaussian,1) // Gaussian PDF
70};
71
72#endif
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
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
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.
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.
TObject * clone(const char *newname) const override
Definition RooGaussian.h:34
RooGaussian(const char *name, const char *title, RooAbsReal &_x, RooAbsReal &_mean, RooAbsReal &_sigma)
Definition RooGaussian.h:28
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.
RooAbsReal const & getX() const
Get the x variable.
Definition RooGaussian.h:45
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.
RooAbsReal const & getMean() const
Get the mean parameter.
Definition RooGaussian.h:48
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.
bool canComputeBatchWithCuda() const override
Definition RooGaussian.h:65
RooRealProxy x
Definition RooGaussian.h:59
RooAbsReal const & getSigma() const
Get the sigma parameter.
Definition RooGaussian.h:51
const T & arg() const
Return reference to object held in proxy.
Mother of all ROOT objects.
Definition TObject.h:41