Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooCBShape.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17/** \class RooCBShape
18 \ingroup Roofit
19
20PDF implementing the Crystal Ball line shape.
21**/
22
23#include "RooCBShape.h"
24
25#include "RooAbsReal.h"
26#include "RooRealVar.h"
27#include "RooMath.h"
28#include "RooBatchCompute.h"
29
30#include "TMath.h"
31
32#include <cmath>
33
34using namespace std;
35
37
38////////////////////////////////////////////////////////////////////////////////
39
41{
42 static const double erflim = 5.0;
43 if( arg > erflim )
44 return 1.0;
45 if( arg < -erflim )
46 return -1.0;
47
48 return RooMath::erf(arg);
49}
50
51////////////////////////////////////////////////////////////////////////////////
52
53RooCBShape::RooCBShape(const char *name, const char *title,
54 RooAbsReal& _m, RooAbsReal& _m0, RooAbsReal& _sigma,
55 RooAbsReal& _alpha, RooAbsReal& _n) :
56 RooAbsPdf(name, title),
57 m("m", "Dependent", this, _m),
58 m0("m0", "M0", this, _m0),
59 sigma("sigma", "Sigma", this, _sigma),
60 alpha("alpha", "Alpha", this, _alpha),
61 n("n", "Order", this, _n)
62{
63}
64
65////////////////////////////////////////////////////////////////////////////////
66
67RooCBShape::RooCBShape(const RooCBShape& other, const char* name) :
68 RooAbsPdf(other, name), m("m", this, other.m), m0("m0", this, other.m0),
69 sigma("sigma", this, other.sigma), alpha("alpha", this, other.alpha),
70 n("n", this, other.n)
71{
72}
73
74////////////////////////////////////////////////////////////////////////////////
75
77 Double_t t = (m-m0)/sigma;
78 if (alpha < 0) t = -t;
79
80 Double_t absAlpha = fabs((Double_t)alpha);
81
82 if (t >= -absAlpha) {
83 return exp(-0.5*t*t);
84 }
85 else {
86 Double_t a = TMath::Power(n/absAlpha,n)*exp(-0.5*absAlpha*absAlpha);
87 Double_t b= n/absAlpha - absAlpha;
88
89 return a/TMath::Power(b - t, n);
90 }
91}
92
93////////////////////////////////////////////////////////////////////////////////
94/// Compute multiple values of Crystal ball Shape distribution.
96 return RooBatchCompute::dispatch->computeCBShape(this, evalData, m->getValues(evalData, normSet), m0->getValues(evalData, normSet), sigma->getValues(evalData, normSet), alpha->getValues(evalData, normSet), n->getValues(evalData, normSet));
97}
98
99////////////////////////////////////////////////////////////////////////////////
100
101Int_t RooCBShape::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
102{
103 if( matchArgs(allVars,analVars,m) )
104 return 1 ;
105
106 return 0;
107}
108
109////////////////////////////////////////////////////////////////////////////////
110
111Double_t RooCBShape::analyticalIntegral(Int_t code, const char* rangeName) const
112{
113 static const double sqrtPiOver2 = 1.2533141373;
114 static const double sqrt2 = 1.4142135624;
115
116 R__ASSERT(code==1);
117 double result = 0.0;
118 bool useLog = false;
119
120 if( fabs(n-1.0) < 1.0e-05 )
121 useLog = true;
122
123 double sig = fabs((Double_t)sigma);
124
125 double tmin = (m.min(rangeName)-m0)/sig;
126 double tmax = (m.max(rangeName)-m0)/sig;
127
128 if(alpha < 0) {
129 double tmp = tmin;
130 tmin = -tmax;
131 tmax = -tmp;
132 }
133
134 double absAlpha = fabs((Double_t)alpha);
135
136 if( tmin >= -absAlpha ) {
137 result += sig*sqrtPiOver2*( ApproxErf(tmax/sqrt2)
138 - ApproxErf(tmin/sqrt2) );
139 }
140 else if( tmax <= -absAlpha ) {
141 double a = TMath::Power(n/absAlpha,n)*exp(-0.5*absAlpha*absAlpha);
142 double b = n/absAlpha - absAlpha;
143
144 if(useLog) {
145 result += a*sig*( log(b-tmin) - log(b-tmax) );
146 }
147 else {
148 result += a*sig/(1.0-n)*( 1.0/(TMath::Power(b-tmin,n-1.0))
149 - 1.0/(TMath::Power(b-tmax,n-1.0)) );
150 }
151 }
152 else {
153 double a = TMath::Power(n/absAlpha,n)*exp(-0.5*absAlpha*absAlpha);
154 double b = n/absAlpha - absAlpha;
155
156 double term1 = 0.0;
157 if(useLog) {
158 term1 = a*sig*( log(b-tmin) - log(n/absAlpha));
159 }
160 else {
161 term1 = a*sig/(1.0-n)*( 1.0/(TMath::Power(b-tmin,n-1.0))
162 - 1.0/(TMath::Power(n/absAlpha,n-1.0)) );
163 }
164
165 double term2 = sig*sqrtPiOver2*( ApproxErf(tmax/sqrt2)
166 - ApproxErf(-absAlpha/sqrt2) );
167
168
169 result += term1 + term2;
170 }
171
172 return result != 0. ? result : 1.E-300;
173}
174
175////////////////////////////////////////////////////////////////////////////////
176/// Advertise that we know the maximum of self for given (m0,alpha,n,sigma)
177
179{
180 RooArgSet dummy ;
181
182 if (matchArgs(vars,dummy,m)) {
183 return 1 ;
184 }
185 return 0 ;
186}
187
188////////////////////////////////////////////////////////////////////////////////
189
191{
192 R__ASSERT(code==1) ;
193
194 // The maximum value for given (m0,alpha,n,sigma)
195 // is 1./ Integral in the variable range
196 return 1.0/analyticalIntegral(1) ;
197}
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
#define ClassImp(name)
Definition Rtypes.h:364
#define R__ASSERT(e)
Definition TError.h:120
char name[80]
Definition TGX11.cxx:110
double exp(double)
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 > computeCBShape(const RooAbsReal *, RunContext &, RooSpan< const double > m, RooSpan< const double > m0, RooSpan< const double > sigma, RooSpan< const double > alpha, RooSpan< const double > n)=0
PDF implementing the Crystal Ball line shape.
Definition RooCBShape.h:24
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
virtual Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooRealProxy n
Definition RooCBShape.h:51
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Double_t ApproxErf(Double_t arg) const
RooRealProxy m0
Definition RooCBShape.h:48
RooRealProxy m
Definition RooCBShape.h:47
virtual Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise that we know the maximum of self for given (m0,alpha,n,sigma)
RooSpan< double > evaluateSpan(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet) const
Compute multiple values of Crystal ball Shape distribution.
RooRealProxy sigma
Definition RooCBShape.h:49
RooRealProxy alpha
Definition RooCBShape.h:50
static std::complex< double > erf(const std::complex< double > z)
complex erf function
Definition RooMath.cxx:580
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.
const Double_t sigma
const Int_t n
Definition legend1.C:16
R__EXTERN RooBatchComputeInterface * dispatch
This dispatch pointer points to an implementation of the compute library, provided one has been loade...
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition TMath.h:735
This struct enables passing computation data around between elements of a computation graph.
Definition RunContext.h:31
auto * m
Definition textangle.C:8