RooBukinPdf.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitModels                                                     *
00004  * @(#)root/roofit:$Id: RooBukinPdf.cxx 24286 2008-06-16 15:47:04Z wouter $
00005  * Authors:                                                                  *
00006  *   RW, Ruddick William  UC Colorado        wor@slac.stanford.edu           *
00007  *                                                                           *
00008  * Copyright (c) 2000-2005, Regents of the University of California          *
00009  *                          and Stanford University. All rights reserved.    *
00010  *                                                                           *
00011  * Redistribution and use in source and binary forms,                        *
00012  * with or without modification, are permitted according to the terms        *
00013  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00014  *****************************************************************************/
00015 
00016 //////////////////////////////////////////////////////////////////////////////
00017 //
00018 // BEGIN_HTML
00019 // RooBukinPdf implements the NovosibirskA function 
00020 // END_HTML
00021 //
00022 
00023 // Original Fortran Header below
00024 /*****************************************************************************
00025  * Fitting function for asymmetric peaks with 6 free parameters:             *
00026  *     Ap   - peak value                                                     *
00027  *     Xp   - peak position                                                  *
00028  *     sigp - FWHM divided by 2*sqrt(2*log(2))=2.35                          *
00029  *     xi   - peak asymmetry parameter                                       *
00030  *     rho1 - parameter of the "left tail"                                   *
00031  *     rho2 - parameter of the "right tail"                                  *
00032  *   ---------------------------------------------                           *
00033  *       May 26, 2003                                                        *
00034  *       A.Bukin, Budker INP, Novosibirsk                                    *
00035  *       Documentation:                                                      *
00036  *       http://www.slac.stanford.edu/BFROOT/www/Organization/CollabMtgs/2003/detJuly2003/Tues3a/bukin.ps 
00037  *   -------------------------------------------                             *
00038  *****************************************************************************/
00039 
00040 #include "RooFit.h"
00041 
00042 #include <math.h>
00043 
00044 
00045 #include "RooBukinPdf.h"
00046 #include "RooRealVar.h"
00047 #include "TMath.h"
00048 
00049 ClassImp(RooBukinPdf)
00050 
00051 
00052 
00053 //_____________________________________________________________________________
00054 RooBukinPdf::RooBukinPdf(const char *name, const char *title,
00055                          RooAbsReal& _x,    RooAbsReal& _Xp,
00056                          RooAbsReal& _sigp, RooAbsReal& _xi,
00057                          RooAbsReal& _rho1, RooAbsReal& _rho2) :
00058   // The two addresses refer to our first dependent variable and
00059   // parameter, respectively, as declared in the rdl file
00060   RooAbsPdf(name, title),
00061   x("x","x",this,_x),
00062   Xp("Xp","Xp",this,_Xp),
00063   sigp("sigp","sigp",this,_sigp),
00064   xi("xi","xi",this,_xi),
00065   rho1("rho1","rho1",this,_rho1),
00066   rho2("rho2","rho2",this,_rho2)
00067 {
00068   // Constructor
00069   consts = 2*sqrt(2*log(2.));
00070 }
00071 
00072 
00073 
00074 //_____________________________________________________________________________
00075 RooBukinPdf::RooBukinPdf(const RooBukinPdf& other, const char *name):
00076   RooAbsPdf(other,name),
00077   x("x",this,other.x),
00078   Xp("Xp",this,other.Xp),
00079   sigp("sigp",this,other.sigp),
00080   xi("xi",this,other.xi),
00081   rho1("rho1",this,other.rho1),
00082   rho2("rho2",this,other.rho2)
00083 
00084 {
00085   // Copy constructor
00086   consts = 2*sqrt(2*log(2.));
00087 }
00088 
00089 
00090 
00091 //_____________________________________________________________________________
00092 Double_t RooBukinPdf::evaluate() const 
00093 {
00094   // Implementation 
00095 
00096   double r1=0,r2=0,r3=0,r4=0,r5=0,hp=0;
00097   double x1 = 0,x2 = 0;
00098   double fit_result = 0;
00099   
00100   hp=sigp*consts;
00101   r3=log(2.);
00102   r4=sqrt(TMath::Power(xi,2)+1);
00103   r1=xi/r4;  
00104 
00105   if(TMath::Abs(xi) > exp(-6.)){
00106     r5=xi/log(r4+xi);
00107   }
00108   else
00109     r5=1;
00110     
00111   x1 = Xp + (hp / 2) * (r1-1);
00112   x2 = Xp + (hp / 2) * (r1+1);
00113   
00114   //--- Left Side
00115   if(x < x1){
00116     r2=rho1*TMath::Power((x-x1)/(Xp-x1),2)-r3 + 4 * r3 * (x-x1)/hp * r5 * r4/TMath::Power((r4-xi),2);
00117   }
00118 
00119 
00120   //--- Center
00121   else if(x < x2) {
00122     if(TMath::Abs(xi) > exp(-6.)) {
00123       r2=log(1 + 4 * xi * r4 * (x-Xp)/hp)/log(1+2*xi*(xi-r4));
00124       r2=-r3*(TMath::Power(r2,2));
00125     }
00126     else{
00127       r2=-4*r3*TMath::Power(((x-Xp)/hp),2);  
00128     }
00129   }
00130   
00131 
00132   //--- Right Side
00133   else {
00134     r2=rho2*TMath::Power((x-x2)/(Xp-x2),2)-r3 - 4 * r3 * (x-x2)/hp * r5 * r4/TMath::Power((r4+xi),2);
00135   }
00136 
00137   if(TMath::Abs(r2) > 100){
00138     fit_result = 0;  
00139   }
00140   else{
00141     //---- Normalize the result
00142     fit_result = exp(r2);
00143   }
00144   
00145   return fit_result;
00146   
00147 }

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