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

#include <G4LCapture.hh>

+ Inheritance diagram for G4LCapture:

Public Member Functions

 G4LCapture (const G4String &name="G4LCapture")
 
 ~G4LCapture ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
const G4HadronicInteractionGetMyPointer () const
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 67 of file G4LCapture.hh.

Constructor & Destructor Documentation

◆ G4LCapture()

G4LCapture::G4LCapture ( const G4String name = "G4LCapture")

Definition at line 51 of file G4LCapture.cc.

53{
54 SetMinEnergy(0.0*GeV);
56 // Description();
57}
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
#define DBL_MAX
Definition: templates.hh:83

◆ ~G4LCapture()

G4LCapture::~G4LCapture ( )

Definition at line 60 of file G4LCapture.cc.

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4LCapture::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
virtual

Implements G4HadronicInteraction.

Definition at line 80 of file G4LCapture.cc.

81{
84
85 G4double N = targetNucleus.GetA_asInt();
86 G4double Z = targetNucleus.GetZ_asInt();
87
88 const G4LorentzVector theMom = aTrack.Get4Momentum();
89 G4double P = theMom.vect().mag()/GeV;
90 G4double Px = theMom.vect().x()/GeV;
91 G4double Py = theMom.vect().y()/GeV;
92 G4double Pz = theMom.vect().z()/GeV;
93 G4double E = theMom.e()/GeV;
94 G4double E0 = aTrack.GetDefinition()->GetPDGMass()/GeV;
95 G4double Q = aTrack.GetDefinition()->GetPDGCharge();
96 if (verboseLevel > 1) {
97 G4cout << "G4LCapture:ApplyYourself: incident particle:" << G4endl;
98 G4cout << "P " << P << " GeV/c" << G4endl;
99 G4cout << "Px " << Px << " GeV/c" << G4endl;
100 G4cout << "Py " << Py << " GeV/c" << G4endl;
101 G4cout << "Pz " << Pz << " GeV/c" << G4endl;
102 G4cout << "E " << E << " GeV" << G4endl;
103 G4cout << "mass " << E0 << " GeV" << G4endl;
104 G4cout << "charge " << Q << G4endl;
105 }
106
107 // GHEISHA ADD operation to get total energy, mass, charge:
108
109 if (verboseLevel > 1) {
110 G4cout << "G4LCapture:ApplyYourself: material:" << G4endl;
111 G4cout << "A " << N << G4endl;
112 G4cout << "Z " << Z << G4endl;
113 G4cout << "atomic mass " <<
114 Atomas(N, Z) << "GeV" << G4endl;
115 }
116 E = E + Atomas(N, Z);
117 G4double E02 = E*E - P*P;
118 E0 = std::sqrt(std::abs(E02));
119 if (E02 < 0) E0 = -E0;
120 Q = Q + Z;
121 if (verboseLevel > 1) {
122 G4cout << "G4LCapture:ApplyYourself: total:" << G4endl;
123 G4cout << "E " << E << " GeV" << G4endl;
124 G4cout << "mass " << E0 << " GeV" << G4endl;
125 G4cout << "charge " << Q << G4endl;
126 }
127 Px = -Px;
128 Py = -Py;
129 Pz = -Pz;
130
131 // Make a gamma...
132
133 G4double p;
134 if (Z == 1 && N == 1) { // special case for hydrogen
135 p = 0.0022;
136 } else {
137 G4double ran = G4RandGauss::shoot();
138 p = 0.0065 + ran*0.0010;
139 }
140
141 G4double ran1 = G4UniformRand();
142 G4double ran2 = G4UniformRand();
143 G4double cost = -1. + 2.*ran1;
144 G4double sint = std::sqrt(std::abs(1. - cost*cost));
145 G4double phi = ran2*twopi;
146
147 G4double px = p*sint*std::sin(phi);
148 G4double py = p*sint*std::cos(phi);
149 G4double pz = p*cost;
150 G4double e = p;
151
152 G4double a = px*Px + py*Py + pz*Pz;
153 a = (a/(E + E0) - e)/E0;
154
155 px = px + a*Px;
156 py = py + a*Py;
157 pz = pz + a*Pz;
158
159 G4DynamicParticle* aGamma;
161 G4ThreeVector(px*GeV, py*GeV, pz*GeV));
163
164 // Make another gamma if there is sufficient energy left over...
165
166 G4double xp = 0.008 - p;
167 if (xp > 0.) {
168 if (Z > 1 || N > 1) {
169 ran1 = G4UniformRand();
170 ran2 = G4UniformRand();
171 cost = -1. + 2.*ran1;
172 sint = std::sqrt(std::abs(1. - cost*cost));
173 phi = ran2*twopi;
174
175 px = xp*sint*std::sin(phi);
176 py = xp*sint*std::cos(phi);
177 pz = xp*cost;
178 e = xp;
179
180 a = px*Px + py*Py + pz*Pz;
181 a = (a/(E + E0) - e)/E0;
182
183 px = px + a*Px;
184 py = py + a*Py;
185 pz = pz + a*Pz;
186
188 G4ThreeVector(px*GeV, py*GeV, pz*GeV));
190 }
191 }
192 return &theParticleChange;
193}
@ stopAndKill
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
double z() const
double x() const
double y() const
double mag() const
Hep3Vector vect() const
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP)
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double GetPDGCharge() const

Referenced by G4NeutronHPorLCaptureModel::ApplyYourself().

◆ GetFatalEnergyCheckLevels()

const std::pair< G4double, G4double > G4LCapture::GetFatalEnergyCheckLevels ( ) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 195 of file G4LCapture.cc.

196{
197 // max energy non-conservation is mass of heavy nucleus
198 return std::pair<G4double, G4double>(5*perCent,250*GeV);
199}

◆ ModelDescription()

void G4LCapture::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 66 of file G4LCapture.cc.

67{
68 outFile << "G4LCapture is one of the Low Energy Parameterized\n"
69 << "(LEP) models used to implement neutron capture on nuclei.\n"
70 << "It is a re-engineered version of the GHEISHA code of\n"
71 << "H. Fesefeldt which simply adds the neutron mass and energy\n"
72 << "to the target nucleus, and emits gammas isotropically as\n"
73 << "long as there is sufficient excitation energy in the\n"
74 << "daughter nucleus. The model is applicable to all incident\n"
75 << "neutron energies.\n";
76}

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