ROOT
git-r3/HEAD
Reference Guide
Loading...
Searching...
No Matches
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
/** \class RooDstD0BG
20
\ingroup Roofit
21
22
Special p.d.f shape that can be used to model the background of
23
D*-D0 mass difference distributions. It computes
24
25
\f[
26
\mathrm{RooDSTD0}(m \, | \, m_0, A, B, C) =
27
\left(1 - \exp\left(-\frac{m - m_0}{C}\right) \right)
28
\cdot \left(\frac{m}{m_0}\right)^A + B
29
\cdot \left(\frac{m}{m_0} - 1 \right)
30
\f]
31
**/
32
33
#include <
RooDstD0BG.h
>
34
35
#include <
RooBatchCompute.h
>
36
37
#include <
TMath.h
>
38
39
#include <cmath>
40
41
42
////////////////////////////////////////////////////////////////////////////////
43
44
RooDstD0BG::RooDstD0BG
(
const
char
*
name
,
const
char
*title,
RooAbsReal
&_dm,
RooAbsReal
&_dm0,
RooAbsReal
&_c,
45
RooAbsReal
&_a,
RooAbsReal
&_b)
46
:
RooAbsPdf
(
name
, title),
47
dm
(
"dm"
,
"Dstar-D0 Mass Diff"
, this, _dm),
48
dm0
(
"dm0"
,
"Threshold"
, this, _dm0),
49
C
(
"C"
,
"Shape Parameter"
, this, _c),
50
A
(
"A"
,
"Shape Parameter 2"
, this, _a),
51
B
(
"B"
,
"Shape Parameter 3"
, this, _b)
52
{
53
}
54
55
////////////////////////////////////////////////////////////////////////////////
56
57
RooDstD0BG::RooDstD0BG
(
const
RooDstD0BG
&other,
const
char
*
name
)
58
:
RooAbsPdf
(other,
name
),
59
dm
(
"dm"
, this, other.
dm
),
60
dm0
(
"dm0"
, this, other.
dm0
),
61
C
(
"C"
, this, other.
C
),
62
A
(
"A"
, this, other.
A
),
63
B
(
"B"
, this, other.
B
)
64
{
65
}
66
67
////////////////////////////////////////////////////////////////////////////////
68
69
double
RooDstD0BG::evaluate
()
const
70
{
71
double
arg =
dm
-
dm0
;
72
if
(arg <= 0)
73
return
0;
74
double
ratio =
dm
/
dm0
;
75
double
val = (1 - std::exp(-arg /
C
)) * std::pow(ratio,
A
) +
B
* (ratio - 1);
76
77
return
(val > 0 ? val : 0);
78
}
79
80
////////////////////////////////////////////////////////////////////////////////
81
/// Compute multiple values of D*-D0 mass difference distribution.
82
void
RooDstD0BG::doEval
(
RooFit::EvalContext
& ctx)
const
83
{
84
RooBatchCompute::compute
(ctx.
config
(
this
),
RooBatchCompute::DstD0BG
, ctx.
output
(),
85
{ctx.at(dm), ctx.at(dm0), ctx.at(C), ctx.at(A), ctx.at(B)});
86
}
87
88
////////////////////////////////////////////////////////////////////////////////
89
/// if (matchArgs(allVars,analVars,dm)) return 1 ;
90
91
Int_t
RooDstD0BG::getAnalyticalIntegral
(
RooArgSet
&
/*allVars*/
,
RooArgSet
&
/*analVars*/
,
92
const
char
*
/*rangeName*/
)
const
93
{
94
return
0;
95
}
96
97
////////////////////////////////////////////////////////////////////////////////
98
99
double
RooDstD0BG::analyticalIntegral
(
Int_t
/*code*/
,
const
char
*
/*rangeName*/
)
const
100
{
101
// switch (code) {
102
// case 1: {
103
// double min = dm.min(rangeName);
104
// double max = dm.max(rangeName);
105
// if (max <= dm0)
106
// return 0;
107
// else if (min < dm0)
108
// min = dm0;
109
//
110
// bool doNumerical = false;
111
// if (A != 0)
112
// doNumerical = true;
113
// else if (B < 0) {
114
// // If b<0, pdf can be negative at large dm, the integral should
115
// // only up to where pdf hits zero. Better solution should be
116
// // solve the zero and use it as max.
117
// // Here we check this whether pdf(max) < 0. If true, let numerical
118
// // integral take care of. ( kind of ugly!)
119
// if (1 - exp(-(max - dm0) / C) + B * (max / dm0 - 1) < 0)
120
// doNumerical = true;
121
// }
122
// if (!doNumerical) {
123
// return (max - min) + C * exp(dm0 / C) * (exp(-max / C) - exp(-min / C)) +
124
// B * (0.5 * (max * max - min * min) / dm0 - (max - min));
125
// } else {
126
// // In principle the integral for a!=0 can be done analytically.
127
// // It involves incomplete Gamma function, TMath::Gamma(a+1,m/c),
128
// // which is not defined for a < -1. And the whole expression is
129
// // not stable for m/c >> 1.
130
// // Do numerical integral
131
// RooArgSet vset(dm.arg(),"vset");
132
// std::unique_ptr<RooAbsFunc> func{bindVars(vset)};
133
// RooRombergIntegrator integrator(*func,min,max);
134
// return integrator.integral();
135
// }
136
// }
137
// }
138
139
assert(0);
140
return
0.0;
141
}
RooBatchCompute.h
RooDstD0BG.h
Int_t
int Int_t
Signed integer 4 bytes (int).
Definition
RtypesCore.h:59
name
char name[80]
Definition
TGX11.cxx:148
TMath.h
RooAbsPdf::RooAbsPdf
RooAbsPdf()
Default constructor.
Definition
RooAbsPdf.cxx:237
RooAbsPdf::RooAbsReal
friend class RooAbsReal
Definition
RooAbsPdf.h:342
RooArgSet
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition
RooArgSet.h:24
RooDstD0BG::B
RooRealProxy B
Definition
RooDstD0BG.h:41
RooDstD0BG::dm
RooRealProxy dm
Definition
RooDstD0BG.h:39
RooDstD0BG::A
RooRealProxy A
Definition
RooDstD0BG.h:41
RooDstD0BG::dm0
RooRealProxy dm0
Definition
RooDstD0BG.h:40
RooDstD0BG::evaluate
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition
RooDstD0BG.cxx:69
RooDstD0BG::RooDstD0BG
RooDstD0BG()
Definition
RooDstD0BG.h:28
RooDstD0BG::getAnalyticalIntegral
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
if (matchArgs(allVars,analVars,dm)) return 1 ;
Definition
RooDstD0BG.cxx:91
RooDstD0BG::C
RooRealProxy C
Definition
RooDstD0BG.h:41
RooDstD0BG::doEval
void doEval(RooFit::EvalContext &) const override
Compute multiple values of D*-D0 mass difference distribution.
Definition
RooDstD0BG.cxx:82
RooDstD0BG::analyticalIntegral
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Definition
RooDstD0BG.cxx:99
RooFit::EvalContext
Definition
EvalContext.h:84
RooFit::EvalContext::output
std::span< double > output()
Definition
EvalContext.h:112
RooFit::EvalContext::config
RooBatchCompute::Config config(RooAbsArg const *arg) const
Definition
EvalContext.cxx:73
RooBatchCompute::compute
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})
Definition
RooBatchCompute.h:212
RooBatchCompute::DstD0BG
@ DstD0BG
Definition
RooBatchCompute.h:86
roofit
roofit
src
RooDstD0BG.cxx
ROOTgit-r3/HEAD - Reference Guide Generated on
(GVA Time) using Doxygen 1.16.1