00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "RooFit.h"
00027
00028 #include "Riostream.h"
00029 #include "Riostream.h"
00030 #include "RooRealVar.h"
00031 #include "RooRandom.h"
00032 #include "RooBCPEffDecay.h"
00033 #include "RooRealIntegral.h"
00034
00035 ClassImp(RooBCPEffDecay)
00036 ;
00037
00038
00039
00040
00041 RooBCPEffDecay::RooBCPEffDecay(const char *name, const char *title,
00042 RooRealVar& t, RooAbsCategory& tag,
00043 RooAbsReal& tau, RooAbsReal& dm,
00044 RooAbsReal& avgMistag, RooAbsReal& CPeigenval,
00045 RooAbsReal& a, RooAbsReal& b,
00046 RooAbsReal& effRatio, RooAbsReal& delMistag,
00047 const RooResolutionModel& model, DecayType type) :
00048 RooAbsAnaConvPdf(name,title,model,t),
00049 _absLambda("absLambda","Absolute value of lambda",this,a),
00050 _argLambda("argLambda","Arg(Lambda)",this,b),
00051 _effRatio("effRatio","B0/B0bar efficiency ratio",this,effRatio),
00052 _CPeigenval("CPeigenval","CP eigen value",this,CPeigenval),
00053 _avgMistag("avgMistag","Average mistag rate",this,avgMistag),
00054 _delMistag("delMistag","Delta mistag rate",this,delMistag),
00055 _t("t","time",this,t),
00056 _tau("tau","decay time",this,tau),
00057 _dm("dm","mixing frequency",this,dm),
00058 _tag("tag","CP state",this,tag),
00059 _genB0Frac(0),
00060 _type(type)
00061 {
00062
00063 switch(type) {
00064 case SingleSided:
00065 _basisExp = declareBasis("exp(-@0/@1)",RooArgList(tau,dm)) ;
00066 _basisSin = declareBasis("exp(-@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
00067 _basisCos = declareBasis("exp(-@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
00068 break ;
00069 case Flipped:
00070 _basisExp = declareBasis("exp(@0)/@1)",RooArgList(tau,dm)) ;
00071 _basisSin = declareBasis("exp(@0/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
00072 _basisCos = declareBasis("exp(@0/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
00073 break ;
00074 case DoubleSided:
00075 _basisExp = declareBasis("exp(-abs(@0)/@1)",RooArgList(tau,dm)) ;
00076 _basisSin = declareBasis("exp(-abs(@0)/@1)*sin(@0*@2)",RooArgList(tau,dm)) ;
00077 _basisCos = declareBasis("exp(-abs(@0)/@1)*cos(@0*@2)",RooArgList(tau,dm)) ;
00078 break ;
00079 }
00080 }
00081
00082
00083
00084
00085 RooBCPEffDecay::RooBCPEffDecay(const RooBCPEffDecay& other, const char* name) :
00086 RooAbsAnaConvPdf(other,name),
00087 _absLambda("absLambda",this,other._absLambda),
00088 _argLambda("argLambda",this,other._argLambda),
00089 _effRatio("effRatio",this,other._effRatio),
00090 _CPeigenval("CPeigenval",this,other._CPeigenval),
00091 _avgMistag("avgMistag",this,other._avgMistag),
00092 _delMistag("delMistag",this,other._delMistag),
00093 _t("t",this,other._t),
00094 _tau("tau",this,other._tau),
00095 _dm("dm",this,other._dm),
00096 _tag("tag",this,other._tag),
00097 _genB0Frac(other._genB0Frac),
00098 _type(other._type),
00099 _basisExp(other._basisExp),
00100 _basisSin(other._basisSin),
00101 _basisCos(other._basisCos)
00102 {
00103
00104 }
00105
00106
00107
00108
00109 RooBCPEffDecay::~RooBCPEffDecay()
00110 {
00111
00112 }
00113
00114
00115
00116
00117 Double_t RooBCPEffDecay::coefficient(Int_t basisIndex) const
00118 {
00119
00120
00121
00122 if (basisIndex==_basisExp) {
00123
00124 return (1 - _tag*_delMistag)*(1+_absLambda*_absLambda)/2 ;
00125
00126 }
00127
00128 if (basisIndex==_basisSin) {
00129
00130 return -1*_tag*(1-2*_avgMistag)*_CPeigenval*_absLambda*_argLambda ;
00131
00132 }
00133
00134 if (basisIndex==_basisCos) {
00135
00136 return -1*_tag*(1-2*_avgMistag)*(1-_absLambda*_absLambda)/2 ;
00137
00138 }
00139
00140 return 0 ;
00141 }
00142
00143
00144
00145
00146 Int_t RooBCPEffDecay::getCoefAnalyticalIntegral(Int_t , RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
00147 {
00148 if (rangeName) return 0 ;
00149
00150 if (matchArgs(allVars,analVars,_tag)) return 1 ;
00151 return 0 ;
00152 }
00153
00154
00155
00156
00157 Double_t RooBCPEffDecay::coefAnalyticalIntegral(Int_t basisIndex, Int_t code, const char* ) const
00158 {
00159 switch(code) {
00160
00161 case 0: return coefficient(basisIndex) ;
00162
00163
00164 case 1:
00165 if (basisIndex==_basisExp) {
00166 return (1+_absLambda*_absLambda) ;
00167 }
00168
00169 if (basisIndex==_basisSin || basisIndex==_basisCos) {
00170 return 0 ;
00171 }
00172 break ;
00173
00174 default:
00175 assert(0) ;
00176 }
00177
00178 return 0 ;
00179 }
00180
00181
00182
00183
00184 Int_t RooBCPEffDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK) const
00185 {
00186 if (staticInitOK) {
00187 if (matchArgs(directVars,generateVars,_t,_tag)) return 2 ;
00188 }
00189 if (matchArgs(directVars,generateVars,_t)) return 1 ;
00190 return 0 ;
00191 }
00192
00193
00194
00195
00196 void RooBCPEffDecay::initGenerator(Int_t code)
00197 {
00198 if (code==2) {
00199
00200 Double_t sumInt = RooRealIntegral("sumInt","sum integral",*this,RooArgSet(_t.arg(),_tag.arg())).getVal() ;
00201 _tag = 1 ;
00202 Double_t b0Int = RooRealIntegral("mixInt","mix integral",*this,RooArgSet(_t.arg())).getVal() ;
00203 _genB0Frac = b0Int/sumInt ;
00204 }
00205 }
00206
00207
00208
00209
00210 void RooBCPEffDecay::generateEvent(Int_t code)
00211 {
00212
00213 if (code==2) {
00214 Double_t rand = RooRandom::uniform() ;
00215 _tag = (rand<=_genB0Frac) ? 1 : -1 ;
00216 }
00217
00218
00219 while(1) {
00220 Double_t rand = RooRandom::uniform() ;
00221 Double_t tval(0) ;
00222
00223 switch(_type) {
00224 case SingleSided:
00225 tval = -_tau*log(rand);
00226 break ;
00227 case Flipped:
00228 tval= +_tau*log(rand);
00229 break ;
00230 case DoubleSided:
00231 tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
00232 break ;
00233 }
00234
00235
00236 Double_t maxDil = 1.0 ;
00237 Double_t al2 = _absLambda*_absLambda ;
00238 Double_t maxAcceptProb = (1+al2) + fabs(maxDil*_CPeigenval*_absLambda*_argLambda) + fabs(maxDil*(1-al2)/2);
00239 Double_t acceptProb = (1+al2)/2*(1-_tag*_delMistag)
00240 - (_tag*(1-2*_avgMistag))*(_CPeigenval*_absLambda*_argLambda)*sin(_dm*tval)
00241 - (_tag*(1-2*_avgMistag))*(1-al2)/2*cos(_dm*tval);
00242
00243 Bool_t accept = maxAcceptProb*RooRandom::uniform() < acceptProb ? kTRUE : kFALSE ;
00244
00245 if (tval<_t.max() && tval>_t.min() && accept) {
00246 _t = tval ;
00247 break ;
00248 }
00249 }
00250
00251 }
00252