RooVoigtian.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitModels                                                     *
00004  * @(#)root/roofit:$Id: RooVoigtian.cxx 24286 2008-06-16 15:47:04Z wouter $
00005  * Authors:                                                                  *
00006  *   TS, Thomas Schietinger, SLAC,           schieti@slac.stanford.edu       *
00007  *                                                                           *
00008  * Copyright (c) 2000-2005, Regents of the University of California          *
00009  *                          and Stanford University. All rights reserved.    *
00010  *                                                                           *
00011  * Redistribution and use in source and binary forms,                        *
00012  * with or without modification, are permitted according to the terms        *
00013  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00014  *****************************************************************************/
00015 
00016 //////////////////////////////////////////////////////////////////////////////
00017 //
00018 // BEGIN_HTML
00019 // RooVoigtian is an efficient implementation of the convolution of a 
00020 // Breit-Wigner with a Gaussian, making use of the complex error function.
00021 // RooFitCore provides two algorithms for the evaluation of the complex error 
00022 // function (the default CERNlib C335 algorithm, and a faster, look-up-table 
00023 // based method). By default, RooVoigtian employs the default (CERNlib) 
00024 // algorithm. Select the faster algorithm either in the constructor, or with
00025 // the selectFastAlgorithm() method.
00026 // END_HTML
00027 //
00028 
00029 
00030 #include "RooFit.h"
00031 
00032 #include "Riostream.h"
00033 #include "Riostream.h"
00034 #include <math.h>
00035 
00036 #include "RooVoigtian.h"
00037 #include "RooAbsReal.h"
00038 #include "RooRealVar.h"
00039 #include "RooComplex.h"
00040 #include "RooMath.h"
00041 
00042 ClassImp(RooVoigtian)
00043 
00044 
00045 //_____________________________________________________________________________
00046 RooVoigtian::RooVoigtian(const char *name, const char *title,
00047                          RooAbsReal& _x, RooAbsReal& _mean,
00048                          RooAbsReal& _width, RooAbsReal& _sigma,
00049                          Bool_t doFast) :
00050   RooAbsPdf(name,title),
00051   x("x","Dependent",this,_x),
00052   mean("mean","Mean",this,_mean),
00053   width("width","Breit-Wigner Width",this,_width),
00054   sigma("sigma","Gauss Width",this,_sigma),
00055   _doFast(doFast)
00056 {
00057   _invRootPi= 1./sqrt(atan2(0.,-1.));
00058 }
00059 
00060 
00061 
00062 //_____________________________________________________________________________
00063 RooVoigtian::RooVoigtian(const RooVoigtian& other, const char* name) : 
00064   RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
00065   width("width",this,other.width),sigma("sigma",this,other.sigma),
00066   _doFast(other._doFast)
00067 {
00068   _invRootPi= 1./sqrt(atan2(0.,-1.));
00069 }
00070 
00071 
00072 
00073 //_____________________________________________________________________________
00074 Double_t RooVoigtian::evaluate() const
00075 {
00076   Double_t s = (sigma>0) ? sigma : -sigma ;
00077   Double_t w = (width>0) ? width : -width ;
00078 
00079   Double_t coef= -0.5/(s*s);
00080   Double_t arg = x - mean;
00081 
00082   // return constant for zero width and sigma
00083   if (s==0. && w==0.) return 1.;
00084 
00085   // Breit-Wigner for zero sigma
00086   if (s==0.) return (1./(arg*arg+0.25*w*w));
00087 
00088   // Gauss for zero width
00089   if (w==0.) return exp(coef*arg*arg);
00090 
00091   // actual Voigtian for non-trivial width and sigma
00092   Double_t c = 1./(sqrt(2.)*s);
00093   Double_t a = 0.5*c*w;
00094   Double_t u = c*arg;
00095   RooComplex z(u,a) ;
00096   RooComplex v(0.) ;
00097 
00098   if (_doFast) {
00099     v = RooMath::FastComplexErrFunc(z);
00100   } else {
00101     v = RooMath::ComplexErrFunc(z);
00102   }
00103   return c*_invRootPi*v.re();
00104 
00105 }

Generated on Tue Jul 5 14:55:27 2011 for ROOT_528-00b_version by  doxygen 1.5.1