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