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
33
34/// \brief Constructor for RooBCPGenDecay.
35///
36/// Creates an instance of RooBCPGenDecay 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] tag The CP state category.
42/// \param[in] tau The decay time parameter.
43/// \param[in] dm The mixing frequency parameter.
44/// \param[in] avgMistag The average mistag rate parameter.
45/// \param[in] avgC Coefficient of cos term.
46/// \param[in] avgS Coefficient of sin term.
47/// \param[in] delMistag Delta mistag rate parameter.
48/// \param[in] mu Tag efficiency difference parameter.
49/// \param[in] model The resolution model.
50/// \param[in] type The decay type.
51
52RooBCPGenDecay::RooBCPGenDecay(const char *name, const char *title,
54 RooAbsReal& tau, RooAbsReal& dm,
55 RooAbsReal& avgMistag,
56 RooAbsReal& avgC, RooAbsReal& avgS,
57 RooAbsReal& delMistag,
58 RooAbsReal& mu,
59 const RooResolutionModel& model, DecayType type) :
60 RooAbsAnaConvPdf(name,title,model,t),
61 _avgC("C","Coefficient of cos term",this,avgC),
62 _avgS("S","Coefficient of sin term",this,avgS),
63 _avgMistag("avgMistag","Average mistag rate",this,avgMistag),
64 _delMistag("delMistag","Delta mistag rate",this,delMistag),
65 _mu("mu","Tag efficiency difference",this,mu),
66 _t("t","time",this,t),
67 _tau("tau","decay time",this,tau),
68 _dm("dm","mixing frequency",this,dm),
69 _tag("tag","CP state",this,tag),
70 _genB0Frac(0),
71 _type(type)
72{
73 switch(type) {
74 case SingleSided:
75 _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau,dm)) ;
76 _basisSin = declareBasis("exp(-@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
77 _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
78 break ;
79 case Flipped:
80 _basisExp = declareBasis("exp(@0)/@1)",RooArgList(tau,dm)) ;
81 _basisSin = declareBasis("exp(@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
82 _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
83 break ;
84 case DoubleSided:
85 _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau,dm)) ;
86 _basisSin = declareBasis("exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
87 _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
88 break ;
89 }
90}
91
92////////////////////////////////////////////////////////////////////////////////
93/// Copy constructor
94
97 _avgC("C",this,other._avgC),
98 _avgS("S",this,other._avgS),
99 _avgMistag("avgMistag",this,other._avgMistag),
100 _delMistag("delMistag",this,other._delMistag),
101 _mu("mu",this,other._mu),
102 _t("t",this,other._t),
103 _tau("tau",this,other._tau),
104 _dm("dm",this,other._dm),
105 _tag("tag",this,other._tag),
106 _genB0Frac(other._genB0Frac),
107 _type(other._type),
108 _basisExp(other._basisExp),
109 _basisSin(other._basisSin),
110 _basisCos(other._basisCos)
111{
112}
113
114////////////////////////////////////////////////////////////////////////////////
115/// B0 : _tag = +1
116/// B0bar : _tag = -1
117
118double RooBCPGenDecay::coefficient(Int_t basisIndex) const
119{
120 if (basisIndex==_basisExp) {
121 //exp term: (1 -/+ dw + mu*_tag*w)
122 return (1 - _tag*_delMistag + _mu*_tag*(1. - 2.*_avgMistag)) ;
123 // = 1 + _tag*deltaDil/2 + _mu*avgDil
124 }
125
126 if (basisIndex==_basisSin) {
127 //sin term: (+/- (1-2w) + _mu*(1 -/+ delw))*S
128 return (_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgS ;
129 // = (_tag*avgDil + _mu*(1 + tag*deltaDil/2)) * S
130 }
131
132 if (basisIndex==_basisCos) {
133 //cos term: -(+/- (1-2w) + _mu*(1 -/+ delw))*C
134 return -1.*(_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgC ;
135 // = -(_tag*avgDil + _mu*(1 + _tag*deltaDil/2) )* C
136 }
137
138 return 0 ;
139}
140
141////////////////////////////////////////////////////////////////////////////////
142
143Int_t RooBCPGenDecay::getCoefAnalyticalIntegral(Int_t /*code*/, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
144{
145 if (rangeName) return 0 ;
146 if (matchArgs(allVars,analVars,_tag)) return 1 ;
147 return 0 ;
148}
149
150////////////////////////////////////////////////////////////////////////////////
151
152double RooBCPGenDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* /*rangeName*/) const
153{
154 switch(code) {
155 // No integration
156 case 0: return coefficient(basisIndex) ;
157
158 // Integration over 'tag'
159 case 1:
160 if (basisIndex==_basisExp) {
161 return 2 ;
162 }
163
164 if (basisIndex==_basisSin) {
165 return 2*_mu*_avgS ;
166 }
167 if (basisIndex==_basisCos) {
168 return -2*_mu*_avgC ;
169 }
170 break ;
171
172 default:
173 assert(0) ;
174 }
175
176 return 0 ;
177}
178
179////////////////////////////////////////////////////////////////////////////////
180
181Int_t RooBCPGenDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool staticInitOK) const
182{
183 if (staticInitOK) {
184 if (matchArgs(directVars,generateVars,_t,_tag)) return 2 ;
185 }
186 if (matchArgs(directVars,generateVars,_t)) return 1 ;
187 return 0 ;
188}
189
190////////////////////////////////////////////////////////////////////////////////
191
193{
194 if (code==2) {
195 // Calculate the fraction of mixed events to generate
196 double sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tag.arg())).getVal() ;
197 _tag = 1 ;
198 double b0Int = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
199 _genB0Frac = b0Int/sumInt ;
200 }
201}
202
203////////////////////////////////////////////////////////////////////////////////
204/// Generate mix-state dependent
205
207{
208 if (code==2) {
209 double rand = RooRandom::uniform() ;
210 _tag = (rand<=_genB0Frac) ? 1 : -1 ;
211 }
212
213 // Generate delta-t dependent
214 while(true) {
215 double rand = RooRandom::uniform() ;
216 double tval(0) ;
217
218 switch(_type) {
219 case SingleSided:
220 tval = -_tau*log(rand);
221 break ;
222 case Flipped:
223 tval= +_tau*log(rand);
224 break ;
225 case DoubleSided:
226 tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
227 break ;
228 }
229
230 // Accept event if T is in generated range
231 double maxDil = 1.0 ;
232// 2 in next line is conservative and inefficient - allows for delMistag=1!
233 double maxAcceptProb = 2 + std::abs(maxDil*_avgS) + std::abs(maxDil*_avgC);
234 double acceptProb = (1-_tag*_delMistag + _mu*_tag*(1. - 2.*_avgMistag))
235 + (_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgS*sin(_dm*tval)
236 - (_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgC*cos(_dm*tval);
237
238 bool accept = maxAcceptProb*RooRandom::uniform() < acceptProb ? true : false ;
239
240 if (tval<_t.max() && tval>_t.min() && accept) {
241 _t = tval ;
242 break ;
243 }
244 }
245
246}
#define ClassImp(name)
Definition Rtypes.h:382
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: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
Interface for one-time initialization to setup the generator for the specified code.
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
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
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.