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 "RooFit.h"
35#include "RooAbsReal.h"
36#include "RooRealVar.h"
37#include "RooIntegrator1D.h"
38#include "RooAbsFunc.h"
39#include "RooBatchCompute.h"
40
41#include "TMath.h"
42
43#include <cmath>
44using namespace std;
45
47
48////////////////////////////////////////////////////////////////////////////////
49
50RooDstD0BG::RooDstD0BG(const char *name, const char *title,
51 RooAbsReal& _dm, RooAbsReal& _dm0,
52 RooAbsReal& _c, RooAbsReal& _a, RooAbsReal& _b) :
53 RooAbsPdf(name,title),
54 dm("dm","Dstar-D0 Mass Diff",this, _dm),
55 dm0("dm0","Threshold",this, _dm0),
56 C("C","Shape Parameter",this, _c),
57 A("A","Shape Parameter 2",this, _a),
58 B("B","Shape Parameter 3",this, _b)
59{
60}
61
62////////////////////////////////////////////////////////////////////////////////
63
64RooDstD0BG::RooDstD0BG(const RooDstD0BG& other, const char *name) :
65 RooAbsPdf(other,name), dm("dm",this,other.dm), dm0("dm0",this,other.dm0),
66 C("C",this,other.C), A("A",this,other.A), B("B",this,other.B)
67{
68}
69
70////////////////////////////////////////////////////////////////////////////////
71
73{
74 Double_t arg= dm- dm0;
75 if (arg <= 0 ) return 0;
76 Double_t ratio= dm/dm0;
77 Double_t val= (1- exp(-arg/C))* TMath::Power(ratio, A) + B*(ratio-1);
78
79 return (val > 0 ? val : 0) ;
80}
81
82////////////////////////////////////////////////////////////////////////////////
83/// Compute multiple values of D*-D0 mass difference distribution.
85 return RooBatchCompute::dispatch->computeDstD0BG(this, evalData, dm->getValues(evalData, normSet), dm0->getValues(evalData, normSet), C->getValues(evalData, normSet), A->getValues(evalData, normSet), B->getValues(evalData, normSet));
86}
87
88
89////////////////////////////////////////////////////////////////////////////////
90/// if (matchArgs(allVars,analVars,dm)) return 1 ;
91
92Int_t RooDstD0BG::getAnalyticalIntegral(RooArgSet& /*allVars*/, RooArgSet& /*analVars*/, const char* /*rangeName*/) const
93{
94 return 0 ;
95}
96
97////////////////////////////////////////////////////////////////////////////////
98
99Double_t RooDstD0BG::analyticalIntegral(Int_t code, const char* rangeName) const
100{
101 switch(code) {
102 case 1:
103 {
104 Double_t min= dm.min(rangeName);
105 Double_t max= dm.max(rangeName);
106 if (max <= dm0 ) return 0;
107 else if (min < dm0) min = dm0;
108
109 Bool_t doNumerical= kFALSE;
110 if ( A != 0 ) doNumerical= kTRUE;
111 else if (B < 0) {
112 // If b<0, pdf can be negative at large dm, the integral should
113 // only up to where pdf hits zero. Better solution should be
114 // solve the zero and use it as max.
115 // Here we check this whether pdf(max) < 0. If true, let numerical
116 // integral take care of. ( kind of ugly!)
117 if ( 1- exp(-(max-dm0)/C) + B*(max/dm0 -1) < 0) doNumerical= kTRUE;
118 }
119 if ( ! doNumerical ) {
120 return (max-min)+ C* exp(dm0/C)* (exp(-max/C)- exp(-min/C)) +
121 B * (0.5* (max*max - min*min)/dm0 - (max- min));
122 } else {
123 // In principle the integral for a!=0 can be done analytically.
124 // It involves incomplete Gamma function, TMath::Gamma(a+1,m/c),
125 // which is not defined for a < -1. And the whole expression is
126 // not stable for m/c >> 1.
127 // Do numerical integral
128 RooArgSet vset(dm.arg(),"vset");
129 RooAbsFunc *func= bindVars(vset);
130 RooIntegrator1D integrator(*func,min,max);
131 return integrator.integral();
132 }
133 }
134 }
135
136 assert(0) ;
137 return 0 ;
138}
const Bool_t kFALSE
Definition RtypesCore.h:92
const Bool_t kTRUE
Definition RtypesCore.h:91
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
double exp(double)
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:61
virtual RooSpan< const double > getValues(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet=nullptr) const
by this change, please consult the release notes for ROOT 6.24 for guidance on how to make this trans...
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).
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
virtual RooSpan< double > computeDstD0BG(const RooAbsReal *, RunContext &, RooSpan< const double > dm, RooSpan< const double > dm0, RooSpan< const double > C, RooSpan< const double > A, RooSpan< const double > B)=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
RooRealProxy A
Definition RooDstD0BG.h:45
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooRealProxy dm0
Definition RooDstD0BG.h:44
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 ;
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooSpan< double > evaluateSpan(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet) const
Compute multiple values of D*-D0 mass difference distribution.
RooIntegrator1D implements an adaptive one-dimensional numerical integration algorithm.
virtual Double_t integral(const Double_t *yvec=0)
Calculate numeric integral at given set of function binding parameters.
A simple container to hold a batch of data values.
Definition RooSpan.h:34
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.
R__EXTERN RooBatchComputeInterface * dispatch
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)
Definition TMath.h:735
This struct enables passing computation data around between elements of a computation graph.
Definition RunContext.h:31