Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooBMixDecay.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 RooBMixDecay
18 \ingroup Roofit
19
20Class RooBMixDecay is a RooAbsAnaConvPdf implementation that describes
21the decay of B mesons with the effects of B0/B0bar mixing.
22This function can be analytically convolved with any RooResolutionModel implementation
23**/
24
25#include "Riostream.h"
26#include "TMath.h"
27#include "RooRealVar.h"
28#include "RooBMixDecay.h"
29#include "RooRealIntegral.h"
30#include "RooRandom.h"
31#include "RooBatchCompute.h"
32
34
35/// \brief Constructor for RooBMixDecay.
36///
37/// Creates an instance of RooBMixDecay 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] mixState The mixing state category.
43/// \param[in] tagFlav The flavour of tagged B0 category.
44/// \param[in] tau The mixing life time parameter.
45/// \param[in] dm The mixing frequency parameter.
46/// \param[in] mistag The mistag rate parameter.
47/// \param[in] delMistag The delta mistag rate parameter.
48/// \param[in] model The resolution model.
49/// \param[in] type The decay type.
50
51RooBMixDecay::RooBMixDecay(const char *name, const char *title,
52 RooRealVar& t, RooAbsCategory& mixState,
53 RooAbsCategory& tagFlav,
54 RooAbsReal& tau, RooAbsReal& dm,
55 RooAbsReal& mistag, RooAbsReal& delMistag,
56 const RooResolutionModel& model,
58 RooAbsAnaConvPdf(name,title,model,t),
59 _type(type),
60 _mistag("mistag","Mistag rate",this,mistag),
61 _delMistag("delMistag","Delta mistag rate",this,delMistag),
62 _mixState("mixState","Mixing state",this,mixState),
63 _tagFlav("tagFlav","Flavour of tagged B0",this,tagFlav),
64 _tau("tau","Mixing life time",this,tau),
65 _dm("dm","Mixing frequency",this,dm),
66 _t("_t","time",this,t), _genMixFrac(0)
67{
68 switch(type) {
69 case SingleSided:
70 _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau)) ;
71 _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
72 break ;
73 case Flipped:
74 _basisExp = declareBasis("exp(@0/@1)",RooArgList(tau)) ;
75 _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
76 break ;
77 case DoubleSided:
78 _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau)) ;
79 _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
80 break ;
81 }
82}
83
84////////////////////////////////////////////////////////////////////////////////
85/// Copy constructor
86
87RooBMixDecay::RooBMixDecay(const RooBMixDecay& other, const char* name) :
89 _type(other._type),
90 _mistag("mistag",this,other._mistag),
91 _delMistag("delMistag",this,other._delMistag),
92 _mixState("mixState",this,other._mixState),
93 _tagFlav("tagFlav",this,other._tagFlav),
94 _tau("tau",this,other._tau),
95 _dm("dm",this,other._dm),
96 _t("t",this,other._t),
97 _basisExp(other._basisExp),
98 _basisCos(other._basisCos),
99 _genMixFrac(other._genMixFrac),
100 _genFlavFrac(other._genFlavFrac),
101 _genFlavFracMix(other._genFlavFracMix),
102 _genFlavFracUnmix(other._genFlavFracUnmix)
103{
104}
105
106////////////////////////////////////////////////////////////////////////////////
107/// Comp with tFit MC: must be (1 - tagFlav*...)
108
109double RooBMixDecay::coefficient(Int_t basisIndex) const
110{
111 if (basisIndex==_basisExp) {
112 return (1 - _tagFlav*_delMistag) ;
113 }
114
115 if (basisIndex==_basisCos) {
116 return _mixState*(1-2*_mistag) ;
117 }
118
119 return 0 ;
120}
121
123{
125 {ctx.at(&_convSet[0]), ctx.at(&_convSet[1]), ctx.at(_tagFlav),
126 ctx.at(_delMistag), ctx.at(_mixState), ctx.at(_mistag)});
127}
128
129////////////////////////////////////////////////////////////////////////////////
130/// cout << "RooBMixDecay::getCoefAI " ; allVars.Print("1") ;
131
132Int_t RooBMixDecay::getCoefAnalyticalIntegral(Int_t /*code*/, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
133{
134 if (rangeName) {
135 return 0 ;
136 }
137
138 if (matchArgs(allVars,analVars,_mixState,_tagFlav)) return 3 ;
139 if (matchArgs(allVars,analVars,_mixState)) return 2 ;
140 if (matchArgs(allVars,analVars,_tagFlav)) return 1 ;
141 return 0 ;
142}
143
144////////////////////////////////////////////////////////////////////////////////
145
146double RooBMixDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* /*rangeName*/) const
147{
148 switch(code) {
149 // No integration
150 case 0: return coefficient(basisIndex) ;
151
152 // Integration over 'mixState' and 'tagFlav'
153 case 3:
154 if (basisIndex==_basisExp) {
155 return 4.0 ;
156 }
157 if (basisIndex==_basisCos) {
158 return 0.0 ;
159 }
160 break ;
161
162 // Integration over 'mixState'
163 case 2:
164 if (basisIndex==_basisExp) {
165 return 2.0*coefficient(basisIndex) ;
166 }
167 if (basisIndex==_basisCos) {
168 return 0.0 ;
169 }
170 break ;
171
172 // Integration over 'tagFlav'
173 case 1:
174 if (basisIndex==_basisExp) {
175 return 2.0 ;
176 }
177 if (basisIndex==_basisCos) {
178 return 2.0*coefficient(basisIndex) ;
179 }
180 break ;
181
182 default:
183 assert(0) ;
184 }
185
186 return 0 ;
187}
188
189////////////////////////////////////////////////////////////////////////////////
190
191Int_t RooBMixDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool staticInitOK) const
192{
193 if (staticInitOK) {
194 if (matchArgs(directVars,generateVars,_t,_mixState,_tagFlav)) return 4 ;
195 if (matchArgs(directVars,generateVars,_t,_mixState)) return 3 ;
196 if (matchArgs(directVars,generateVars,_t,_tagFlav)) return 2 ;
197 }
198
199 if (matchArgs(directVars,generateVars,_t)) return 1 ;
200 return 0 ;
201}
202
203////////////////////////////////////////////////////////////////////////////////
204
206{
207 switch (code) {
208 case 2:
209 {
210 // Calculate the fraction of B0bar events to generate
211 double sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tagFlav.arg())).getVal() ;
212 _tagFlav = 1 ; // B0
213 double flavInt = RooRealIntegral("flavInt","flav integral",*this,RooArgSet(_t.arg())).getVal() ;
214 _genFlavFrac = flavInt/sumInt ;
215 break ;
216 }
217 case 3:
218 {
219 // Calculate the fraction of mixed events to generate
220 double sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_mixState.arg())).getVal() ;
221 _mixState = -1 ; // mixed
222 double mixInt = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
223 _genMixFrac = mixInt/sumInt ;
224 break ;
225 }
226 case 4:
227 {
228 // Calculate the fraction of mixed events to generate
229 double sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_mixState.arg(),_tagFlav.arg())).getVal() ;
230 _mixState = -1 ; // mixed
231 double mixInt = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg(),_tagFlav.arg())).getVal() ;
232 _genMixFrac = mixInt/sumInt ;
233
234 // Calculate the fraction of B0bar tags for mixed and unmixed
235 RooRealIntegral dtInt("mixInt","mix integral",*this,RooArgSet(_t.arg())) ;
236 _mixState = -1 ; // Mixed
237 _tagFlav = 1 ; // B0
238 _genFlavFracMix = dtInt.getVal() / mixInt ;
239 _mixState = 1 ; // Unmixed
240 _tagFlav = 1 ; // B0
241 _genFlavFracUnmix = dtInt.getVal() / (sumInt - mixInt) ;
242 break ;
243 }
244 }
245}
246
247////////////////////////////////////////////////////////////////////////////////
248/// Generate mix-state dependent
249
251{
252 switch(code) {
253 case 2:
254 {
255 double rand = RooRandom::uniform() ;
256 _tagFlav = (Int_t) ((rand<=_genFlavFrac) ? 1 : -1) ;
257 break ;
258 }
259 case 3:
260 {
261 double rand = RooRandom::uniform() ;
262 _mixState = (Int_t) ((rand<=_genMixFrac) ? -1 : 1) ;
263 break ;
264 }
265 case 4:
266 {
267 double rand = RooRandom::uniform() ;
268 _mixState = (Int_t) ((rand<=_genMixFrac) ? -1 : 1) ;
269
270 rand = RooRandom::uniform() ;
271 double genFlavFrac = (_mixState==-1) ? _genFlavFracMix : _genFlavFracUnmix ;
272 _tagFlav = (Int_t) ((rand<=genFlavFrac) ? 1 : -1) ;
273 break ;
274 }
275 }
276
277 // Generate delta-t dependent
278 while(true) {
279 double rand = RooRandom::uniform() ;
280 double tval(0) ;
281
282 switch(_type) {
283 case SingleSided:
284 tval = -_tau*log(rand);
285 break ;
286 case Flipped:
287 tval= +_tau*log(rand);
288 break ;
289 case DoubleSided:
290 tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
291 break ;
292 }
293
294 // Accept event if T is in generated range
295 double dil = 1-2.*_mistag ;
296 double maxAcceptProb = 1 + std::abs(_delMistag) + std::abs(dil) ;
297 double acceptProb = (1-_tagFlav*_delMistag) + _mixState*dil*cos(_dm*tval);
298 bool mixAccept = maxAcceptProb*RooRandom::uniform() < acceptProb ? true : false ;
299
300 if (tval<_t.max() && tval>_t.min() && mixAccept) {
301 _t = tval ;
302 break ;
303 }
304 }
305}
int Int_t
Definition RtypesCore.h:45
#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
Class RooBMixDecay is a RooAbsAnaConvPdf implementation that describes the decay of B mesons with the...
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...
RooRealProxy _tau
double _genMixFrac
double coefAnalyticalIntegral(Int_t coef, Int_t code, const char *rangeName=nullptr) const override
Default implementation of function implementing advertised integrals.
DecayType _type
double _genFlavFracMix
RooRealProxy _mistag
void doEval(RooFit::EvalContext &) const override
Base function for computing multiple values of a RooAbsReal.
RooRealProxy _t
void initGenerator(Int_t code) override
Interface for one-time initialization to setup the generator for the specified code.
void generateEvent(Int_t code) override
Generate mix-state dependent.
RooCategoryProxy _mixState
Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
cout << "RooBMixDecay::getCoefAI " ; allVars.Print("1") ;
double coefficient(Int_t basisIndex) const override
Comp with tFit MC: must be (1 - tagFlav*...)
double _genFlavFracUnmix
RooRealProxy _dm
double _genFlavFrac
do not persist
RooCategoryProxy _tagFlav
RooRealProxy _delMistag
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
static double uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition RooRandom.cxx:78
Performs hybrid numerical/analytical integrals of RooAbsReal objects.
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.
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})