Logo ROOT  
Reference Guide
RooBCPEffDecay.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17/** \class RooBCPEffDecay
18 \ingroup Roofit
19
20PDF describing decay time distribution of B meson including effects of standard model CP violation.
21This function can be analytically convolved with any RooResolutionModel implementation.
22*/
23
24
25#include "RooFit.h"
26
27#include "Riostream.h"
28#include "Riostream.h"
29#include "RooRealVar.h"
30#include "RooRandom.h"
31#include "RooBCPEffDecay.h"
32#include "RooRealIntegral.h"
33
34using namespace std;
35
37
38////////////////////////////////////////////////////////////////////////////////
39/// Constructor.
40
41RooBCPEffDecay::RooBCPEffDecay(const char *name, const char *title,
43 RooAbsReal& tau, RooAbsReal& dm,
44 RooAbsReal& avgMistag, RooAbsReal& CPeigenval,
46 RooAbsReal& effRatio, RooAbsReal& delMistag,
47 const RooResolutionModel& model, DecayType type) :
48 RooAbsAnaConvPdf(name,title,model,t),
49 _absLambda("absLambda","Absolute value of lambda",this,a),
50 _argLambda("argLambda","Arg(Lambda)",this,b),
51 _effRatio("effRatio","B0/B0bar efficiency ratio",this,effRatio),
52 _CPeigenval("CPeigenval","CP eigen value",this,CPeigenval),
53 _avgMistag("avgMistag","Average mistag rate",this,avgMistag),
54 _delMistag("delMistag","Delta mistag rate",this,delMistag),
55 _t("t","time",this,t),
56 _tau("tau","decay time",this,tau),
57 _dm("dm","mixing frequency",this,dm),
58 _tag("tag","CP state",this,tag),
59 _genB0Frac(0),
60 _type(type)
61{
62 switch(type) {
63 case SingleSided:
64 _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau,dm)) ;
65 _basisSin = declareBasis("exp(-@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
66 _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
67 break ;
68 case Flipped:
69 _basisExp = declareBasis("exp(@0)/@1)",RooArgList(tau,dm)) ;
70 _basisSin = declareBasis("exp(@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
71 _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
72 break ;
73 case DoubleSided:
74 _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau,dm)) ;
75 _basisSin = declareBasis("exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
76 _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
77 break ;
78 }
79}
80
81////////////////////////////////////////////////////////////////////////////////
82/// Copy constructor.
83
86 _absLambda("absLambda",this,other._absLambda),
87 _argLambda("argLambda",this,other._argLambda),
88 _effRatio("effRatio",this,other._effRatio),
89 _CPeigenval("CPeigenval",this,other._CPeigenval),
90 _avgMistag("avgMistag",this,other._avgMistag),
91 _delMistag("delMistag",this,other._delMistag),
92 _t("t",this,other._t),
93 _tau("tau",this,other._tau),
94 _dm("dm",this,other._dm),
95 _tag("tag",this,other._tag),
96 _genB0Frac(other._genB0Frac),
97 _type(other._type),
98 _basisExp(other._basisExp),
99 _basisSin(other._basisSin),
100 _basisCos(other._basisCos)
101{
102}
103
104////////////////////////////////////////////////////////////////////////////////
105/// Destructor.
106
108{
109}
110
111////////////////////////////////////////////////////////////////////////////////
112/// B0 : _tag = +1
113///
114/// B0bar : _tag = -1
115/// \param[in] basisIndex
116
118{
119 if (basisIndex==_basisExp) {
120 //exp term: (1 -/+ dw)(1+a^2)/2
121 return (1 - _tag*_delMistag)*(1+_absLambda*_absLambda)/2 ;
122 // = 1 + _tag*deltaDil/2
123 }
124
125 if (basisIndex==_basisSin) {
126 //sin term: +/- (1-2w)*ImLambda(= -etaCP * |l| * sin2b)
128 // = _tag*avgDil * ...
129 }
130
131 if (basisIndex==_basisCos) {
132 //cos term: +/- (1-2w)*(1-a^2)/2
133 return -1*_tag*(1-2*_avgMistag)*(1-_absLambda*_absLambda)/2 ;
134 // = -_tag*avgDil * ...
135 }
136
137 return 0 ;
138}
139
140////////////////////////////////////////////////////////////////////////////////
141
142Int_t RooBCPEffDecay::getCoefAnalyticalIntegral(Int_t /*code*/, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
143{
144 if (rangeName) return 0 ;
145
146 if (matchArgs(allVars,analVars,_tag)) return 1 ;
147 return 0 ;
148}
149
150////////////////////////////////////////////////////////////////////////////////
151
152Double_t RooBCPEffDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* /*rangeName*/) const
153{
154 switch(code) {
155 // No integration
156 case 0: return coefficient(basisIndex) ;
157
158 // Integration over 'tag'
159 case 1:
160 if (basisIndex==_basisExp) {
161 return (1+_absLambda*_absLambda) ;
162 }
163
164 if (basisIndex==_basisSin || basisIndex==_basisCos) {
165 return 0 ;
166 }
167 break ;
168
169 default:
170 assert(0) ;
171 }
172
173 return 0 ;
174}
175
176////////////////////////////////////////////////////////////////////////////////
177
178Int_t RooBCPEffDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK) const
179{
180 if (staticInitOK) {
181 if (matchArgs(directVars,generateVars,_t,_tag)) return 2 ;
182 }
183 if (matchArgs(directVars,generateVars,_t)) return 1 ;
184 return 0 ;
185}
186
187////////////////////////////////////////////////////////////////////////////////
188
190{
191 if (code==2) {
192 // Calculate the fraction of mixed events to generate
193 Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tag.arg())).getVal() ;
194 _tag = 1 ;
195 Double_t b0Int = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
196 _genB0Frac = b0Int/sumInt ;
197 }
198}
199
200////////////////////////////////////////////////////////////////////////////////
201/// Generates mix-state dependent.
202/// \param[in] code
203
205{
206 if (code==2) {
208 _tag = (rand<=_genB0Frac) ? 1 : -1 ;
209 }
210
211 // Generate delta-t dependent
212 while(1) {
214 Double_t tval(0) ;
215
216 switch(_type) {
217 case SingleSided:
218 tval = -_tau*log(rand);
219 break ;
220 case Flipped:
221 tval= +_tau*log(rand);
222 break ;
223 case DoubleSided:
224 tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
225 break ;
226 }
227
228 // Accept event if T is in generated range
229 Double_t maxDil = 1.0 ;
231 Double_t maxAcceptProb = (1+al2) + fabs(maxDil*_CPeigenval*_absLambda*_argLambda) + fabs(maxDil*(1-al2)/2);
232 Double_t acceptProb = (1+al2)/2*(1-_tag*_delMistag)
234 - (_tag*(1-2*_avgMistag))*(1-al2)/2*cos(_dm*tval);
235
236 Bool_t accept = maxAcceptProb*RooRandom::uniform() < acceptProb ? kTRUE : kFALSE ;
237
238 if (tval<_t.max() && tval>_t.min() && accept) {
239 _t = tval ;
240 break ;
241 }
242 }
243
244}
#define b(i)
Definition: RSha256.hxx:100
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
int type
Definition: TGX11.cxx:120
double cos(double)
double sin(double)
double log(double)
RooAbsAnaConvPdf is the base class for PDFs that represent a physics model that can be analytically c...
Int_t declareBasis(const char *expression, const RooArgList &params)
Declare a basis function for use in this physics model.
friend class RooArgSet
Definition: RooAbsArg.h:551
RooAbsCategory is the common abstract base class for objects that represent a discrete value with a f...
friend class RooRealIntegral
Definition: RooAbsPdf.h:314
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:59
Bool_t 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:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
PDF describing decay time distribution of B meson including effects of standard model CP violation.
DecayType _type
Double_t _genB0Frac
virtual Double_t coefAnalyticalIntegral(Int_t coef, Int_t code, const char *rangeName=0) const
Default implementation of function implementing advertised integrals.
RooRealProxy _t
RooRealProxy _avgMistag
virtual Double_t coefficient(Int_t basisIndex) const
B0 : _tag = +1.
virtual Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Default implementation of function advertising integration capabilities.
RooRealProxy _delMistag
void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
RooRealProxy _dm
RooRealProxy _absLambda
void generateEvent(Int_t code)
Generates mix-state dependent.
RooRealProxy _argLambda
virtual ~RooBCPEffDecay()
Destructor.
RooCategoryProxy _tag
RooRealProxy _CPeigenval
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
RooRealProxy _tau
const RooAbsCategory & arg() const
static Double_t uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition: RooRandom.cxx:84
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...
Double_t min(const char *rname=0) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
const T & arg() const
Return reference to object held in proxy.
Double_t max(const char *rname=0) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
auto * a
Definition: textangle.C:12