Logo ROOT   6.18/05
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
24**/
25
26#include "RooFit.h"
27
28#include "Riostream.h"
29#include "Riostream.h"
30#include <math.h>
31#include "TMath.h"
32
33#include "RooDstD0BG.h"
34#include "RooAbsReal.h"
35#include "RooRealVar.h"
36#include "RooIntegrator1D.h"
37#include "RooAbsFunc.h"
38
39using namespace std;
40
42
43////////////////////////////////////////////////////////////////////////////////
44
45RooDstD0BG::RooDstD0BG(const char *name, const char *title,
46 RooAbsReal& _dm, RooAbsReal& _dm0,
47 RooAbsReal& _c, RooAbsReal& _a, RooAbsReal& _b) :
48 RooAbsPdf(name,title),
49 dm("dm","Dstar-D0 Mass Diff",this, _dm),
50 dm0("dm0","Threshold",this, _dm0),
51 C("C","Shape Parameter",this, _c),
52 A("A","Shape Parameter 2",this, _a),
53 B("B","Shape Parameter 3",this, _b)
54{
55}
56
57////////////////////////////////////////////////////////////////////////////////
58
59RooDstD0BG::RooDstD0BG(const RooDstD0BG& other, const char *name) :
60 RooAbsPdf(other,name), dm("dm",this,other.dm), dm0("dm0",this,other.dm0),
61 C("C",this,other.C), A("A",this,other.A), B("B",this,other.B)
62{
63}
64
65////////////////////////////////////////////////////////////////////////////////
66
68{
69 Double_t arg= dm- dm0;
70 if (arg <= 0 ) return 0;
71 Double_t ratio= dm/dm0;
72 Double_t val= (1- exp(-arg/C))* TMath::Power(ratio, A) + B*(ratio-1);
73
74 return (val > 0 ? val : 0) ;
75}
76
77
78////////////////////////////////////////////////////////////////////////////////
79/// if (matchArgs(allVars,analVars,dm)) return 1 ;
80
81Int_t RooDstD0BG::getAnalyticalIntegral(RooArgSet& /*allVars*/, RooArgSet& /*analVars*/, const char* /*rangeName*/) const
82{
83 return 0 ;
84}
85
86////////////////////////////////////////////////////////////////////////////////
87
88Double_t RooDstD0BG::analyticalIntegral(Int_t code, const char* rangeName) const
89{
90 switch(code) {
91 case 1:
92 {
93 Double_t min= dm.min(rangeName);
94 Double_t max= dm.max(rangeName);
95 if (max <= dm0 ) return 0;
96 else if (min < dm0) min = dm0;
97
98 Bool_t doNumerical= kFALSE;
99 if ( A != 0 ) doNumerical= kTRUE;
100 else if (B < 0) {
101 // If b<0, pdf can be negative at large dm, the integral should
102 // only up to where pdf hits zero. Better solution should be
103 // solve the zero and use it as max.
104 // Here we check this whether pdf(max) < 0. If true, let numerical
105 // integral take care of. ( kind of ugly!)
106 if ( 1- exp(-(max-dm0)/C) + B*(max/dm0 -1) < 0) doNumerical= kTRUE;
107 }
108 if ( ! doNumerical ) {
109 return (max-min)+ C* exp(dm0/C)* (exp(-max/C)- exp(-min/C)) +
110 B * (0.5* (max*max - min*min)/dm0 - (max- min));
111 } else {
112 // In principle the integral for a!=0 can be done analytically.
113 // It involves incomplete Gamma function, TMath::Gamma(a+1,m/c),
114 // which is not defined for a < -1. And the whole expression is
115 // not stable for m/c >> 1.
116 // Do numerical integral
117 RooArgSet vset(dm.arg(),"vset");
118 RooAbsFunc *func= bindVars(vset);
119 RooIntegrator1D integrator(*func,min,max);
120 return integrator.integral();
121 }
122 }
123 }
124
125 assert(0) ;
126 return 0 ;
127}
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:365
char name[80]
Definition: TGX11.cxx:109
double exp(double)
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
Definition: RooAbsFunc.h:23
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
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:28
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.
Definition: RooDstD0BG.cxx:88
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 ;
Definition: RooDstD0BG.cxx:81
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition: RooDstD0BG.cxx:67
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.
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
static double B[]
static double A[]
static double C[]
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:723