Logo ROOT  
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
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 "RooFit.h"
30#include "RooAbsReal.h"
31#include "RooRealVar.h"
32#include "RooMath.h"
33#include "RooBatchCompute.h"
34
35#include <cmath>
36#include <complex>
37using namespace std;
38
40
41////////////////////////////////////////////////////////////////////////////////
42
43RooVoigtian::RooVoigtian(const char *name, const char *title,
44 RooAbsReal& _x, RooAbsReal& _mean,
45 RooAbsReal& _width, RooAbsReal& _sigma,
46 Bool_t doFast) :
47 RooAbsPdf(name,title),
48 x("x","Dependent",this,_x),
49 mean("mean","Mean",this,_mean),
50 width("width","Breit-Wigner Width",this,_width),
51 sigma("sigma","Gauss Width",this,_sigma),
52 _doFast(doFast)
53{
54
55}
56
57////////////////////////////////////////////////////////////////////////////////
58
59RooVoigtian::RooVoigtian(const RooVoigtian& other, const char* name) :
60 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
61 width("width",this,other.width),sigma("sigma",this,other.sigma),
62 _doFast(other._doFast)
63{
64
65}
66
67////////////////////////////////////////////////////////////////////////////////
68
70{
71 Double_t s = (sigma>0) ? sigma : -sigma ;
72 Double_t w = (width>0) ? width : -width ;
73
74 Double_t coef= -0.5/(s*s);
75 Double_t arg = x - mean;
76
77 // return constant for zero width and sigma
78 if (s==0. && w==0.) return 1.;
79
80 // Breit-Wigner for zero sigma
81 if (s==0.) return (1./(arg*arg+0.25*w*w));
82
83 // Gauss for zero width
84 if (w==0.) return exp(coef*arg*arg);
85
86 // actual Voigtian for non-trivial width and sigma
87 Double_t c = 1./(sqrt(2.)*s);
88 Double_t a = 0.5*c*w;
89 Double_t u = c*arg;
90 std::complex<Double_t> z(u,a) ;
91 std::complex<Double_t> v(0.) ;
92
93 if (_doFast) {
95 } else {
97 }
98 return c * v.real();
99}
100
101////////////////////////////////////////////////////////////////////////////////
102/// Compute multiple values of Voigtian distribution.
104 return RooBatchCompute::dispatch->computeVoigtian(this, evalData, x->getValues(evalData, normSet), mean->getValues(evalData, normSet), width->getValues(evalData, normSet), sigma->getValues(evalData, normSet));
105}
106
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
#define ClassImp(name)
Definition Rtypes.h:364
include TDocParser_001 C image html pict1_TDocParser_001 png width
char name[80]
Definition TGX11.cxx:110
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:61
virtual RooSpan< const double > getValues(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet=nullptr) const
by this change, please consult the release notes for ROOT 6.24 for guidance on how to make this trans...
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
virtual RooSpan< double > computeVoigtian(const RooAbsReal *, RunContext &, RooSpan< const double > x, RooSpan< const double > mean, RooSpan< const double > width, RooSpan< const double > sigma)=0
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
A simple container to hold a batch of data values.
Definition RooSpan.h:34
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
RooSpan< double > evaluateSpan(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet) const
Compute multiple values of Voigtian distribution.
Bool_t _doFast
Definition RooVoigtian.h:54
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooRealProxy sigma
Definition RooVoigtian.h:47
RooRealProxy mean
Definition RooVoigtian.h:45
RooRealProxy width
Definition RooVoigtian.h:46
const Double_t sigma
Double_t x[n]
Definition legend1.C:17
R__EXTERN RooBatchComputeInterface * dispatch
This dispatch pointer points to an implementation of the compute library, provided one has been loade...
This struct enables passing computation data around between elements of a computation graph.
Definition RunContext.h:31