Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4INCL::KinematicsUtils Namespace Reference

Functions

void transformToLocalEnergyFrame (Nucleus const *const n, Particle *const p)
 
G4double getLocalEnergy (Nucleus const *const n, Particle *const p)
 
ThreeVector makeBoostVector (Particle const *const p1, Particle const *const p2)
 
G4double totalEnergyInCM (Particle const *const p1, Particle const *const p2)
 
G4double squareTotalEnergyInCM (Particle const *const p1, Particle const *const p2)
 
G4double momentumInCM (Particle const *const p1, Particle const *const p2)
 gives the momentum in the CM frame of two particles.
 
G4double momentumInCM (const G4double E, const G4double M1, const G4double M2)
 
G4double momentumInLab (Particle const *const p1, Particle const *const p2)
 gives the momentum in the lab frame of two particles.
 
G4double momentumInLab (const G4double s, const G4double m1, const G4double m2)
 
G4double sumTotalEnergies (const ParticleList &)
 
ThreeVector sumMomenta (const ParticleList &)
 
G4double energy (const ThreeVector &p, const G4double m)
 
G4double invariantMass (const G4double E, const ThreeVector &p)
 
G4double squareInvariantMass (const G4double E, const ThreeVector &p)
 
G4double gammaFromKineticEnergy (const ParticleSpecies &p, const G4double EKin)
 

Function Documentation

◆ energy()

G4double G4INCL::KinematicsUtils::energy ( const ThreeVector p,
const G4double  m 
)

Definition at line 158 of file G4INCLKinematicsUtils.cc.

158 {
159 return std::sqrt(p.mag2() + m*m);
160 }
G4double mag2() const

◆ gammaFromKineticEnergy()

G4double G4INCL::KinematicsUtils::gammaFromKineticEnergy ( const ParticleSpecies p,
const G4double  EKin 
)

Definition at line 170 of file G4INCLKinematicsUtils.cc.

170 {
171 G4double mass;
172 if(p.theType==Composite)
173 mass = ParticleTable::getTableMass(p.theA, p.theZ, p.theS);
174 else
175 mass = ParticleTable::getTableParticleMass(p.theType);
176 return (1.+EKin/mass);
177 }
double G4double
Definition: G4Types.hh:83

◆ getLocalEnergy()

G4double G4INCL::KinematicsUtils::getLocalEnergy ( Nucleus const *const  n,
Particle *const  p 
)

Definition at line 53 of file G4INCLKinematicsUtils.cc.

53 {
54// assert(!p->isMeson()); // No local energy for mesons
55
56 G4double vloc = 0.0;
57 const G4double r = p->getPosition().mag();
58 const G4double mass = p->getMass();
59
60 // Local energy is constant outside the surface
61 if(r > n->getUniverseRadius()) {
62 INCL_WARN("Tried to evaluate local energy for a particle outside the maximum radius."
63 << '\n' << p->print() << '\n'
64 << "Maximum radius = " << n->getDensity()->getMaximumRadius() << '\n'
65 << "Universe radius = " << n->getUniverseRadius() << '\n');
66 return 0.0;
67 }
68
69 G4double pfl0 = 0.0;
70 const ParticleType t = p->getType();
71 const G4double kinE = p->getKineticEnergy();
72 if(kinE <= n->getPotential()->getFermiEnergy(t)) {
73 pfl0 = n->getPotential()->getFermiMomentum(p);
74 } else {
75 const G4double tf0 = p->getPotentialEnergy() - n->getPotential()->getSeparationEnergy(p);
76 if(tf0<0.0) return 0.0;
77 pfl0 = std::sqrt(tf0*(tf0 + 2.0*mass));
78 }
79 const G4double pReflection = p->getReflectionMomentum()/pfl0;
80 const G4double reflectionRadius = n->getDensity()->getMaxRFromP(p->getType(), pReflection);
81 const G4double pNominal = p->getMomentum().mag()/pfl0;
82 const G4double nominalReflectionRadius = n->getDensity()->getMaxRFromP(p->getType(), pNominal);
83 const G4double pl = pfl0*n->getDensity()->getMinPFromR(t, r*nominalReflectionRadius/reflectionRadius);
84 vloc = std::sqrt(pl*pl + mass*mass) - mass;
85
86 return vloc;
87 }
#define INCL_WARN(x)
G4double getPotentialEnergy() const
Get the particle potential energy.
const G4INCL::ThreeVector & getPosition() const
G4double getKineticEnergy() const
Get the particle kinetic energy.
G4double getReflectionMomentum() const
Return the reflection momentum.
const G4INCL::ThreeVector & getMomentum() const
G4INCL::ParticleType getType() const
std::string print() const
G4double getMass() const
Get the cached particle mass.
G4double mag() const

Referenced by transformToLocalEnergyFrame().

◆ invariantMass()

G4double G4INCL::KinematicsUtils::invariantMass ( const G4double  E,
const ThreeVector p 
)

Definition at line 162 of file G4INCLKinematicsUtils.cc.

162 {
163 return std::sqrt(squareInvariantMass(E, p));
164 }
G4double squareInvariantMass(const G4double E, const ThreeVector &p)

◆ makeBoostVector()

ThreeVector G4INCL::KinematicsUtils::makeBoostVector ( Particle const *const  p1,
Particle const *const  p2 
)

Definition at line 89 of file G4INCLKinematicsUtils.cc.

89 {
90 const G4double totalEnergy = p1->getEnergy() + p2->getEnergy();
91 return ((p1->getMomentum() + p2->getMomentum())/totalEnergy);
92 }

Referenced by G4INCL::InteractionAvatar::preInteraction(), and squareTotalEnergyInCM().

◆ momentumInCM() [1/2]

G4double G4INCL::KinematicsUtils::momentumInCM ( const G4double  E,
const G4double  M1,
const G4double  M2 
)

Definition at line 119 of file G4INCLKinematicsUtils.cc.

119 {
120 return 0.5*std::sqrt((E*E - std::pow(M1 + M2, 2))
121 *(E*E - std::pow(M1 - M2, 2)))/E;
122 }

◆ momentumInCM() [2/2]

G4double G4INCL::KinematicsUtils::momentumInCM ( Particle const *const  p1,
Particle const *const  p2 
)

gives the momentum in the CM frame of two particles.

The formula is the following:

\[ p_{CM}^2 = \frac{z^2 - m_1^2 m_2^2}{2 z + m_1^2 + m_2^2} \]

where $z$ is the scalar product of the momentum four-vectors:

\[ z = E_1 E_2 - \vec{p}_1\cdot\vec{p}_2 \]

Parameters
p1pointer to particle 1
p2pointer to particle 2
Returns
the absolute value of the momentum of any of the two particles in the CM frame, in MeV/c.

Definition at line 107 of file G4INCLKinematicsUtils.cc.

107 {
108 const G4double m1sq = std::pow(p1->getMass(),2);
109 const G4double m2sq = std::pow(p2->getMass(),2);
110 const G4double z = p1->getEnergy()*p2->getEnergy() - p1->getMomentum().dot(p2->getMomentum());
111 G4double pcm2 = (z*z-m1sq*m2sq)/(2*z+m1sq+m2sq);
112 if(pcm2 < 0.0) {
113 INCL_ERROR("momentumInCM: pcm2 == " << pcm2 << " < 0.0" << '\n');
114 pcm2 = 0.0;
115 }
116 return std::sqrt(pcm2);
117 }
#define INCL_ERROR(x)

Referenced by G4INCL::DeltaDecayChannel::computeDecayTime(), G4INCL::Nucleus::decayOutgoingDeltas(), G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), G4INCL::CrossSectionsMultiPionsAndResonances::etaNToPiN(), G4INCL::DeltaDecayChannel::fillFinalState(), G4INCL::DeltaProductionChannel::fillFinalState(), G4INCL::NKbToLpiChannel::fillFinalState(), G4INCL::NKbToNKbChannel::fillFinalState(), G4INCL::NKbToSpiChannel::fillFinalState(), G4INCL::NKElasticChannel::fillFinalState(), G4INCL::NKToNKChannel::fillFinalState(), G4INCL::NpiToLKChannel::fillFinalState(), G4INCL::NpiToSKChannel::fillFinalState(), G4INCL::NYElasticChannel::fillFinalState(), G4INCL::PionResonanceDecayChannel::fillFinalState(), G4INCL::RecombinationChannel::fillFinalState(), G4INCL::SigmaZeroDecayChannel::fillFinalState(), G4INCL::StrangeAbsorbtionChannel::fillFinalState(), G4INCL::PhaseSpaceKopylov::generate(), G4INCL::NKbElasticChannel::KaonMomentum(), and G4INCL::CrossSectionsMultiPionsAndResonances::omegaNToPiN().

◆ momentumInLab() [1/2]

G4double G4INCL::KinematicsUtils::momentumInLab ( const G4double  s,
const G4double  m1,
const G4double  m2 
)

Definition at line 124 of file G4INCLKinematicsUtils.cc.

124 {
125 const G4double m1sq = m1*m1;
126 const G4double m2sq = m2*m2;
127 G4double plab2 = (s*s-2*s*(m1sq+m2sq)+(m1sq-m2sq)*(m1sq-m2sq))/(4*m2sq);
128 if(plab2 < 0.0) {
129 INCL_ERROR("momentumInLab: plab2 == " << plab2 << " < 0.0; m1sq == " << m1sq << "; m2sq == " << m2sq << "; s == " << s << '\n');
130 plab2 = 0.0;
131 }
132 return std::sqrt(plab2);
133 }

◆ momentumInLab() [2/2]

G4double G4INCL::KinematicsUtils::momentumInLab ( Particle const *const  p1,
Particle const *const  p2 
)

gives the momentum in the lab frame of two particles.

Assumes particle 1 carries all the momentum and particle 2 is at rest.

The formula is the following:

\[ p_{lab}^2 = \frac{s^2 - 2 s (m_1^2 + m_2^2) + {(m_1^2 - m_2^2)}^2}{4 m_2^2} \]

Parameters
p1pointer to particle 1
p2pointer to particle 2
Returns
the absolute value of the momentum of particle 1 in the lab frame, in MeV/c

Definition at line 135 of file G4INCLKinematicsUtils.cc.

135 {
136 const G4double m1 = p1->getMass();
137 const G4double m2 = p2->getMass();
138 const G4double s = squareTotalEnergyInCM(p1, p2);
139 return momentumInLab(s, m1, m2);
140 }
G4double squareTotalEnergyInCM(Particle const *const p1, Particle const *const p2)
G4double momentumInLab(Particle const *const p1, Particle const *const p2)
gives the momentum in the lab frame of two particles.

Referenced by G4INCL::CrossSectionsINCL46::elasticNNLegacy(), G4INCL::CrossSectionsMultiPionsAndResonances::etaNElastic(), G4INCL::CrossSectionsMultiPionsAndResonances::etaNToPiN(), G4INCL::CrossSectionsMultiPionsAndResonances::etaNToPiPiN(), G4INCL::DeltaProductionChannel::fillFinalState(), G4INCL::ElasticChannel::fillFinalState(), G4INCL::EtaNElasticChannel::fillFinalState(), G4INCL::EtaNToPiNChannel::fillFinalState(), G4INCL::NNToMissingStrangenessChannel::fillFinalState(), G4INCL::NpiToMissingStrangenessChannel::fillFinalState(), G4INCL::NKbElasticChannel::KaonMomentum(), G4INCL::NKbToLpiChannel::KaonMomentum(), G4INCL::NKbToNKbChannel::KaonMomentum(), G4INCL::NKbToSpiChannel::KaonMomentum(), G4INCL::NpiToLKChannel::KaonMomentum(), G4INCL::NpiToSKChannel::KaonMomentum(), momentumInLab(), G4INCL::CrossSectionsINCL46::NDeltaToNN(), G4INCL::CrossSectionsStrangeness::NKbelastic(), G4INCL::CrossSectionsStrangeness::NKbToNKb(), G4INCL::CrossSectionsStrangeness::NKbToNKb2pi(), G4INCL::CrossSectionsStrangeness::NKbToNKbpi(), G4INCL::CrossSectionsStrangeness::NKbToS2pi(), G4INCL::CrossSectionsStrangeness::NKbToSpi(), G4INCL::CrossSectionsStrangeness::NKelastic(), G4INCL::CrossSectionsStrangeness::NKToNK(), G4INCL::CrossSectionsStrangeness::NKToNK2pi(), G4INCL::CrossSectionsStrangeness::NKToNKpi(), G4INCL::CrossSectionsStrangeness::NLToNS(), G4INCL::CrossSectionsMultiPions::NNElastic(), G4INCL::CrossSectionsMultiPions::NNElasticFixed(), G4INCL::CrossSectionsMultiPions::NNOnePiOrDelta(), G4INCL::CrossSectionsMultiPions::NNThreePi(), G4INCL::CrossSectionsStrangeness::NNToMissingStrangeness(), G4INCL::CrossSectionsINCL46::NNToNDelta(), G4INCL::CrossSectionsStrangeness::NNToNLK(), G4INCL::CrossSectionsStrangeness::NNToNSK(), G4INCL::CrossSectionsMultiPions::NNTotFixed(), G4INCL::CrossSectionsMultiPions::NNTwoPi(), G4INCL::CrossSectionsStrangeness::NpiToLK2pi(), G4INCL::CrossSectionsStrangeness::NpiToLKpi(), G4INCL::CrossSectionsStrangeness::NpiToMissingStrangeness(), G4INCL::CrossSectionsStrangeness::NpiToNKKb(), G4INCL::CrossSectionsStrangeness::NpiToSK2pi(), G4INCL::CrossSectionsStrangeness::NpiToSKpi(), G4INCL::CrossSectionsStrangeness::NSToNL(), G4INCL::CrossSectionsStrangeness::NSToNS(), G4INCL::CrossSectionsStrangeness::NYelastic(), G4INCL::CrossSectionsMultiPionsAndResonances::omegaNElastic(), G4INCL::CrossSectionsMultiPionsAndResonances::omegaNInelastic(), G4INCL::CrossSectionsMultiPionsAndResonances::omegaNToPiN(), G4INCL::CrossSectionsStrangeness::p_kmToL_pp_pm(), G4INCL::CrossSectionsStrangeness::p_kmToL_pz(), G4INCL::CrossSectionsStrangeness::p_pimToLK0(), G4INCL::CrossSectionsStrangeness::p_pimToSmKp(), G4INCL::CrossSectionsStrangeness::p_pimToSzKz(), G4INCL::CrossSectionsStrangeness::p_pipToSpKp(), G4INCL::CrossSectionsStrangeness::p_pizToSzKp(), G4INCL::CrossSectionsMultiPions::piMinuspIne(), G4INCL::CrossSectionsMultiPions::piMinuspOnePi(), G4INCL::CrossSectionsMultiPionsAndResonances::piMinuspToEtaN(), G4INCL::CrossSectionsMultiPionsAndResonances::piMinuspToOmegaN(), G4INCL::CrossSectionsMultiPions::piMinuspTwoPi(), G4INCL::CrossSectionsMultiPions::piNIne(), G4INCL::CrossSectionsMultiPions::piNOnePi(), G4INCL::CrossSectionsMultiPions::piNToxPiN(), G4INCL::CrossSectionsMultiPions::piNTwoPi(), G4INCL::CrossSectionsMultiPions::piPluspIne(), G4INCL::CrossSectionsMultiPions::piPluspOnePi(), and G4INCL::CrossSectionsMultiPions::piPluspTwoPi().

◆ squareInvariantMass()

G4double G4INCL::KinematicsUtils::squareInvariantMass ( const G4double  E,
const ThreeVector p 
)

Definition at line 166 of file G4INCLKinematicsUtils.cc.

166 {
167 return E*E - p.mag2();
168 }

Referenced by invariantMass().

◆ squareTotalEnergyInCM()

◆ sumMomenta()

ThreeVector G4INCL::KinematicsUtils::sumMomenta ( const ParticleList pl)

Definition at line 150 of file G4INCLKinematicsUtils.cc.

150 {
151 ThreeVector p(0.0, 0.0, 0.0);
152 for(ParticleIter i=pl.begin(), e=pl.end(); i!=e; ++i) {
153 p += (*i)->getMomentum();
154 }
155 return p;
156 }
ParticleList::const_iterator ParticleIter

◆ sumTotalEnergies()

G4double G4INCL::KinematicsUtils::sumTotalEnergies ( const ParticleList pl)

Definition at line 142 of file G4INCLKinematicsUtils.cc.

142 {
143 G4double E = 0.0;
144 for(ParticleIter i=pl.begin(), e=pl.end(); i!=e; ++i) {
145 E += (*i)->getEnergy();
146 }
147 return E;
148 }

◆ totalEnergyInCM()

G4double G4INCL::KinematicsUtils::totalEnergyInCM ( Particle const *const  p1,
Particle const *const  p2 
)

Definition at line 94 of file G4INCLKinematicsUtils.cc.

94 {
95 return std::sqrt(squareTotalEnergyInCM(p1,p2));
96 }

Referenced by G4INCL::CrossSectionsMultiPionsAndResonances::etaNToPiN(), G4INCL::DeltaProductionChannel::fillFinalState(), G4INCL::EtaNToPiPiNChannel::fillFinalState(), G4INCL::NDeltaEtaProductionChannel::fillFinalState(), G4INCL::NDeltaOmegaProductionChannel::fillFinalState(), G4INCL::NDeltaToDeltaLKChannel::fillFinalState(), G4INCL::NDeltaToDeltaSKChannel::fillFinalState(), G4INCL::NDeltaToNLKChannel::fillFinalState(), G4INCL::NDeltaToNNKKbChannel::fillFinalState(), G4INCL::NDeltaToNSKChannel::fillFinalState(), G4INCL::NKbToL2piChannel::fillFinalState(), G4INCL::NKbToNKb2piChannel::fillFinalState(), G4INCL::NKbToNKbpiChannel::fillFinalState(), G4INCL::NKbToS2piChannel::fillFinalState(), G4INCL::NKToNK2piChannel::fillFinalState(), G4INCL::NKToNKpiChannel::fillFinalState(), G4INCL::NLToNSChannel::fillFinalState(), G4INCL::NNEtaToMultiPionsChannel::fillFinalState(), G4INCL::NNOmegaToMultiPionsChannel::fillFinalState(), G4INCL::NNToMissingStrangenessChannel::fillFinalState(), G4INCL::NNToMultiPionsChannel::fillFinalState(), G4INCL::NNToNLK2piChannel::fillFinalState(), G4INCL::NNToNLKChannel::fillFinalState(), G4INCL::NNToNLKpiChannel::fillFinalState(), G4INCL::NNToNNEtaChannel::fillFinalState(), G4INCL::NNToNNKKbChannel::fillFinalState(), G4INCL::NNToNNOmegaChannel::fillFinalState(), G4INCL::NNToNSK2piChannel::fillFinalState(), G4INCL::NNToNSKChannel::fillFinalState(), G4INCL::NNToNSKpiChannel::fillFinalState(), G4INCL::NpiToLK2piChannel::fillFinalState(), G4INCL::NpiToLKpiChannel::fillFinalState(), G4INCL::NpiToMissingStrangenessChannel::fillFinalState(), G4INCL::NpiToNKKbChannel::fillFinalState(), G4INCL::NpiToSK2piChannel::fillFinalState(), G4INCL::NpiToSKpiChannel::fillFinalState(), G4INCL::NSToNLChannel::fillFinalState(), G4INCL::NSToNSChannel::fillFinalState(), G4INCL::OmegaNToPiPiNChannel::fillFinalState(), G4INCL::PiNToEtaChannel::fillFinalState(), G4INCL::PiNToMultiPionsChannel::fillFinalState(), G4INCL::RecombinationChannel::fillFinalState(), G4INCL::StrangeAbsorbtionChannel::fillFinalState(), G4INCL::CrossSectionsStrangeness::NDeltaToNNKKb(), G4INCL::CrossSectionsMultiPions::NNOnePi(), G4INCL::CrossSectionsMultiPions::NNOnePiOrDelta(), G4INCL::CrossSectionsMultiPions::NNThreePi(), G4INCL::CrossSectionsINCL46::NNToNDelta(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNDeltaEta(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNDeltaOmega(), G4INCL::CrossSectionsStrangeness::NNToNLK2pi(), G4INCL::CrossSectionsStrangeness::NNToNLKpi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEta(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaExclu(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaFourPi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaOnePi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaOnePiOrDelta(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaThreePi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaTwoPi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaxPi(), G4INCL::CrossSectionsStrangeness::NNToNNKKb(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmega(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaExclu(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaFourPi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaOnePi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaOnePiOrDelta(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaThreePi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaTwoPi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaxPi(), G4INCL::CrossSectionsStrangeness::NNToNSK2pi(), G4INCL::CrossSectionsStrangeness::NNToNSKpi(), G4INCL::CrossSectionsMultiPions::NNTwoPi(), G4INCL::CrossSectionsMultiPionsAndResonances::omegaNToPiN(), G4INCL::CrossSectionsMultiPionsAndResonances::piMinuspToEtaN(), G4INCL::CrossSectionsMultiPionsAndResonances::piMinuspToOmegaN(), G4INCL::CrossSectionsINCL46::piNToDelta(), G4INCL::CrossSectionsMultiPions::piNToDelta(), and G4INCL::CrossSectionsMultiPions::piNTot().

◆ transformToLocalEnergyFrame()

void G4INCL::KinematicsUtils::transformToLocalEnergyFrame ( Nucleus const *const  n,
Particle *const  p 
)

Definition at line 45 of file G4INCLKinematicsUtils.cc.

45 {
46// assert(!p->isMeson()); // No local energy for mesons
47 const G4double localEnergy = getLocalEnergy(n, p);
48 const G4double localTotalEnergy = p->getEnergy() - localEnergy;
49 p->setEnergy(localTotalEnergy);
51 }
G4double getEnergy() const
const ThreeVector & adjustMomentumFromEnergy()
Rescale the momentum to match the total energy.
void setEnergy(G4double energy)
G4double getLocalEnergy(Nucleus const *const n, Particle *const p)

Referenced by G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), and G4INCL::InteractionAvatar::preInteractionLocalEnergy().