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 "RooFit.h"
30#include "RooAbsReal.h"
31#include "RooRealVar.h"
32#include "RooMath.h"
33#include "BatchHelpers.h"
34#include "RooVDTHeaders.h"
35
36#include <cmath>
37#include <complex>
38using namespace std;
39
41
42////////////////////////////////////////////////////////////////////////////////
43
44RooVoigtian::RooVoigtian(const char *name, const char *title,
45 RooAbsReal& _x, RooAbsReal& _mean,
46 RooAbsReal& _width, RooAbsReal& _sigma,
47 Bool_t doFast) :
48 RooAbsPdf(name,title),
49 x("x","Dependent",this,_x),
50 mean("mean","Mean",this,_mean),
51 width("width","Breit-Wigner Width",this,_width),
52 sigma("sigma","Gauss Width",this,_sigma),
53 _doFast(doFast)
54{
55
56}
57
58////////////////////////////////////////////////////////////////////////////////
59
60RooVoigtian::RooVoigtian(const RooVoigtian& other, const char* name) :
61 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
62 width("width",this,other.width),sigma("sigma",this,other.sigma),
63 _doFast(other._doFast)
64{
65
66}
67
68////////////////////////////////////////////////////////////////////////////////
69
71{
72 Double_t s = (sigma>0) ? sigma : -sigma ;
73 Double_t w = (width>0) ? width : -width ;
74
75 Double_t coef= -0.5/(s*s);
76 Double_t arg = x - mean;
77
78 // return constant for zero width and sigma
79 if (s==0. && w==0.) return 1.;
80
81 // Breit-Wigner for zero sigma
82 if (s==0.) return (1./(arg*arg+0.25*w*w));
83
84 // Gauss for zero width
85 if (w==0.) return exp(coef*arg*arg);
86
87 // actual Voigtian for non-trivial width and sigma
88 Double_t c = 1./(sqrt(2.)*s);
89 Double_t a = 0.5*c*w;
90 Double_t u = c*arg;
91 std::complex<Double_t> z(u,a) ;
92 std::complex<Double_t> v(0.) ;
93
94 if (_doFast) {
96 } else {
98 }
99 return c * v.real();
100}
101
102////////////////////////////////////////////////////////////////////////////////
103
104namespace {
105//Author: Emmanouil Michalainas, CERN 11 September 2019
106
107template<class Tx, class Tmean, class Twidth, class Tsigma>
108void compute( size_t batchSize, double * __restrict output,
109 Tx X, Tmean M, Twidth W, Tsigma S)
110{
111 constexpr double invSqrt2 = 0.707106781186547524400844362105;
112 for (size_t i=0; i<batchSize; i++) {
113 const double arg = (X[i]-M[i])*(X[i]-M[i]);
114 if (S[i]==0.0 && W[i]==0.0) {
115 output[i] = 1.0;
116 } else if (S[i]==0.0) {
117 output[i] = 1/(arg+0.25*W[i]*W[i]);
118 } else if (W[i]==0.0) {
119 output[i] = _rf_fast_exp(-0.5*arg/(S[i]*S[i]));
120 } else {
121 output[i] = invSqrt2/S[i];
122 }
123 }
124
125 for (size_t i=0; i<batchSize; i++) {
126 if (S[i]!=0.0 && W[i]!=0.0) {
127 if (output[i] < 0) output[i] = -output[i];
128 const double factor = W[i]>0.0 ? 0.5 : -0.5;
129 std::complex<Double_t> z( output[i]*(X[i]-M[i]) , factor*output[i]*W[i] );
130 output[i] *= RooMath::faddeeva(z).real();
131 }
132 }
133}
134};
135
136RooSpan<double> RooVoigtian::evaluateBatch(std::size_t begin, std::size_t batchSize) const {
137 using namespace BatchHelpers;
138
139 EvaluateInfo info = getInfo( {&x, &mean, &width, &sigma}, begin, batchSize );
140 if (info.nBatches == 0) {
141 return {};
142 }
143 auto output = _batchData.makeWritableBatchUnInit(begin, batchSize);
144 auto xData = x.getValBatch(begin, info.size);
145
146 if (info.nBatches==1 && !xData.empty()) {
147 compute(info.size, output.data(), xData.data(),
151 }
152 else {
153 compute(info.size, output.data(),
158 }
159 return output;
160}
161
#define c(i)
Definition: RSha256.hxx:101
double _rf_fast_exp(double x)
VDT headers for RooFit.
Definition: RooVDTHeaders.h:47
#define ClassImp(name)
Definition: Rtypes.h:361
include TDocParser_001 C image html pict1_TDocParser_001 png width
Definition: TDocParser.cxx:121
char name[80]
Definition: TGX11.cxx:109
double sqrt(double)
double exp(double)
RooSpan< double > makeWritableBatchUnInit(std::size_t begin, std::size_t batchSize, const RooArgSet *const normSet=nullptr, Tag_t ownerTag=kUnspecified)
Make a batch and return a span pointing to the pdf-local memory.
Definition: BatchData.cxx:118
Little adapter that gives a bracket operator to types that don't have one.
Definition: BatchHelpers.h:58
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:60
BatchHelpers::BatchData _batchData
Definition: RooAbsReal.h:450
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:32
RooSpan< const typename T::value_type > getValBatch(std::size_t begin, std::size_t batchSize) const
Retrieve a batch of real or category data.
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 > evaluateBatch(std::size_t begin, std::size_t batchSize) const
Evaluate function for a batch of input data points.
Bool_t _doFast
Definition: RooVoigtian.h:55
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition: RooVoigtian.cxx:70
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
EvaluateInfo getInfo(std::vector< const RooRealProxy * > parameters, size_t begin, size_t batchSize)
RooArgSet S(const RooAbsArg &v1)
static constexpr double s
auto * a
Definition: textangle.C:12
static void output(int code)
Definition: gifencode.c:226