Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooBifurGauss.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * Abi Soffer, Colorado State University, abi@slac.stanford.edu *
7 * *
8 * Copyright (c) 2000-2005, Regents of the University of California, *
9 * Colorado State University *
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 RooBifurGauss
18 \ingroup Roofit
19
20Bifurcated Gaussian p.d.f with different widths on left and right
21side of maximum value.
22**/
23
24#include <RooBifurGauss.h>
25
26#include "RooBatchCompute.h"
27
29
31
32////////////////////////////////////////////////////////////////////////////////
33
34RooBifurGauss::RooBifurGauss(const char *name, const char *title, RooAbsReal &_x, RooAbsReal &_mean,
35 RooAbsReal &_sigmaL, RooAbsReal &_sigmaR)
36 : RooAbsPdf(name, title),
37 x("x", "Dependent", this, _x),
38 mean("mean", "Mean", this, _mean),
39 sigmaL("sigmaL", "Left Sigma", this, _sigmaL),
40 sigmaR("sigmaR", "Right Sigma", this, _sigmaR)
41
42{
43}
44
45////////////////////////////////////////////////////////////////////////////////
46
48 : RooAbsPdf(other, name),
49 x("x", this, other.x),
50 mean("mean", this, other.mean),
51 sigmaL("sigmaL", this, other.sigmaL),
52 sigmaR("sigmaR", this, other.sigmaR)
53{
54}
55
56////////////////////////////////////////////////////////////////////////////////
57
59{
61}
62
63////////////////////////////////////////////////////////////////////////////////
64
66{
67 ctx.addResult(this, ctx.buildCall("RooFit::Detail::MathFuncs::bifurGauss", x, mean, sigmaL, sigmaR));
68}
69
70////////////////////////////////////////////////////////////////////////////////
71/// Compute multiple values of BifurGauss distribution.
73{
75 {ctx.at(x),ctx.at(mean),ctx.at(sigmaL),ctx.at(sigmaR)});
76}
77
78////////////////////////////////////////////////////////////////////////////////
79
80Int_t RooBifurGauss::getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char * /*rangeName*/) const
81{
82 if (matchArgs(allVars, analVars, x))
83 return 1;
84 if (matchArgs(allVars, analVars, mean))
85 return 2;
86 return 0;
87}
88
89////////////////////////////////////////////////////////////////////////////////
90
91double RooBifurGauss::analyticalIntegral(Int_t code, const char *rangeName) const
92{
93 auto &constant = code == 1 ? mean : x;
94 auto &integrand = code == 1 ? x : mean;
95
96 return RooFit::Detail::MathFuncs::bifurGaussIntegral(integrand.min(rangeName), integrand.max(rangeName),
97 constant, sigmaL, sigmaR);
98}
99
100////////////////////////////////////////////////////////////////////////////////
101
102std::string RooBifurGauss::buildCallToAnalyticIntegral(Int_t code, const char *rangeName,
104{
105 auto &constant = code == 1 ? mean : x;
106 auto &integrand = code == 1 ? x : mean;
107
108 return ctx.buildCall("RooFit::Detail::MathFuncs::bifurGaussIntegral", integrand.min(rangeName),
109 integrand.max(rangeName), constant, sigmaL, sigmaR);
110}
#define ClassImp(name)
Definition Rtypes.h:377
char name[80]
Definition TGX11.cxx:110
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
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:55
Bifurcated Gaussian p.d.f with different widths on left and right side of maximum value.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
void doEval(RooFit::EvalContext &) const override
Compute multiple values of BifurGauss distribution.
RooRealProxy mean
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 sigmaL
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooRealProxy sigmaR
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 x
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
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
Double_t x[n]
Definition legend1.C:17
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})
double bifurGaussIntegral(double xMin, double xMax, double mean, double sigmaL, double sigmaR)
Definition MathFuncs.h:450
double bifurGauss(double x, double mean, double sigmaL, double sigmaR)
Definition MathFuncs.h:108