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////////////////////////////////////////////////////////////////////////////////
36/// Constructor.
37
38RooBCPEffDecay::RooBCPEffDecay(const char *name, const char *title,
40 RooAbsReal& tau, RooAbsReal& dm,
41 RooAbsReal& avgMistag, RooAbsReal& CPeigenval,
43 RooAbsReal& effRatio, RooAbsReal& delMistag,
44 const RooResolutionModel& model, DecayType type) :
45 RooAbsAnaConvPdf(name,title,model,t),
46 _absLambda("absLambda","Absolute value of lambda",this,a),
47 _argLambda("argLambda","Arg(Lambda)",this,b),
48 _effRatio("effRatio","B0/B0bar efficiency ratio",this,effRatio),
49 _CPeigenval("CPeigenval","CP eigen value",this,CPeigenval),
50 _avgMistag("avgMistag","Average mistag rate",this,avgMistag),
51 _delMistag("delMistag","Delta mistag rate",this,delMistag),
52 _t("t","time",this,t),
53 _tau("tau","decay time",this,tau),
54 _dm("dm","mixing frequency",this,dm),
55 _tag("tag","CP state",this,tag),
56 _genB0Frac(0),
57 _type(type)
58{
59 switch(type) {
60 case SingleSided:
61 _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau,dm)) ;
62 _basisSin = declareBasis("exp(-@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
63 _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
64 break ;
65 case Flipped:
66 _basisExp = declareBasis("exp(@0)/@1)",RooArgList(tau,dm)) ;
67 _basisSin = declareBasis("exp(@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
68 _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
69 break ;
70 case DoubleSided:
71 _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau,dm)) ;
72 _basisSin = declareBasis("exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
73 _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
74 break ;
75 }
76}
77
78////////////////////////////////////////////////////////////////////////////////
79/// Copy constructor.
80
83 _absLambda("absLambda",this,other._absLambda),
84 _argLambda("argLambda",this,other._argLambda),
85 _effRatio("effRatio",this,other._effRatio),
86 _CPeigenval("CPeigenval",this,other._CPeigenval),
87 _avgMistag("avgMistag",this,other._avgMistag),
88 _delMistag("delMistag",this,other._delMistag),
89 _t("t",this,other._t),
90 _tau("tau",this,other._tau),
91 _dm("dm",this,other._dm),
92 _tag("tag",this,other._tag),
93 _genB0Frac(other._genB0Frac),
94 _type(other._type),
95 _basisExp(other._basisExp),
96 _basisSin(other._basisSin),
97 _basisCos(other._basisCos)
98{
99}
100
101////////////////////////////////////////////////////////////////////////////////
102/// Destructor.
103
105{
106}
107
108////////////////////////////////////////////////////////////////////////////////
109/// B0 : _tag = +1
110///
111/// B0bar : _tag = -1
112/// \param[in] basisIndex
113
114double RooBCPEffDecay::coefficient(Int_t basisIndex) const
115{
116 if (basisIndex==_basisExp) {
117 //exp term: (1 -/+ dw)(1+a^2)/2
118 return (1 - _tag*_delMistag)*(1+_absLambda*_absLambda)/2 ;
119 // = 1 + _tag*deltaDil/2
120 }
121
122 if (basisIndex==_basisSin) {
123 //sin term: +/- (1-2w)*ImLambda(= -etaCP * |l| * sin2b)
125 // = _tag*avgDil * ...
126 }
127
128 if (basisIndex==_basisCos) {
129 //cos term: +/- (1-2w)*(1-a^2)/2
130 return -1*_tag*(1-2*_avgMistag)*(1-_absLambda*_absLambda)/2 ;
131 // = -_tag*avgDil * ...
132 }
133
134 return 0 ;
135}
136
137////////////////////////////////////////////////////////////////////////////////
138
139Int_t RooBCPEffDecay::getCoefAnalyticalIntegral(Int_t /*code*/, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
140{
141 if (rangeName) return 0 ;
142
143 if (matchArgs(allVars,analVars,_tag)) return 1 ;
144 return 0 ;
145}
146
147////////////////////////////////////////////////////////////////////////////////
148
149double RooBCPEffDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* /*rangeName*/) const
150{
151 switch(code) {
152 // No integration
153 case 0: return coefficient(basisIndex) ;
154
155 // Integration over 'tag'
156 case 1:
157 if (basisIndex==_basisExp) {
158 return (1+_absLambda*_absLambda) ;
159 }
160
161 if (basisIndex==_basisSin || basisIndex==_basisCos) {
162 return 0 ;
163 }
164 break ;
165
166 default:
167 assert(0) ;
168 }
169
170 return 0 ;
171}
172
173////////////////////////////////////////////////////////////////////////////////
174
175Int_t RooBCPEffDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool staticInitOK) const
176{
177 if (staticInitOK) {
178 if (matchArgs(directVars,generateVars,_t,_tag)) return 2 ;
179 }
180 if (matchArgs(directVars,generateVars,_t)) return 1 ;
181 return 0 ;
182}
183
184////////////////////////////////////////////////////////////////////////////////
185
187{
188 if (code==2) {
189 // Calculate the fraction of mixed events to generate
190 double sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tag.arg())).getVal() ;
191 _tag = 1 ;
192 double b0Int = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
193 _genB0Frac = b0Int/sumInt ;
194 }
195}
196
197////////////////////////////////////////////////////////////////////////////////
198/// Generates mix-state dependent.
199/// \param[in] code
200
202{
203 if (code==2) {
204 double rand = RooRandom::uniform() ;
205 _tag = (rand<=_genB0Frac) ? 1 : -1 ;
206 }
207
208 // Generate delta-t dependent
209 while(1) {
210 double rand = RooRandom::uniform() ;
211 double tval(0) ;
212
213 switch(_type) {
214 case SingleSided:
215 tval = -_tau*log(rand);
216 break ;
217 case Flipped:
218 tval= +_tau*log(rand);
219 break ;
220 case DoubleSided:
221 tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
222 break ;
223 }
224
225 // Accept event if T is in generated range
226 double maxDil = 1.0 ;
227 double al2 = _absLambda*_absLambda ;
228 double maxAcceptProb = (1+al2) + std::abs(maxDil*_CPeigenval*_absLambda*_argLambda) + std::abs(maxDil*(1-al2)/2);
229 double acceptProb = (1+al2)/2*(1-_tag*_delMistag)
230 - (_tag*(1-2*_avgMistag))*(_CPeigenval*_absLambda*_argLambda)*sin(_dm*tval)
231 - (_tag*(1-2*_avgMistag))*(1-al2)/2*cos(_dm*tval);
232
233 bool accept = maxAcceptProb*RooRandom::uniform() < acceptProb ? true : false ;
234
235 if (tval<_t.max() && tval>_t.min() && accept) {
236 _t = tval ;
237 break ;
238 }
239 }
240
241}
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
#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
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.
A space to attach TBranches.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
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
~RooBCPEffDecay() override
Destructor.
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
RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
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.