RooLognormal.cxx

Go to the documentation of this file.
00001  /***************************************************************************** 
00002   * Project: RooFit                                                           * 
00003   * @(#)root/roofit:$Id: RooLognormal.cxx 34064 2010-06-22 15:05:19Z wouter $ *
00004   *                                                                           * 
00005   * RooFit Lognormal PDF                                                      *
00006   *                                                                           * 
00007   * Author: Gregory Schott and Stefan Schmitz                                 *
00008   *                                                                           * 
00009   *****************************************************************************/ 
00010 
00011 //////////////////////////////////////////////////////////////////////////////
00012 //
00013 // BEGIN_HTML
00014 // RooFit Lognormal PDF. The two parameters are:
00015 //  - m0: the median of the distribution
00016 //  - k=exp(sigma): sigma is called the shape parameter in the TMath parametrization
00017 //
00018 //   Lognormal(x,m0,k) = 1/(sqrt(2*pi)*ln(k)*x)*exp(-ln^2(x/m0)/(2*ln^2(k)))
00019 //
00020 // The parametrization here is physics driven and differs from the ROOT::Math::lognormal_pdf(x,m,s,x0) with:
00021 //  - m = log(m0)
00022 //  - s = log(k)
00023 //  - x0 = 0
00024 // END_HTML
00025 
00026 
00027 #include "RooFit.h"
00028 
00029 #include "Riostream.h"
00030 #include "Riostream.h"
00031 #include <math.h>
00032 
00033 #include "RooLognormal.h"
00034 #include "RooAbsReal.h"
00035 #include "RooRealVar.h"
00036 #include "RooRandom.h"
00037 #include "RooMath.h"
00038 #include "TMath.h"
00039 
00040 #include <Math/SpecFuncMathCore.h>
00041 #include <Math/PdfFuncMathCore.h>
00042 #include <Math/ProbFuncMathCore.h>
00043 
00044 ClassImp(RooLognormal)
00045 
00046 
00047 //_____________________________________________________________________________
00048 RooLognormal::RooLognormal(const char *name, const char *title,
00049                          RooAbsReal& _x, RooAbsReal& _m0,
00050                          RooAbsReal& _k) :
00051   RooAbsPdf(name,title),
00052   x("x","Observable",this,_x),
00053   m0("m0","m0",this,_m0),
00054   k("k","k",this,_k)
00055 {
00056 }
00057 
00058 
00059 
00060 //_____________________________________________________________________________
00061 RooLognormal::RooLognormal(const RooLognormal& other, const char* name) : 
00062   RooAbsPdf(other,name), x("x",this,other.x), m0("m0",this,other.m0),
00063   k("k",this,other.k)
00064 {
00065 }
00066 
00067 
00068 
00069 //_____________________________________________________________________________
00070 Double_t RooLognormal::evaluate() const
00071 {
00072   // ln(k)<1 would correspond to sigma < 0 in the parametrization
00073   // resulting by transforming a normal random variable in its
00074   // standard parametrization to a lognormal random variable
00075   // => treat ln(k) as -ln(k) for k<1
00076 
00077   Double_t xv = x;
00078   Double_t ln_k = TMath::Abs(TMath::Log(k));
00079   Double_t ln_m0 = TMath::Log(m0);
00080   Double_t x0 = 0;
00081 
00082   Double_t ret = ROOT::Math::lognormal_pdf(xv,ln_m0,ln_k,x0);
00083   return ret ;
00084 }
00085 
00086 
00087 
00088 //_____________________________________________________________________________
00089 Int_t RooLognormal::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const 
00090 {
00091   if (matchArgs(allVars,analVars,x)) return 1 ;
00092   return 0 ;
00093 }
00094 
00095 
00096 
00097 //_____________________________________________________________________________
00098 Double_t RooLognormal::analyticalIntegral(Int_t code, const char* rangeName) const 
00099 {
00100   assert(code==1) ;
00101 
00102   static const Double_t root2 = sqrt(2.) ;
00103 
00104   Double_t ln_k = TMath::Abs(TMath::Log(k));
00105   Double_t ret = 0.5*( RooMath::erf( TMath::Log(x.max(rangeName)/m0)/(root2*ln_k) ) - RooMath::erf( TMath::Log(x.min(rangeName)/m0)/(root2*ln_k) ) ) ;
00106 
00107   return ret ;
00108 }
00109 
00110 
00111 
00112 
00113 //_____________________________________________________________________________
00114 Int_t RooLognormal::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
00115 {
00116   if (matchArgs(directVars,generateVars,x)) return 1 ;  
00117   return 0 ;
00118 }
00119 
00120 
00121 
00122 //_____________________________________________________________________________
00123 void RooLognormal::generateEvent(Int_t code)
00124 {
00125   assert(code==1) ;
00126 
00127   Double_t xgen ;
00128   while(1) {    
00129     xgen = TMath::Exp(RooRandom::randomGenerator()->Gaus(TMath::Log(m0),TMath::Log(k)));
00130     if (xgen<=x.max() && xgen>=x.min()) {
00131       x = xgen ;
00132       break;
00133     }
00134   }
00135 
00136   return;
00137 }
00138 
00139 

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