ROOT  6.06/09
Reference Guide
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 //////////////////////////////////////////////////////////////////////////////
21 //
22 // BEGIN_HTML
23 // Implement standard CP physics model with S and C (no mention of lambda)
24 // Suitably stolen and modified from RooBCPEffDecay
25 // END_HTML
26 //
27 
28 #include "RooRealVar.h"
29 #include "RooRandom.h"
30 #include "RooBCPGenDecay.h"
31 #include "RooRealIntegral.h"
32 
33 using namespace std;
34 
36 ;
37 
38 
39 
40 ////////////////////////////////////////////////////////////////////////////////
41 /// Constructor
42 
43 RooBCPGenDecay::RooBCPGenDecay(const char *name, const char *title,
44  RooRealVar& t, RooAbsCategory& tag,
45  RooAbsReal& tau, RooAbsReal& dm,
46  RooAbsReal& avgMistag,
47  RooAbsReal& a, RooAbsReal& b,
48  RooAbsReal& delMistag,
49  RooAbsReal& mu,
50  const RooResolutionModel& model, DecayType type) :
51  RooAbsAnaConvPdf(name,title,model,t),
52  _avgC("C","Coefficient of cos term",this,a),
53  _avgS("S","Coefficient of sin term",this,b),
54  _avgMistag("avgMistag","Average mistag rate",this,avgMistag),
55  _delMistag("delMistag","Delta mistag rate",this,delMistag),
56  _mu("mu","Tagg efficiency difference",this,mu),
57  _t("t","time",this,t),
58  _tau("tau","decay time",this,tau),
59  _dm("dm","mixing frequency",this,dm),
60  _tag("tag","CP state",this,tag),
61  _genB0Frac(0),
62  _type(type)
63 {
64  switch(type) {
65  case SingleSided:
66  _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau,dm)) ;
67  _basisSin = declareBasis("exp(-@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
68  _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
69  break ;
70  case Flipped:
71  _basisExp = declareBasis("exp(@0)/@1)",RooArgList(tau,dm)) ;
72  _basisSin = declareBasis("exp(@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
73  _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
74  break ;
75  case DoubleSided:
76  _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau,dm)) ;
77  _basisSin = declareBasis("exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
78  _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
79  break ;
80  }
81 }
82 
83 
84 
85 ////////////////////////////////////////////////////////////////////////////////
86 /// Copy constructor
87 
88 RooBCPGenDecay::RooBCPGenDecay(const RooBCPGenDecay& other, const char* name) :
89  RooAbsAnaConvPdf(other,name),
90  _avgC("C",this,other._avgC),
91  _avgS("S",this,other._avgS),
92  _avgMistag("avgMistag",this,other._avgMistag),
93  _delMistag("delMistag",this,other._delMistag),
94  _mu("mu",this,other._mu),
95  _t("t",this,other._t),
96  _tau("tau",this,other._tau),
97  _dm("dm",this,other._dm),
98  _tag("tag",this,other._tag),
99  _genB0Frac(other._genB0Frac),
100  _type(other._type),
101  _basisExp(other._basisExp),
102  _basisSin(other._basisSin),
103  _basisCos(other._basisCos)
104 {
105 }
106 
107 
108 
109 ////////////////////////////////////////////////////////////////////////////////
110 /// Destructor
111 
113 {
114 }
115 
116 
117 
118 ////////////////////////////////////////////////////////////////////////////////
119 /// B0 : _tag = +1
120 /// B0bar : _tag = -1
121 
123 {
124  if (basisIndex==_basisExp) {
125  //exp term: (1 -/+ dw + mu*_tag*w)
126  return (1 - _tag*_delMistag + _mu*_tag*(1. - 2.*_avgMistag)) ;
127  // = 1 + _tag*deltaDil/2 + _mu*avgDil
128  }
129 
130  if (basisIndex==_basisSin) {
131  //sin term: (+/- (1-2w) + _mu*(1 -/+ delw))*S
132  return (_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgS ;
133  // = (_tag*avgDil + _mu*(1 + tag*deltaDil/2)) * S
134  }
135 
136  if (basisIndex==_basisCos) {
137  //cos term: -(+/- (1-2w) + _mu*(1 -/+ delw))*C
138  return -1.*(_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgC ;
139  // = -(_tag*avgDil + _mu*(1 + _tag*deltaDil/2) )* C
140  }
141 
142  return 0 ;
143 }
144 
145 
146 
147 ////////////////////////////////////////////////////////////////////////////////
148 
149 Int_t RooBCPGenDecay::getCoefAnalyticalIntegral(Int_t /*code*/, RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
150 {
151  if (rangeName) return 0 ;
152  if (matchArgs(allVars,analVars,_tag)) return 1 ;
153  return 0 ;
154 }
155 
156 
157 
158 ////////////////////////////////////////////////////////////////////////////////
159 
160 Double_t RooBCPGenDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* /*rangeName*/) const
161 {
162  switch(code) {
163  // No integration
164  case 0: return coefficient(basisIndex) ;
165 
166  // Integration over 'tag'
167  case 1:
168  if (basisIndex==_basisExp) {
169  return 2 ;
170  }
171 
172  if (basisIndex==_basisSin) {
173  return 2*_mu*_avgS ;
174  }
175  if (basisIndex==_basisCos) {
176  return -2*_mu*_avgC ;
177  }
178  break ;
179 
180  default:
181  assert(0) ;
182  }
183 
184  return 0 ;
185 }
186 
187 
188 
189 ////////////////////////////////////////////////////////////////////////////////
190 
191 Int_t RooBCPGenDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK) const
192 {
193  if (staticInitOK) {
194  if (matchArgs(directVars,generateVars,_t,_tag)) return 2 ;
195  }
196  if (matchArgs(directVars,generateVars,_t)) return 1 ;
197  return 0 ;
198 }
199 
200 
201 
202 ////////////////////////////////////////////////////////////////////////////////
203 
205 {
206  if (code==2) {
207  // Calculate the fraction of mixed events to generate
208  Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tag.arg())).getVal() ;
209  _tag = 1 ;
210  Double_t b0Int = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
211  _genB0Frac = b0Int/sumInt ;
212  }
213 }
214 
215 
216 
217 ////////////////////////////////////////////////////////////////////////////////
218 /// Generate mix-state dependent
219 
221 {
222  if (code==2) {
223  Double_t rand = RooRandom::uniform() ;
224  _tag = (rand<=_genB0Frac) ? 1 : -1 ;
225  }
226 
227  // Generate delta-t dependent
228  while(1) {
229  Double_t rand = RooRandom::uniform() ;
230  Double_t tval(0) ;
231 
232  switch(_type) {
233  case SingleSided:
234  tval = -_tau*log(rand);
235  break ;
236  case Flipped:
237  tval= +_tau*log(rand);
238  break ;
239  case DoubleSided:
240  tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
241  break ;
242  }
243 
244  // Accept event if T is in generated range
245  Double_t maxDil = 1.0 ;
246 // 2 in next line is conservative and inefficient - allows for delMistag=1!
247  Double_t maxAcceptProb = 2 + fabs(maxDil*_avgS) + fabs(maxDil*_avgC);
248  Double_t acceptProb = (1-_tag*_delMistag + _mu*_tag*(1. - 2.*_avgMistag))
249  + (_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgS*sin(_dm*tval)
250  - (_tag*(1-2*_avgMistag) + _mu*(1. - _tag*_delMistag))*_avgC*cos(_dm*tval);
251 
252  Bool_t accept = maxAcceptProb*RooRandom::uniform() < acceptProb ? kTRUE : kFALSE ;
253 
254  if (tval<_t.max() && tval>_t.min() && accept) {
255  _t = tval ;
256  break ;
257  }
258  }
259 
260 }
261 
RooRealProxy _mu
virtual Int_t getCoefAnalyticalIntegral(Int_t coef, RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Default implementation of function advertising integration capabilities.
ClassImp(RooBCPGenDecay)
RooRealProxy _avgC
#define assert(cond)
Definition: unittest.h:542
RooRealProxy _dm
virtual ~RooBCPGenDecay()
Destructor.
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
STL namespace.
virtual Double_t coefAnalyticalIntegral(Int_t coef, Int_t code, const char *rangeName=0) const
Default implementation of function implementing advertised integrals.
double cos(double)
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
Double_t _genB0Frac
virtual Double_t coefficient(Int_t basisIndex) const
B0 : _tag = +1 B0bar : _tag = -1.
RooRealProxy _tau
RooCategoryProxy _tag
friend class RooArgSet
Definition: RooAbsArg.h:469
double sin(double)
DecayType _type
RooRealProxy _delMistag
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
void initGenerator(Int_t code)
Interface for one-time initialization to setup the generator for the specified code.
static Double_t uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition: RooRandom.cxx:83
double Double_t
Definition: RtypesCore.h:55
RooRealProxy _avgMistag
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
void generateEvent(Int_t code)
Generate mix-state dependent.
int type
Definition: TGX11.cxx:120
Int_t declareBasis(const char *expression, const RooArgList &params)
Declare a basis function for use in this physics model.
#define name(a, b)
Definition: linkTestLib0.cpp:5
friend class RooRealIntegral
Definition: RooAbsPdf.h:294
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
const RooAbsCategory & arg() const
RooRealProxy _t
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
const Bool_t kTRUE
Definition: Rtypes.h:91
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
RooRealProxy _avgS
double log(double)