Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LightTargetCollider Class Reference

#include <G4LightTargetCollider.hh>

+ Inheritance diagram for G4LightTargetCollider:

Public Member Functions

 G4LightTargetCollider ()
 
virtual ~G4LightTargetCollider ()
 
void collide (G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)
 
void setVerboseLevel (G4int verbose=0)
 
- Public Member Functions inherited from G4CascadeColliderBase
 G4CascadeColliderBase (const G4String &name, G4int verbose=0)
 
virtual ~G4CascadeColliderBase ()
 
virtual void rescatter (G4InuclParticle *, G4KineticTrackVector *, G4V3DNucleus *, G4CollisionOutput &)
 
virtual void setVerboseLevel (G4int verbose=0)
 
- Public Member Functions inherited from G4VCascadeCollider
 G4VCascadeCollider (const G4String &name, G4int verbose=0)
 
virtual ~G4VCascadeCollider ()
 
virtual void collide (G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)=0
 
virtual void setVerboseLevel (G4int verbose=0)
 

Additional Inherited Members

- Protected Member Functions inherited from G4CascadeColliderBase
virtual G4bool useEPCollider (G4InuclParticle *bullet, G4InuclParticle *target) const
 
virtual G4bool inelasticInteractionPossible (G4InuclParticle *bullet, G4InuclParticle *target, G4double ekin) const
 
virtual G4bool validateOutput (G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
 
virtual G4bool validateOutput (const G4Fragment &fragment, G4CollisionOutput &output)
 
virtual G4bool validateOutput (G4InuclParticle *bullet, G4InuclParticle *target, const std::vector< G4InuclElementaryParticle > &particles)
 
- Protected Member Functions inherited from G4VCascadeCollider
virtual void setName (const G4String &name)
 
- Protected Attributes inherited from G4CascadeColliderBase
G4InteractionCase interCase
 
G4CascadeCheckBalancebalance
 
- Protected Attributes inherited from G4VCascadeCollider
G4String theName
 
G4int verboseLevel
 

Detailed Description

Definition at line 52 of file G4LightTargetCollider.hh.

Constructor & Destructor Documentation

◆ G4LightTargetCollider()

G4LightTargetCollider::G4LightTargetCollider ( )

Definition at line 55 of file G4LightTargetCollider.cc.

56 : G4CascadeColliderBase("G4LightTargetCollider"),
57 theElementaryParticleCollider(new G4ElementaryParticleCollider)
58{
59 mP = G4Proton::Proton()->GetPDGMass()/CLHEP::GeV;
60 mN = G4Neutron::Neutron()->GetPDGMass()/CLHEP::GeV;
61 mD = G4Deuteron::Deuteron()->GetPDGMass()/CLHEP::GeV;
62 pFermiD = 0.045; // Fermi momentum of nucleon in deuteron Hulthen potential
63}
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4Proton * Proton()
Definition: G4Proton.cc:92

◆ ~G4LightTargetCollider()

G4LightTargetCollider::~G4LightTargetCollider ( )
virtual

Definition at line 65 of file G4LightTargetCollider.cc.

65 {
66 delete theElementaryParticleCollider;
67}

Member Function Documentation

◆ collide()

void G4LightTargetCollider::collide ( G4InuclParticle bullet,
G4InuclParticle target,
G4CollisionOutput globalOutput 
)
virtual

Implements G4VCascadeCollider.

Definition at line 78 of file G4LightTargetCollider.cc.

81{
82 if (verboseLevel) {
83 G4cout << " >>> G4LightTargetCollider::collide" << G4endl;
84 G4cout << " Projectile: " << bullet->getDefinition()->GetParticleName() << G4endl;
85 G4cout << " Target: " << target->getDefinition()->GetParticleName() << G4endl;
86 }
87
88 G4double ke = bullet->getKineticEnergy();
89
90 if (target->getDefinition() == G4Proton::Proton() ) {
91 if (ke < 0.1447) {
92 // Below threshold lab energy for pi0 creation
93 globalOutput.trivialise(bullet, target);
94 } else {
95 // Need inelastic cross section in this class if we want to replace ElementaryParticleCollider
96 // with SingleNucleonScattering
97 theElementaryParticleCollider->collide(bullet, target, globalOutput);
98 if (globalOutput.numberOfOutgoingParticles() == 0) globalOutput.trivialise(bullet, target);
99 }
100
101 } else if (target->getDefinition() == G4Deuteron::Deuteron() ) {
102
103 if (ke < mP + mN - mD) {
104 // Should not happen as long as inelastic cross section is zero
105 G4Exception("G4LightTargetCollider::collide()","HAD_BERT_201",
106 JustWarning, "Projectile energy below reaction threshold");
107 globalOutput.trivialise(bullet, target);
108
109 } else {
110 // Get p, n and deuteron cross sections; use lab energy to access
113 G4double gammaDXS = GammaDCrossSection(ke);
114
115 G4double probP = 0.0;
116 G4double probN = 0.0;
117 // Highest threshold is 0.152 (for gamma p -> n pi+)
118 // Because of Fermi momentum in deuteron, raise this to 0.159
119 if (ke > 0.159) {
120 G4double totalDXS = gammaPXS + gammaNXS + gammaDXS;
121 probP = gammaPXS/totalDXS;
122 probN = (gammaPXS+gammaNXS)/totalDXS;
123 }
124
125 G4double rndm = G4UniformRand();
126 if (rndm < probP) {
127 // Generate Fermi momenta of bullet and target
128 G4ThreeVector fermiMomentum = pFermiD*G4RandomDirection();
129 G4LorentzVector protonMomentum(fermiMomentum, std::sqrt(mP*mP + pFermiD*pFermiD) );
130 G4LorentzVector neutronMomentum(-fermiMomentum, std::sqrt(mN*mN + pFermiD*pFermiD) );
131
132 G4LorentzVector bulletMomentum = bullet->getMomentum();
133 G4ThreeVector betacm = bulletMomentum.findBoostToCM(protonMomentum);
134
135 // First boost bullet and target so that target is at rest
136 G4ThreeVector toProtonRest = -protonMomentum.boostVector();
137 protonMomentum.boost(toProtonRest);
138 bulletMomentum.boost(toProtonRest);
139
140 G4InuclElementaryParticle projectile(bulletMomentum, bullet->getDefinition() );
141 G4InuclElementaryParticle targetNucleon(protonMomentum, G4Proton::Proton() );
142 G4InuclElementaryParticle spectatorNucleon(neutronMomentum, G4Neutron::Neutron() );
143 ScatteringProducts products = SingleNucleonScattering(projectile, targetNucleon);
144
145 // Particles from SingleNucleonScattering are in CM frame of projectile
146 // and moving proton. Transform back to lab frame with -betacm, then
147 // add them to outgoing list.
148 globalOutput.reset();
149 G4LorentzVector temp;
150 for (G4int i = 0; i < G4int(products.size()); i++) {
151 temp = products[i].getMomentum();
152 temp.boost(-betacm);
153 products[i].setMomentum(temp);
154 globalOutput.addOutgoingParticle(products[i]);
155 }
156
157 // Add the recoil nucleon unmodified
158 globalOutput.addOutgoingParticle(spectatorNucleon);
159
160 } else if (rndm < probN) {
161 G4ThreeVector fermiMomentum = pFermiD*G4RandomDirection();
162 G4LorentzVector protonMomentum(fermiMomentum, std::sqrt(mP*mP + pFermiD*pFermiD) );
163 G4LorentzVector neutronMomentum(-fermiMomentum, std::sqrt(mN*mN + pFermiD*pFermiD) );
164
165 G4LorentzVector bulletMomentum = bullet->getMomentum();
166 G4ThreeVector betacm = bulletMomentum.findBoostToCM(neutronMomentum);
167
168 // First boost bullet and target so that target is at rest
169 G4ThreeVector toNeutronRest = -neutronMomentum.boostVector();
170 neutronMomentum.boost(toNeutronRest);
171 bulletMomentum.boost(toNeutronRest);
172
173 G4InuclElementaryParticle projectile(bulletMomentum, bullet->getDefinition() );
174 G4InuclElementaryParticle targetNucleon(neutronMomentum, G4Neutron::Neutron() );
175 G4InuclElementaryParticle spectatorNucleon(protonMomentum, G4Proton::Proton() );
176
177 ScatteringProducts products = SingleNucleonScattering(projectile, targetNucleon);
178
179 // Particles from SingleNucleonScattering are in CM frame of projectile
180 // and moving neutron. Transform back to lab frame with -betacm, then add
181 // them to outgoing list
182 globalOutput.reset();
183 G4LorentzVector temp;
184 for (G4int i = 0; i < G4int(products.size()); i++) {
185 temp = products[i].getMomentum();
186 temp.boost(-betacm);
187 products[i].setMomentum(temp);
188 globalOutput.addOutgoingParticle(products[i]);
189 }
190
191 // Add the recoil nucleon unmodified
192 globalOutput.addOutgoingParticle(spectatorNucleon);
193
194 } else {
195 NucleonPair products = AbsorptionOnDeuteron(bullet);
196 globalOutput.reset();
197 globalOutput.addOutgoingParticle(products.first);
198 globalOutput.addOutgoingParticle(products.second);
199 }
200 } // Energy above threshold ?
201
202 // Test code
203 // G4int numPart = globalOutput.numberOfOutgoingParticles();
204 // std::vector<G4InuclElementaryParticle> testList = globalOutput.getOutgoingParticles();
205 // G4LorentzVector sumP;
206 // G4cout << " Global output " << G4endl;
207 // for (G4int i = 0; i < numPart; i++) {
208 // sumP += testList[i].getMomentum();
209 // G4cout << testList[i] << G4endl;
210 // }
211 // G4cout << " Global 4-momentum sum = " << sumP << G4endl;
212 // G4cout << " Initial lab energy = " << mD + bullet->getEnergy() << G4endl;
213
214 } else {
215 G4Exception("G4LightTargetCollider::collide()","HAD_BERT_203",
216 FatalException, "Scattering from this target not implemented");
217 }
218
219 return;
220}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::vector< G4InuclElementaryParticle > ScatteringProducts
std::pair< G4InuclElementaryParticle, G4InuclElementaryParticle > NucleonPair
G4ThreeVector G4RandomDirection()
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
HepLorentzVector & boost(double, double, double)
Hep3Vector findBoostToCM() const
static const G4CascadeChannel * GetTable(G4int initialState)
virtual G4double getCrossSection(double ke) const =0
G4int numberOfOutgoingParticles() const
void addOutgoingParticle(const G4InuclElementaryParticle &particle)
void trivialise(G4InuclParticle *bullet, G4InuclParticle *target)
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
const G4ParticleDefinition * getDefinition() const
G4double getKineticEnergy() const
G4LorentzVector getMomentum() const
const G4String & GetParticleName() const

Referenced by G4CascadeInterface::ApplyYourself().

◆ setVerboseLevel()

void G4LightTargetCollider::setVerboseLevel ( G4int  verbose = 0)
virtual

Reimplemented from G4CascadeColliderBase.

Definition at line 71 of file G4LightTargetCollider.cc.

71 {
73 theElementaryParticleCollider->setVerboseLevel(verboseLevel);
75}
virtual void setVerboseLevel(G4int verbose=0)
void setVerboseLevel(G4int verbose)

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