Logo ROOT   6.18/05
Reference Guide
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
19RooVoigtian is an efficient implementation of the convolution of a
20Breit-Wigner with a Gaussian, making use of the complex error function.
21RooFitCore provides two algorithms for the evaluation of the complex error
22function (the default CERNlib C335 algorithm, and a faster, look-up-table
23based method). By default, RooVoigtian employs the default (CERNlib)
24algorithm. Select the faster algorithm either in the constructor, or with
25the selectFastAlgorithm() method.
26**/
27
28#include <cmath>
29#include <complex>
30
31#include "RooFit.h"
32
33#include "Riostream.h"
34
35#include "RooVoigtian.h"
36#include "RooAbsReal.h"
37#include "RooRealVar.h"
38#include "RooMath.h"
39
40using namespace std;
41
43
44////////////////////////////////////////////////////////////////////////////////
45
46RooVoigtian::RooVoigtian(const char *name, const char *title,
47 RooAbsReal& _x, RooAbsReal& _mean,
48 RooAbsReal& _width, RooAbsReal& _sigma,
49 Bool_t doFast) :
50 RooAbsPdf(name,title),
51 x("x","Dependent",this,_x),
52 mean("mean","Mean",this,_mean),
53 width("width","Breit-Wigner Width",this,_width),
54 sigma("sigma","Gauss Width",this,_sigma),
55 _doFast(doFast)
56{
57 _invRootPi= 1./sqrt(atan2(0.,-1.));
58}
59
60////////////////////////////////////////////////////////////////////////////////
61
62RooVoigtian::RooVoigtian(const RooVoigtian& other, const char* name) :
63 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
64 width("width",this,other.width),sigma("sigma",this,other.sigma),
65 _doFast(other._doFast)
66{
67 _invRootPi= 1./sqrt(atan2(0.,-1.));
68}
69
70////////////////////////////////////////////////////////////////////////////////
71
73{
74 Double_t s = (sigma>0) ? sigma : -sigma ;
75 Double_t w = (width>0) ? width : -width ;
76
77 Double_t coef= -0.5/(s*s);
78 Double_t arg = x - mean;
79
80 // return constant for zero width and sigma
81 if (s==0. && w==0.) return 1.;
82
83 // Breit-Wigner for zero sigma
84 if (s==0.) return (1./(arg*arg+0.25*w*w));
85
86 // Gauss for zero width
87 if (w==0.) return exp(coef*arg*arg);
88
89 // actual Voigtian for non-trivial width and sigma
90 Double_t c = 1./(sqrt(2.)*s);
91 Double_t a = 0.5*c*w;
92 Double_t u = c*arg;
93 std::complex<Double_t> z(u,a) ;
94 std::complex<Double_t> v(0.) ;
95
96 if (_doFast) {
98 } else {
100 }
101 return c*_invRootPi*v.real();
102
103}
SVector< double, 2 > v
Definition: Dict.h:5
#define c(i)
Definition: RSha256.hxx:101
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
#define ClassImp(name)
Definition: Rtypes.h:365
include TDocParser_001 C image html pict1_TDocParser_001 png width
Definition: TDocParser.cxx:121
char name[80]
Definition: TGX11.cxx:109
double atan2(double, double)
double sqrt(double)
double exp(double)
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
static std::complex< double > faddeeva(std::complex< double > z)
evaluate Faddeeva function for complex argument
Definition: RooMath.cxx:542
static std::complex< double > faddeeva_fast(std::complex< double > z)
evaluate Faddeeva function for complex argument (fast version)
Definition: RooMath.cxx:549
RooVoigtian is an efficient implementation of the convolution of a Breit-Wigner with a Gaussian,...
Definition: RooVoigtian.h:24
RooRealProxy x
Definition: RooVoigtian.h:44
Bool_t _doFast
Definition: RooVoigtian.h:54
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition: RooVoigtian.cxx:72
RooRealProxy sigma
Definition: RooVoigtian.h:47
RooRealProxy mean
Definition: RooVoigtian.h:45
RooRealProxy width
Definition: RooVoigtian.h:46
Double_t _invRootPi
Definition: RooVoigtian.h:53
const Double_t sigma
Double_t x[n]
Definition: legend1.C:17
static constexpr double s
auto * a
Definition: textangle.C:12