ROOT  6.06/09
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 //
18 // BEGIN_HTML
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 // END_HTML
27 //
28 
29 #include <cmath>
30 #include <complex>
31 
32 #include "RooFit.h"
33 
34 #include "Riostream.h"
35 
36 #include "RooVoigtian.h"
37 #include "RooAbsReal.h"
38 #include "RooRealVar.h"
39 #include "RooMath.h"
40 
41 using namespace std;
42 
44 
45 
46 ////////////////////////////////////////////////////////////////////////////////
47 
48 RooVoigtian::RooVoigtian(const char *name, const char *title,
49  RooAbsReal& _x, RooAbsReal& _mean,
50  RooAbsReal& _width, RooAbsReal& _sigma,
51  Bool_t doFast) :
52  RooAbsPdf(name,title),
53  x("x","Dependent",this,_x),
54  mean("mean","Mean",this,_mean),
55  width("width","Breit-Wigner Width",this,_width),
56  sigma("sigma","Gauss Width",this,_sigma),
57  _doFast(doFast)
58 {
59  _invRootPi= 1./sqrt(atan2(0.,-1.));
60 }
61 
62 
63 
64 ////////////////////////////////////////////////////////////////////////////////
65 
66 RooVoigtian::RooVoigtian(const RooVoigtian& other, const char* name) :
67  RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
68  width("width",this,other.width),sigma("sigma",this,other.sigma),
69  _doFast(other._doFast)
70 {
71  _invRootPi= 1./sqrt(atan2(0.,-1.));
72 }
73 
74 
75 
76 ////////////////////////////////////////////////////////////////////////////////
77 
79 {
80  Double_t s = (sigma>0) ? sigma : -sigma ;
81  Double_t w = (width>0) ? width : -width ;
82 
83  Double_t coef= -0.5/(s*s);
84  Double_t arg = x - mean;
85 
86  // return constant for zero width and sigma
87  if (s==0. && w==0.) return 1.;
88 
89  // Breit-Wigner for zero sigma
90  if (s==0.) return (1./(arg*arg+0.25*w*w));
91 
92  // Gauss for zero width
93  if (w==0.) return exp(coef*arg*arg);
94 
95  // actual Voigtian for non-trivial width and sigma
96  Double_t c = 1./(sqrt(2.)*s);
97  Double_t a = 0.5*c*w;
98  Double_t u = c*arg;
99  std::complex<Double_t> z(u,a) ;
100  std::complex<Double_t> v(0.) ;
101 
102  if (_doFast) {
103  v = RooMath::faddeeva_fast(z);
104  } else {
105  v = RooMath::faddeeva(z);
106  }
107  return c*_invRootPi*v.real();
108 
109 }
RooRealProxy mean
Definition: RooVoigtian.h:45
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
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:824
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)
#define name(a, b)
Definition: linkTestLib0.cpp:5
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
ClassImp(RooVoigtian) RooVoigtian
Definition: RooVoigtian.cxx:43
Double_t evaluate() const
Definition: RooVoigtian.cxx:78
double exp(double)
RooRealProxy x
Definition: RooVoigtian.h:44