ROOT
v6-32
Reference Guide
Loading...
Searching...
No Matches
RooVoigtian.cxx
Go to the documentation of this file.
1
/*****************************************************************************
2
* Project: RooFit *
3
* Package: RooFitModels *
4
* @(#)root/roofit:$Id$
5
* Authors: *
6
* TS, Thomas Schietinger, SLAC, schieti@slac.stanford.edu *
7
* *
8
* Copyright (c) 2000-2005, Regents of the University of California *
9
* and Stanford University. All rights reserved. *
10
* *
11
* Redistribution and use in source and binary forms, *
12
* with or without modification, are permitted according to the terms *
13
* listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
14
*****************************************************************************/
15
16
/** \class RooVoigtian
17
\ingroup Roofit
18
19
RooVoigtian is an efficient implementation of the convolution of a
20
Breit-Wigner with a Gaussian, making use of the complex error function.
21
RooFitCore provides two algorithms for the evaluation of the complex error
22
function (the default CERNlib C335 algorithm, and a faster, look-up-table
23
based method). By default, RooVoigtian employs the default (CERNlib)
24
algorithm. Select the faster algorithm either in the constructor, or with
25
the selectFastAlgorithm() method.
26
27
\note The "width" parameter that determines the Breit-Wigner shape
28
represents the **full width at half maximum (FWHM)** of the
29
Breit-Wigner (often referred to as \f$\Gamma\f$ or \f$2\gamma\f$).
30
**/
31
32
#include <
RooVoigtian.h
>
33
34
#include <
RooMath.h
>
35
#include <
RooBatchCompute.h
>
36
37
#include <cmath>
38
#include <complex>
39
40
ClassImp
(
RooVoigtian
);
41
42
////////////////////////////////////////////////////////////////////////////////
43
/// Construct a RooVoigtian PDF, which represents the convolution of a
44
/// Breit-Wigner with a Gaussian.
45
/// \param name Name that identifies the PDF in computations.
46
/// \param title Title for plotting.
47
/// \param _x The observable for the PDF.
48
/// \param _mean The mean of the distribution.
49
/// \param _width The **full width at half maximum (FWHM)** of the Breit-Wigner
50
/// (often referred to as \f$\Gamma\f$ or \f$2\gamma\f$).
51
/// \param _sigma The width of the Gaussian distribution.
52
/// \param doFast Use the faster look-up-table-based method for the evaluation
53
/// of the complex error function.
54
55
RooVoigtian::RooVoigtian
(
const
char
*
name
,
const
char
*title,
56
RooAbsReal
& _x,
RooAbsReal
& _mean,
57
RooAbsReal
&
_width
,
RooAbsReal
& _sigma,
58
bool
doFast
) :
59
RooAbsPdf
(
name
,title),
60
x
(
"x"
,
"Dependent"
,
this
,_x),
61
mean(
"mean"
,
"Mean"
,
this
,_mean),
62
width
(
"width"
,
"Breit-Wigner Width"
,
this
,
_width
),
63
sigma
(
"sigma"
,
"Gauss Width"
,
this
,_sigma),
64
_doFast(
doFast
)
65
{
66
67
}
68
69
////////////////////////////////////////////////////////////////////////////////
70
71
RooVoigtian::RooVoigtian
(
const
RooVoigtian
&
other
,
const
char
*
name
) :
72
RooAbsPdf
(
other
,
name
),
x
(
"x"
,
this
,
other
.
x
), mean(
"mean"
,
this
,
other
.mean),
73
width
(
"width"
,
this
,
other
.
width
),
sigma
(
"sigma"
,
this
,
other
.
sigma
),
74
_doFast(
other
._doFast)
75
{
76
77
}
78
79
////////////////////////////////////////////////////////////////////////////////
80
81
double
RooVoigtian::evaluate
()
const
82
{
83
double
s = (
sigma
>0) ?
sigma
: -
sigma
;
84
double
w
= (
width
>0) ?
width
: -
width
;
85
86
double
coef= -0.5/(s*s);
87
double
arg =
x
-
mean
;
88
89
// return constant for zero width and sigma
90
if
(s==0. &&
w
==0.)
return
1.;
91
92
// Breit-Wigner for zero sigma
93
if
(s==0.)
return
(1./(arg*arg+0.25*
w
*
w
));
94
95
// Gauss for zero width
96
if
(
w
==0.)
return
std::exp(coef*arg*arg);
97
98
// actual Voigtian for non-trivial width and sigma
99
double
c
= 1./(sqrt(2.)*s);
100
double
a
= 0.5*
c
*
w
;
101
double
u
=
c
*arg;
102
std::complex<double> z(
u
,
a
) ;
103
std::complex<double>
v
(0.) ;
104
105
if
(
_doFast
) {
106
v
=
RooMath::faddeeva_fast
(z);
107
}
else
{
108
v
=
RooMath::faddeeva
(z);
109
}
110
return
c
*
v
.real();
111
}
112
113
////////////////////////////////////////////////////////////////////////////////
114
/// Compute multiple values of Voigtian distribution.
115
void
RooVoigtian::doEval
(
RooFit::EvalContext
&ctx)
const
116
{
117
RooBatchCompute::compute
(ctx.
config
(
this
),
RooBatchCompute::Voigtian
, ctx.
output
(),
118
{ctx.at(x), ctx.at(mean), ctx.at(width), ctx.at(sigma)});
119
}
c
#define c(i)
Definition
RSha256.hxx:101
a
#define a(i)
Definition
RSha256.hxx:99
RooBatchCompute.h
RooMath.h
RooVoigtian.h
ClassImp
#define ClassImp(name)
Definition
Rtypes.h:377
w
winID w
Definition
TGWin32VirtualGLProxy.cxx:39
width
Option_t Option_t width
Definition
TGWin32VirtualXProxy.cxx:56
name
char name[80]
Definition
TGX11.cxx:110
ROOT::Detail::TRangeCast
Definition
TCollection.h:311
RooAbsPdf
Abstract interface for all probability density functions.
Definition
RooAbsPdf.h:40
RooAbsReal
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition
RooAbsReal.h:59
RooFit::EvalContext
Definition
EvalContext.h:84
RooFit::EvalContext::output
std::span< double > output()
Definition
EvalContext.h:112
RooFit::EvalContext::config
RooBatchCompute::Config config(RooAbsArg const *arg) const
Definition
EvalContext.cxx:73
RooMath::faddeeva
static std::complex< double > faddeeva(std::complex< double > z)
evaluate Faddeeva function for complex argument
Definition
RooMath.cxx:30
RooMath::faddeeva_fast
static std::complex< double > faddeeva_fast(std::complex< double > z)
evaluate Faddeeva function for complex argument (fast version)
Definition
RooMath.cxx:35
RooVoigtian
RooVoigtian is an efficient implementation of the convolution of a Breit-Wigner with a Gaussian,...
Definition
RooVoigtian.h:22
RooVoigtian::x
RooRealProxy x
Definition
RooVoigtian.h:42
RooVoigtian::doEval
void doEval(RooFit::EvalContext &) const override
Compute multiple values of Voigtian distribution.
Definition
RooVoigtian.cxx:115
RooVoigtian::evaluate
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition
RooVoigtian.cxx:81
RooVoigtian::_doFast
bool _doFast
Definition
RooVoigtian.h:53
RooVoigtian::sigma
RooRealProxy sigma
Definition
RooVoigtian.h:45
RooVoigtian::RooVoigtian
RooVoigtian()
Definition
RooVoigtian.h:24
RooVoigtian::mean
RooRealProxy mean
Definition
RooVoigtian.h:43
RooVoigtian::width
RooRealProxy width
Definition
RooVoigtian.h:44
sigma
const Double_t sigma
Definition
h1analysisProxy.h:11
x
Double_t x[n]
Definition
legend1.C:17
RooBatchCompute::compute
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})
Definition
RooBatchCompute.h:205
RooBatchCompute::Voigtian
@ Voigtian
Definition
RooBatchCompute.h:106
v
@ v
Definition
rootcling_impl.cxx:3687
roofit
roofit
src
RooVoigtian.cxx
ROOT v6-32 - Reference Guide Generated on Thu Mar 6 2025 14:28:50 (GVA Time) using Doxygen 1.10.0