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
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047 #include "RooFit.h"
00048
00049 #include "Riostream.h"
00050 #include "Riostream.h"
00051 #include "RooRealVar.h"
00052 #include "RooRandom.h"
00053 #include "RooNonCPEigenDecay.h"
00054 #include "TMath.h"
00055 #include "RooRealIntegral.h"
00056
00057 ClassImp(RooNonCPEigenDecay);
00058
00059 #define Debug_RooNonCPEigenDecay 1
00060
00061
00062
00063 RooNonCPEigenDecay::RooNonCPEigenDecay( const char *name, const char *title,
00064 RooRealVar& t,
00065 RooAbsCategory& tag,
00066 RooAbsReal& tau,
00067 RooAbsReal& dm,
00068 RooAbsReal& avgW,
00069 RooAbsReal& delW,
00070 RooAbsCategory& rhoQ,
00071 RooAbsReal& correctQ,
00072 RooAbsReal& wQ,
00073 RooAbsReal& acp,
00074 RooAbsReal& C,
00075 RooAbsReal& delC,
00076 RooAbsReal& S,
00077 RooAbsReal& delS,
00078 const RooResolutionModel& model,
00079 DecayType type )
00080 : RooAbsAnaConvPdf( name, title, model, t ),
00081 _acp ( "acp", "acp", this, acp ),
00082 _avgC ( "C", "C", this, C ),
00083 _delC ( "delC", "delC", this, delC ),
00084 _avgS ( "S", "S", this, S ),
00085 _delS ( "delS", "delS", this, delS ),
00086 _avgW ( "avgW", "Average mistag rate",this, avgW ),
00087 _delW ( "delW", "Shift mistag rate", this, delW ),
00088 _t ( "t", "time", this, t ),
00089 _tau ( "tau", "decay time", this, tau ),
00090 _dm ( "dm", "mixing frequency", this, dm ),
00091 _tag ( "tag", "CP state", this, tag ),
00092 _rhoQ ( "rhoQ", "Charge of the rho", this, rhoQ ),
00093 _correctQ ( "correctQ", "correction of rhoQ", this, correctQ ),
00094 _wQ ( "wQ", "mischarge", this, wQ ),
00095 _genB0Frac ( 0 ),
00096 _genRhoPlusFrac( 0 ),
00097 _type ( type )
00098 {
00099
00100
00101 switch(type) {
00102 case SingleSided:
00103 _basisExp = declareBasis( "exp(-@0/@1)", RooArgList( tau ) );
00104 _basisSin = declareBasis( "exp(-@0/@1)*sin(@0*@2)", RooArgList( tau, dm ) );
00105 _basisCos = declareBasis( "exp(-@0/@1)*cos(@0*@2)", RooArgList( tau, dm ) );
00106 break;
00107 case Flipped:
00108 _basisExp = declareBasis( "exp(@0)/@1)", RooArgList( tau ) );
00109 _basisSin = declareBasis( "exp(@0/@1)*sin(@0*@2)", RooArgList( tau, dm ) );
00110 _basisCos = declareBasis( "exp(@0/@1)*cos(@0*@2)", RooArgList( tau, dm ) );
00111 break;
00112 case DoubleSided:
00113 _basisExp = declareBasis( "exp(-abs(@0)/@1)", RooArgList( tau ) );
00114 _basisSin = declareBasis( "exp(-abs(@0)/@1)*sin(@0*@2)", RooArgList( tau, dm ) );
00115 _basisCos = declareBasis( "exp(-abs(@0)/@1)*cos(@0*@2)", RooArgList( tau, dm ) );
00116 break;
00117 }
00118 }
00119
00120
00121
00122 RooNonCPEigenDecay::RooNonCPEigenDecay( const char *name, const char *title,
00123 RooRealVar& t,
00124 RooAbsCategory& tag,
00125 RooAbsReal& tau,
00126 RooAbsReal& dm,
00127 RooAbsReal& avgW,
00128 RooAbsReal& delW,
00129 RooAbsCategory& rhoQ,
00130 RooAbsReal& correctQ,
00131 RooAbsReal& acp,
00132 RooAbsReal& C,
00133 RooAbsReal& delC,
00134 RooAbsReal& S,
00135 RooAbsReal& delS,
00136 const RooResolutionModel& model,
00137 DecayType type )
00138 : RooAbsAnaConvPdf( name, title, model, t ),
00139 _acp ( "acp", "acp", this, acp ),
00140 _avgC ( "C", "C", this, C ),
00141 _delC ( "delC", "delC", this, delC ),
00142 _avgS ( "S", "S", this, S ),
00143 _delS ( "delS", "delS", this, delS ),
00144 _avgW ( "avgW", "Average mistag rate",this, avgW ),
00145 _delW ( "delW", "Shift mistag rate", this, delW ),
00146 _t ( "t", "time", this, t ),
00147 _tau ( "tau", "decay time", this, tau ),
00148 _dm ( "dm", "mixing frequency", this, dm ),
00149 _tag ( "tag", "CP state", this, tag ),
00150 _rhoQ ( "rhoQ", "Charge of the rho", this, rhoQ ),
00151 _correctQ ( "correctQ", "correction of rhoQ", this, correctQ ),
00152 _genB0Frac ( 0 ),
00153 _genRhoPlusFrac( 0 ),
00154 _type ( type )
00155 {
00156
00157
00158 _wQ = RooRealProxy( "wQ", "mischarge", this, *(new RooRealVar( "wQ", "wQ", 0 )) );
00159
00160 switch(type) {
00161 case SingleSided:
00162 _basisExp = declareBasis( "exp(-@0/@1)", RooArgList( tau ) );
00163 _basisSin = declareBasis( "exp(-@0/@1)*sin(@0*@2)", RooArgList( tau, dm ) );
00164 _basisCos = declareBasis( "exp(-@0/@1)*cos(@0*@2)", RooArgList( tau, dm ) );
00165 break;
00166 case Flipped:
00167 _basisExp = declareBasis( "exp(@0)/@1)", RooArgList( tau ) );
00168 _basisSin = declareBasis( "exp(@0/@1)*sin(@0*@2)", RooArgList( tau, dm ) );
00169 _basisCos = declareBasis( "exp(@0/@1)*cos(@0*@2)", RooArgList( tau, dm ) );
00170 break;
00171 case DoubleSided:
00172 _basisExp = declareBasis( "exp(-abs(@0)/@1)", RooArgList( tau ) );
00173 _basisSin = declareBasis( "exp(-abs(@0)/@1)*sin(@0*@2)", RooArgList( tau, dm ) );
00174 _basisCos = declareBasis( "exp(-abs(@0)/@1)*cos(@0*@2)", RooArgList( tau, dm ) );
00175 break;
00176 }
00177 }
00178
00179
00180
00181 RooNonCPEigenDecay::RooNonCPEigenDecay( const RooNonCPEigenDecay& other, const char* name )
00182 : RooAbsAnaConvPdf( other, name ),
00183 _acp ( "acp", this, other._acp ),
00184 _avgC ( "C", this, other._avgC ),
00185 _delC ( "delC", this, other._delC ),
00186 _avgS ( "S", this, other._avgS ),
00187 _delS ( "delS", this, other._delS ),
00188 _avgW ( "avgW", this, other._avgW ),
00189 _delW ( "delW", this, other._delW ),
00190 _t ( "t", this, other._t ),
00191 _tau ( "tau", this, other._tau ),
00192 _dm ( "dm", this, other._dm ),
00193 _tag ( "tag", this, other._tag ),
00194 _rhoQ ( "rhoQ", this, other._rhoQ ),
00195 _correctQ ( "correctQ", this, other._correctQ ),
00196 _wQ ( "wQ", this, other._wQ ),
00197 _genB0Frac ( other._genB0Frac ),
00198 _genRhoPlusFrac( other._genRhoPlusFrac ),
00199 _type ( other._type ),
00200 _basisExp ( other._basisExp ),
00201 _basisSin ( other._basisSin ),
00202 _basisCos ( other._basisCos )
00203 {
00204
00205 }
00206
00207
00208
00209 RooNonCPEigenDecay::~RooNonCPEigenDecay( void )
00210 {
00211
00212 }
00213
00214
00215
00216 Double_t RooNonCPEigenDecay::coefficient( Int_t basisIndex ) const
00217 {
00218
00219
00220
00221
00222
00223
00224 Int_t rhoQc = _rhoQ * int(_correctQ);
00225 assert( rhoQc == 1 || rhoQc == -1 );
00226
00227 Double_t a_sin_p = _avgS + _delS;
00228 Double_t a_sin_m = _avgS - _delS;
00229 Double_t a_cos_p = _avgC + _delC;
00230 Double_t a_cos_m = _avgC - _delC;
00231
00232 if (basisIndex == _basisExp) {
00233 if (rhoQc == -1 || rhoQc == +1)
00234 return (1 + rhoQc*_acp*(1 - 2*_wQ))*(1 + 0.5*_tag*(2*_delW));
00235 else
00236 return 1;
00237 }
00238
00239 if (basisIndex == _basisSin) {
00240
00241 if (rhoQc == -1)
00242
00243 return - ((1 - _acp)*a_sin_m*(1 - _wQ) + (1 + _acp)*a_sin_p*_wQ)*(1 - 2*_avgW)*_tag;
00244
00245 else if (rhoQc == +1)
00246
00247 return - ((1 + _acp)*a_sin_p*(1 - _wQ) + (1 - _acp)*a_sin_m*_wQ)*(1 - 2*_avgW)*_tag;
00248
00249 else
00250 return - _tag*((a_sin_p + a_sin_m)/2)*(1 - 2*_avgW);
00251 }
00252
00253 if (basisIndex == _basisCos) {
00254
00255 if ( rhoQc == -1)
00256
00257 return + ((1 - _acp)*a_cos_m*(1 - _wQ) + (1 + _acp)*a_cos_p*_wQ)*(1 - 2*_avgW)*_tag;
00258
00259 else if (rhoQc == +1)
00260
00261 return + ((1 + _acp)*a_cos_p*(1 - _wQ) + (1 - _acp)*a_cos_m*_wQ)*(1 - 2*_avgW)*_tag;
00262
00263 else
00264 return _tag*((a_cos_p + a_cos_m)/2)*(1 - 2*_avgW);
00265 }
00266
00267 return 0;
00268 }
00269
00270
00271
00272
00273 Int_t RooNonCPEigenDecay::getCoefAnalyticalIntegral( Int_t , RooArgSet& allVars,
00274 RooArgSet& analVars, const char* rangeName ) const
00275 {
00276 if (rangeName) return 0 ;
00277
00278 if (matchArgs( allVars, analVars, _tag, _rhoQ )) return 3;
00279 if (matchArgs( allVars, analVars, _rhoQ )) return 2;
00280 if (matchArgs( allVars, analVars, _tag )) return 1;
00281
00282 return 0;
00283 }
00284
00285
00286
00287 Double_t RooNonCPEigenDecay::coefAnalyticalIntegral( Int_t basisIndex,
00288 Int_t code, const char* ) const
00289 {
00290
00291 Int_t rhoQc = _rhoQ*int(_correctQ);
00292
00293 Double_t a_sin_p = _avgS + _delS;
00294 Double_t a_sin_m = _avgS - _delS;
00295 Double_t a_cos_p = _avgC + _delC;
00296 Double_t a_cos_m = _avgC - _delC;
00297
00298 switch(code) {
00299
00300
00301 case 0: return coefficient(basisIndex);
00302
00303
00304 case 1:
00305 if (basisIndex == _basisExp) return 2*(1 + rhoQc*_acp*(1 - 2*_wQ));
00306 if (basisIndex == _basisSin || basisIndex==_basisCos) return 0;
00307 assert( kFALSE );
00308
00309
00310 case 2:
00311 if (basisIndex == _basisExp) return 2*(1 + 0.5*_tag*(2.*_delW));
00312
00313 if (basisIndex == _basisSin)
00314
00315 return - ( (1 - _acp)*a_sin_m + (1 + _acp)*a_sin_p )*(1 - 2*_avgW)*_tag;
00316
00317 if (basisIndex == _basisCos)
00318
00319 return + ( (1 - _acp)*a_cos_m + (1 + _acp)*a_cos_p )*(1 - 2*_avgW)*_tag;
00320
00321 assert( kFALSE );
00322
00323
00324 case 3:
00325 if (basisIndex == _basisExp) return 2*2;
00326 if (basisIndex == _basisSin || basisIndex==_basisCos) return 0;
00327 assert( kFALSE );
00328
00329 default:
00330 assert( kFALSE );
00331 }
00332
00333 return 0;
00334 }
00335
00336
00337
00338 Int_t RooNonCPEigenDecay::getGenerator( const RooArgSet& directVars,
00339 RooArgSet& generateVars, Bool_t staticInitOK ) const
00340 {
00341 if (staticInitOK) {
00342 if (matchArgs( directVars, generateVars, _t, _tag, _rhoQ )) return 4;
00343 if (matchArgs( directVars, generateVars, _t, _rhoQ )) return 3;
00344 if (matchArgs( directVars, generateVars, _t, _tag )) return 2;
00345 }
00346 if (matchArgs( directVars, generateVars, _t )) return 1;
00347 return 0;
00348 }
00349
00350
00351
00352 void RooNonCPEigenDecay::initGenerator( Int_t code )
00353 {
00354
00355 if (code == 2 || code == 4) {
00356
00357 Double_t sumInt1 = RooRealIntegral( "sumInt1", "sum integral1", *this,
00358 RooArgSet( _t.arg(), _tag.arg(), _rhoQ.arg() )
00359 ).getVal();
00360 _tag = -1;
00361 Double_t b0Int1 = RooRealIntegral( "mixInt1", "mix integral1", *this,
00362 RooArgSet( _t.arg(), _rhoQ.arg() )
00363 ).getVal();
00364 _genB0Frac = b0Int1/sumInt1;
00365
00366 if (Debug_RooNonCPEigenDecay == 1)
00367 cout << " o RooNonCPEigenDecay::initgenerator: genB0Frac : "
00368 << _genB0Frac
00369 << ", tag dilution: " << (1 - 2*_avgW)
00370 << endl;
00371 }
00372
00373 if (code == 3 || code == 4) {
00374
00375 Double_t sumInt2 = RooRealIntegral( "sumInt2", "sum integral2", *this,
00376 RooArgSet( _t.arg(), _tag.arg(), _rhoQ.arg() )
00377 ).getVal();
00378 _rhoQ = 1;
00379 Double_t b0Int2 = RooRealIntegral( "mixInt2", "mix integral2", *this,
00380 RooArgSet( _t.arg(), _tag.arg() )
00381 ).getVal();
00382 _genRhoPlusFrac = b0Int2/sumInt2;
00383
00384 if (Debug_RooNonCPEigenDecay == 1)
00385 cout << " o RooNonCPEigenDecay::initgenerator: genRhoPlusFrac: "
00386 << _genRhoPlusFrac << endl;
00387 }
00388 }
00389
00390
00391
00392
00393 void RooNonCPEigenDecay::generateEvent( Int_t code )
00394 {
00395
00396
00397 while (kTRUE) {
00398
00399
00400 if (code != 1) {
00401 if (code != 3) _tag = (RooRandom::uniform()<=0.5) ? -1 : +1;
00402 if (code != 2) _rhoQ = (RooRandom::uniform()<=0.5) ? 1 : -1;
00403 }
00404
00405
00406
00407
00408 Double_t a_sin_p = _avgS + _delS;
00409 Double_t a_sin_m = _avgS - _delS;
00410 Double_t a_cos_p = _avgC + _delC;
00411 Double_t a_cos_m = _avgC - _delC;
00412
00413
00414 double a1 = 1 + sqrt(TMath::Power(a_cos_m, 2) + TMath::Power(a_sin_m, 2));
00415 double a2 = 1 + sqrt(TMath::Power(a_cos_p, 2) + TMath::Power(a_sin_p, 2));
00416
00417 Double_t maxAcceptProb = (1.10 + TMath::Abs(_acp)) * (a1 > a2 ? a1 : a2);
00418
00419
00420 Double_t rand = RooRandom::uniform();
00421 Double_t tval(0);
00422
00423 switch(_type) {
00424
00425 case SingleSided:
00426 tval = -_tau*log(rand);
00427 break;
00428
00429 case Flipped:
00430 tval = +_tau*log(rand);
00431 break;
00432
00433 case DoubleSided:
00434 tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5));
00435 break;
00436 }
00437
00438
00439 Double_t expC = coefficient( _basisExp );
00440 Double_t sinC = coefficient( _basisSin );
00441 Double_t cosC = coefficient( _basisCos );
00442
00443
00444 Double_t acceptProb = expC + sinC*sin(_dm*tval) + cosC*cos(_dm*tval);
00445
00446
00447 assert( acceptProb <= maxAcceptProb );
00448
00449
00450 Bool_t accept = maxAcceptProb*RooRandom::uniform() < acceptProb ? kTRUE : kFALSE;
00451
00452 if (accept && tval<_t.max() && tval>_t.min()) {
00453 _t = tval;
00454 break;
00455 }
00456 }
00457 }
00458