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#include "RooRealVar.h"
36#include "RooAbsFunc.h"
37#include "RooBatchCompute.h"
38
39#include "TMath.h"
40
41#include <cmath>
42using namespace std;
43
45
46////////////////////////////////////////////////////////////////////////////////
47
48RooDstD0BG::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
62RooDstD0BG::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
71{
72 double arg= dm- dm0;
73 if (arg <= 0 ) return 0;
74 double ratio= dm/dm0;
75 double val= (1- exp(-arg/C))* TMath::Power(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.
82void RooDstD0BG::computeBatch(double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
83{
85 {dataMap.at(dm), dataMap.at(dm0), dataMap.at(C), dataMap.at(A), dataMap.at(B)});
86}
87
88////////////////////////////////////////////////////////////////////////////////
89/// if (matchArgs(allVars,analVars,dm)) return 1 ;
90
91Int_t RooDstD0BG::getAnalyticalIntegral(RooArgSet& /*allVars*/, RooArgSet& /*analVars*/, const char* /*rangeName*/) const
92{
93 return 0 ;
94}
95
96////////////////////////////////////////////////////////////////////////////////
97
98double RooDstD0BG::analyticalIntegral(Int_t code, const char* rangeName) const
99{
100 switch(code) {
101 case 1:
102 {
103 double min= dm.min(rangeName);
104 double max= dm.max(rangeName);
105 if (max <= dm0 ) return 0;
106 else if (min < dm0) min = dm0;
107
108 bool doNumerical= false;
109 if ( A != 0 ) doNumerical= true;
110 else if (B < 0) {
111 // If b<0, pdf can be negative at large dm, the integral should
112 // only up to where pdf hits zero. Better solution should be
113 // solve the zero and use it as max.
114 // Here we check this whether pdf(max) < 0. If true, let numerical
115 // integral take care of. ( kind of ugly!)
116 if ( 1- exp(-(max-dm0)/C) + B*(max/dm0 -1) < 0) doNumerical= true;
117 }
118 if ( ! doNumerical ) {
119 return (max-min)+ C* exp(dm0/C)* (exp(-max/C)- exp(-min/C)) +
120 B * (0.5* (max*max - min*min)/dm0 - (max- min));
121 } else {
122 // In principle the integral for a!=0 can be done analytically.
123 // It involves incomplete Gamma function, TMath::Gamma(a+1,m/c),
124 // which is not defined for a < -1. And the whole expression is
125 // not stable for m/c >> 1.
126 // Do numerical integral
127 RooArgSet vset(dm.arg(),"vset");
128 std::unique_ptr<RooAbsFunc> func{bindVars(vset)};
129 RooRombergIntegrator integrator(*func,min,max);
130 return integrator.integral();
131 }
132 }
133 }
134
135 assert(0) ;
136 return 0 ;
137}
#define ClassImp(name)
Definition Rtypes.h:377
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
RooFit::OwningPtr< RooAbsFunc > bindVars(const RooArgSet &vars, const RooArgSet *nset=nullptr, 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:55
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:44
RooRealProxy dm
Definition RooDstD0BG.h:42
RooRealProxy A
Definition RooDstD0BG.h:44
RooRealProxy dm0
Definition RooDstD0BG.h:43
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
void computeBatch(double *output, size_t nEvents, RooFit::Detail::DataMap const &) const override
Compute multiple values of D*-D0 mass difference distribution.
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:44
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooBatchCompute::Config config(RooAbsArg const *arg) const
Definition DataMap.cxx:40
RooRombergIntegrator implements an adaptive numerical integration algorithm.
double integral(const double *yvec=nullptr) override
double max(const char *rname=nullptr) 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.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
void compute(Config cfg, Computer comp, RestrictArr output, size_t size, const VarVector &vars, ArgVector &extraArgs)
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:721
static void output()