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 "RooRealVar.h"
26#include "RooMath.h"
27#include "RooBatchCompute.h"
28
29#include "TMath.h"
30
31#include <cmath>
32
34
35////////////////////////////////////////////////////////////////////////////////
36
37double RooCBShape::ApproxErf(double arg) const
38{
39 static const double erflim = 5.0;
40 if( arg > erflim )
41 return 1.0;
42 if( arg < -erflim )
43 return -1.0;
44
45 return RooMath::erf(arg);
46}
47
48////////////////////////////////////////////////////////////////////////////////
49
50RooCBShape::RooCBShape(const char *name, const char *title,
51 RooAbsReal& _m, RooAbsReal& _m0, RooAbsReal& _sigma,
52 RooAbsReal& _alpha, RooAbsReal& _n) :
53 RooAbsPdf(name, title),
54 m("m", "Dependent", this, _m),
55 m0("m0", "M0", this, _m0),
56 sigma("sigma", "Sigma", this, _sigma),
57 alpha("alpha", "Alpha", this, _alpha),
58 n("n", "Order", this, _n)
59{
60}
61
62////////////////////////////////////////////////////////////////////////////////
63
64RooCBShape::RooCBShape(const RooCBShape& other, const char* name) :
65 RooAbsPdf(other, name), m("m", this, other.m), m0("m0", this, other.m0),
66 sigma("sigma", this, other.sigma), alpha("alpha", this, other.alpha),
67 n("n", this, other.n)
68{
69}
70
71////////////////////////////////////////////////////////////////////////////////
72
73double RooCBShape::evaluate() const {
74 double t = (m-m0)/sigma;
75 if (alpha < 0) t = -t;
76
77 double absAlpha = std::abs((double)alpha);
78
79 if (t >= -absAlpha) {
80 return exp(-0.5*t*t);
81 }
82 else {
83 double a = TMath::Power(n/absAlpha,n)*exp(-0.5*absAlpha*absAlpha);
84 double b= n/absAlpha - absAlpha;
85
86 return a/TMath::Power(b - t, n);
87 }
88}
89
90////////////////////////////////////////////////////////////////////////////////
91/// Compute multiple values of Crystal ball Shape distribution.
93{
95 {ctx.at(m), ctx.at(m0), ctx.at(sigma), ctx.at(alpha), ctx.at(n)});
96}
97
98////////////////////////////////////////////////////////////////////////////////
99
100Int_t RooCBShape::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
101{
102 if( matchArgs(allVars,analVars,m) )
103 return 1 ;
104
105 return 0;
106}
107
108////////////////////////////////////////////////////////////////////////////////
109
110double RooCBShape::analyticalIntegral(Int_t code, const char* rangeName) const
111{
112 static const double sqrtPiOver2 = 1.2533141373;
113 static const double sqrt2 = 1.4142135624;
114
115 R__ASSERT(code==1);
116 double result = 0.0;
117 bool useLog = false;
118
119 if( std::abs(n-1.0) < 1.0e-05 )
120 useLog = true;
121
122 double sig = std::abs((double)sigma);
123
124 double tmin = (m.min(rangeName)-m0)/sig;
125 double tmax = (m.max(rangeName)-m0)/sig;
126
127 if(alpha < 0) {
128 double tmp = tmin;
129 tmin = -tmax;
130 tmax = -tmp;
131 }
132
133 double absAlpha = std::abs((double)alpha);
134
135 if( tmin >= -absAlpha ) {
136 result += sig*sqrtPiOver2*( ApproxErf(tmax/sqrt2)
137 - ApproxErf(tmin/sqrt2) );
138 }
139 else if( tmax <= -absAlpha ) {
140 double a = TMath::Power(n/absAlpha,n)*exp(-0.5*absAlpha*absAlpha);
141 double b = n/absAlpha - absAlpha;
142
143 if(useLog) {
144 result += a*sig*( log(b-tmin) - log(b-tmax) );
145 }
146 else {
147 result += a*sig/(1.0-n)*( 1.0/(TMath::Power(b-tmin,n-1.0))
148 - 1.0/(TMath::Power(b-tmax,n-1.0)) );
149 }
150 }
151 else {
152 double a = TMath::Power(n/absAlpha,n)*exp(-0.5*absAlpha*absAlpha);
153 double b = n/absAlpha - absAlpha;
154
155 double term1 = 0.0;
156 if(useLog) {
157 term1 = a*sig*( log(b-tmin) - log(n/absAlpha));
158 }
159 else {
160 term1 = a*sig/(1.0-n)*( 1.0/(TMath::Power(b-tmin,n-1.0))
161 - 1.0/(TMath::Power(n/absAlpha,n-1.0)) );
162 }
163
164 double term2 = sig*sqrtPiOver2*( ApproxErf(tmax/sqrt2)
165 - ApproxErf(-absAlpha/sqrt2) );
166
167
168 result += term1 + term2;
169 }
170
171 return result != 0. ? result : 1.E-300;
172}
173
174////////////////////////////////////////////////////////////////////////////////
175/// Advertise that we know the maximum of self for given (m0,alpha,n,sigma)
176
178{
179 RooArgSet dummy ;
180
181 if (matchArgs(vars,dummy,m)) {
182 return 1 ;
183 }
184 return 0 ;
185}
186
187////////////////////////////////////////////////////////////////////////////////
188
189double RooCBShape::maxVal(Int_t code) const
190{
191 R__ASSERT(code==1) ;
192
193 // The maximum value for given (m0,alpha,n,sigma)
194 // is 1./ Integral in the variable range
195 return 1.0/analyticalIntegral(1) ;
196}
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
#define ClassImp(name)
Definition Rtypes.h:377
#define R__ASSERT(e)
Definition TError.h:118
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
char name[80]
Definition TGX11.cxx:110
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
bool 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:55
PDF implementing the Crystal Ball line shape.
Definition RooCBShape.h:24
double ApproxErf(double arg) const
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
void doEval(RooFit::EvalContext &) const override
Compute multiple values of Crystal ball Shape distribution.
RooRealProxy n
Definition RooCBShape.h:49
RooRealProxy m0
Definition RooCBShape.h:46
RooRealProxy m
Definition RooCBShape.h:45
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooRealProxy sigma
Definition RooCBShape.h:47
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooRealProxy alpha
Definition RooCBShape.h:48
Int_t getMaxVal(const RooArgSet &vars) const override
Advertise that we know the maximum of self for given (m0,alpha,n,sigma)
double maxVal(Int_t code) const override
Return maximum value for set of observables identified by code assigned in getMaxVal.
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
static std::complex< double > erf(const std::complex< double > z)
complex erf function
Definition RooMath.cxx:59
double max(const char *rname=nullptr) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
const Double_t sigma
const Int_t n
Definition legend1.C:16
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:721
TMarker m
Definition textangle.C:8