ROOT  6.06/09
Reference Guide
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 //////////////////////////////////////////////////////////////////////////////
18 //
19 // BEGIN_HTML
20 // P.d.f implementing the Crystall Ball line shape
21 // END_HTML
22 //
23 
24 #include "RooFit.h"
25 
26 #include "Riostream.h"
27 #include "Riostream.h"
28 #include <math.h>
29 
30 #include "RooCBShape.h"
31 #include "RooAbsReal.h"
32 #include "RooRealVar.h"
33 #include "RooMath.h"
34 #include "TMath.h"
35 
36 #include "TError.h"
37 
38 using namespace std;
39 
41 ;
42 
43 
44 
45 ////////////////////////////////////////////////////////////////////////////////
46 
48 {
49  static const double erflim = 5.0;
50  if( arg > erflim )
51  return 1.0;
52  if( arg < -erflim )
53  return -1.0;
54 
55  return RooMath::erf(arg);
56 }
57 
58 
59 
60 ////////////////////////////////////////////////////////////////////////////////
61 
62 RooCBShape::RooCBShape(const char *name, const char *title,
63  RooAbsReal& _m, RooAbsReal& _m0, RooAbsReal& _sigma,
64  RooAbsReal& _alpha, RooAbsReal& _n) :
65  RooAbsPdf(name, title),
66  m("m", "Dependent", this, _m),
67  m0("m0", "M0", this, _m0),
68  sigma("sigma", "Sigma", this, _sigma),
69  alpha("alpha", "Alpha", this, _alpha),
70  n("n", "Order", this, _n)
71 {
72 }
73 
74 
75 ////////////////////////////////////////////////////////////////////////////////
76 
77 RooCBShape::RooCBShape(const RooCBShape& other, const char* name) :
78  RooAbsPdf(other, name), m("m", this, other.m), m0("m0", this, other.m0),
79  sigma("sigma", this, other.sigma), alpha("alpha", this, other.alpha),
80  n("n", this, other.n)
81 {
82 }
83 
84 
85 ////////////////////////////////////////////////////////////////////////////////
86 
88  Double_t t = (m-m0)/sigma;
89  if (alpha < 0) t = -t;
90 
91  Double_t absAlpha = fabs((Double_t)alpha);
92 
93  if (t >= -absAlpha) {
94  return exp(-0.5*t*t);
95  }
96  else {
97  Double_t a = TMath::Power(n/absAlpha,n)*exp(-0.5*absAlpha*absAlpha);
98  Double_t b= n/absAlpha - absAlpha;
99 
100  return a/TMath::Power(b - t, n);
101  }
102 }
103 
104 
105 ////////////////////////////////////////////////////////////////////////////////
106 
107 Int_t RooCBShape::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
108 {
109  if( matchArgs(allVars,analVars,m) )
110  return 1 ;
111 
112  return 0;
113 }
114 
115 
116 
117 ////////////////////////////////////////////////////////////////////////////////
118 
119 Double_t RooCBShape::analyticalIntegral(Int_t code, const char* rangeName) const
120 {
121  static const double sqrtPiOver2 = 1.2533141373;
122  static const double sqrt2 = 1.4142135624;
123 
124  R__ASSERT(code==1);
125  double result = 0.0;
126  bool useLog = false;
127 
128  if( fabs(n-1.0) < 1.0e-05 )
129  useLog = true;
130 
131  double sig = fabs((Double_t)sigma);
132 
133  double tmin = (m.min(rangeName)-m0)/sig;
134  double tmax = (m.max(rangeName)-m0)/sig;
135 
136  if(alpha < 0) {
137  double tmp = tmin;
138  tmin = -tmax;
139  tmax = -tmp;
140  }
141 
142  double absAlpha = fabs((Double_t)alpha);
143 
144  if( tmin >= -absAlpha ) {
145  result += sig*sqrtPiOver2*( ApproxErf(tmax/sqrt2)
146  - ApproxErf(tmin/sqrt2) );
147  }
148  else if( tmax <= -absAlpha ) {
149  double a = TMath::Power(n/absAlpha,n)*exp(-0.5*absAlpha*absAlpha);
150  double b = n/absAlpha - absAlpha;
151 
152  if(useLog) {
153  result += a*sig*( log(b-tmin) - log(b-tmax) );
154  }
155  else {
156  result += a*sig/(1.0-n)*( 1.0/(TMath::Power(b-tmin,n-1.0))
157  - 1.0/(TMath::Power(b-tmax,n-1.0)) );
158  }
159  }
160  else {
161  double a = TMath::Power(n/absAlpha,n)*exp(-0.5*absAlpha*absAlpha);
162  double b = n/absAlpha - absAlpha;
163 
164  double term1 = 0.0;
165  if(useLog) {
166  term1 = a*sig*( log(b-tmin) - log(n/absAlpha));
167  }
168  else {
169  term1 = a*sig/(1.0-n)*( 1.0/(TMath::Power(b-tmin,n-1.0))
170  - 1.0/(TMath::Power(n/absAlpha,n-1.0)) );
171  }
172 
173  double term2 = sig*sqrtPiOver2*( ApproxErf(tmax/sqrt2)
174  - ApproxErf(-absAlpha/sqrt2) );
175 
176 
177  result += term1 + term2;
178  }
179 
180  return result;
181 }
182 
183 
184 
185 ////////////////////////////////////////////////////////////////////////////////
186 /// Advertise that we know the maximum of self for given (m0,alpha,n,sigma)
187 
189 {
190  RooArgSet dummy ;
191 
192  if (matchArgs(vars,dummy,m)) {
193  return 1 ;
194  }
195  return 0 ;
196 }
197 
198 
199 
200 ////////////////////////////////////////////////////////////////////////////////
201 
203 {
204  R__ASSERT(code==1) ;
205 
206  // The maximum value for given (m0,alpha,n,sigma)
207  return 1.0 ;
208 }
209 
210 
RooRealProxy m
Definition: RooCBShape.h:47
RooRealProxy m0
Definition: RooCBShape.h:48
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise that we know the maximum of self for given (m0,alpha,n,sigma)
Definition: RooCBShape.cxx:188
#define R__ASSERT(e)
Definition: TError.h:98
int Int_t
Definition: RtypesCore.h:41
TArc * a
Definition: textangle.C:12
STL namespace.
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
virtual Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
Definition: RooCBShape.cxx:119
Double_t evaluate() const
Definition: RooCBShape.cxx:87
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
Definition: RooCBShape.cxx:202
TMarker * m
Definition: textangle.C:8
virtual Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
Definition: RooCBShape.cxx:107
ClassImp(RooCBShape)
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooRealProxy sigma
Definition: RooCBShape.h:49
static RooMathCoreReg dummy
#define name(a, b)
Definition: linkTestLib0.cpp:5
RooRealProxy alpha
Definition: RooCBShape.h:50
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
RooRealProxy n
Definition: RooCBShape.h:51
static std::complex< double > erf(const std::complex< double > z)
complex erf function
Definition: RooMath.cxx:584
double result[121]
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
double exp(double)
const Int_t n
Definition: legend1.C:16
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
double log(double)
Double_t ApproxErf(Double_t arg) const
Definition: RooCBShape.cxx:47