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