Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPFinalState Class Referenceabstract

#include <G4ParticleHPFinalState.hh>

+ Inheritance diagram for G4ParticleHPFinalState:

Public Member Functions

 G4ParticleHPFinalState ()
 
virtual ~G4ParticleHPFinalState ()
 
void Init (G4double A, G4double Z, const G4String &dirName, const G4String &aFSType, G4ParticleDefinition *p)
 
virtual void Init (G4double A, G4double Z, G4int M, const G4String &dirName, const G4String &aFSType, G4ParticleDefinition *)=0
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &)
 
virtual G4ParticleHPFinalStateNew ()=0
 
G4bool HasXsec () const
 
G4bool HasFSData () const
 
G4bool HasAnyData () const
 
virtual G4double GetXsec (G4double) const
 
virtual G4ParticleHPVectorGetXsec () const
 
void SetA_Z (G4double anA, G4double aZ, G4int aM=0)
 
G4double GetZ () const
 
G4double GetN () const
 
G4double GetA () const
 
G4int GetM () const
 
void SetAZMs (const G4ParticleHPDataUsed &used)
 
void SetAZMs (G4double anA, G4double aZ, G4int aM, const G4ParticleHPDataUsed &used)
 
void SetProjectile (G4ParticleDefinition *projectile)
 
G4ParticleHPFinalStateoperator= (const G4ParticleHPFinalState &right)=delete
 
 G4ParticleHPFinalState (const G4ParticleHPFinalState &)=delete
 

Protected Member Functions

void adjust_final_state (G4LorentzVector)
 

Protected Attributes

G4ParticleDefinitiontheProjectile {nullptr}
 
G4ParticleHPManagerfManager
 
G4IonTableionTable
 
G4int theBaseA {0}
 
G4int theBaseZ {0}
 
G4int theBaseM {0}
 
G4int theNDLDataZ {0}
 
G4int theNDLDataA {0}
 
G4int theNDLDataM {0}
 
G4int secID {-1}
 
G4bool hasXsec {true}
 
G4bool hasFSData {true}
 
G4bool hasAnyData {true}
 
G4ParticleHPNames theNames
 
G4Cache< G4HadFinalState * > theResult
 

Detailed Description

Definition at line 49 of file G4ParticleHPFinalState.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPFinalState() [1/2]

G4ParticleHPFinalState::G4ParticleHPFinalState ( )

Definition at line 47 of file G4ParticleHPFinalState.cc.

48{
50 theResult.Put(nullptr);
53}
static G4IonTable * GetIonTable()
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
G4ParticleHPManager * fManager
G4ParticleDefinition * theProjectile
G4Cache< G4HadFinalState * > theResult
static G4ParticleHPManager * GetInstance()

Referenced by G4FissionLibrary::G4FissionLibrary(), G4ParticleHPFinalState(), G4ParticleHPInelasticCompFS::G4ParticleHPInelasticCompFS(), G4FissionLibrary::New(), G4NeutronHPCaptureFS::New(), G4ParticleHP2AInelasticFS::New(), G4ParticleHP2N2AInelasticFS::New(), G4ParticleHP2NAInelasticFS::New(), G4ParticleHP2NDInelasticFS::New(), G4ParticleHP2NInelasticFS::New(), G4ParticleHP2NPInelasticFS::New(), G4ParticleHP2PInelasticFS::New(), G4ParticleHP3AInelasticFS::New(), G4ParticleHP3NAInelasticFS::New(), G4ParticleHP3NInelasticFS::New(), G4ParticleHP3NPInelasticFS::New(), G4ParticleHP4NInelasticFS::New(), G4ParticleHPAInelasticFS::New(), G4ParticleHPD2AInelasticFS::New(), G4ParticleHPDAInelasticFS::New(), G4ParticleHPDInelasticFS::New(), G4ParticleHPElasticFS::New(), G4ParticleHPFCFissionFS::New(), G4ParticleHPFFFissionFS::New(), New(), G4ParticleHPFissionFS::New(), G4ParticleHPFSFissionFS::New(), G4ParticleHPHe3InelasticFS::New(), G4ParticleHPInelasticBaseFS::New(), G4ParticleHPInelasticCompFS::New(), G4ParticleHPLCFissionFS::New(), G4ParticleHPN2AInelasticFS::New(), G4ParticleHPN2PInelasticFS::New(), G4ParticleHPN3AInelasticFS::New(), G4ParticleHPNAInelasticFS::New(), G4ParticleHPND2AInelasticFS::New(), G4ParticleHPNDInelasticFS::New(), G4ParticleHPNHe3InelasticFS::New(), G4ParticleHPNInelasticFS::New(), G4ParticleHPNPAInelasticFS::New(), G4ParticleHPNPInelasticFS::New(), G4ParticleHPNT2AInelasticFS::New(), G4ParticleHPNTInelasticFS::New(), G4ParticleHPNXInelasticFS::New(), G4ParticleHPPAInelasticFS::New(), G4ParticleHPPDInelasticFS::New(), G4ParticleHPPInelasticFS::New(), G4ParticleHPPTInelasticFS::New(), G4ParticleHPSCFissionFS::New(), G4ParticleHPT2AInelasticFS::New(), G4ParticleHPTCFissionFS::New(), G4ParticleHPTInelasticFS::New(), and operator=().

◆ ~G4ParticleHPFinalState()

G4ParticleHPFinalState::~G4ParticleHPFinalState ( )
virtual

Definition at line 55 of file G4ParticleHPFinalState.cc.

56{
57 delete theResult.Get();
58}

◆ G4ParticleHPFinalState() [2/2]

G4ParticleHPFinalState::G4ParticleHPFinalState ( const G4ParticleHPFinalState & )
delete

Member Function Documentation

◆ adjust_final_state()

void G4ParticleHPFinalState::adjust_final_state ( G4LorentzVector init_4p_lab)
protected

Definition at line 60 of file G4ParticleHPFinalState.cc.

61{
62 G4double minimum_energy = 1 * keV;
63
64 if (fManager->GetDoNotAdjustFinalState()) return;
65
66 auto nSecondaries = (G4int)theResult.Get()->GetNumberOfSecondaries();
67
68 G4int sum_Z = 0;
69 G4int sum_A = 0;
70 G4int max_SecZ = 0;
71 G4int max_SecA = 0;
72 G4int imaxA = -1;
73 for (G4int i = 0; i < nSecondaries; ++i) {
74 auto ptr = theResult.Get()->GetSecondary(i)->GetParticle()->GetDefinition();
75 sum_Z += ptr->GetAtomicNumber();
76 max_SecZ = std::max(max_SecZ, ptr->GetAtomicNumber());
77 sum_A += ptr->GetAtomicMass();
78 max_SecA = std::max(max_SecA, ptr->GetAtomicMass());
79 if (ptr->GetAtomicMass() == max_SecA)
80 imaxA = i;
81#ifdef G4PHPDEBUG
82 if (fManager->GetDEBUG())
83 G4cout << "G4ParticleHPFinalState::adjust_final_stat SECO " << i << " "
84 << ptr->GetParticleName() << G4endl;
85#endif
86 }
87
88 G4ParticleDefinition* resi_pd = nullptr;
89
90 G4int baseZNew = theBaseZ;
91 G4int baseANew = theBaseA;
93 baseANew++;
94 }
95 else if (theProjectile == G4Proton::Proton()) {
96 baseZNew++;
97 baseANew++;
98 }
99 else if (theProjectile == G4Deuteron::Deuteron()) {
100 baseZNew++;
101 baseANew += 2;
102 }
103 else if (theProjectile == G4Triton::Triton()) {
104 baseZNew++;
105 baseANew += 3;
106 }
107 else if (theProjectile == G4He3::He3()) {
108 baseZNew += 2;
109 baseANew += 3;
110 }
111 else if (theProjectile == G4Alpha::Alpha()) {
112 baseZNew += 2;
113 baseANew += 4;
114 }
115
116#ifdef G4PHPDEBUG
117 if (fManager->GetDEBUG())
118 G4cout << "G4ParticleHPFinalState::adjust_final_stat BaseZ " << baseZNew << " BaseA "
119 << baseANew << " sum_Z " << sum_Z << " sum_A " << sum_A << G4endl;
120#endif
121
122 G4bool needOneMoreSec = false;
123 G4ParticleDefinition* oneMoreSec_pd = nullptr;
124 if (baseZNew == sum_Z && baseANew == sum_A) {
125 // All secondaries are already created;
126 resi_pd = theResult.Get()->GetSecondary(imaxA)->GetParticle()->GetDefinition();
127 }
128 else {
129 if (max_SecA >= baseANew - sum_A) {
130 // Most heavy secondary is interpreted as residual
131 resi_pd = theResult.Get()->GetSecondary(imaxA)->GetParticle()->GetDefinition();
132 needOneMoreSec = true;
133 }
134 else {
135 // creation of residual is required
136 resi_pd = ionTable->GetIon(baseZNew - sum_Z, baseANew - sum_A, 0.0);
137 }
138
139 if (needOneMoreSec) {
140 if (baseZNew == sum_Z && baseANew == sum_A) {
141 // In this case, one neutron is added to secondaries
142 if (baseANew - sum_A > 1)
143 G4cout << "More than one neutron is required for the balance of baryon number!"
144 << G4endl;
145 oneMoreSec_pd = G4Neutron::Neutron();
146 }
147 else {
148#ifdef G4PHPDEBUG
149 if (fManager->GetDEBUG())
150 G4cout << this << "G4ParticleHPFinalState oneMoreSec_pd Z "
151 << baseZNew << " - " << sum_Z
152 << " A " << baseANew << " - " << sum_A << " projectile "
153 << theProjectile->GetParticleName() << G4endl;
154#endif
155 oneMoreSec_pd =
156 G4IonTable::GetIonTable()->GetIon(baseZNew - sum_Z, baseANew - sum_A, 0.0);
157 if (oneMoreSec_pd == nullptr) {
159 ed << "G4ParticleHPFinalState oneMoreSec_pd Z=" << baseZNew << " - " << sum_Z
160 << " A=" << baseANew << " - " << sum_A << " projectile "
161 << theProjectile->GetParticleName();
162 G4Exception("G4ParticleHPFinalState:adjust_final_state", "hadr01", JustWarning,
163 ed, "No adjustment will be done!");
164 return;
165 }
166 }
167 }
168
169 if (resi_pd == nullptr) {
170 // theNDLDataZ,A has the Z and A of used NDL file
171 G4int ndlZNew = theNDLDataZ;
172 G4int ndlANew = theNDLDataA;
174 ndlANew++;
175 }
176 else if (theProjectile == G4Proton::Proton()) {
177 ndlZNew++;
178 ndlANew++;
179 }
180 else if (theProjectile == G4Deuteron::Deuteron()) {
181 ndlZNew++;
182 ndlANew += 2;
183 }
184 else if (theProjectile == G4Triton::Triton()) {
185 ndlZNew++;
186 ndlANew += 3;
187 }
188 else if (theProjectile == G4He3::He3()) {
189 ndlZNew += 2;
190 ndlANew += 3;
191 }
192 else if (theProjectile == G4Alpha::Alpha()) {
193 ndlZNew += 2;
194 ndlANew += 4;
195 }
196 // theNDLDataZ,A has the Z and A of used NDL file
197 if (ndlZNew == sum_Z && ndlANew == sum_A) {
198 G4int dif_Z = theNDLDataZ - theBaseZ;
199 G4int dif_A = theNDLDataA - theBaseA;
200 resi_pd = ionTable->GetIon(max_SecZ - dif_Z, max_SecA - dif_A, 0.0);
201 if (resi_pd == nullptr) {
203 ed << "resi_pd Z=" << max_SecZ << " - " << dif_Z << " A="
204 << max_SecA << " - " << dif_A << " projectile "
205 << theProjectile->GetParticleName();
206 G4Exception("G4ParticleHPFinalState:adjust_final_state", "hadr02", JustWarning,
207 ed, "No adjustment will be done!");
208 return;
209 }
210
211 for (G4int i = 0; i < nSecondaries; ++i) {
212 auto ptr = theResult.Get()->GetSecondary(i)->GetParticle();
213 if (ptr->GetDefinition()->GetAtomicNumber() == max_SecZ &&
214 ptr->GetDefinition()->GetAtomicMass() == max_SecA)
215 {
216 G4ThreeVector p = ptr->GetMomentum() * resi_pd->GetPDGMass()
217 / ionTable->GetIon(max_SecZ, max_SecA, 0.0)->GetPDGMass();
218 ptr->SetDefinition(resi_pd);
219 ptr->SetMomentum(p);
220 }
221 }
222 }
223 }
224 }
225
226 G4LorentzVector secs_4p_lab(0.0);
227
228 auto n_sec = (G4int)theResult.Get()->GetNumberOfSecondaries();
229 G4double fast = 0;
230 G4double slow = 1;
231 G4int ifast = 0;
232 G4int islow = 0;
233 G4int ires = -1;
234
235 for (G4int i = 0; i < n_sec; ++i) {
236 auto ptr = theResult.Get()->GetSecondary(i)->GetParticle();
237 secs_4p_lab += ptr->Get4Momentum();
238
239 G4double beta = 0;
240 if (ptr->GetDefinition() != G4Gamma::Gamma()) {
241 beta = ptr->Get4Momentum().beta();
242 }
243 else {
244 beta = 1.0;
245 }
246
247 if (ptr->GetDefinition() == resi_pd) ires = i;
248
249 if (slow > beta && beta != 0) {
250 slow = beta;
251 islow = i;
252 }
253
254 if (fast <= beta) {
255 if (fast != 1) {
256 fast = beta;
257 ifast = i;
258 }
259 else {
260 // fast is already photon then check E
261 G4double e = ptr->Get4Momentum().e();
262 if (e > theResult.Get()->GetSecondary(ifast)->GetParticle()->Get4Momentum().e()) {
263 // among photons, the highest E becomes the fastest
264 ifast = i;
265 }
266 }
267 }
268 }
269
270 G4LorentzVector dif_4p = init_4p_lab - secs_4p_lab;
271
272 G4LorentzVector p4(0);
273 if (ires == -1) {
274 // Create and Add Residual Nucleus
275 ires = nSecondaries;
276 nSecondaries += 1;
277
278 auto res = new G4DynamicParticle(resi_pd, dif_4p.v());
279 theResult.Get()->AddSecondary(res, secID);
280
281 p4 = res->Get4Momentum();
282 if (slow > p4.beta()) {
283 slow = p4.beta();
284 islow = ires;
285 }
286 dif_4p = init_4p_lab - (secs_4p_lab + p4);
287 }
288
289 if (needOneMoreSec && oneMoreSec_pd != nullptr)
290 //
291 // fca: this is not a fix, this is a crash avoidance...
292 // fca: the baryon number is still wrong, most probably because it
293 // fca: should have been decreased, but since we could not create a particle
294 // fca: we just do not add it
295 //
296 {
297 nSecondaries += 1;
298 auto one = new G4DynamicParticle(oneMoreSec_pd, dif_4p.v());
299 theResult.Get()->AddSecondary(one, secID);
300 p4 = one->Get4Momentum();
301 if (slow > p4.beta()) {
302 slow = p4.beta();
303 islow = nSecondaries - 1; // Because the first is 0th, so the last becomes "nSecondaries-1"
304 }
305 dif_4p = init_4p_lab - (secs_4p_lab + p4);
306 }
307
308 // Which is bigger dif_p or dif_e
309
310 if (dif_4p.v().mag() < std::abs(dif_4p.e())) {
311 // Adjust p
312 if (minimum_energy < dif_4p.v().mag() && dif_4p.v().mag() < 1 * MeV) {
313 nSecondaries += 1;
314 theResult.Get()->AddSecondary(new G4DynamicParticle(G4Gamma::Gamma(), dif_4p.v()), secID);
315 }
316 }
317 else {
318 // dif_p > dif_e
319 // at first momentum
320 // Move residual momentum
321 p4 = theResult.Get()->GetSecondary(ires)->GetParticle()->Get4Momentum();
322 theResult.Get()->GetSecondary(ires)->GetParticle()->SetMomentum(p4.v() + dif_4p.v());
323 dif_4p = init_4p_lab -
324 (secs_4p_lab - p4 + theResult.Get()->GetSecondary(ires)->GetParticle()->Get4Momentum());
325 }
326
327 G4double dif_e = dif_4p.e() - (dif_4p.v()).mag();
328
329 if (dif_e > 0) {
330 // create 2 gamma
331
332 nSecondaries += 2;
333 G4double e1 = (dif_4p.e() - dif_4p.v().mag()) / 2;
334
335 if (minimum_energy < e1) {
337 theResult.Get()->AddSecondary(new G4DynamicParticle(G4Gamma::Gamma(), e1 * dir), secID);
338 theResult.Get()->AddSecondary(new G4DynamicParticle(G4Gamma::Gamma(), -e1 * dir), secID);
339 }
340 }
341 else // dif_e < 0
342 {
343 // At first reduce KE of the fastest secondary;
344 auto ptr = theResult.Get()->GetSecondary(ifast)->GetParticle();
345 G4double ke0 = ptr->GetKineticEnergy();
346 G4ThreeVector p0 = ptr->GetMomentum();
347 G4ThreeVector dir = p0.unit();
348
349 if (ke0 + dif_e > 0) {
350 ptr->SetKineticEnergy(ke0 + dif_e);
351 G4ThreeVector dp = p0 - theResult.Get()->GetSecondary(ifast)->GetParticle()->GetMomentum();
352 G4ThreeVector p = theResult.Get()->GetSecondary(islow)->GetParticle()->GetMomentum();
353 theResult.Get()->GetSecondary(islow)->GetParticle()->SetMomentum(p + dp);
354 }
355 }
356}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
CLHEP::HepLorentzVector G4LorentzVector
G4ThreeVector G4RandomDirection()
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
double mag() const
Hep3Vector v() const
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4He3 * He3()
Definition G4He3.cc:90
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
static G4Proton * Proton()
Definition G4Proton.cc:90
static G4Triton * Triton()
Definition G4Triton.cc:90

Referenced by G4ParticleHPInelasticBaseFS::BaseApply(), and G4ParticleHPInelasticCompFS::CompositeApply().

◆ ApplyYourself()

◆ GetA()

G4double G4ParticleHPFinalState::GetA ( ) const
inline

Definition at line 92 of file G4ParticleHPFinalState.hh.

92{ return theBaseA; }

◆ GetM()

G4int G4ParticleHPFinalState::GetM ( ) const
inline

Definition at line 93 of file G4ParticleHPFinalState.hh.

◆ GetN()

G4double G4ParticleHPFinalState::GetN ( ) const
inline

Definition at line 91 of file G4ParticleHPFinalState.hh.

91{ return theBaseA; }

◆ GetXsec() [1/2]

virtual G4ParticleHPVector * G4ParticleHPFinalState::GetXsec ( ) const
inlinevirtual

Reimplemented in G4ParticleHPFissionBaseFS, G4ParticleHPInelasticBaseFS, and G4ParticleHPInelasticCompFS.

Definition at line 82 of file G4ParticleHPFinalState.hh.

82{ return nullptr; }

◆ GetXsec() [2/2]

virtual G4double G4ParticleHPFinalState::GetXsec ( G4double ) const
inlinevirtual

Reimplemented in G4ParticleHPFissionBaseFS, G4ParticleHPInelasticBaseFS, and G4ParticleHPInelasticCompFS.

Definition at line 81 of file G4ParticleHPFinalState.hh.

81{ return 0.; }

◆ GetZ()

G4double G4ParticleHPFinalState::GetZ ( ) const
inline

Definition at line 90 of file G4ParticleHPFinalState.hh.

90{ return theBaseZ; }

◆ HasAnyData()

G4bool G4ParticleHPFinalState::HasAnyData ( ) const
inline

Definition at line 79 of file G4ParticleHPFinalState.hh.

◆ HasFSData()

G4bool G4ParticleHPFinalState::HasFSData ( ) const
inline

◆ HasXsec()

G4bool G4ParticleHPFinalState::HasXsec ( ) const
inline

Definition at line 77 of file G4ParticleHPFinalState.hh.

◆ Init() [1/2]

void G4ParticleHPFinalState::Init ( G4double A,
G4double Z,
const G4String & dirName,
const G4String & aFSType,
G4ParticleDefinition * p )
inline

Definition at line 56 of file G4ParticleHPFinalState.hh.

58 {
59 theProjectile = p;
60 Init(A, Z, 0, dirName, aFSType, p);
61 }
const G4double A[17]
void Init(G4double A, G4double Z, const G4String &dirName, const G4String &aFSType, G4ParticleDefinition *p)

Referenced by Init().

◆ Init() [2/2]

◆ New()

◆ operator=()

G4ParticleHPFinalState & G4ParticleHPFinalState::operator= ( const G4ParticleHPFinalState & right)
delete

◆ SetA_Z()

void G4ParticleHPFinalState::SetA_Z ( G4double anA,
G4double aZ,
G4int aM = 0 )
inline

Definition at line 84 of file G4ParticleHPFinalState.hh.

85 {
86 theBaseA = G4lrint(anA);
87 theBaseZ = G4lrint(aZ);
88 theBaseM = aM;
89 }
int G4lrint(double ad)
Definition templates.hh:134

Referenced by G4ParticleHPElasticFS::Init(), G4ParticleHPInelasticBaseFS::Init(), and G4ParticleHPInelasticCompFS::Init().

◆ SetAZMs() [1/2]

◆ SetAZMs() [2/2]

void G4ParticleHPFinalState::SetAZMs ( G4double anA,
G4double aZ,
G4int aM,
const G4ParticleHPDataUsed & used )
inline

Definition at line 102 of file G4ParticleHPFinalState.hh.

103 {
104 theBaseA = G4lrint(anA);
105 theBaseZ = G4lrint(aZ);
106 theBaseM = aM;
107 theNDLDataA = G4lrint(used.GetA());
108 theNDLDataZ = G4lrint(used.GetZ());
109 theNDLDataM = used.GetM();
110 }

◆ SetProjectile()

void G4ParticleHPFinalState::SetProjectile ( G4ParticleDefinition * projectile)
inline

Definition at line 112 of file G4ParticleHPFinalState.hh.

113 {
114 theProjectile = projectile;
115 }

Referenced by G4ParticleHPChannel::Register().

Member Data Documentation

◆ fManager

◆ hasAnyData

◆ hasFSData

◆ hasXsec

◆ ionTable

◆ secID

G4int G4ParticleHPFinalState::secID {-1}
protected

Definition at line 135 of file G4ParticleHPFinalState.hh.

135{-1};

Referenced by adjust_final_state(), G4NeutronHPCaptureFS::ApplyYourself(), G4ParticleHPElasticFS::ApplyYourself(), G4ParticleHPFissionFS::ApplyYourself(), G4ParticleHPInelasticBaseFS::BaseApply(), G4ParticleHPInelasticCompFS::CompositeApply(), G4NeutronHPCaptureFS::G4NeutronHPCaptureFS(), G4ParticleHP2AInelasticFS::G4ParticleHP2AInelasticFS(), G4ParticleHP2N2AInelasticFS::G4ParticleHP2N2AInelasticFS(), G4ParticleHP2NAInelasticFS::G4ParticleHP2NAInelasticFS(), G4ParticleHP2NDInelasticFS::G4ParticleHP2NDInelasticFS(), G4ParticleHP2NInelasticFS::G4ParticleHP2NInelasticFS(), G4ParticleHP2NPInelasticFS::G4ParticleHP2NPInelasticFS(), G4ParticleHP2PInelasticFS::G4ParticleHP2PInelasticFS(), G4ParticleHP3AInelasticFS::G4ParticleHP3AInelasticFS(), G4ParticleHP3NAInelasticFS::G4ParticleHP3NAInelasticFS(), G4ParticleHP3NInelasticFS::G4ParticleHP3NInelasticFS(), G4ParticleHP3NPInelasticFS::G4ParticleHP3NPInelasticFS(), G4ParticleHP4NInelasticFS::G4ParticleHP4NInelasticFS(), G4ParticleHPAInelasticFS::G4ParticleHPAInelasticFS(), G4ParticleHPD2AInelasticFS::G4ParticleHPD2AInelasticFS(), G4ParticleHPDAInelasticFS::G4ParticleHPDAInelasticFS(), G4ParticleHPDInelasticFS::G4ParticleHPDInelasticFS(), G4ParticleHPElasticFS::G4ParticleHPElasticFS(), G4ParticleHPFissionFS::G4ParticleHPFissionFS(), G4ParticleHPHe3InelasticFS::G4ParticleHPHe3InelasticFS(), G4ParticleHPN2AInelasticFS::G4ParticleHPN2AInelasticFS(), G4ParticleHPN2PInelasticFS::G4ParticleHPN2PInelasticFS(), G4ParticleHPN3AInelasticFS::G4ParticleHPN3AInelasticFS(), G4ParticleHPNAInelasticFS::G4ParticleHPNAInelasticFS(), G4ParticleHPND2AInelasticFS::G4ParticleHPND2AInelasticFS(), G4ParticleHPNDInelasticFS::G4ParticleHPNDInelasticFS(), G4ParticleHPNHe3InelasticFS::G4ParticleHPNHe3InelasticFS(), G4ParticleHPNInelasticFS::G4ParticleHPNInelasticFS(), G4ParticleHPNPAInelasticFS::G4ParticleHPNPAInelasticFS(), G4ParticleHPNPInelasticFS::G4ParticleHPNPInelasticFS(), G4ParticleHPNT2AInelasticFS::G4ParticleHPNT2AInelasticFS(), G4ParticleHPNTInelasticFS::G4ParticleHPNTInelasticFS(), G4ParticleHPNXInelasticFS::G4ParticleHPNXInelasticFS(), G4ParticleHPPAInelasticFS::G4ParticleHPPAInelasticFS(), G4ParticleHPPDInelasticFS::G4ParticleHPPDInelasticFS(), G4ParticleHPPInelasticFS::G4ParticleHPPInelasticFS(), G4ParticleHPPTInelasticFS::G4ParticleHPPTInelasticFS(), G4ParticleHPT2AInelasticFS::G4ParticleHPT2AInelasticFS(), and G4ParticleHPTInelasticFS::G4ParticleHPTInelasticFS().

◆ theBaseA

◆ theBaseM

G4int G4ParticleHPFinalState::theBaseM {0}
protected

Definition at line 130 of file G4ParticleHPFinalState.hh.

130{0};

Referenced by GetM(), SetA_Z(), and SetAZMs().

◆ theBaseZ

◆ theNames

◆ theNDLDataA

◆ theNDLDataM

G4int G4ParticleHPFinalState::theNDLDataM {0}
protected

Definition at line 133 of file G4ParticleHPFinalState.hh.

133{0};

Referenced by G4ParticleHPElasticFS::Init(), SetAZMs(), and SetAZMs().

◆ theNDLDataZ

◆ theProjectile

◆ theResult

G4Cache<G4HadFinalState*> G4ParticleHPFinalState::theResult
protected

Definition at line 143 of file G4ParticleHPFinalState.hh.

Referenced by adjust_final_state(), G4FissionLibrary::ApplyYourself(), G4NeutronHPCaptureFS::ApplyYourself(), G4ParticleHP2AInelasticFS::ApplyYourself(), G4ParticleHP2N2AInelasticFS::ApplyYourself(), G4ParticleHP2NAInelasticFS::ApplyYourself(), G4ParticleHP2NDInelasticFS::ApplyYourself(), G4ParticleHP2NInelasticFS::ApplyYourself(), G4ParticleHP2NPInelasticFS::ApplyYourself(), G4ParticleHP2PInelasticFS::ApplyYourself(), G4ParticleHP3AInelasticFS::ApplyYourself(), G4ParticleHP3NAInelasticFS::ApplyYourself(), G4ParticleHP3NInelasticFS::ApplyYourself(), G4ParticleHP3NPInelasticFS::ApplyYourself(), G4ParticleHP4NInelasticFS::ApplyYourself(), G4ParticleHPAInelasticFS::ApplyYourself(), G4ParticleHPD2AInelasticFS::ApplyYourself(), G4ParticleHPDAInelasticFS::ApplyYourself(), G4ParticleHPDInelasticFS::ApplyYourself(), G4ParticleHPElasticFS::ApplyYourself(), G4ParticleHPFissionFS::ApplyYourself(), G4ParticleHPHe3InelasticFS::ApplyYourself(), G4ParticleHPN2AInelasticFS::ApplyYourself(), G4ParticleHPN2PInelasticFS::ApplyYourself(), G4ParticleHPN3AInelasticFS::ApplyYourself(), G4ParticleHPNAInelasticFS::ApplyYourself(), G4ParticleHPND2AInelasticFS::ApplyYourself(), G4ParticleHPNDInelasticFS::ApplyYourself(), G4ParticleHPNHe3InelasticFS::ApplyYourself(), G4ParticleHPNInelasticFS::ApplyYourself(), G4ParticleHPNPAInelasticFS::ApplyYourself(), G4ParticleHPNPInelasticFS::ApplyYourself(), G4ParticleHPNT2AInelasticFS::ApplyYourself(), G4ParticleHPNTInelasticFS::ApplyYourself(), G4ParticleHPNXInelasticFS::ApplyYourself(), G4ParticleHPPAInelasticFS::ApplyYourself(), G4ParticleHPPDInelasticFS::ApplyYourself(), G4ParticleHPPInelasticFS::ApplyYourself(), G4ParticleHPPTInelasticFS::ApplyYourself(), G4ParticleHPT2AInelasticFS::ApplyYourself(), G4ParticleHPTInelasticFS::ApplyYourself(), G4ParticleHPInelasticBaseFS::BaseApply(), G4ParticleHPInelasticCompFS::CompositeApply(), G4ParticleHPFinalState(), and ~G4ParticleHPFinalState().


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