Logo ROOT   6.07/09
Reference Guide
RooDstD0BG.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * @(#)root/roofit:$Id$
5  * Authors: *
6  * UE, Ulrik Egede, RAL, U.Egede@rl.ac.uk *
7  * MT, Max Turri, UC Santa Cruz turri@slac.stanford.edu *
8  * CC, Chih-hsiang Cheng, Stanford chcheng@slac.stanford.edu *
9  * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
10  * *
11  * Copyright (c) 2000-2005, Regents of the University of California *
12  * RAL and Stanford University. All rights reserved.*
13  * *
14  * Redistribution and use in source and binary forms, *
15  * with or without modification, are permitted according to the terms *
16  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
17  *****************************************************************************/
18 
19 /**
20 \file RooDstD0BG.cxx
21 \class RooDstD0BG
22 \ingroup Roofit
23 
24 Special p.d.f shape that can be used to model the background of
25 D*-D0 mass difference distributions
26 **/
27 
28 #include "RooFit.h"
29 
30 #include "Riostream.h"
31 #include "Riostream.h"
32 #include <math.h>
33 #include "TMath.h"
34 
35 #include "RooDstD0BG.h"
36 #include "RooAbsReal.h"
37 #include "RooRealVar.h"
38 #include "RooIntegrator1D.h"
39 #include "RooAbsFunc.h"
40 
41 using namespace std;
42 
44 
45 
46 ////////////////////////////////////////////////////////////////////////////////
47 
48 RooDstD0BG::RooDstD0BG(const char *name, const char *title,
49  RooAbsReal& _dm, RooAbsReal& _dm0,
50  RooAbsReal& _c, RooAbsReal& _a, RooAbsReal& _b) :
51  RooAbsPdf(name,title),
52  dm("dm","Dstar-D0 Mass Diff",this, _dm),
53  dm0("dm0","Threshold",this, _dm0),
54  C("C","Shape Parameter",this, _c),
55  A("A","Shape Parameter 2",this, _a),
56  B("B","Shape Parameter 3",this, _b)
57 {
58 }
59 
60 
61 ////////////////////////////////////////////////////////////////////////////////
62 
63 RooDstD0BG::RooDstD0BG(const RooDstD0BG& other, const char *name) :
64  RooAbsPdf(other,name), dm("dm",this,other.dm), dm0("dm0",this,other.dm0),
65  C("C",this,other.C), A("A",this,other.A), B("B",this,other.B)
66 {
67 }
68 
69 
70 ////////////////////////////////////////////////////////////////////////////////
71 
73 {
74  Double_t arg= dm- dm0;
75  if (arg <= 0 ) return 0;
76  Double_t ratio= dm/dm0;
77  Double_t val= (1- exp(-arg/C))* TMath::Power(ratio, A) + B*(ratio-1);
78 
79  return (val > 0 ? val : 0) ;
80 }
81 
82 
83 ////////////////////////////////////////////////////////////////////////////////
84 /// if (matchArgs(allVars,analVars,dm)) return 1 ;
85 
86 Int_t RooDstD0BG::getAnalyticalIntegral(RooArgSet& /*allVars*/, RooArgSet& /*analVars*/, const char* /*rangeName*/) const
87 {
88  return 0 ;
89 }
90 
91 
92 ////////////////////////////////////////////////////////////////////////////////
93 
94 Double_t RooDstD0BG::analyticalIntegral(Int_t code, const char* rangeName) const
95 {
96  switch(code) {
97  case 1:
98  {
99  Double_t min= dm.min(rangeName);
100  Double_t max= dm.max(rangeName);
101  if (max <= dm0 ) return 0;
102  else if (min < dm0) min = dm0;
103 
104  Bool_t doNumerical= kFALSE;
105  if ( A != 0 ) doNumerical= kTRUE;
106  else if (B < 0) {
107  // If b<0, pdf can be negative at large dm, the integral should
108  // only up to where pdf hits zero. Better solution should be
109  // solve the zero and use it as max.
110  // Here we check this whether pdf(max) < 0. If true, let numerical
111  // integral take care of. ( kind of ugly!)
112  if ( 1- exp(-(max-dm0)/C) + B*(max/dm0 -1) < 0) doNumerical= kTRUE;
113  }
114  if ( ! doNumerical ) {
115  return (max-min)+ C* exp(dm0/C)* (exp(-max/C)- exp(-min/C)) +
116  B * (0.5* (max*max - min*min)/dm0 - (max- min));
117  } else {
118  // In principle the integral for a!=0 can be done analytically.
119  // It involves incomplete Gamma function, TMath::Gamma(a+1,m/c),
120  // which is not defined for a < -1. And the whole expression is
121  // not stable for m/c >> 1.
122  // Do numerical integral
123  RooArgSet vset(dm.arg(),"vset");
124  RooAbsFunc *func= bindVars(vset);
125  RooIntegrator1D integrator(*func,min,max);
126  return integrator.integral();
127  }
128  }
129  }
130 
131  assert(0) ;
132  return 0 ;
133 }
static double B[]
RooRealProxy dm
Definition: RooDstD0BG.h:43
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
Special p.d.f shape that can be used to model the background of D*-D0 mass difference distributions...
Definition: RooDstD0BG.h:26
STL namespace.
static double A[]
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
RooRealProxy C
Definition: RooDstD0BG.h:45
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
if (matchArgs(allVars,analVars,dm)) return 1 ;
Definition: RooDstD0BG.cxx:86
RooRealProxy dm0
Definition: RooDstD0BG.h:44
RooRealProxy A
Definition: RooDstD0BG.h:45
static double C[]
virtual Double_t integral(const Double_t *yvec=0)
Calculate numeric integral at given set of function binding parameters.
RooIntegrator1D implements an adaptive one-dimensional numerical integration algorithm.
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooRealProxy B
Definition: RooDstD0BG.h:45
Double_t evaluate() const
Definition: RooDstD0BG.cxx:72
double func(double *x, double *p)
Definition: stressTF1.cxx:213
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
double exp(double)
const Bool_t kTRUE
Definition: Rtypes.h:91
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
Definition: RooAbsFunc.h:23
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
RooAbsFunc * bindVars(const RooArgSet &vars, const RooArgSet *nset=0, Bool_t clipInvalid=kFALSE) const
Create an interface adaptor f(vars) that binds us to the specified variables (in arbitrary order)...
char name[80]
Definition: TGX11.cxx:109
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
Definition: RooDstD0BG.cxx:94