FairRoot/PandaRoot
Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
chigen::PythiaChiGen Class Reference

#include <PythiaChiGen.h>

Public Member Functions

 PythiaChiGen (PartonicModel &partonicModel, int maxTries=30, bool doHadronization=true)
 
PartonicModelgetPartonicModel () const
 
virtual ~PythiaChiGen ()
 
Event * next ()
 
 ClassDef (PythiaChiGen, 1)
 

Private Member Functions

void calculateKinematics ()
 
void calculateHadronRemnants (double x, double hadronPz, Vec4 &u, Vec4 &d)
 
void decay (Event &event)
 
void fillColorEvent (Event &event)
 

Private Attributes

double eCM
 
PartonicModelpartonicModel
 
int maxTries
 
double protonPz
 
bool doHadronization
 
Vec4 uMomentum
 
Vec4 dMomentum
 
Vec4 ubarMomentum
 
Vec4 dbarMomentum
 
Vec4 mesonMomentum
 
Vec4 gluonMomentum
 
Vec4 boost
 
EvtStdHep evtstdhep
 

Static Private Attributes

static int COLORS [] = {101, 102, 103}
 

Detailed Description

Definition at line 19 of file PythiaChiGen.h.

Constructor & Destructor Documentation

chigen::PythiaChiGen::PythiaChiGen ( PartonicModel partonicModel,
int  maxTries = 30,
bool  doHadronization = true 
)
Parameters
partonicModelpartonic model
maxTriesmax errors number
doHadronizationperform hadronization

Definition at line 16 of file PythiaChiGen.cxx.

References doHadronization, chigen::PartonicModel::eCM, eCM, chigen::ensure_chigen_is_initialized(), maxTries, protonPz, and sqrt().

18 : partonicModel(partonicModelRef) {
20  eCM = partonicModelRef.eCM;
21  protonPz = sqrt(pow2(eCM) / 4 - 1);
22  maxTries = tries;
23  doHadronization = hadronize;
24 }
void ensure_chigen_is_initialized()
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
PartonicModel & partonicModel
Definition: PythiaChiGen.h:52
chigen::PythiaChiGen::~PythiaChiGen ( )
virtual

Definition at line 26 of file PythiaChiGen.cxx.

26  {
27 }

Member Function Documentation

void chigen::PythiaChiGen::calculateHadronRemnants ( double  x,
double  hadronPz,
Vec4 &  u,
Vec4 &  d 
)
private

Calculates remnants kinematics

Parameters
xx-fraction
hadronPzhadron z-momentum
uu (anti) quark
dd (anti) quark

Definition at line 63 of file PythiaChiGen.cxx.

References cos(), chigen::ChiGenRandomEngine::flat(), phi, chigen::random::random_engine, sin(), and theta.

63  {
64 
65  Vec4 totalRemnant(0, 0, (1 - x) * hadronPz, eCM / 2 - x * abs(hadronPz));
66  double halfRemnantInvMass = totalRemnant.mCalc() / 2;
67  //make remnants pT
68  double theta = chigen::random::random_engine->flat() * M_PI;
69  double sine = halfRemnantInvMass * sin(theta), cosine = halfRemnantInvMass * cos(theta);
70  u.p(0, sine, cosine, halfRemnantInvMass);
71  d.p(0, -sine, -cosine, halfRemnantInvMass);
72 
73  //rotate vectors in XY plane
74  double phi = chigen::random::random_engine->flat()* 2 * M_PI;
75  u.rot(0, phi);
76  d.rot(0, phi);
77 
78  //boost to lab frame
79  u.bst(totalRemnant);
80  d.bst(totalRemnant);
81 }
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
TObjArray * d
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
Double_t x
ChiGenRandomEngine * random_engine
void chigen::PythiaChiGen::calculateKinematics ( )
private

Calculates hadron-level kinematics

Definition at line 31 of file PythiaChiGen.cxx.

References chigen::ChiGenRandomEngine::flat(), chigen::random::random_engine, and sqrt().

31  {
32 
34  double pL = (partonicModel.tH - partonicModel.uH) / (2 * sqrt(partonicModel.sH));
35  double mesonE = sqrt(partonicModel.mesonMass2 + pow2(pT) + pow2(pL));
36  double gluonE = sqrt(pow2(pT) + pow2(pL));
37 
38 
39  //quarkonia
40  mesonMomentum = Vec4(0., pT, pL, mesonE);
41  //gluon
42  gluonMomentum = Vec4(0., -pT, -pL, gluonE);
43 
44  //rotate vectors in XY plane
45  double hardPhiAngle = chigen::random::random_engine->flat()* 2 * M_PI;
46  mesonMomentum.rot(0, hardPhiAngle);
47  gluonMomentum.rot(0, hardPhiAngle);
48 
49  //boost from hard c.m. to hadronic c.m.
50  boost = Vec4(0., 0., (partonicModel.x1 - partonicModel.x2) * protonPz,
52  //do not boost meson momentum until EvtGen perform decay!!!
53  //mesonMomentum.bst(boost) <-- this will destroy polarization
54  gluonMomentum.bst(boost);
55 
56  //remnants kinematics
57  //proton remnants
59  //antiproton remnants
61 }
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
void calculateHadronRemnants(double x, double hadronPz, Vec4 &u, Vec4 &d)
ChiGenRandomEngine * random_engine
PartonicModel & partonicModel
Definition: PythiaChiGen.h:52
chigen::PythiaChiGen::ClassDef ( PythiaChiGen  ,
 
)
void chigen::PythiaChiGen::decay ( Event &  event)
private

Perform decay of polarized chi-meson.

Parameters
eventevent reference

Definition at line 83 of file PythiaChiGen.cxx.

References chigen::evtgen::evt_gen, chigen::isCharmonia(), and p.

83  {
84  //lets decay polarized meson via EvtGen
85  EvtVector4R chiMomentum(
86  mesonMomentum.e(),
87  mesonMomentum.px(),
88  mesonMomentum.py(),
89  mesonMomentum.pz());
90 
91 
92  EvtParticle *chiMeson =
93  EvtParticleFactory::particleFactory(partonicModel.evtId, chiMomentum);
94  chiMeson->setSpinDensityForward(partonicModel.spinDensity);
95 
96  //decay
97  chigen::evtgen::evt_gen->generateDecay(chiMeson);
98  //write info
99  evtstdhep.init();
100  chiMeson->makeStdHep(evtstdhep);
101 
102  assert(evtstdhep.getNPart() == 5);
103  int id;
104  EvtVector4R vertexPoint, evtMomentum;
105  for (int particleNumber = 0, size = evtstdhep.getNPart();
106  particleNumber < size; ++particleNumber) {
107  vertexPoint = evtstdhep.getX4(particleNumber);
108  evtMomentum = evtstdhep.getP4(particleNumber);
109 
110  Pythia8::Particle p;
111  p.id(id = evtstdhep.getStdHepID(particleNumber));
112  p.px(evtMomentum.get(1));
113  p.py(evtMomentum.get(2));
114  p.pz(evtMomentum.get(3));
115  p.e(evtMomentum.get(0));
116  //boost from hard c.m. frame to hadronic c.m. frame
117  p.bst(boost);
118 
119  p.daughter1(evtstdhep.getFirstDaughter(particleNumber) + 1);
120  p.daughter2(evtstdhep.getLastDaughter(particleNumber) + 1);
121  p.mother1(evtstdhep.getFirstMother(particleNumber) + 1);
122  p.mother2(evtstdhep.getLastMother(particleNumber) + 1);
123  if (chigen::isCharmonia(id))
124  p.status(-33);
125  else
126  p.status(91);
127 
128  event.append(p);
129  }
130 }
bool isCharmonia(int pdgCode)
Double_t p
Definition: anasim.C:58
EvtSpinDensity spinDensity
Definition: PartonicModel.h:52
PartonicModel & partonicModel
Definition: PythiaChiGen.h:52
void chigen::PythiaChiGen::fillColorEvent ( Event &  event)
private

Fills the event by hadronic remnants & gluon

Parameters
eventevent reference

Definition at line 132 of file PythiaChiGen.cxx.

132  {
133  // first setup colors and anti-colors
134  // the following color mask is adopted:
135  //
136  // u (bar) d (bar) g
137  // color : 101 102 103
138  // anti-color: COLORS[0] COLORS[1] COLORS[2]
139 
140  // randomly shuffle colors
141  std::random_shuffle(COLORS, COLORS + 3); //todo use same seed as in rndm
142 
143  //already added chi-meson in decay()
144  // //filling the event record:
145  // //quarkonia
146  // event.append(partonicModel->pdgId/*id*/, 21 /*status code*/,
147  // 0, 0/*colors*/,
148  // mesonMomentum /*4-momentum*/, partonicModel->mesonMass /*mass*/);
149 
150  //gluon
151  event.append(21/*id*/, 23 /*status code (should be hadronized)*/,
152  103, COLORS[2]/*colors*/,
153  gluonMomentum /*4-momentum*/, 0.0 /*mass*/);
154 
155 
156  //u-quark
157  event.append(2/*id*/, 23 /*status code (should be hadronized)*/,
158  101, 0/*colors*/,
159  uMomentum/*4-momentum*/, 0.0/* mass*/);
160 
161  //d-quark
162  event.append(1/*id*/, 23 /*status code (should be hadronized)*/,
163  102, 0/*colors*/,
164  dMomentum/*4-momentum*/, 0.0/* mass*/);
165 
166  //u-bar
167  event.append(-2/*id*/, 23 /*status code (should be hadronized)*/,
168  0, COLORS[0]/*colors*/,
169  ubarMomentum/*4-momentum*/, 0.0/* mass*/);
170 
171  //d-bar
172  event.append(-1/*id*/, 23 /*status code (should be hadronized)*/,
173  0, COLORS[1]/*colors*/,
174  dbarMomentum/*4-momentum*/, 0.0 /*mass*/);
175 }
static int COLORS[]
Definition: PythiaChiGen.h:92
PartonicModel& chigen::PythiaChiGen::getPartonicModel ( ) const
inline

Definition at line 30 of file PythiaChiGen.h.

30  {
31  return partonicModel;
32  }
PartonicModel & partonicModel
Definition: PythiaChiGen.h:52
Event * chigen::PythiaChiGen::next ( )

Calculates and returns the next event

Returns
event

Definition at line 177 of file PythiaChiGen.cxx.

References next, and chigen::pythia::pythia.

177  {
178  Event& event = chigen::pythia::pythia->event;
179  //reset event is necessary
180  event.reset();
181 
182  //number of attempts
183  int attempt = 0;
184  while (!partonicModel.next() && attempt++ < maxTries);
185  //parton level failed
186  if (attempt >= maxTries)
187  Fatal("PythiaHadronizationModel::next", "Can't set kinematics at parton level.");
188 
189  attempt = 0;
191  decay(event);
192  fillColorEvent(event);
193 
194  if (doHadronization) {
195  while (!chigen::pythia::pythia->next() && attempt++ < maxTries) {
196  event.reset();
198  decay(event);
199  fillColorEvent(event);
200  }
201  //hadron level failed
202  if (attempt >= maxTries)
203  Fatal("PythiaHadronizationModel::next", "Failed to hadronize remnants.");
204  }
205  return &event;
206 }
void fillColorEvent(Event &event)
void decay(Event &event)
Pythia8::Pythia * pythia
PartonicModel & partonicModel
Definition: PythiaChiGen.h:52

Member Data Documentation

Vec4 chigen::PythiaChiGen::boost
private

Boost from hard c.m. frame to hadronic c.m. frame

Definition at line 80 of file PythiaChiGen.h.

int chigen::PythiaChiGen::COLORS = {101, 102, 103}
staticprivate

Colors of particles:

u-quark : 101 d-quark : 102 gluon (color) : 103

Definition at line 92 of file PythiaChiGen.h.

Vec4 chigen::PythiaChiGen::dbarMomentum
private

Definition at line 72 of file PythiaChiGen.h.

Vec4 chigen::PythiaChiGen::dMomentum
private

Definition at line 68 of file PythiaChiGen.h.

bool chigen::PythiaChiGen::doHadronization
private

Do hadronization or not

Definition at line 64 of file PythiaChiGen.h.

Referenced by PythiaChiGen().

double chigen::PythiaChiGen::eCM
private

Hadronic c.m. energy

Definition at line 48 of file PythiaChiGen.h.

Referenced by PythiaChiGen().

EvtStdHep chigen::PythiaChiGen::evtstdhep
private

EvtGen data on chi_c decay

Definition at line 84 of file PythiaChiGen.h.

Vec4 chigen::PythiaChiGen::gluonMomentum
private

Definition at line 76 of file PythiaChiGen.h.

int chigen::PythiaChiGen::maxTries
private

Maximum number of tries to hadronize with diff kinematics

Definition at line 56 of file PythiaChiGen.h.

Referenced by PythiaChiGen().

Vec4 chigen::PythiaChiGen::mesonMomentum
private

Hard process particles

Definition at line 76 of file PythiaChiGen.h.

PartonicModel& chigen::PythiaChiGen::partonicModel
private

Underlying parton model

Definition at line 52 of file PythiaChiGen.h.

double chigen::PythiaChiGen::protonPz
private

Proton momentum

Definition at line 60 of file PythiaChiGen.h.

Referenced by PythiaChiGen().

Vec4 chigen::PythiaChiGen::ubarMomentum
private

Antiproton remnants

Definition at line 72 of file PythiaChiGen.h.

Vec4 chigen::PythiaChiGen::uMomentum
private

Proton remnants

Definition at line 68 of file PythiaChiGen.h.


The documentation for this class was generated from the following files: