Logo ROOT  
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
22Special p.d.f shape that can be used to model the background of
23D*-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
42
43////////////////////////////////////////////////////////////////////////////////
44
45RooDstD0BG::RooDstD0BG(const char *name, const char *title, RooAbsReal &_dm, RooAbsReal &_dm0, RooAbsReal &_c,
46 RooAbsReal &_a, RooAbsReal &_b)
47 : RooAbsPdf(name, title),
48 dm("dm", "Dstar-D0 Mass Diff", this, _dm),
49 dm0("dm0", "Threshold", this, _dm0),
50 C("C", "Shape Parameter", this, _c),
51 A("A", "Shape Parameter 2", this, _a),
52 B("B", "Shape Parameter 3", this, _b)
53{
54}
55
56////////////////////////////////////////////////////////////////////////////////
57
58RooDstD0BG::RooDstD0BG(const RooDstD0BG &other, const char *name)
59 : RooAbsPdf(other, name),
60 dm("dm", this, other.dm),
61 dm0("dm0", this, other.dm0),
62 C("C", this, other.C),
63 A("A", this, other.A),
64 B("B", this, other.B)
65{
66}
67
68////////////////////////////////////////////////////////////////////////////////
69
71{
72 double arg = dm - dm0;
73 if (arg <= 0)
74 return 0;
75 double ratio = dm / dm0;
76 double val = (1 - std::exp(-arg / C)) * std::pow(ratio, A) + B * (ratio - 1);
77
78 return (val > 0 ? val : 0);
79}
80
81////////////////////////////////////////////////////////////////////////////////
82/// Compute multiple values of D*-D0 mass difference distribution.
84{
86 {ctx.at(dm), ctx.at(dm0), ctx.at(C), ctx.at(A), ctx.at(B)});
87}
88
89////////////////////////////////////////////////////////////////////////////////
90/// if (matchArgs(allVars,analVars,dm)) return 1 ;
91
93 const char * /*rangeName*/) const
94{
95 return 0;
96}
97
98////////////////////////////////////////////////////////////////////////////////
99
100double RooDstD0BG::analyticalIntegral(Int_t /*code*/, const char * /*rangeName*/) const
101{
102 // switch (code) {
103 // case 1: {
104 // double min = dm.min(rangeName);
105 // double max = dm.max(rangeName);
106 // if (max <= dm0)
107 // return 0;
108 // else if (min < dm0)
109 // min = dm0;
110 //
111 // bool doNumerical = false;
112 // if (A != 0)
113 // doNumerical = true;
114 // else if (B < 0) {
115 // // If b<0, pdf can be negative at large dm, the integral should
116 // // only up to where pdf hits zero. Better solution should be
117 // // solve the zero and use it as max.
118 // // Here we check this whether pdf(max) < 0. If true, let numerical
119 // // integral take care of. ( kind of ugly!)
120 // if (1 - exp(-(max - dm0) / C) + B * (max / dm0 - 1) < 0)
121 // doNumerical = true;
122 // }
123 // if (!doNumerical) {
124 // return (max - min) + C * exp(dm0 / C) * (exp(-max / C) - exp(-min / C)) +
125 // B * (0.5 * (max * max - min * min) / dm0 - (max - min));
126 // } else {
127 // // In principle the integral for a!=0 can be done analytically.
128 // // It involves incomplete Gamma function, TMath::Gamma(a+1,m/c),
129 // // which is not defined for a < -1. And the whole expression is
130 // // not stable for m/c >> 1.
131 // // Do numerical integral
132 // RooArgSet vset(dm.arg(),"vset");
133 // std::unique_ptr<RooAbsFunc> func{bindVars(vset)};
134 // RooRombergIntegrator integrator(*func,min,max);
135 // return integrator.integral();
136 // }
137 // }
138 // }
139
140 assert(0);
141 return 0.0;
142}
#define ClassImp(name)
Definition Rtypes.h:382
char name[80]
Definition TGX11.cxx:110
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Special p.d.f shape that can be used to model the background of D*-D0 mass difference distributions.
Definition RooDstD0BG.h:26
RooRealProxy B
Definition RooDstD0BG.h:41
RooRealProxy dm
Definition RooDstD0BG.h:39
RooRealProxy A
Definition RooDstD0BG.h:41
RooRealProxy dm0
Definition RooDstD0BG.h:40
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
if (matchArgs(allVars,analVars,dm)) return 1 ;
RooRealProxy C
Definition RooDstD0BG.h:41
void doEval(RooFit::EvalContext &) const override
Compute multiple values of D*-D0 mass difference distribution.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})