Logo ROOT  
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/** \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#include "RooAbsReal.h"
35#include "RooRealVar.h"
36#include "RooIntegrator1D.h"
37#include "RooAbsFunc.h"
38#include "RooBatchCompute.h"
39
40#include "TMath.h"
41
42#include <cmath>
43using namespace std;
44
46
47////////////////////////////////////////////////////////////////////////////////
48
49RooDstD0BG::RooDstD0BG(const char *name, const char *title,
50 RooAbsReal& _dm, RooAbsReal& _dm0,
51 RooAbsReal& _c, RooAbsReal& _a, RooAbsReal& _b) :
52 RooAbsPdf(name,title),
53 dm("dm","Dstar-D0 Mass Diff",this, _dm),
54 dm0("dm0","Threshold",this, _dm0),
55 C("C","Shape Parameter",this, _c),
56 A("A","Shape Parameter 2",this, _a),
57 B("B","Shape Parameter 3",this, _b)
58{
59}
60
61////////////////////////////////////////////////////////////////////////////////
62
63RooDstD0BG::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
72{
73 double arg= dm- dm0;
74 if (arg <= 0 ) return 0;
75 double ratio= dm/dm0;
76 double val= (1- exp(-arg/C))* TMath::Power(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.
83void RooDstD0BG::computeBatch(cudaStream_t* stream, double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
84{
86 dispatch->compute(stream, RooBatchCompute::DstD0BG, output, nEvents,
87 {dataMap.at(dm), dataMap.at(dm0), dataMap.at(C), dataMap.at(A), dataMap.at(B)});
88}
89
90////////////////////////////////////////////////////////////////////////////////
91/// if (matchArgs(allVars,analVars,dm)) return 1 ;
92
93Int_t RooDstD0BG::getAnalyticalIntegral(RooArgSet& /*allVars*/, RooArgSet& /*analVars*/, 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 {
105 double min= dm.min(rangeName);
106 double max= dm.max(rangeName);
107 if (max <= dm0 ) return 0;
108 else if (min < dm0) min = dm0;
109
110 bool doNumerical= false;
111 if ( A != 0 ) doNumerical= true;
112 else if (B < 0) {
113 // If b<0, pdf can be negative at large dm, the integral should
114 // only up to where pdf hits zero. Better solution should be
115 // solve the zero and use it as max.
116 // Here we check this whether pdf(max) < 0. If true, let numerical
117 // integral take care of. ( kind of ugly!)
118 if ( 1- exp(-(max-dm0)/C) + B*(max/dm0 -1) < 0) doNumerical= true;
119 }
120 if ( ! doNumerical ) {
121 return (max-min)+ C* exp(dm0/C)* (exp(-max/C)- exp(-min/C)) +
122 B * (0.5* (max*max - min*min)/dm0 - (max- min));
123 } else {
124 // In principle the integral for a!=0 can be done analytically.
125 // It involves incomplete Gamma function, TMath::Gamma(a+1,m/c),
126 // which is not defined for a < -1. And the whole expression is
127 // not stable for m/c >> 1.
128 // Do numerical integral
129 RooArgSet vset(dm.arg(),"vset");
130 RooAbsFunc *func= bindVars(vset);
131 RooIntegrator1D integrator(*func,min,max);
132 return integrator.integral();
133 }
134 }
135 }
136
137 assert(0) ;
138 return 0 ;
139}
#define ClassImp(name)
Definition: Rtypes.h:375
char name[80]
Definition: TGX11.cxx:110
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
Definition: RooAbsFunc.h:27
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:64
RooAbsFunc * bindVars(const RooArgSet &vars, const RooArgSet *nset=0, bool clipInvalid=false) const
Create an interface adaptor f(vars) that binds us to the specified variables (in arbitrary order).
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:57
virtual void compute(cudaStream_t *, Computer, RestrictArr, size_t, const VarVector &, const ArgVector &={})=0
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:45
RooRealProxy dm
Definition: RooDstD0BG.h:43
void computeBatch(cudaStream_t *, double *output, size_t nEvents, RooFit::Detail::DataMap const &) const override
Compute multiple values of D*-D0 mass difference distribution.
Definition: RooDstD0BG.cxx:83
RooRealProxy A
Definition: RooDstD0BG.h:45
RooRealProxy dm0
Definition: RooDstD0BG.h:44
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition: RooDstD0BG.cxx:71
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const override
if (matchArgs(allVars,analVars,dm)) return 1 ;
Definition: RooDstD0BG.cxx:93
RooRealProxy C
Definition: RooDstD0BG.h:45
double analyticalIntegral(Int_t code, const char *rangeName=0) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Definition: RooDstD0BG.cxx:100
auto & at(RooAbsArg const *arg, RooAbsArg const *=nullptr)
Definition: DataMap.h:88
RooIntegrator1D implements an adaptive one-dimensional numerical integration algorithm.
double integral(const double *yvec=0) override
Calculate numeric integral at given set of function binding parameters.
double min(const char *rname=0) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
double max(const char *rname=0) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
const T & arg() const
Return reference to object held in proxy.
RVec< PromoteType< T > > exp(const RVec< T > &v)
Definition: RVec.hxx:1744
static double B[]
static double A[]
static double C[]
R__EXTERN RooBatchComputeInterface * dispatchCUDA
R__EXTERN RooBatchComputeInterface * dispatchCPU
This dispatch pointer points to an implementation of the compute library, provided one has been loade...
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition: TMath.h:718
static void output()