RooArgusBG.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitModels                                                     *
00004  * @(#)root/roofit:$Id: RooArgusBG.cxx 24286 2008-06-16 15:47:04Z 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 // RooArgusBg is a RooAbsPdf implementation describing the ARGUS background 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 "RooArgusBG.h"
00031 #include "RooRealVar.h"
00032 #include "RooRealConstant.h"
00033 #include "RooMath.h"
00034 #include "TMath.h"
00035 
00036 ClassImp(RooArgusBG)
00037 
00038 
00039 //_____________________________________________________________________________
00040 RooArgusBG::RooArgusBG(const char *name, const char *title,
00041                        RooAbsReal& _m, RooAbsReal& _m0, RooAbsReal& _c) :
00042   RooAbsPdf(name, title), 
00043   m("m","Mass",this,_m),
00044   m0("m0","Resonance mass",this,_m0),
00045   c("c","Slope parameter",this,_c),
00046   p("p","Power",this,(RooRealVar&)RooRealConstant::value(0.5))
00047 {
00048 }
00049 
00050 
00051 //_____________________________________________________________________________
00052 RooArgusBG::RooArgusBG(const char *name, const char *title,
00053                        RooAbsReal& _m, RooAbsReal& _m0, RooAbsReal& _c, RooAbsReal& _p) :
00054   RooAbsPdf(name, title), 
00055   m("m","Mass",this,_m),
00056   m0("m0","Resonance mass",this,_m0),
00057   c("c","Slope parameter",this,_c),
00058   p("p","Power",this,_p)
00059 {
00060 }
00061 
00062 
00063 //_____________________________________________________________________________
00064 RooArgusBG::RooArgusBG(const RooArgusBG& other, const char* name) :
00065   RooAbsPdf(other,name), 
00066   m("m",this,other.m), 
00067   m0("m0",this,other.m0), 
00068   c("c",this,other.c),
00069   p("p",this,other.p)
00070 {
00071 }
00072 
00073 
00074 
00075 //_____________________________________________________________________________
00076 Double_t RooArgusBG::evaluate() const {
00077   Double_t t= m/m0;
00078   if(t >= 1) return 0;
00079 
00080   Double_t u= 1 - t*t;
00081   //cout << "c = " << c << " result = " << m*TMath::Power(u,p)*exp(c*u) << endl ; 
00082   return m*TMath::Power(u,p)*exp(c*u) ;
00083 }
00084 
00085 
00086 
00087 //_____________________________________________________________________________
00088 Int_t RooArgusBG::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
00089 {
00090   if (p.arg().isConstant()) {
00091     // We can integrate over m if power = 0.5
00092     if (matchArgs(allVars,analVars,m) && p == 0.5) return 1;
00093   }
00094   return 0;
00095 
00096 }
00097 
00098 
00099 
00100 //_____________________________________________________________________________
00101 Double_t RooArgusBG::analyticalIntegral(Int_t code, const char* rangeName) const
00102 {
00103   assert(code==1);
00104   // Formula for integration over m when p=0.5
00105   static const Double_t pi = atan2(0.0,-1.0);
00106   Double_t min = (m.min(rangeName) < m0) ? m.min(rangeName) : m0;
00107   Double_t max = (m.max(rangeName) < m0) ? m.max(rangeName) : m0;
00108   Double_t f1 = (1.-TMath::Power(min/m0,2));
00109   Double_t f2 = (1.-TMath::Power(max/m0,2));
00110   Double_t aLow, aHigh ;
00111   aLow  = -0.5*m0*m0*(exp(c*f1)*sqrt(f1)/c + 0.5/TMath::Power(-c,1.5)*sqrt(pi)*RooMath::erf(sqrt(-c*f1)));
00112   aHigh = -0.5*m0*m0*(exp(c*f2)*sqrt(f2)/c + 0.5/TMath::Power(-c,1.5)*sqrt(pi)*RooMath::erf(sqrt(-c*f2)));
00113   Double_t area = aHigh - aLow;
00114   //cout << "c = " << c << "aHigh = " << aHigh << " aLow = " << aLow << " area = " << area << endl ;
00115   return area;
00116 
00117 }
00118 
00119 

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