Logo ROOT   6.08/07
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 /**
17 \file RooVoigtian.cxx
18 \class RooVoigtian
19 \ingroup Roofit
20 
21 RooVoigtian is an efficient implementation of the convolution of a
22 Breit-Wigner with a Gaussian, making use of the complex error function.
23 RooFitCore provides two algorithms for the evaluation of the complex error
24 function (the default CERNlib C335 algorithm, and a faster, look-up-table
25 based method). By default, RooVoigtian employs the default (CERNlib)
26 algorithm. Select the faster algorithm either in the constructor, or with
27 the selectFastAlgorithm() method.
28 **/
29 
30 #include <cmath>
31 #include <complex>
32 
33 #include "RooFit.h"
34 
35 #include "Riostream.h"
36 
37 #include "RooVoigtian.h"
38 #include "RooAbsReal.h"
39 #include "RooRealVar.h"
40 #include "RooMath.h"
41 
42 using namespace std;
43 
45 
46 
47 ////////////////////////////////////////////////////////////////////////////////
48 
49 RooVoigtian::RooVoigtian(const char *name, const char *title,
50  RooAbsReal& _x, RooAbsReal& _mean,
51  RooAbsReal& _width, RooAbsReal& _sigma,
52  Bool_t doFast) :
53  RooAbsPdf(name,title),
54  x("x","Dependent",this,_x),
55  mean("mean","Mean",this,_mean),
56  width("width","Breit-Wigner Width",this,_width),
57  sigma("sigma","Gauss Width",this,_sigma),
58  _doFast(doFast)
59 {
60  _invRootPi= 1./sqrt(atan2(0.,-1.));
61 }
62 
63 
64 
65 ////////////////////////////////////////////////////////////////////////////////
66 
67 RooVoigtian::RooVoigtian(const RooVoigtian& other, const char* name) :
68  RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
69  width("width",this,other.width),sigma("sigma",this,other.sigma),
70  _doFast(other._doFast)
71 {
72  _invRootPi= 1./sqrt(atan2(0.,-1.));
73 }
74 
75 
76 
77 ////////////////////////////////////////////////////////////////////////////////
78 
80 {
81  Double_t s = (sigma>0) ? sigma : -sigma ;
82  Double_t w = (width>0) ? width : -width ;
83 
84  Double_t coef= -0.5/(s*s);
85  Double_t arg = x - mean;
86 
87  // return constant for zero width and sigma
88  if (s==0. && w==0.) return 1.;
89 
90  // Breit-Wigner for zero sigma
91  if (s==0.) return (1./(arg*arg+0.25*w*w));
92 
93  // Gauss for zero width
94  if (w==0.) return exp(coef*arg*arg);
95 
96  // actual Voigtian for non-trivial width and sigma
97  Double_t c = 1./(sqrt(2.)*s);
98  Double_t a = 0.5*c*w;
99  Double_t u = c*arg;
100  std::complex<Double_t> z(u,a) ;
101  std::complex<Double_t> v(0.) ;
102 
103  if (_doFast) {
104  v = RooMath::faddeeva_fast(z);
105  } else {
106  v = RooMath::faddeeva(z);
107  }
108  return c*_invRootPi*v.real();
109 
110 }
RooRealProxy mean
Definition: RooVoigtian.h:45
return c
Double_t _invRootPi
Definition: RooVoigtian.h:53
Bool_t _doFast
Definition: RooVoigtian.h:54
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
RooRealProxy width
Definition: RooVoigtian.h:46
STL namespace.
double sqrt(double)
Double_t x[n]
Definition: legend1.C:17
static std::complex< double > faddeeva(std::complex< double > z)
evaluate Faddeeva function for complex argument
Definition: RooMath.cxx:546
const Double_t sigma
static std::complex< double > faddeeva_fast(std::complex< double > z)
evaluate Faddeeva function for complex argument (fast version)
Definition: RooMath.cxx:553
SVector< double, 2 > v
Definition: Dict.h:5
Double_t Mean(Long64_t n, const T *a, const Double_t *w=0)
Definition: TMath.h:811
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
double atan2(double, double)
Double_t evaluate() const
Definition: RooVoigtian.cxx:79
static std::shared_ptr< std::function< double(double)> > Gauss
Definition: NeuralNet.icc:71
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
RooRealProxy sigma
Definition: RooVoigtian.h:47
double exp(double)
RooRealProxy x
Definition: RooVoigtian.h:44
RooVoigtian is an efficient implementation of the convolution of a Breit-Wigner with a Gaussian...
Definition: RooVoigtian.h:24
char name[80]
Definition: TGX11.cxx:109