Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
RooBDecay.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * PL, Parker C Lund, UC Irvine *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
9 * *
10 * Copyright (c) 2000-2005, Regents of the University of California *
11 * and Stanford University. All rights reserved. *
12 * *
13 * Redistribution and use in source and binary forms, *
14 * with or without modification, are permitted according to the terms *
15 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
16 *****************************************************************************/
17
18
19/** \class RooBDecay
20 \ingroup Roofit
21
22Most general description of B decay time distribution with effects
23of CP violation, mixing and life time differences. This function can
24be analytically convolved with any RooResolutionModel implementation.
25**/
26
27#include "Riostream.h"
28#include "TMath.h"
29#include "RooBDecay.h"
30#include "RooRealVar.h"
31#include "RooRandom.h"
32
33#include "TError.h"
34
35
36/// \brief Constructor for RooBDecay.
37///
38/// Creates an instance of RooBDecay with the specified parameters.
39///
40/// \param[in] name The name of the PDF.
41/// \param[in] title The title of the PDF.
42/// \param[in] t The time variable.
43/// \param[in] tau The average decay time parameter.
44/// \param[in] dgamma The Delta Gamma parameter.
45/// \param[in] f0 The Cosh Coefficient.
46/// \param[in] f1 The Sinh Coefficient.
47/// \param[in] f2 The Cos Coefficient.
48/// \param[in] f3 The Sin Coefficient.
49/// \param[in] dm The Delta Mass parameter.
50/// \param[in] model The resolution model.
51/// \param[in] type The decay type.
52
53RooBDecay::RooBDecay(const char *name, const char* title,
56 RooAbsReal& dm, const RooResolutionModel& model, DecayType type) :
57 RooAbsAnaConvPdf(name, title, model, t),
58 _t("t", "time", this, t),
59 _tau("tau", "Average Decay Time", this, tau),
60 _dgamma("dgamma", "Delta Gamma", this, dgamma),
61 _f0("f0", "Cosh Coefficient", this, f0),
62 _f1("f1", "Sinh Coefficient", this, f1),
63 _f2("f2", "Cos Coefficient", this, f2),
64 _f3("f3", "Sin Coefficient", this, f3),
65 _dm("dm", "Delta Mass", this, dm),
66 _type(type)
67
68{
69 //Constructor
70 switch(type)
71 {
72 case SingleSided:
73 _basisCosh = declareBasis("exp(-@0/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
74 _basisSinh = declareBasis("exp(-@0/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
75 _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau, dm));
76 _basisSin = declareBasis("exp(-@0/@1)*sin(@0*@2)",RooArgList(tau, dm));
77 break;
78 case Flipped:
79 _basisCosh = declareBasis("exp(@0/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
80 _basisSinh = declareBasis("exp(@0/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
81 _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau, dm));
82 _basisSin = declareBasis("exp(@0/@1)*sin(@0*@2)",RooArgList(tau, dm));
83 break;
84 case DoubleSided:
85 _basisCosh = declareBasis("exp(-abs(@0)/@1)*cosh(@0*@2/2)", RooArgList(tau,dgamma));
86 _basisSinh = declareBasis("exp(-abs(@0)/@1)*sinh(@0*@2/2)", RooArgList(tau,dgamma));
87 _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau, dm));
88 _basisSin = declareBasis("exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau, dm));
89 break;
90 }
91}
92
93////////////////////////////////////////////////////////////////////////////////
94///Copy constructor
95
98 _t("t", this, other._t),
99 _tau("tau", this, other._tau),
100 _dgamma("dgamma", this, other._dgamma),
101 _f0("f0", this, other._f0),
102 _f1("f1", this, other._f1),
103 _f2("f2", this, other._f2),
104 _f3("f3", this, other._f3),
105 _dm("dm", this, other._dm),
106 _basisCosh(other._basisCosh),
107 _basisSinh(other._basisSinh),
108 _basisCos(other._basisCos),
109 _basisSin(other._basisSin),
110 _type(other._type)
111{
112}
113
114////////////////////////////////////////////////////////////////////////////////
115
117{
119 {
120 return _f0;
121 }
123 {
124 return _f1;
125 }
126 if(basisIndex == _basisCos)
127 {
128 return _f2;
129 }
130 if(basisIndex == _basisSin)
131 {
132 return _f3;
133 }
134
135 return 0 ;
136}
137
138////////////////////////////////////////////////////////////////////////////////
139
141{
143 {
144 return _f0.arg().getVariables();
145 }
147 {
148 return _f1.arg().getVariables();
149 }
150 if(basisIndex == _basisCos)
151 {
152 return _f2.arg().getVariables();
153 }
154 if(basisIndex == _basisSin)
155 {
156 return _f3.arg().getVariables();
157 }
158
159 return nullptr;
160}
161
162////////////////////////////////////////////////////////////////////////////////
163
165{
166 if(coef == _basisCosh)
167 {
168 return _f0.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
169 }
170 if(coef == _basisSinh)
171 {
172 return _f1.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
173 }
174 if(coef == _basisCos)
175 {
176 return _f2.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
177 }
178 if(coef == _basisSin)
179 {
180 return _f3.arg().getAnalyticalIntegral(allVars,analVars,rangeName) ;
181 }
182
183 return 0 ;
184}
185
186////////////////////////////////////////////////////////////////////////////////
187
188double RooBDecay::coefAnalyticalIntegral(Int_t coef, Int_t code, const char* rangeName) const
189{
190 if(coef == _basisCosh)
191 {
192 return _f0.arg().analyticalIntegral(code,rangeName) ;
193 }
194 if(coef == _basisSinh)
195 {
196 return _f1.arg().analyticalIntegral(code,rangeName) ;
197 }
198 if(coef == _basisCos)
199 {
200 return _f2.arg().analyticalIntegral(code,rangeName) ;
201 }
202 if(coef == _basisSin)
203 {
204 return _f3.arg().analyticalIntegral(code,rangeName) ;
205 }
206
207 return 0 ;
208}
209
210////////////////////////////////////////////////////////////////////////////////
211
213{
214 if (matchArgs(directVars, generateVars, _t)) return 1;
215 return 0;
216}
217
218////////////////////////////////////////////////////////////////////////////////
219
221{
222 R__ASSERT(code==1);
223 double gammamin = 1/_tau-std::abs(_dgamma)/2;
224 while(true) {
225 double t = -log(RooRandom::uniform())/gammamin;
226 if (_type == Flipped || (_type == DoubleSided && RooRandom::uniform() <0.5) ) t *= -1;
227 if ( t<_t.min() || t>_t.max() ) continue;
228
229 double dgt = _dgamma*t/2;
230 double dmt = _dm*t;
231 double ft = std::abs(t);
232 double f = exp(-ft/_tau)*(_f0*cosh(dgt)+_f1*sinh(dgt)+_f2*cos(dmt)+_f3*sin(dmt));
233 if(f < 0) {
234 std::cout << "RooBDecay::generateEvent(" << GetName() << ") ERROR: PDF value less than zero" << std::endl;
235 ::abort();
236 }
237 double w = 1.001*exp(-ft*gammamin)*(std::abs(_f0)+std::abs(_f1)+sqrt(_f2*_f2+_f3*_f3));
238 if(w < f) {
239 std::cout << "RooBDecay::generateEvent(" << GetName() << ") ERROR: Envelope function less than p.d.f. " << std::endl;
240 std::cout << _f0 << std::endl;
241 std::cout << _f1 << std::endl;
242 std::cout << _f2 << std::endl;
243 std::cout << _f3 << std::endl;
244 ::abort();
245 }
246 if(w*RooRandom::uniform() > f) continue;
247 _t = t;
248 break;
249 }
250}
#define f(i)
Definition RSha256.hxx:104
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:110
Base class for PDFs that represent a physics model that can be analytically convolved with a resoluti...
Int_t declareBasis(const char *expression, const RooArgList &params)
Declare a basis function for use in this physics model.
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
bool matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Most general description of B decay time distribution with effects of CP violation,...
Definition RooBDecay.h:25
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const override
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
double coefAnalyticalIntegral(Int_t coef, Int_t code, const char *rangeName=nullptr) const override
Default implementation of function implementing advertised integrals.
RooRealProxy _t
Definition RooBDecay.h:57
RooRealProxy _f2
Definition RooBDecay.h:62
Int_t _basisCosh
Definition RooBDecay.h:65
DecayType _type
Definition RooBDecay.h:70
Int_t _basisSin
Definition RooBDecay.h:68
RooRealProxy _tau
Definition RooBDecay.h:58
void generateEvent(Int_t code) override
Interface for generation of an event using the algorithm corresponding to the specified code.
RooRealProxy _f0
Definition RooBDecay.h:60
@ DoubleSided
Definition RooBDecay.h:29
@ SingleSided
Definition RooBDecay.h:29
double coefficient(Int_t basisIndex) const override
RooRealProxy _dgamma
Definition RooBDecay.h:59
RooRealProxy _f3
Definition RooBDecay.h:63
Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Default implementation of function advertising integration capabilities.
RooFit::OwningPtr< RooArgSet > coefVars(Int_t coefIdx) const override
Return set of parameters with are used exclusively by the coefficient functions.
Int_t _basisSinh
Definition RooBDecay.h:66
RooRealProxy _dm
Definition RooBDecay.h:64
Int_t _basisCos
Definition RooBDecay.h:67
RooRealProxy _f1
Definition RooBDecay.h:61
static double uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition RooRandom.cxx:77
Variable that can be changed from the outside.
Definition RooRealVar.h:37
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...
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.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
TF1 * f1
Definition legend1.C:11
T * OwningPtr
An alias for raw pointers for indicating that the return type of a RooFit function is an owning point...
Definition Config.h:35