BOSS
7.0.1
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtRexc.cc
Go to the documentation of this file.
1
//--------------------------------------------------------------------------
2
// pingrg-2014-10-16
3
// Model REXC : used to decay vgam to the final state as in ConExc model
4
// Name: REXC: R-scan exclusive decay model
5
// Decay cards:
6
// Decay vgam
7
// 1 REXC;
8
// Enddecay
9
//------------------------------------------------------------------------
10
11
#include "
EvtGenBase/EvtPatches.hh
"
12
#include <stdlib.h>
13
#include "
EvtGenBase/EvtParticle.hh
"
14
#include "
EvtGenBase/EvtGenKine.hh
"
15
#include "
EvtGenBase/EvtPDL.hh
"
16
#include "
EvtGenBase/EvtReport.hh
"
17
#include "
EvtGenModels/EvtRexc.hh
"
18
#include "
EvtGenBase/EvtRandom.hh
"
19
#include <string>
20
21
EvtRexc::~EvtRexc
() {}
22
23
void
EvtRexc::getName
(std::string& model_name){
24
25
model_name=
"REXC"
;
//R-scan exclusive decay model
26
27
}
28
29
EvtDecayBase
*
EvtRexc::clone
(){
30
31
return
new
EvtRexc
;
32
33
}
34
35
36
void
EvtRexc::init
(){
37
38
// check that there are 0 arguments
39
checkNArg
(0);
40
myconexc =
new
EvtConExc
();
41
}
42
43
void
EvtRexc::initProbMax
(){
44
45
noProbMax
();
46
47
}
48
49
void
EvtRexc::decay
(
EvtParticle
*p ){
50
double
mhds = p->
mass
();
51
int
mymode =
EvtConExc::conexcmode
;
52
myconexc->
init_mode
(mymode);
53
//std::cout<<"EvtRexc:: A selected mode is "<<mymode<<" with Mhds= "<<mhds<<std::endl; //debugging
54
int
_ndaugs = myconexc->
getNdaugs
();
55
EvtId
*daugs=myconexc->
getDaugId
();
56
//debugging
57
/*
58
std::cout<<"Ndaugs= "<<_ndaugs<<std::endl;
59
for(int i=0;i<_ndaugs;i++){
60
std::cout<<i<<" "<<EvtPDL::getStdHep(daugs[i])<<std::endl;
61
}
62
*/
63
Loop:
64
double
totmass=0;
65
p->
makeDaughters
(_ndaugs,daugs);
66
for
(
int
i=0;i< _ndaugs;i++){
67
EvtParticle
* di=p->
getDaug
(i);
68
double
xmi=
EvtPDL::getMass
(di->
getId
());
69
di->
setMass
(xmi);
70
totmass += xmi;
71
}
72
if
(totmass > p->
mass
())
goto
Loop;
73
74
double
weight
= p->
initializePhaseSpace
( _ndaugs , daugs);
75
// std::cout<<"weight= "<<weight<<std::endl;
76
if
( (2.5< mhds && mhds<3.092 || mhds>3.1012) && !
angularSampling
(p))
goto
Loop;
77
return ;
78
}
79
80
bool
EvtRexc::angularSampling
(
EvtParticle
* par){
81
bool
tagp,tagk;
82
tagk=0;
83
tagp=0;
84
int
nds = par->
getNDaug
();
85
for
(
int
i=0;i<par->
getNDaug
();i++){
86
EvtId
di=par->
getDaug
(i)->
getId
();
87
EvtVector4R
p4i = par->
getDaug
(i)->
getP4Lab
();
88
int
pdgcode =
EvtPDL::getStdHep
( di );
89
double
alpha
=1;
90
double
angmax = 1+
alpha
;
91
double
costheta =
cos
(p4i.
theta
());
92
double
ang=1+
alpha
*costheta*costheta;
93
double
xratio = ang/angmax;
94
double
xi=
EvtRandom::Flat
(0.,1);
95
//std::cout<<"pdgcode "<<pdgcode<<std::endl;
96
//std::cout<<ang<<" "<<angmax<<" "<<xi<<" "<<xratio<<std::endl;
97
if
(xi>xratio)
return
false
;
98
}
//loop over duaghters
99
return
true
;
100
}
101
102
double
EvtRexc::baryonAng
(
double
mx){
103
double
mp
=0.938;
104
double
u = 0.938/mx;
105
u = u*u;
106
double
u2 = (1+u)*(1+u);
107
double
uu = u*(1+6*u);
108
double
alpha
= (u2-uu)/(u2+uu);
109
return
alpha
;
110
}
111
EvtGenKine.hh
EvtPDL.hh
EvtParticle.hh
EvtPatches.hh
EvtRandom.hh
EvtReport.hh
EvtRexc.hh
alpha
const double alpha
Definition:
FastVertexFit.cxx:4
cos
double cos(const BesAngle a)
Definition:
InstallArea/include/MdcGeom/MdcGeom/BesAngle.h:213
weight
*********Class see also m_nmax DOUBLE PRECISION m_MasPhot DOUBLE PRECISION m_phsu DOUBLE PRECISION m_Xenph DOUBLE PRECISION m_r2 DOUBLE PRECISION m_WtMass INTEGER m_nmax INTEGER m_Nevgen INTEGER m_IsFSR INTEGER m_MarTot *COMMON c_KarFin $ !Output file $ !Event serial number $ !alpha QED at Thomson limit $ !minimum energy at CMS for remooval $ !infrared dimensionless $ !dummy photon IR regulator $ !crude photon multiplicity enhancement factor *EVENT $ !MC crude volume of PhhSpace *Sfactors $ !YFS formfactor IR part only $ !YFS formfactor non IR finite part $ !mass weight
Definition:
KarFin.h:34
EvtConExc
Definition:
EvtConExc.hh:87
EvtConExc::getNdaugs
int getNdaugs()
Definition:
EvtConExc.hh:154
EvtConExc::conexcmode
static int conexcmode
Definition:
EvtConExc.hh:157
EvtConExc::init_mode
void init_mode(int mode)
Definition:
EvtConExc.cc:576
EvtConExc::getDaugId
EvtId * getDaugId()
Definition:
EvtConExc.hh:155
EvtDecayBase
Definition:
EvtDecayBase.hh:33
EvtDecayBase::noProbMax
void noProbMax()
Definition:
EvtDecayBase.cc:304
EvtDecayBase::checkNArg
void checkNArg(int a1, int a2=-1, int a3=-1, int a4=-1)
Definition:
EvtDecayBase.cc:482
EvtId
Definition:
EvtId.hh:27
EvtPDL::getStdHep
static int getStdHep(EvtId id)
Definition:
EvtPDL.hh:56
EvtPDL::getMass
static double getMass(EvtId i)
Definition:
EvtPDL.hh:46
EvtParticle
Definition:
EvtParticle.hh:42
EvtParticle::setMass
void setMass(double m)
Definition:
EvtParticle.hh:372
EvtParticle::makeDaughters
void makeDaughters(int ndaug, EvtId *id)
Definition:
EvtParticle.cc:1177
EvtParticle::getP4Lab
EvtVector4R getP4Lab()
Definition:
EvtParticle.cc:685
EvtParticle::getId
EvtId getId() const
Definition:
EvtParticle.cc:113
EvtParticle::getNDaug
int getNDaug() const
Definition:
EvtParticle.cc:125
EvtParticle::getDaug
EvtParticle * getDaug(int i)
Definition:
EvtParticle.cc:85
EvtParticle::mass
double mass() const
Definition:
EvtParticle.cc:127
EvtParticle::initializePhaseSpace
double initializePhaseSpace(int numdaughter, EvtId *daughters, double poleSize=-1., int whichTwo1=0, int whichTwo2=1)
Definition:
EvtParticle.cc:1071
EvtRandom::Flat
static double Flat()
Definition:
EvtRandom.cc:73
EvtRexc::~EvtRexc
virtual ~EvtRexc()
Definition:
EvtRexc.cc:21
EvtRexc::init
void init()
Definition:
EvtRexc.cc:36
EvtRexc::decay
void decay(EvtParticle *p)
Definition:
EvtRexc.cc:49
EvtRexc::clone
EvtDecayBase * clone()
Definition:
EvtRexc.cc:29
EvtRexc::angularSampling
bool angularSampling(EvtParticle *par)
Definition:
EvtRexc.cc:80
EvtRexc::initProbMax
void initProbMax()
Definition:
EvtRexc.cc:43
EvtRexc::EvtRexc
EvtRexc()
Definition:
EvtRexc.hh:23
EvtRexc::baryonAng
double baryonAng(double mx)
Definition:
EvtRexc.cc:102
EvtRexc::getName
void getName(std::string &name)
Definition:
EvtRexc.cc:23
EvtVector4R
Definition:
EvtVector4R.hh:29
EvtVector4R::theta
double theta()
Definition:
EvtVector4R.cc:249
mp
const double mp
Definition:
incllambda.cxx:45
source
Generator
BesEvtGen
BesEvtGen-00-03-36
src
EvtGen
EvtGenModels
EvtRexc.cc
Generated by
1.9.6