Loading [MathJax]/extensions/tex2jax.js
ROOT  6.06/09
Reference Guide
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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 //
21 // BEGIN_HTML
22 // Special p.d.f shape that can be used to model the background of
23 // D*-D0 mass difference distributions
24 // END_HTML
25 //
26 
27 #include "RooFit.h"
28 
29 #include "Riostream.h"
30 #include "Riostream.h"
31 #include <math.h>
32 #include "TMath.h"
33 
34 #include "RooDstD0BG.h"
35 #include "RooAbsReal.h"
36 #include "RooRealVar.h"
37 #include "RooIntegrator1D.h"
38 #include "RooAbsFunc.h"
39 
40 using namespace std;
41 
43 
44 
45 ////////////////////////////////////////////////////////////////////////////////
46 
47 RooDstD0BG::RooDstD0BG(const char *name, const char *title,
48  RooAbsReal& _dm, RooAbsReal& _dm0,
49  RooAbsReal& _c, RooAbsReal& _a, RooAbsReal& _b) :
50  RooAbsPdf(name,title),
51  dm("dm","Dstar-D0 Mass Diff",this, _dm),
52  dm0("dm0","Threshold",this, _dm0),
53  C("C","Shape Parameter",this, _c),
54  A("A","Shape Parameter 2",this, _a),
55  B("B","Shape Parameter 3",this, _b)
56 {
57 }
58 
59 
60 ////////////////////////////////////////////////////////////////////////////////
61 
62 RooDstD0BG::RooDstD0BG(const RooDstD0BG& other, const char *name) :
63  RooAbsPdf(other,name), dm("dm",this,other.dm), dm0("dm0",this,other.dm0),
64  C("C",this,other.C), A("A",this,other.A), B("B",this,other.B)
65 {
66 }
67 
68 
69 ////////////////////////////////////////////////////////////////////////////////
70 
72 {
73  Double_t arg= dm- dm0;
74  if (arg <= 0 ) return 0;
75  Double_t ratio= dm/dm0;
76  Double_t val= (1- exp(-arg/C))* TMath::Power(ratio, A) + B*(ratio-1);
77 
78  return (val > 0 ? val : 0) ;
79 }
80 
81 
82 ////////////////////////////////////////////////////////////////////////////////
83 /// if (matchArgs(allVars,analVars,dm)) return 1 ;
84 
85 Int_t RooDstD0BG::getAnalyticalIntegral(RooArgSet& /*allVars*/, RooArgSet& /*analVars*/, const char* /*rangeName*/) const
86 {
87  return 0 ;
88 }
89 
90 
91 ////////////////////////////////////////////////////////////////////////////////
92 
93 Double_t RooDstD0BG::analyticalIntegral(Int_t code, const char* rangeName) const
94 {
95  switch(code) {
96  case 1:
97  {
98  Double_t min= dm.min(rangeName);
99  Double_t max= dm.max(rangeName);
100  if (max <= dm0 ) return 0;
101  else if (min < dm0) min = dm0;
102 
103  Bool_t doNumerical= kFALSE;
104  if ( A != 0 ) doNumerical= kTRUE;
105  else if (B < 0) {
106  // If b<0, pdf can be negative at large dm, the integral should
107  // only up to where pdf hits zero. Better solution should be
108  // solve the zero and use it as max.
109  // Here we check this whether pdf(max) < 0. If true, let numerical
110  // integral take care of. ( kind of ugly!)
111  if ( 1- exp(-(max-dm0)/C) + B*(max/dm0 -1) < 0) doNumerical= kTRUE;
112  }
113  if ( ! doNumerical ) {
114  return (max-min)+ C* exp(dm0/C)* (exp(-max/C)- exp(-min/C)) +
115  B * (0.5* (max*max - min*min)/dm0 - (max- min));
116  } else {
117  // In principle the integral for a!=0 can be done analytically.
118  // It involves incomplete Gamma function, TMath::Gamma(a+1,m/c),
119  // which is not defined for a < -1. And the whole expression is
120  // not stable for m/c >> 1.
121  // Do numerical integral
122  RooArgSet vset(dm.arg(),"vset");
123  RooAbsFunc *func= bindVars(vset);
124  RooIntegrator1D integrator(*func,min,max);
125  return integrator.integral();
126  }
127  }
128  }
129 
130  assert(0) ;
131  return 0 ;
132 }
static double B[]
static Vc_ALWAYS_INLINE int_v min(const int_v &x, const int_v &y)
Definition: vector.h:433
#define assert(cond)
Definition: unittest.h:542
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
STL namespace.
ClassImp(RooDstD0BG) RooDstD0BG
Definition: RooDstD0BG.cxx:42
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:85
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.
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:71
double func(double *x, double *p)
Definition: stressTF1.cxx:213
static Vc_ALWAYS_INLINE int_v max(const int_v &x, const int_v &y)
Definition: vector.h:440
#define name(a, b)
Definition: linkTestLib0.cpp:5
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
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)...
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
Definition: RooDstD0BG.cxx:93