Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooBCPGenDecay.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * JGS, Jim Smith , University of Colorado, jgsmith@pizero.colorado.edu *
7 * History:
8 * 15-Aug-2002 JGS Created initial version
9 * 11-Sep-2002 JGS Mods to introduce mu (Mirna van Hoek, JGS, Nick Danielson)
10 * *
11 * Copyright (c) 2000-2005, Regents of the University of California, *
12 * University of Colorado *
13 * and Stanford University. All rights reserved. *
14 * *
15 * Redistribution and use in source and binary forms, *
16 * with or without modification, are permitted according to the terms *
17 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
18 *****************************************************************************/
19
20/** \class RooBCPGenDecay
21 \ingroup Roofit
22
23Implement standard CP physics model with S and C (no mention of lambda)
24Suitably stolen and modified from RooBCPEffDecay
25**/
26
27#include "RooRealVar.h"
28#include "RooRandom.h"
29#include "RooBCPGenDecay.h"
30#include "RooRealIntegral.h"
31
32
33/// \brief Constructor for RooBCPGenDecay.
34///
35/// Creates an instance of RooBCPGenDecay with the specified parameters.
36///
37/// \param[in] name The name of the PDF.
38/// \param[in] title The title of the PDF.
39/// \param[in] t The time variable.
40/// \param[in] tag The CP state category.
41/// \param[in] tau The decay time parameter.
42/// \param[in] dm The mixing frequency parameter.
43/// \param[in] avgMistag The average mistag rate parameter.
44/// \param[in] avgC Coefficient of cos term.
45/// \param[in] avgS Coefficient of sin term.
46/// \param[in] delMistag Delta mistag rate parameter.
47/// \param[in] mu Tag efficiency difference parameter.
48/// \param[in] model The resolution model.
49/// \param[in] type The decay type.
50
51RooBCPGenDecay::RooBCPGenDecay(const char *name, const char *title,
53 RooAbsReal& tau, RooAbsReal& dm,
57 RooAbsReal& mu,
58 const RooResolutionModel& model, DecayType type) :
59 RooAbsAnaConvPdf(name,title,model,t),
60 _avgC("C","Coefficient of cos term",this,avgC),
61 _avgS("S","Coefficient of sin term",this,avgS),
62 _avgMistag("avgMistag","Average mistag rate",this,avgMistag),
63 _delMistag("delMistag","Delta mistag rate",this,delMistag),
64 _mu("mu","Tag efficiency difference",this,mu),
65 _t("t","time",this,t),
66 _tau("tau","decay time",this,tau),
67 _dm("dm","mixing frequency",this,dm),
68 _tag("tag","CP state",this,tag),
69 _genB0Frac(0),
70 _type(type)
71{
72 switch(type) {
73 case SingleSided:
74 _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau,dm)) ;
75 _basisSin = declareBasis("exp(-@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
76 _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
77 break ;
78 case Flipped:
79 _basisExp = declareBasis("exp(@0)/@1)",RooArgList(tau,dm)) ;
80 _basisSin = declareBasis("exp(@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
81 _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
82 break ;
83 case DoubleSided:
84 _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau,dm)) ;
85 _basisSin = declareBasis("exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
86 _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
87 break ;
88 }
89}
90
91////////////////////////////////////////////////////////////////////////////////
92/// Copy constructor
93
96 _avgC("C",this,other._avgC),
97 _avgS("S",this,other._avgS),
98 _avgMistag("avgMistag",this,other._avgMistag),
99 _delMistag("delMistag",this,other._delMistag),
100 _mu("mu",this,other._mu),
101 _t("t",this,other._t),
102 _tau("tau",this,other._tau),
103 _dm("dm",this,other._dm),
104 _tag("tag",this,other._tag),
105 _genB0Frac(other._genB0Frac),
106 _type(other._type),
107 _basisExp(other._basisExp),
108 _basisSin(other._basisSin),
109 _basisCos(other._basisCos)
110{
111}
112
113////////////////////////////////////////////////////////////////////////////////
114/// B0 : _tag = +1
115/// B0bar : _tag = -1
116
118{
119 if (basisIndex==_basisExp) {
120 //exp term: (1 -/+ dw + mu*_tag*w)
121 return (1 - _tag*_delMistag + _mu*_tag*(1. - 2.*_avgMistag)) ;
122 // = 1 + _tag*deltaDil/2 + _mu*avgDil
123 }
124
125 if (basisIndex==_basisSin) {
126 //sin term: (+/- (1-2w) + _mu*(1 -/+ delw))*S
127 return (_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgS ;
128 // = (_tag*avgDil + _mu*(1 + tag*deltaDil/2)) * S
129 }
130
131 if (basisIndex==_basisCos) {
132 //cos term: -(+/- (1-2w) + _mu*(1 -/+ delw))*C
133 return -1.*(_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgC ;
134 // = -(_tag*avgDil + _mu*(1 + _tag*deltaDil/2) )* C
135 }
136
137 return 0 ;
138}
139
140////////////////////////////////////////////////////////////////////////////////
141
142Int_t RooBCPGenDecay::getCoefAnalyticalIntegral(Int_t /*code*/, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
143{
144 if (rangeName) return 0 ;
145 if (matchArgs(allVars,analVars,_tag)) return 1 ;
146 return 0 ;
147}
148
149////////////////////////////////////////////////////////////////////////////////
150
151double RooBCPGenDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* /*rangeName*/) const
152{
153 switch(code) {
154 // No integration
155 case 0: return coefficient(basisIndex) ;
156
157 // Integration over 'tag'
158 case 1:
159 if (basisIndex==_basisExp) {
160 return 2 ;
161 }
162
163 if (basisIndex==_basisSin) {
164 return 2*_mu*_avgS ;
165 }
166 if (basisIndex==_basisCos) {
167 return -2*_mu*_avgC ;
168 }
169 break ;
170
171 default:
172 assert(0) ;
173 }
174
175 return 0 ;
176}
177
178////////////////////////////////////////////////////////////////////////////////
179
181{
182 if (staticInitOK) {
183 if (matchArgs(directVars,generateVars,_t,_tag)) return 2 ;
184 }
185 if (matchArgs(directVars,generateVars,_t)) return 1 ;
186 return 0 ;
187}
188
189////////////////////////////////////////////////////////////////////////////////
190
192{
193 if (code==2) {
194 // Calculate the fraction of mixed events to generate
195 double sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tag.arg())).getVal() ;
196 _tag = 1 ;
197 double b0Int = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
199 }
200}
201
202////////////////////////////////////////////////////////////////////////////////
203/// Generate mix-state dependent
204
206{
207 if (code==2) {
208 double rand = RooRandom::uniform() ;
209 _tag = (rand<=_genB0Frac) ? 1 : -1 ;
210 }
211
212 // Generate delta-t dependent
213 while(true) {
214 double rand = RooRandom::uniform() ;
215 double tval(0) ;
216
217 switch(_type) {
218 case SingleSided:
219 tval = -_tau*log(rand);
220 break ;
221 case Flipped:
222 tval= +_tau*log(rand);
223 break ;
224 case DoubleSided:
225 tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
226 break ;
227 }
228
229 // Accept event if T is in generated range
230 double maxDil = 1.0 ;
231// 2 in next line is conservative and inefficient - allows for delMistag=1!
232 double maxAcceptProb = 2 + std::abs(maxDil*_avgS) + std::abs(maxDil*_avgC);
233 double acceptProb = (1-_tag*_delMistag + _mu*_tag*(1. - 2.*_avgMistag))
234 + (_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgS*sin(_dm*tval)
235 - (_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgC*cos(_dm*tval);
236
238
239 if (tval<_t.max() && tval>_t.min() && accept) {
240 _t = tval ;
241 break ;
242 }
243 }
244
245}
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
Implement standard CP physics model with S and C (no mention of lambda) Suitably stolen and modified ...
RooRealProxy _mu
RooRealProxy _avgC
RooRealProxy _avgS
void initGenerator(Int_t code) override
RooRealProxy _dm
RooRealProxy _tau
void generateEvent(Int_t code) override
Generate mix-state dependent.
RooCategoryProxy _tag
RooRealProxy _t
double coefficient(Int_t basisIndex) const override
B0 : _tag = +1 B0bar : _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 _avgMistag
RooRealProxy _delMistag
double coefAnalyticalIntegral(Int_t coef, Int_t code, const char *rangeName=nullptr) const override
Default implementation of function implementing advertised integrals.
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const override
static double uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition RooRandom.cxx:77
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.