RooDstD0BG.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitModels                                                     *
00004  * @(#)root/roofit:$Id: RooDstD0BG.cxx 24286 2008-06-16 15:47:04Z wouter $
00005  * Authors:                                                                  *
00006  *   UE, Ulrik Egede,     RAL,               U.Egede@rl.ac.uk                *
00007  *   MT, Max Turri,       UC Santa Cruz      turri@slac.stanford.edu         *
00008  *   CC, Chih-hsiang Cheng, Stanford         chcheng@slac.stanford.edu       *
00009  *   DK, David Kirkby,    UC Irvine,         dkirkby@uci.edu                 *
00010  *                                                                           *
00011  * Copyright (c) 2000-2005, Regents of the University of California          *
00012  *                          RAL and Stanford University. All rights reserved.*
00013  *                                                                           *
00014  * Redistribution and use in source and binary forms,                        *
00015  * with or without modification, are permitted according to the terms        *
00016  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00017  *****************************************************************************/
00018 
00019 //////////////////////////////////////////////////////////////////////////////
00020 //
00021 // BEGIN_HTML
00022 // Special p.d.f shape that can be used to model the background of
00023 // D*-D0 mass difference distributions
00024 // END_HTML
00025 //
00026 
00027 #include "RooFit.h"
00028 
00029 #include "Riostream.h"
00030 #include "Riostream.h"
00031 #include <math.h>
00032 #include "TMath.h"
00033 
00034 #include "RooDstD0BG.h"
00035 #include "RooAbsReal.h"
00036 #include "RooRealVar.h"
00037 #include "RooIntegrator1D.h"
00038 #include "RooAbsFunc.h"
00039 
00040 ClassImp(RooDstD0BG) 
00041 
00042 static const char rcsid[] =
00043 "$Id: RooDstD0BG.cxx 24286 2008-06-16 15:47:04Z wouter $";
00044 
00045 
00046 //_____________________________________________________________________________
00047 RooDstD0BG::RooDstD0BG(const char *name, const char *title,
00048                        RooAbsReal& _dm, RooAbsReal& _dm0,
00049                        RooAbsReal& _c, RooAbsReal& _a, RooAbsReal& _b) :
00050   RooAbsPdf(name,title),
00051   dm("dm","Dstar-D0 Mass Diff",this, _dm),
00052   dm0("dm0","Threshold",this, _dm0),
00053   C("C","Shape Parameter",this, _c),
00054   A("A","Shape Parameter 2",this, _a),
00055   B("B","Shape Parameter 3",this, _b)
00056 {
00057 }
00058 
00059 
00060 //_____________________________________________________________________________
00061 RooDstD0BG::RooDstD0BG(const RooDstD0BG& other, const char *name) :
00062   RooAbsPdf(other,name), dm("dm",this,other.dm), dm0("dm0",this,other.dm0),
00063   C("C",this,other.C), A("A",this,other.A), B("B",this,other.B)
00064 {
00065 }
00066 
00067 
00068 //_____________________________________________________________________________
00069 Double_t RooDstD0BG::evaluate() const
00070 {
00071   Double_t arg= dm- dm0;
00072   if (arg <= 0 ) return 0;
00073   Double_t ratio= dm/dm0;
00074   Double_t val= (1- exp(-arg/C))* TMath::Power(ratio, A) + B*(ratio-1);
00075 
00076   return (val > 0 ? val : 0) ;
00077 }
00078 
00079 
00080 //_____________________________________________________________________________
00081 Int_t RooDstD0BG::getAnalyticalIntegral(RooArgSet& /*allVars*/, RooArgSet& /*analVars*/, const char* /*rangeName*/) const 
00082 {
00083   // if (matchArgs(allVars,analVars,dm)) return 1 ;
00084   return 0 ;
00085 }
00086 
00087 
00088 //_____________________________________________________________________________
00089 Double_t RooDstD0BG::analyticalIntegral(Int_t code, const char* rangeName) const 
00090 {
00091   switch(code) {
00092   case 1: 
00093     {
00094       Double_t min= dm.min(rangeName);
00095       Double_t max= dm.max(rangeName);
00096       if (max <= dm0 ) return 0;
00097       else if (min < dm0) min = dm0;
00098 
00099       Bool_t doNumerical= kFALSE;
00100       if ( A != 0 ) doNumerical= kTRUE;
00101       else if (B < 0) {
00102         // If b<0, pdf can be negative at large dm, the integral should
00103         // only up to where pdf hits zero. Better solution should be
00104         // solve the zero and use it as max. 
00105         // Here we check this whether pdf(max) < 0. If true, let numerical
00106         // integral take care of. ( kind of ugly!)
00107         if ( 1- exp(-(max-dm0)/C) + B*(max/dm0 -1) < 0) doNumerical= kTRUE;
00108       }
00109       if ( ! doNumerical ) {
00110         return (max-min)+ C* exp(dm0/C)* (exp(-max/C)- exp(-min/C)) +
00111           B * (0.5* (max*max - min*min)/dm0 - (max- min));
00112       } else {
00113         // In principle the integral for a!=0  can be done analytically. 
00114         // It involves incomplete Gamma function, TMath::Gamma(a+1,m/c), 
00115         // which is not defined for a < -1. And the whole expression is 
00116         // not stable for m/c >> 1.
00117         // Do numerical integral
00118         RooArgSet vset(dm.arg(),"vset");
00119         RooAbsFunc *func= bindVars(vset);
00120         RooIntegrator1D integrator(*func,min,max);
00121         return integrator.integral();
00122       }
00123     }
00124   }
00125   
00126   assert(0) ;
00127   return 0 ;
00128 }

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