Logo ROOT  
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 "RooVoigtian.h"
29#include "RooRealVar.h"
30#include "RooMath.h"
31#include "RooBatchCompute.h"
32
33#include <cmath>
34#include <complex>
35using namespace std;
36
38
39////////////////////////////////////////////////////////////////////////////////
40
41RooVoigtian::RooVoigtian(const char *name, const char *title,
42 RooAbsReal& _x, RooAbsReal& _mean,
43 RooAbsReal& _width, RooAbsReal& _sigma,
44 bool doFast) :
45 RooAbsPdf(name,title),
46 x("x","Dependent",this,_x),
47 mean("mean","Mean",this,_mean),
48 width("width","Breit-Wigner Width",this,_width),
49 sigma("sigma","Gauss Width",this,_sigma),
50 _doFast(doFast)
51{
52
53}
54
55////////////////////////////////////////////////////////////////////////////////
56
57RooVoigtian::RooVoigtian(const RooVoigtian& other, const char* name) :
58 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
59 width("width",this,other.width),sigma("sigma",this,other.sigma),
60 _doFast(other._doFast)
61{
62
63}
64
65////////////////////////////////////////////////////////////////////////////////
66
68{
69 double s = (sigma>0) ? sigma : -sigma ;
70 double w = (width>0) ? width : -width ;
71
72 double coef= -0.5/(s*s);
73 double arg = x - mean;
74
75 // return constant for zero width and sigma
76 if (s==0. && w==0.) return 1.;
77
78 // Breit-Wigner for zero sigma
79 if (s==0.) return (1./(arg*arg+0.25*w*w));
80
81 // Gauss for zero width
82 if (w==0.) return exp(coef*arg*arg);
83
84 // actual Voigtian for non-trivial width and sigma
85 double c = 1./(sqrt(2.)*s);
86 double a = 0.5*c*w;
87 double u = c*arg;
88 std::complex<double> z(u,a) ;
89 std::complex<double> v(0.) ;
90
91 if (_doFast) {
93 } else {
95 }
96 return c * v.real();
97}
98
99////////////////////////////////////////////////////////////////////////////////
100/// Compute multiple values of Voigtian distribution.
101void RooVoigtian::computeBatch(cudaStream_t* stream, double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
102{
104 dispatch->compute(stream, RooBatchCompute::Voigtian, output, nEvents,
105 {dataMap.at(x), dataMap.at(mean), dataMap.at(width), dataMap.at(sigma)});
106}
#define c(i)
Definition: RSha256.hxx:101
#define ClassImp(name)
Definition: Rtypes.h:375
Option_t Option_t width
char name[80]
Definition: TGX11.cxx:110
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:62
virtual void compute(cudaStream_t *, Computer, RestrictArr, size_t, const VarVector &, const ArgVector &={})=0
auto & at(RooAbsArg const *arg, RooAbsArg const *=nullptr)
Definition: DataMap.h:88
static std::complex< double > faddeeva(std::complex< double > z)
evaluate Faddeeva function for complex argument
Definition: RooMath.cxx:31
static std::complex< double > faddeeva_fast(std::complex< double > z)
evaluate Faddeeva function for complex argument (fast version)
Definition: RooMath.cxx:36
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
void computeBatch(cudaStream_t *, double *output, size_t nEvents, RooFit::Detail::DataMap const &) const override
Compute multiple values of Voigtian distribution.
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition: RooVoigtian.cxx:67
bool _doFast
Definition: RooVoigtian.h:55
RooRealProxy sigma
Definition: RooVoigtian.h:47
RooRealProxy mean
Definition: RooVoigtian.h:45
RooRealProxy width
Definition: RooVoigtian.h:46
RVec< PromoteType< T > > exp(const RVec< T > &v)
Definition: RVec.hxx:1761
const Double_t sigma
Double_t x[n]
Definition: legend1.C:17
VecExpr< UnaryOp< Sqrt< T >, VecExpr< A, T, D >, T >, T, D > sqrt(const VecExpr< A, T, D > &rhs)
R__EXTERN RooBatchComputeInterface * dispatchCUDA
R__EXTERN RooBatchComputeInterface * dispatchCPU
This dispatch pointer points to an implementation of the compute library, provided one has been loade...
static constexpr double s
TArc a
Definition: textangle.C:12
static void output()