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