Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 "Riostream.h"
26#include "RooRealVar.h"
27#include "RooRandom.h"
28#include "RooBCPEffDecay.h"
29#include "RooRealIntegral.h"
30
31using namespace std;
32
34
35/// \brief Constructor for RooBCPEffDecay.
36///
37/// Creates an instance of RooBCPEffDecay with the specified parameters.
38///
39/// \param[in] name The name of the PDF.
40/// \param[in] title The title of the PDF.
41/// \param[in] t The time variable.
42/// \param[in] tag The CP state category.
43/// \param[in] tau The decay time parameter.
44/// \param[in] dm The mixing frequency parameter.
45/// \param[in] avgMistag The average mistag rate parameter.
46/// \param[in] CPeigenval The CP eigen value parameter.
47/// \param[in] absLambda The absolute value of the complex lambda parameter.
48/// \param[in] argLambda The argument of the complex lambda parameter.
49/// \param[in] effRatio The B0/B0bar efficiency ratio.
50/// \param[in] delMistag Delta mistag rate parameter.
51/// \param[in] model The resolution model.
52/// \param[in] type The decay type.
53
54RooBCPEffDecay::RooBCPEffDecay(const char *name, const char *title,
56 RooAbsReal& tau, RooAbsReal& dm,
57 RooAbsReal& avgMistag, RooAbsReal& CPeigenval,
58 RooAbsReal& absLambda, RooAbsReal& argLambda,
59 RooAbsReal& effRatio, RooAbsReal& delMistag,
60 const RooResolutionModel& model, DecayType type) :
61 RooAbsAnaConvPdf(name,title,model,t),
62 _absLambda("absLambda","Absolute value of lambda",this,absLambda),
63 _argLambda("argLambda","Arg(Lambda)",this,argLambda),
64 _effRatio("effRatio","B0/B0bar efficiency ratio",this,effRatio),
65 _CPeigenval("CPeigenval","CP eigen value",this,CPeigenval),
66 _avgMistag("avgMistag","Average mistag rate",this,avgMistag),
67 _delMistag("delMistag","Delta mistag rate",this,delMistag),
68 _t("t","time",this,t),
69 _tau("tau","decay time",this,tau),
70 _dm("dm","mixing frequency",this,dm),
71 _tag("tag","CP state",this,tag),
72 _genB0Frac(0),
73 _type(type)
74{
75 switch(type) {
76 case SingleSided:
77 _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau,dm)) ;
78 _basisSin = declareBasis("exp(-@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
79 _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
80 break ;
81 case Flipped:
82 _basisExp = declareBasis("exp(@0)/@1)",RooArgList(tau,dm)) ;
83 _basisSin = declareBasis("exp(@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
84 _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
85 break ;
86 case DoubleSided:
87 _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau,dm)) ;
88 _basisSin = declareBasis("exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
89 _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
90 break ;
91 }
92}
93
94////////////////////////////////////////////////////////////////////////////////
95/// Copy constructor.
96
99 _absLambda("absLambda",this,other._absLambda),
100 _argLambda("argLambda",this,other._argLambda),
101 _effRatio("effRatio",this,other._effRatio),
102 _CPeigenval("CPeigenval",this,other._CPeigenval),
103 _avgMistag("avgMistag",this,other._avgMistag),
104 _delMistag("delMistag",this,other._delMistag),
105 _t("t",this,other._t),
106 _tau("tau",this,other._tau),
107 _dm("dm",this,other._dm),
108 _tag("tag",this,other._tag),
109 _genB0Frac(other._genB0Frac),
110 _type(other._type),
111 _basisExp(other._basisExp),
112 _basisSin(other._basisSin),
113 _basisCos(other._basisCos)
114{
115}
116
117////////////////////////////////////////////////////////////////////////////////
118/// B0 : _tag = +1
119///
120/// B0bar : _tag = -1
121/// \param[in] basisIndex
122
123double RooBCPEffDecay::coefficient(Int_t basisIndex) const
124{
125 if (basisIndex==_basisExp) {
126 //exp term: (1 -/+ dw)(1+a^2)/2
127 return (1 - _tag*_delMistag)*(1+_absLambda*_absLambda)/2 ;
128 // = 1 + _tag*deltaDil/2
129 }
130
131 if (basisIndex==_basisSin) {
132 //sin term: +/- (1-2w)*ImLambda(= -etaCP * |l| * sin2b)
134 // = _tag*avgDil * ...
135 }
136
137 if (basisIndex==_basisCos) {
138 //cos term: +/- (1-2w)*(1-a^2)/2
139 return -1*_tag*(1-2*_avgMistag)*(1-_absLambda*_absLambda)/2 ;
140 // = -_tag*avgDil * ...
141 }
142
143 return 0 ;
144}
145
146////////////////////////////////////////////////////////////////////////////////
147
148Int_t RooBCPEffDecay::getCoefAnalyticalIntegral(Int_t /*code*/, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
149{
150 if (rangeName) return 0 ;
151
152 if (matchArgs(allVars,analVars,_tag)) return 1 ;
153 return 0 ;
154}
155
156////////////////////////////////////////////////////////////////////////////////
157
158double RooBCPEffDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* /*rangeName*/) const
159{
160 switch(code) {
161 // No integration
162 case 0: return coefficient(basisIndex) ;
163
164 // Integration over 'tag'
165 case 1:
166 if (basisIndex==_basisExp) {
167 return (1+_absLambda*_absLambda) ;
168 }
169
170 if (basisIndex==_basisSin || basisIndex==_basisCos) {
171 return 0 ;
172 }
173 break ;
174
175 default:
176 assert(0) ;
177 }
178
179 return 0 ;
180}
181
182////////////////////////////////////////////////////////////////////////////////
183
184Int_t RooBCPEffDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool staticInitOK) const
185{
186 if (staticInitOK) {
187 if (matchArgs(directVars,generateVars,_t,_tag)) return 2 ;
188 }
189 if (matchArgs(directVars,generateVars,_t)) return 1 ;
190 return 0 ;
191}
192
193////////////////////////////////////////////////////////////////////////////////
194
196{
197 if (code==2) {
198 // Calculate the fraction of mixed events to generate
199 double sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tag.arg())).getVal() ;
200 _tag = 1 ;
201 double b0Int = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
202 _genB0Frac = b0Int/sumInt ;
203 }
204}
205
206////////////////////////////////////////////////////////////////////////////////
207/// Generates mix-state dependent.
208/// \param[in] code
209
211{
212 if (code==2) {
213 double rand = RooRandom::uniform() ;
214 _tag = (rand<=_genB0Frac) ? 1 : -1 ;
215 }
216
217 // Generate delta-t dependent
218 while(true) {
219 double rand = RooRandom::uniform() ;
220 double tval(0) ;
221
222 switch(_type) {
223 case SingleSided:
224 tval = -_tau*log(rand);
225 break ;
226 case Flipped:
227 tval= +_tau*log(rand);
228 break ;
229 case DoubleSided:
230 tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
231 break ;
232 }
233
234 // Accept event if T is in generated range
235 double maxDil = 1.0 ;
236 double al2 = _absLambda*_absLambda ;
237 double maxAcceptProb = (1+al2) + std::abs(maxDil*_CPeigenval*_absLambda*_argLambda) + std::abs(maxDil*(1-al2)/2);
238 double acceptProb = (1+al2)/2*(1-_tag*_delMistag)
239 - (_tag*(1-2*_avgMistag))*(_CPeigenval*_absLambda*_argLambda)*sin(_dm*tval)
240 - (_tag*(1-2*_avgMistag))*(1-al2)/2*cos(_dm*tval);
241
242 bool accept = maxAcceptProb*RooRandom::uniform() < acceptProb ? true : false ;
243
244 if (tval<_t.max() && tval>_t.min() && accept) {
245 _t = tval ;
246 break ;
247 }
248 }
249
250}
#define ClassImp(name)
Definition Rtypes.h:377
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.
A space to attach TBranches.
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
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:55
PDF describing decay time distribution of B meson including effects of standard model CP violation.
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.
double coefficient(Int_t basisIndex) const override
B0 : _tag = +1.
Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Default implementation of function advertising integration capabilities.
RooRealProxy _t
RooRealProxy _avgMistag
RooRealProxy _delMistag
void generateEvent(Int_t code) override
Generates mix-state dependent.
RooRealProxy _dm
RooRealProxy _absLambda
RooRealProxy _argLambda
RooCategoryProxy _tag
RooRealProxy _CPeigenval
RooRealProxy _tau
void initGenerator(Int_t code) override
Interface for one-time initialization to setup the generator for the specified code.
static double uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition RooRandom.cxx:81
Performs hybrid numerical/analytical integrals of RooAbsReal objects.
RooRealVar represents a 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.