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

#include <G4PiMinusStopAbsorption.hh>

Public Member Functions

 G4PiMinusStopAbsorption (G4PiMinusStopMaterial *materialAlgo, const G4double Z, const G4double A)
 
 ~G4PiMinusStopAbsorption ()
 
G4DynamicParticleVectorDoAbsorption ()
 
G4double Energy ()
 
G4ThreeVector RecoilMomentum ()
 
G4int NProtons ()
 
G4int NNeutrons ()
 
void SetVerboseLevel (G4int level)
 

Detailed Description

Definition at line 46 of file G4PiMinusStopAbsorption.hh.

Constructor & Destructor Documentation

◆ G4PiMinusStopAbsorption()

G4PiMinusStopAbsorption::G4PiMinusStopAbsorption ( G4PiMinusStopMaterial materialAlgo,
const G4double  Z,
const G4double  A 
)

Definition at line 53 of file G4PiMinusStopAbsorption.cc.

56{
57 G4HadronicDeprecate("G4PiMinusStopAbsorption");
58 _materialAlgo = materialAlgo;
59 _nucleusZ = Z;
60 _nucleusA = A;
61 _level = 0;
62 _absorptionProducts = new G4DynamicParticleVector();
63}
std::vector< G4DynamicParticle * > G4DynamicParticleVector
#define G4HadronicDeprecate(name)

◆ ~G4PiMinusStopAbsorption()

G4PiMinusStopAbsorption::~G4PiMinusStopAbsorption ( )

Definition at line 67 of file G4PiMinusStopAbsorption.cc.

68{
69 // Memory management of materialAlgo needs better thought (MGP)
70 delete _materialAlgo;
71 // Who owns it? Memory management is not clear... (MGP)
72 // _absorptionProducts->clearAndDestroy();
73 delete _absorptionProducts;
74}

Member Function Documentation

◆ DoAbsorption()

G4DynamicParticleVector * G4PiMinusStopAbsorption::DoAbsorption ( )

Definition at line 76 of file G4PiMinusStopAbsorption.cc.

77{
78 std::vector<G4ParticleDefinition*>* defNucleons = _materialAlgo->DefinitionVector();
79
80 G4double newA = _nucleusA;
81 G4double newZ = _nucleusZ;
82
83 if (defNucleons != 0)
84 {
85 for (unsigned int i=0; i<defNucleons->size(); i++)
86 {
87 if ( (*defNucleons)[i] == G4Proton::Proton())
88 {
89 newA = newA - 1;
90 newZ = newZ - 1;
91 }
92 if ((*defNucleons)[i] == G4Neutron::Neutron())
93 { newA = newA - 1; }
94 }
95 }
96
97 G4double binding = G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(_nucleusA) ,static_cast<G4int>(_nucleusZ)) / _nucleusA;
98 G4double mass = G4NucleiProperties::GetNuclearMass(static_cast<G4int>(newA),static_cast<G4int>(newZ));
99
100
101 std::vector<G4LorentzVector*>* p4Nucleons = _materialAlgo->P4Vector(binding,mass);
102
103 if (defNucleons != 0 && p4Nucleons != 0)
104 {
105 unsigned int nNucleons = p4Nucleons->size();
106
107 G4double seen = _materialAlgo->FinalNucleons() / 2.;
108 G4int maxN = nNucleons;
109 if (defNucleons->size() < nNucleons) { maxN = defNucleons->size(); }
110
111 for (G4int i=0; i<maxN; i++)
112 {
113 G4DynamicParticle* product;
114 if ((*defNucleons)[i] == G4Proton::Proton())
115 { product = new G4DynamicParticle(G4Proton::Proton(),*((*p4Nucleons)[i])); }
116 else
117 { product = new G4DynamicParticle(G4Neutron::Neutron(),*((*p4Nucleons)[i])); }
118 G4double ranflat = G4UniformRand();
119
120 if (ranflat < seen)
121 { _absorptionProducts->push_back(product); }
122 else
123 { delete product; }
124 }
125 }
126
127 return _absorptionProducts;
128
129}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4double GetBindingEnergy(const G4int A, const G4int Z)
static G4double GetNuclearMass(const G4double A, const G4double Z)
virtual std::vector< G4ParticleDefinition * > * DefinitionVector()
virtual std::vector< G4LorentzVector * > * P4Vector(const G4double binding, const G4double mass)
virtual G4double FinalNucleons()=0
static G4Proton * Proton()
Definition: G4Proton.cc:93

Referenced by G4PiMinusAbsorptionAtRest::AtRestDoIt().

◆ Energy()

G4double G4PiMinusStopAbsorption::Energy ( )

Definition at line 169 of file G4PiMinusStopAbsorption.cc.

170{
171 G4double energy = 0.;
172 G4double productEnergy = 0.;
173 G4ThreeVector pProducts(0.,0.,0.);
174 G4int nN = 0;
175 G4int nP = 0;
176
177
178 G4int nAbsorptionProducts = _absorptionProducts->size();
179
180 for (int i = 0; i<nAbsorptionProducts; i++)
181 {
182 productEnergy += (*_absorptionProducts)[i]->GetKineticEnergy();
183 pProducts = pProducts + (*_absorptionProducts)[i]->GetMomentum();
184 if ((*_absorptionProducts)[i]->GetDefinition() == G4Neutron::Neutron()) nN++;
185 if ((*_absorptionProducts)[i]->GetDefinition() == G4Proton::Proton()) nP++;
186 }
187
188 G4double productBinding = (G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(_nucleusA),static_cast<G4int>(_nucleusZ)) / _nucleusA) * nAbsorptionProducts;
189 G4double mass = G4NucleiProperties::GetNuclearMass(_nucleusA - (nP + nN),_nucleusZ - nP);
190 G4double pNucleus = pProducts.mag();
191 G4double eNucleus = std::sqrt(pNucleus*pNucleus + mass*mass);
192 G4double tNucleus = eNucleus - mass;
193 G4double temp =
194 G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(_nucleusA - (nP + nN)),static_cast<G4int>(_nucleusZ - nP)) -
195 G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(_nucleusA),static_cast<G4int>(_nucleusZ));
196 energy = productEnergy + productBinding + tNucleus;
197
198 if (_level > 0)
199 {
200 std::cout << "E products " << productEnergy
201 << " Binding " << productBinding << " " << temp << " "
202 << " Tnucleus " << tNucleus
203 << " energy = " << energy << G4endl;
204 }
205
206 return energy;
207}
#define G4endl
Definition: G4ios.hh:52

Referenced by G4PiMinusAbsorptionAtRest::AtRestDoIt().

◆ NNeutrons()

G4int G4PiMinusStopAbsorption::NNeutrons ( )

Definition at line 156 of file G4PiMinusStopAbsorption.cc.

157{
158 G4int n = 0;
159 G4int entries = _absorptionProducts->size();
160 for (int i = 0; i<entries; i++)
161 {
162 if ((*_absorptionProducts)[i]->GetDefinition() == G4Neutron::Neutron())
163 { n = n + 1; }
164 }
165 return n;
166}

Referenced by G4PiMinusAbsorptionAtRest::AtRestDoIt().

◆ NProtons()

G4int G4PiMinusStopAbsorption::NProtons ( )

Definition at line 143 of file G4PiMinusStopAbsorption.cc.

144{
145 G4int n = 0;
146 G4int entries = _absorptionProducts->size();
147 for (int i = 0; i<entries; i++)
148 {
149 if ((*_absorptionProducts)[i]->GetDefinition() == G4Proton::Proton())
150 { n = n + 1; }
151 }
152 return n;
153}

Referenced by G4PiMinusAbsorptionAtRest::AtRestDoIt().

◆ RecoilMomentum()

G4ThreeVector G4PiMinusStopAbsorption::RecoilMomentum ( )

Definition at line 131 of file G4PiMinusStopAbsorption.cc.

132{
133 G4ThreeVector pProducts(0.,0.,0.);
134
135 for (unsigned int i = 0; i< _absorptionProducts->size(); i++)
136 {
137 pProducts = pProducts + (*_absorptionProducts)[i]->GetMomentum();
138 }
139 return pProducts;
140}

Referenced by G4PiMinusAbsorptionAtRest::AtRestDoIt().

◆ SetVerboseLevel()

void G4PiMinusStopAbsorption::SetVerboseLevel ( G4int  level)

Definition at line 209 of file G4PiMinusStopAbsorption.cc.

210{
211 _level = level;
212 return;
213}

Referenced by G4PiMinusAbsorptionAtRest::AtRestDoIt().


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