Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooNovosibirsk.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * DB, Dieter Best, UC Irvine, best@slac.stanford.edu *
7 * HT, Hirohisa Tanaka SLAC tanaka@slac.stanford.edu *
8 * *
9 * Updated version with analytical integral *
10 * MP, Marko Petric, J. Stefan Institute, marko.petric@ijs.si *
11 * *
12 * Copyright (c) 2000-2013, Regents of the University of California *
13 * and Stanford University. All rights reserved. *
14 * *
15 * Redistribution and use in source and binary forms, *
16 * with or without modification, are permitted according to the terms *
17 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
18 *****************************************************************************/
19
20/** \class RooNovosibirsk
21 \ingroup Roofit
22
23RooNovosibirsk implements the Novosibirsk function
24
25Function taken from H. Ikeda et al. NIM A441 (2000), p. 401 (Belle Collaboration)
26
27**/
28#include "RooNovosibirsk.h"
29#include "RooFit.h"
30#include "RooRealVar.h"
31#include "RooBatchCompute.h"
32
33#include "TMath.h"
34
35#include <cmath>
36using namespace std;
37
39
40////////////////////////////////////////////////////////////////////////////////
41
42RooNovosibirsk::RooNovosibirsk(const char *name, const char *title,
43 RooAbsReal& _x, RooAbsReal& _peak,
44 RooAbsReal& _width, RooAbsReal& _tail) :
45 // The two addresses refer to our first dependent variable and
46 // parameter, respectively, as declared in the rdl file
47 RooAbsPdf(name, title),
48 x("x","x",this,_x),
49 width("width","width",this,_width),
50 peak("peak","peak",this,_peak),
51 tail("tail","tail",this,_tail)
52{
53}
54
55////////////////////////////////////////////////////////////////////////////////
56
58 RooAbsPdf(other,name),
59 x("x",this,other.x),
60 width("width",this,other.width),
61 peak("peak",this,other.peak),
62 tail("tail",this,other.tail)
63{
64}
65
66////////////////////////////////////////////////////////////////////////////////
67///If tail=eta=0 the Belle distribution becomes gaussian
68
70{
71 if (TMath::Abs(tail) < 1.e-7) {
72 return TMath::Exp( -0.5 * TMath::Power( ( (x - peak) / width ), 2 ));
73 }
74
75 Double_t arg = 1.0 - ( x - peak ) * tail / width;
76
77 if (arg < 1.e-7) {
78 //Argument of logarithm negative. Real continuation -> function equals zero
79 return 0.0;
80 }
81
82 Double_t log = TMath::Log(arg);
83 static const Double_t xi = 2.3548200450309494; // 2 Sqrt( Ln(4) )
84
85 Double_t width_zero = ( 2.0 / xi ) * TMath::ASinH( tail * xi * 0.5 );
86 Double_t width_zero2 = width_zero * width_zero;
87 Double_t exponent = ( -0.5 / (width_zero2) * log * log ) - ( width_zero2 * 0.5 );
88
89 return TMath::Exp(exponent) ;
90}
91////////////////////////////////////////////////////////////////////////////////
92/// Compute multiple values of Novosibirsk distribution.
94 return RooBatchCompute::dispatch->computeNovosibirsk(this, evalData, x->getValues(evalData, normSet), peak->getValues(evalData, normSet), width->getValues(evalData, normSet), tail->getValues(evalData, normSet));
95}
96
97////////////////////////////////////////////////////////////////////////////////
98
99Int_t RooNovosibirsk::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* ) const
100{
101 if (matchArgs(allVars,analVars,x)) return 1 ;
102 if (matchArgs(allVars,analVars,peak)) return 2 ;
103
104 //The other two integrals over tali and width are not integrable
105
106 return 0 ;
107}
108
109////////////////////////////////////////////////////////////////////////////////
110
111Double_t RooNovosibirsk::analyticalIntegral(Int_t code, const char* rangeName) const
112{
113 assert(code==1 || code==2) ;
114
115 //The range is defined as [A,B]
116
117 //Numerical values need for the evaluation of the integral
118 static const Double_t sqrt2 = 1.4142135623730950; // Sqrt(2)
119 static const Double_t sqlog2 = 0.832554611157697756; //Sqrt( Log(2) )
120 static const Double_t sqlog4 = 1.17741002251547469; //Sqrt( Log(4) )
121 static const Double_t log4 = 1.38629436111989062; //Log(2)
122 static const Double_t rootpiby2 = 1.2533141373155003; // Sqrt(pi/2)
123 static const Double_t sqpibylog2 = 2.12893403886245236; //Sqrt( pi/Log(2) )
124
125 if (code==1) {
126 Double_t A = x.min(rangeName);
127 Double_t B = x.max(rangeName);
128
129 Double_t result = 0;
130
131
132 //If tail==0 the function becomes gaussian, thus we return a Gaussian integral
133 if (TMath::Abs(tail) < 1.e-7) {
134
135 Double_t xscale = sqrt2*width;
136
137 result = rootpiby2*width*(TMath::Erf((B-peak)/xscale)-TMath::Erf((A-peak)/xscale));
138
139 return result;
140
141 }
142
143 // If the range is not defined correctly the function becomes complex
144 Double_t log_argument_A = ( (peak - A)*tail + width ) / width ;
145 Double_t log_argument_B = ( (peak - B)*tail + width ) / width ;
146
147 //lower limit
148 if ( log_argument_A < 1.e-7) {
149 log_argument_A = 1.e-7;
150 }
151
152 //upper limit
153 if ( log_argument_B < 1.e-7) {
154 log_argument_B = 1.e-7;
155 }
156
157 Double_t term1 = TMath::ASinH( tail * sqlog4 );
158 Double_t term1_2 = term1 * term1;
159
160 //Calculate the error function arguments
161 Double_t erf_termA = ( term1_2 - log4 * TMath::Log( log_argument_A ) ) / ( 2 * term1 * sqlog2 );
162 Double_t erf_termB = ( term1_2 - log4 * TMath::Log( log_argument_B ) ) / ( 2 * term1 * sqlog2 );
163
164 result = 0.5 / tail * width * term1 * ( TMath::Erf(erf_termB) - TMath::Erf(erf_termA)) * sqpibylog2;
165
166 return result;
167
168 } else if (code==2) {
169 Double_t A = x.min(rangeName);
170 Double_t B = x.max(rangeName);
171
172 Double_t result = 0;
173
174
175 //If tail==0 the function becomes gaussian, thus we return a Gaussian integral
176 if (TMath::Abs(tail) < 1.e-7) {
177
178 Double_t xscale = sqrt2*width;
179
180 result = rootpiby2*width*(TMath::Erf((B-x)/xscale)-TMath::Erf((A-x)/xscale));
181
182 return result;
183
184 }
185
186 // If the range is not defined correctly the function becomes complex
187 Double_t log_argument_A = ( (A - x)*tail + width ) / width;
188 Double_t log_argument_B = ( (B - x)*tail + width ) / width;
189
190 //lower limit
191 if ( log_argument_A < 1.e-7) {
192 log_argument_A = 1.e-7;
193 }
194
195 //upper limit
196 if ( log_argument_B < 1.e-7) {
197 log_argument_B = 1.e-7;
198 }
199
200 Double_t term1 = TMath::ASinH( tail * sqlog4 );
201 Double_t term1_2 = term1 * term1;
202
203 //Calculate the error function arguments
204 Double_t erf_termA = ( term1_2 - log4 * TMath::Log( log_argument_A ) ) / ( 2 * term1 * sqlog2 );
205 Double_t erf_termB = ( term1_2 - log4 * TMath::Log( log_argument_B ) ) / ( 2 * term1 * sqlog2 );
206
207 result = 0.5 / tail * width * term1 * ( TMath::Erf(erf_termB) - TMath::Erf(erf_termA)) * sqpibylog2;
208
209 return result;
210
211 }
212
213 // Emit fatal error
214 coutF(Eval) << "Error in RooNovosibirsk::analyticalIntegral" << std::endl;
215
216 // Put dummy return here to avoid compiler warnings
217 return 1.0 ;
218}
#define coutF(a)
#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 log(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...
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
virtual RooSpan< double > computeNovosibirsk(const RooAbsReal *, RunContext &, RooSpan< const double > x, RooSpan< const double > peak, RooSpan< const double > width, RooSpan< const double > tail)=0
RooNovosibirsk implements the Novosibirsk function.
RooSpan< double > evaluateSpan(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet) const
Compute multiple values of Novosibirsk distribution.
RooRealProxy width
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
Double_t evaluate() const
If tail=eta=0 the Belle distribution becomes gaussian.
RooRealProxy peak
RooRealProxy tail
RooRealProxy x
A simple container to hold a batch of data values.
Definition RooSpan.h:34
double min(const char *rname=0) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
double max(const char *rname=0) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
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...
Double_t Exp(Double_t x)
Definition TMath.h:727
Double_t Erf(Double_t x)
Computation of the error function erf(x).
Definition TMath.cxx:184
Double_t ASinH(Double_t)
Definition TMath.cxx:64
Double_t Log(Double_t x)
Definition TMath.h:760
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition TMath.h:735
Short_t Abs(Short_t d)
Definition TMathBase.h:120
This struct enables passing computation data around between elements of a computation graph.
Definition RunContext.h:31