RooCBShape.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitModels                                                     *
00004  * @(#)root/roofit:$Id: RooCBShape.cxx 28259 2009-04-16 16:21:16Z wouter $
00005  * Authors:                                                                  *
00006  *   WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu       *
00007  *   DK, David Kirkby,    UC Irvine,         dkirkby@uci.edu                 *
00008  *                                                                           *
00009  * Copyright (c) 2000-2005, Regents of the University of California          *
00010  *                          and Stanford University. All rights reserved.    *
00011  *                                                                           *
00012  * Redistribution and use in source and binary forms,                        *
00013  * with or without modification, are permitted according to the terms        *
00014  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00015  *****************************************************************************/
00016 
00017 //////////////////////////////////////////////////////////////////////////////
00018 //
00019 // BEGIN_HTML
00020 // P.d.f implementing the Crystall Ball line shape
00021 // END_HTML
00022 //
00023 
00024 #include "RooFit.h"
00025 
00026 #include "Riostream.h"
00027 #include "Riostream.h"
00028 #include <math.h>
00029 
00030 #include "RooCBShape.h"
00031 #include "RooAbsReal.h"
00032 #include "RooRealVar.h"
00033 #include "RooMath.h"
00034 #include "TMath.h"
00035 
00036 ClassImp(RooCBShape)
00037 ;
00038 
00039   
00040 
00041 //_____________________________________________________________________________
00042 Double_t RooCBShape::ApproxErf(Double_t arg) const 
00043 {
00044   static const double erflim = 5.0;
00045   if( arg > erflim )
00046     return 1.0;
00047   if( arg < -erflim )
00048     return -1.0;
00049   
00050   return RooMath::erf(arg);
00051 }
00052 
00053 
00054 static const char rcsid[] =
00055 "$Id: RooCBShape.cxx 28259 2009-04-16 16:21:16Z wouter $";
00056 
00057 
00058 //_____________________________________________________________________________
00059 RooCBShape::RooCBShape(const char *name, const char *title,
00060                        RooAbsReal& _m, RooAbsReal& _m0, RooAbsReal& _sigma,
00061                        RooAbsReal& _alpha, RooAbsReal& _n) :
00062   RooAbsPdf(name, title),
00063   m("m", "Dependent", this, _m),
00064   m0("m0", "M0", this, _m0),
00065   sigma("sigma", "Sigma", this, _sigma),
00066   alpha("alpha", "Alpha", this, _alpha),
00067   n("n", "Order", this, _n)
00068 {
00069 }
00070 
00071 
00072 //_____________________________________________________________________________
00073 RooCBShape::RooCBShape(const RooCBShape& other, const char* name) :
00074   RooAbsPdf(other, name), m("m", this, other.m), m0("m0", this, other.m0),
00075   sigma("sigma", this, other.sigma), alpha("alpha", this, other.alpha),
00076   n("n", this, other.n)
00077 {
00078 }
00079 
00080 
00081 //_____________________________________________________________________________
00082 Double_t RooCBShape::evaluate() const {
00083 
00084   Double_t t = (m-m0)/sigma;
00085   if (alpha < 0) t = -t;
00086 
00087   Double_t absAlpha = fabs((Double_t)alpha);
00088 
00089   if (t >= -absAlpha) {
00090     return exp(-0.5*t*t);
00091   }
00092   else {
00093     Double_t a =  TMath::Power(n/absAlpha,n)*exp(-0.5*absAlpha*absAlpha);
00094     Double_t b= n/absAlpha - absAlpha; 
00095 
00096     return a/TMath::Power(b - t, n);
00097   }
00098 }
00099 
00100 
00101 //_____________________________________________________________________________
00102 Int_t RooCBShape::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
00103 {
00104   if( matchArgs(allVars,analVars,m) )
00105     return 1 ;
00106   
00107   return 0;
00108 }
00109 
00110 
00111 
00112 //_____________________________________________________________________________
00113 Double_t RooCBShape::analyticalIntegral(Int_t code, const char* rangeName) const
00114 {
00115   static const double sqrtPiOver2 = 1.2533141373;
00116   static const double sqrt2 = 1.4142135624;
00117 
00118   assert(code==1);
00119   double result = 0.0;
00120   bool useLog = false;
00121   
00122   if( fabs(n-1.0) < 1.0e-05 )
00123     useLog = true;
00124   
00125   double sig = fabs((Double_t)sigma);
00126   
00127   double tmin = (m.min(rangeName)-m0)/sig;
00128   double tmax = (m.max(rangeName)-m0)/sig;
00129   
00130   if(alpha < 0) {
00131     double tmp = tmin;
00132     tmin = -tmax;
00133     tmax = -tmp;
00134   }
00135 
00136   double absAlpha = fabs((Double_t)alpha);
00137   
00138   if( tmin >= -absAlpha ) {
00139     result += sig*sqrtPiOver2*(   ApproxErf(tmax/sqrt2)
00140                                 - ApproxErf(tmin/sqrt2) );
00141   }
00142   else if( tmax <= -absAlpha ) {
00143     double a = TMath::Power(n/absAlpha,n)*exp(-0.5*absAlpha*absAlpha);
00144     double b = n/absAlpha - absAlpha;
00145     
00146     if(useLog) {
00147       result += a*sig*( log(b-tmin) - log(b-tmax) );
00148     }
00149     else {
00150       result += a*sig/(1.0-n)*(   1.0/(TMath::Power(b-tmin,n-1.0))
00151                                 - 1.0/(TMath::Power(b-tmax,n-1.0)) );
00152     }
00153   }
00154   else {
00155     double a = TMath::Power(n/absAlpha,n)*exp(-0.5*absAlpha*absAlpha);
00156     double b = n/absAlpha - absAlpha;
00157     
00158     double term1 = 0.0;
00159     if(useLog) {
00160       term1 = a*sig*(  log(b-tmin) - log(n/absAlpha));
00161     }
00162     else {
00163       term1 = a*sig/(1.0-n)*(   1.0/(TMath::Power(b-tmin,n-1.0))
00164                               - 1.0/(TMath::Power(n/absAlpha,n-1.0)) );
00165     }
00166     
00167     double term2 = sig*sqrtPiOver2*(   ApproxErf(tmax/sqrt2)
00168                                      - ApproxErf(-absAlpha/sqrt2) );
00169     
00170     
00171     result += term1 + term2;
00172   }
00173   
00174   return result;
00175 }
00176 
00177 
00178 
00179 //_____________________________________________________________________________
00180 Int_t RooCBShape::getMaxVal(const RooArgSet& vars) const 
00181 {
00182   // Advertise that we know the maximum of self for given (m0,alpha,n,sigma)
00183   RooArgSet dummy ;
00184 
00185   if (matchArgs(vars,dummy,m)) {
00186     return 1 ;  
00187   }
00188   return 0 ;  
00189 }
00190 
00191 
00192 
00193 //_____________________________________________________________________________
00194 Double_t RooCBShape::maxVal(Int_t code) const
00195 {
00196   assert(code==1) ;
00197 
00198   // The maximum value for given (m0,alpha,n,sigma)
00199   return 1.0 ;
00200 }
00201 
00202 

Generated on Tue Jul 5 14:54:35 2011 for ROOT_528-00b_version by  doxygen 1.5.1