Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4LEDeuteronInelastic.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// Hadronic Process: Deuteron Inelastic Process
27// J.L. Chuma, TRIUMF, 25-Feb-1997
28// J.L. Chuma, 08-May-2001: Update original incident passed back in vec[0]
29// from NuclearReaction
30
32#include "G4SystemOfUnits.hh"
33#include "Randomize.hh"
34#include "G4Electron.hh"
35
36void G4LEDeuteronInelastic::ModelDescription(std::ostream& outFile) const
37{
38 outFile << "G4LEDeuteronInelastic is one of the Low Energy Parameterized\n"
39 << "(LEP) models used to implement inelastic deuteron scattering\n"
40 << "from nuclei. It is a re-engineered version of the GHEISHA\n"
41 << "code of H. Fesefeldt. It divides the initial collision\n"
42 << "products into backward- and forward-going clusters which are\n"
43 << "then decayed into final state hadrons. The model does not\n"
44 << "conserve energy on an event-by-event basis. It may be\n"
45 << "applied to deuterons with initial energies between 0 and 10\n"
46 << "TeV.\n";
47}
48
49
52 G4Nucleus& targetNucleus)
53{
55 const G4HadProjectile* originalIncident = &aTrack;
56
57 if (verboseLevel > 1) {
58 const G4Material *targetMaterial = aTrack.GetMaterial();
59 G4cout << "G4LEDeuteronInelastic::ApplyYourself called" << G4endl;
60 G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
61 G4cout << "target material = " << targetMaterial->GetName() << ", ";
62 }
63
64 // Work-around for lack of model above 100 MeV
65 if (originalIncident->GetKineticEnergy()/MeV > 100. ||
66 originalIncident->GetKineticEnergy() <= 0.1*MeV) {
70 return &theParticleChange;
71 }
72
73 G4double A = targetNucleus.GetA_asInt();
74 G4double Z = targetNucleus.GetZ_asInt();
75 G4double theAtomicMass = targetNucleus.AtomicMass(A, Z);
76 G4double massVec[9];
77 massVec[0] = targetNucleus.AtomicMass( A+2.0, Z+1.0 );
78 massVec[1] = targetNucleus.AtomicMass( A+1.0, Z+1.0 );
79 massVec[2] = targetNucleus.AtomicMass( A+1.0, Z );
80 massVec[3] = theAtomicMass;
81 massVec[4] = 0.;
82 if (A > 1.0 && A-1.0 > Z)
83 massVec[4] = targetNucleus.AtomicMass(A-1.0, Z);
84 massVec[5] = 0.;
85 if (A > 2.0 && Z > 1.0 && A-2.0 > Z-1.0)
86 massVec[5] = targetNucleus.AtomicMass(A-2.0, Z-1.0);
87 massVec[6] = 0.;
88 if (A > Z+1.0)
89 massVec[6] = targetNucleus.AtomicMass(A, Z+1.0);
90 massVec[7] = massVec[3];
91 massVec[8] = 0.;
92 if (Z > 1.0) massVec[8] = targetNucleus.AtomicMass(A,Z-1.0);
93
94 G4FastVector<G4ReactionProduct,4> vec; // vec will contain the secondary particles
95 G4int vecLen = 0;
96 vec.Initialize( 0 );
97
98 theReactionDynamics.NuclearReaction(vec, vecLen, originalIncident,
99 targetNucleus, theAtomicMass, massVec);
100
101 G4double p = vec[0]->GetMomentum().mag();
102 theParticleChange.SetMomentumChange( vec[0]->GetMomentum() * (1.0/p) );
103 theParticleChange.SetEnergyChange( vec[0]->GetKineticEnergy() );
104 delete vec[0];
105
106 if (vecLen <= 1)
107 {
111 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
112 return &theParticleChange;
113 }
114
116 for (G4int i=1; i<vecLen; ++i) {
117 pd = new G4DynamicParticle();
118 pd->SetDefinition( vec[i]->GetDefinition() );
119 pd->SetMomentum( vec[i]->GetMomentum() );
121 delete vec[i];
122 }
123
124 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
125 return &theParticleChange;
126}
127
128 /* end of file */
129
@ isAlive
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
Hep3Vector unit() const
Hep3Vector vect() const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4ThreeVector &momentum)
void Initialize(G4int items)
Definition: G4FastVector.hh:63
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4Material * GetMaterial() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4ReactionDynamics theReactionDynamics
void DoIsotopeCounting(const G4HadProjectile *theProjectile, const G4Nucleus &aNucleus)
virtual void ModelDescription(std::ostream &outFile) const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
const G4String & GetName() const
Definition: G4Material.hh:177
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double AtomicMass(const G4double A, const G4double Z) const
Definition: G4Nucleus.cc:240
void NuclearReaction(G4FastVector< G4ReactionProduct, 4 > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4Nucleus &aNucleus, const G4double theAtomicMass, const G4double *massVec)