Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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
33
34/// \brief Constructor for RooBMixDecay.
35///
36/// Creates an instance of RooBMixDecay with the specified parameters.
37///
38/// \param[in] name The name of the PDF.
39/// \param[in] title The title of the PDF.
40/// \param[in] t The time variable.
41/// \param[in] mixState The mixing state category.
42/// \param[in] tagFlav The flavour of tagged B0 category.
43/// \param[in] tau The mixing life time parameter.
44/// \param[in] dm The mixing frequency parameter.
45/// \param[in] mistag The mistag rate parameter.
46/// \param[in] delMistag The delta mistag rate parameter.
47/// \param[in] model The resolution model.
48/// \param[in] type The decay type.
49
50RooBMixDecay::RooBMixDecay(const char *name, const char *title,
53 RooAbsReal& tau, RooAbsReal& dm,
55 const RooResolutionModel& model,
57 RooAbsAnaConvPdf(name,title,model,t),
58 _type(type),
59 _mistag("mistag","Mistag rate",this,mistag),
60 _delMistag("delMistag","Delta mistag rate",this,delMistag),
61 _mixState("mixState","Mixing state",this,mixState),
62 _tagFlav("tagFlav","Flavour of tagged B0",this,tagFlav),
63 _tau("tau","Mixing life time",this,tau),
64 _dm("dm","Mixing frequency",this,dm),
65 _t("_t","time",this,t), _genMixFrac(0)
66{
67 switch(type) {
68 case SingleSided:
69 _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau)) ;
70 _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
71 break ;
72 case Flipped:
73 _basisExp = declareBasis("exp(@0/@1)",RooArgList(tau)) ;
74 _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
75 break ;
76 case DoubleSided:
77 _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau)) ;
78 _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
79 break ;
80 }
81}
82
83////////////////////////////////////////////////////////////////////////////////
84/// Copy constructor
85
88 _type(other._type),
89 _mistag("mistag",this,other._mistag),
90 _delMistag("delMistag",this,other._delMistag),
91 _mixState("mixState",this,other._mixState),
92 _tagFlav("tagFlav",this,other._tagFlav),
93 _tau("tau",this,other._tau),
94 _dm("dm",this,other._dm),
95 _t("t",this,other._t),
96 _basisExp(other._basisExp),
97 _basisCos(other._basisCos),
98 _genMixFrac(other._genMixFrac),
99 _genFlavFrac(other._genFlavFrac),
100 _genFlavFracMix(other._genFlavFracMix),
101 _genFlavFracUnmix(other._genFlavFracUnmix)
102{
103}
104
105////////////////////////////////////////////////////////////////////////////////
106/// Comp with tFit MC: must be (1 - tagFlav*...)
107
109{
110 if (basisIndex==_basisExp) {
111 return (1 - _tagFlav*_delMistag) ;
112 }
113
114 if (basisIndex==_basisCos) {
115 return _mixState*(1-2*_mistag) ;
116 }
117
118 return 0 ;
119}
120
122{
124 {ctx.at(&_convSet[0]), ctx.at(&_convSet[1]), ctx.at(_tagFlav),
125 ctx.at(_delMistag), ctx.at(_mixState), ctx.at(_mistag)});
126}
127
128////////////////////////////////////////////////////////////////////////////////
129/// std::cout << "RooBMixDecay::getCoefAI " ; allVars.Print("1") ;
130
131Int_t RooBMixDecay::getCoefAnalyticalIntegral(Int_t /*code*/, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
132{
133 if (rangeName) {
134 return 0 ;
135 }
136
137 if (matchArgs(allVars,analVars,_mixState,_tagFlav)) return 3 ;
138 if (matchArgs(allVars,analVars,_mixState)) return 2 ;
139 if (matchArgs(allVars,analVars,_tagFlav)) return 1 ;
140 return 0 ;
141}
142
143////////////////////////////////////////////////////////////////////////////////
144
145double RooBMixDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* /*rangeName*/) const
146{
147 switch(code) {
148 // No integration
149 case 0: return coefficient(basisIndex) ;
150
151 // Integration over 'mixState' and 'tagFlav'
152 case 3:
153 if (basisIndex==_basisExp) {
154 return 4.0 ;
155 }
156 if (basisIndex==_basisCos) {
157 return 0.0 ;
158 }
159 break ;
160
161 // Integration over 'mixState'
162 case 2:
163 if (basisIndex==_basisExp) {
164 return 2.0*coefficient(basisIndex) ;
165 }
166 if (basisIndex==_basisCos) {
167 return 0.0 ;
168 }
169 break ;
170
171 // Integration over 'tagFlav'
172 case 1:
173 if (basisIndex==_basisExp) {
174 return 2.0 ;
175 }
176 if (basisIndex==_basisCos) {
177 return 2.0*coefficient(basisIndex) ;
178 }
179 break ;
180
181 default:
182 assert(0) ;
183 }
184
185 return 0 ;
186}
187
188////////////////////////////////////////////////////////////////////////////////
189
191{
192 if (staticInitOK) {
196 }
197
198 if (matchArgs(directVars,generateVars,_t)) return 1 ;
199 return 0 ;
200}
201
202////////////////////////////////////////////////////////////////////////////////
203
205{
206 switch (code) {
207 case 2:
208 {
209 // Calculate the fraction of B0bar events to generate
210 double sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tagFlav.arg())).getVal() ;
211 _tagFlav = 1 ; // B0
212 double flavInt = RooRealIntegral("flavInt","flav integral",*this,RooArgSet(_t.arg())).getVal() ;
214 break ;
215 }
216 case 3:
217 {
218 // Calculate the fraction of mixed events to generate
219 double sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_mixState.arg())).getVal() ;
220 _mixState = -1 ; // mixed
221 double mixInt = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
223 break ;
224 }
225 case 4:
226 {
227 // Calculate the fraction of mixed events to generate
228 double sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_mixState.arg(),_tagFlav.arg())).getVal() ;
229 _mixState = -1 ; // mixed
230 double mixInt = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg(),_tagFlav.arg())).getVal() ;
232
233 // Calculate the fraction of B0bar tags for mixed and unmixed
234 RooRealIntegral dtInt("mixInt","mix integral",*this,RooArgSet(_t.arg())) ;
235 _mixState = -1 ; // Mixed
236 _tagFlav = 1 ; // B0
237 _genFlavFracMix = dtInt.getVal() / mixInt ;
238 _mixState = 1 ; // Unmixed
239 _tagFlav = 1 ; // B0
240 _genFlavFracUnmix = dtInt.getVal() / (sumInt - mixInt) ;
241 break ;
242 }
243 }
244}
245
246////////////////////////////////////////////////////////////////////////////////
247/// Generate mix-state dependent
248
250{
251 switch(code) {
252 case 2:
253 {
254 double rand = RooRandom::uniform() ;
255 _tagFlav = (Int_t) ((rand<=_genFlavFrac) ? 1 : -1) ;
256 break ;
257 }
258 case 3:
259 {
260 double rand = RooRandom::uniform() ;
261 _mixState = (Int_t) ((rand<=_genMixFrac) ? -1 : 1) ;
262 break ;
263 }
264 case 4:
265 {
266 double rand = RooRandom::uniform() ;
267 _mixState = (Int_t) ((rand<=_genMixFrac) ? -1 : 1) ;
268
271 _tagFlav = (Int_t) ((rand<=genFlavFrac) ? 1 : -1) ;
272 break ;
273 }
274 }
275
276 // Generate delta-t dependent
277 while(true) {
278 double rand = RooRandom::uniform() ;
279 double tval(0) ;
280
281 switch(_type) {
282 case SingleSided:
283 tval = -_tau*log(rand);
284 break ;
285 case Flipped:
286 tval= +_tau*log(rand);
287 break ;
288 case DoubleSided:
289 tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
290 break ;
291 }
292
293 // Accept event if T is in generated range
294 double dil = 1-2.*_mistag ;
295 double maxAcceptProb = 1 + std::abs(_delMistag) + std::abs(dil) ;
296 double acceptProb = (1-_tagFlav*_delMistag) + _mixState*dil*cos(_dm*tval);
298
299 if (tval<_t.max() && tval>_t.min() && mixAccept) {
300 _t = tval ;
301 break ;
302 }
303 }
304}
int Int_t
Definition RtypesCore.h:45
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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.
friend class RooRealIntegral
Definition RooAbsArg.h:572
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
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
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
std::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:77
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={})