Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4InelasticInteraction.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// $Id$
27//
28// Hadronic Process: Inelastic Interaction
29// original by H.P. Wellisch
30// modified by J.L. Chuma, TRIUMF, 22-Nov-1996
31// Last modified: 27-Mar-1997
32// J.P. Wellisch: 23-Apr-97: throw G4HadronicException(__FILE__, __LINE__, removed
33// J.P. Wellisch: 24-Apr-97: correction for SetUpPions
34// Modified by J.L. Chuma, 30-Apr-97: added originalTarget to CalculateMomenta
35// since TwoBody needed to reset the target particle
36// J.L. Chuma, 20-Jun-97: Modified CalculateMomenta to correct the decision process
37// for whether to use GenerateXandPt or TwoCluster
38// J.L. Chuma, 06-Aug-97: added original incident particle, before Fermi motion and
39// evaporation effects are included, needed for calculating
40// self absorption and corrections for single particle spectra
41// HPW removed misunderstanding of LocalEnergyDeposit, 11.04.98.
42// 23-Jan-2009 V.Ivanchenko move constructor and destructor to the body
43
46#include "G4SystemOfUnits.hh"
47#include "Randomize.hh"
49#include "G4IsoResult.hh"
51
52G4IsoParticleChange* G4InelasticInteraction::theIsoResult = 0;
53G4IsoParticleChange* G4InelasticInteraction::theOldIsoResult = 0;
54
56 : G4HadronicInteraction(name), isotopeProduction(false), cache(0.0)
57{}
58
60{}
61
62// Pmltpc used in Cascade functions
65 G4double b, G4double c)
66{
67 const G4double expxu = 82.; // upper bound for arg. of exp
68 const G4double expxl = -expxu; // lower bound for arg. of exp
69 G4double npf = 0.0;
70 G4double nmf = 0.0;
71 G4double nzf = 0.0;
72 G4int i;
73 for (i = 2; i <= npos; i++) npf += std::log((G4double)i);
74 for (i = 2; i <= nneg; i++) nmf += std::log((G4double)i);
75 for (i = 2; i <= nzero; i++) nzf += std::log((G4double)i);
76 G4double r = std::min(expxu, std::max(expxl,
77 -(npos-nneg+nzero+b)*(npos-nneg+nzero+b)/(2*c*c*n*n)
78 - npf - nmf - nzf ) );
79 return std::exp(r);
80}
81
82
84 const G4ReactionProduct& currentParticle,
85 const G4ReactionProduct& targetParticle,
86 G4ReactionProduct& leadParticle)
87{
88 // the following was in GenerateXandPt and TwoCluster
89 // add a parameter to the GenerateXandPt function telling it about the strange particle
90 //
91 // assumes that the original particle was a strange particle
92 //
93 G4bool lead = false;
94 if ((currentParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
95 (currentParticle.GetDefinition() != G4Proton::Proton()) &&
96 (currentParticle.GetDefinition() != G4Neutron::Neutron()) ) {
97 lead = true;
98 leadParticle = currentParticle; // set lead to the incident particle
99
100 } else if ((targetParticle.GetMass() >= G4KaonPlus::KaonPlus()->GetPDGMass()) &&
101 (targetParticle.GetDefinition() != G4Proton::Proton()) &&
102 (targetParticle.GetDefinition() != G4Neutron::Neutron() ) ) {
103 lead = true;
104 leadParticle = targetParticle; // set lead to the target particle
105 }
106
107 return lead;
108}
109
110
111void
113 const G4int nzero,
115 G4int& vecLen)
116{
117 if (npos + nneg + nzero == 0) return;
118 G4int i;
120
121 for (i = 0; i < npos; ++i) {
122 p = new G4ReactionProduct;
124 (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
125 vec.SetElement( vecLen++, p );
126 }
127 for (i = npos; i < npos+nneg; ++i) {
128 p = new G4ReactionProduct;
130 (G4UniformRand() < 0.5) ? p->SetSide(-1) : p->SetSide(1);
131 vec.SetElement( vecLen++, p );
132 }
133 for (i = npos+nneg; i < npos+nneg+nzero; ++i) {
134 p = new G4ReactionProduct;
136 (G4UniformRand() < 0.5) ? p->SetSide(-1) : p->SetSide(1);
137 vec.SetElement( vecLen++, p );
138 }
139}
140
141
142void
144 const G4double energy, // MeV, <0 means annihilation channels
145 G4double &n,
146 G4double &anpn )
147 {
148 const G4double expxu = 82.; // upper bound for arg. of exp
149 const G4double expxl = -expxu; // lower bound for arg. of exp
150 const G4int numSec = 60;
151 //
152 // the only difference between the calculation for annihilation channels
153 // and normal is the starting value, iBegin, for the loop below
154 //
155 G4int iBegin = 1;
156 G4double en = energy;
157 if( energy < 0.0 )
158 {
159 iBegin = 2;
160 en *= -1.0;
161 }
162 //
163 // number of total particles vs. centre of mass Energy - 2*proton mass
164 //
165 G4double aleab = std::log(en/GeV);
166 n = 3.62567 + aleab*(0.665843 + aleab*(0.336514 + aleab*(0.117712 + 0.0136912*aleab)));
167 n -= 2.0;
168 //
169 // normalization constant for kno-distribution
170 //
171 anpn = 0.0;
172 G4double test, temp;
173 for( G4int i=iBegin; i<=numSec; ++i )
174 {
175 temp = pi*i/(2.0*n*n);
176 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(i*i)/(n*n) ) ) );
177 if( temp < 1.0 )
178 {
179 if( test >= 1.0e-10 )anpn += temp*test;
180 }
181 else
182 anpn += temp*test;
183 }
184 }
185
188 G4int& vecLen,
189 const G4HadProjectile* originalIncident,
190 const G4DynamicParticle* originalTarget,
191 G4ReactionProduct& modifiedOriginal, // Fermi motion and evap. effects included
192 G4Nucleus& targetNucleus,
193 G4ReactionProduct& currentParticle,
194 G4ReactionProduct& targetParticle,
195 G4bool& incidentHasChanged,
196 G4bool& targetHasChanged,
197 G4bool quasiElastic)
198{
199 cache = 0;
200 what = originalIncident->Get4Momentum().vect();
201
202 theReactionDynamics.ProduceStrangeParticlePairs(vec, vecLen, modifiedOriginal,
203 originalTarget, currentParticle,
204 targetParticle, incidentHasChanged,
205 targetHasChanged);
206
207 if (quasiElastic) {
208 theReactionDynamics.TwoBody(vec, vecLen,
209 modifiedOriginal, originalTarget,
210 currentParticle, targetParticle,
211 targetNucleus, targetHasChanged);
212 return;
213 }
214 G4ReactionProduct leadingStrangeParticle;
215 G4bool leadFlag = MarkLeadingStrangeParticle(currentParticle,
216 targetParticle,
217 leadingStrangeParticle);
218
219 // Note: the number of secondaries can be reduced in GenerateXandPt and TwoCluster
220 G4bool finishedGenXPt = false;
221 G4bool annihilation = false;
222 if( originalIncident->GetDefinition()->GetPDGEncoding() < 0 &&
223 currentParticle.GetMass() == 0.0 && targetParticle.GetMass() == 0.0 )
224 {
225 // original was an anti-particle and annihilation has taken place
226 annihilation = true;
227 G4double ekcor = 1.0;
228 G4double ek = originalIncident->GetKineticEnergy();
229 G4double ekOrg = ek;
230
231 const G4double tarmas = originalTarget->GetDefinition()->GetPDGMass();
232 if( ek > 1.0*GeV )ekcor = 1./(ek/GeV);
233 const G4double atomicWeight = targetNucleus.GetA_asInt();
234 ek = 2*tarmas + ek*(1.+ekcor/atomicWeight);
235 G4double tkin = targetNucleus.Cinema(ek);
236 ek += tkin;
237 ekOrg += tkin;
238 // modifiedOriginal.SetKineticEnergy( ekOrg );
239 //
240 // evaporation -- re-calculate black track energies
241 // this was Done already just before the cascade
242 //
243 tkin = targetNucleus.AnnihilationEvaporationEffects(ek, ekOrg);
244 ekOrg -= tkin;
245 ekOrg = std::max( 0.0001*GeV, ekOrg );
246 modifiedOriginal.SetKineticEnergy( ekOrg );
247 G4double amas = originalIncident->GetDefinition()->GetPDGMass();
248 G4double et = ekOrg + amas;
249 G4double p = std::sqrt( std::abs(et*et-amas*amas) );
250 G4double pp = modifiedOriginal.GetMomentum().mag();
251 if( pp > 0.0 )
252 {
253 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
254 modifiedOriginal.SetMomentum( momentum * (p/pp) );
255 }
256 if( ekOrg <= 0.0001 )
257 {
258 modifiedOriginal.SetKineticEnergy( 0.0 );
259 modifiedOriginal.SetMomentum( 0.0, 0.0, 0.0 );
260 }
261 }
262 const G4double twsup[] = { 1.0, 0.7, 0.5, 0.3, 0.2, 0.1 };
263 G4double rand1 = G4UniformRand();
264 G4double rand2 = G4UniformRand();
265
266 // Cache current, target, and secondaries
267 G4ReactionProduct saveCurrent = currentParticle;
268 G4ReactionProduct saveTarget = targetParticle;
269 std::vector<G4ReactionProduct> savevec;
270 for (G4int i = 0; i < vecLen; i++) savevec.push_back(*vec[i]);
271
272 if (annihilation ||
273 vecLen >= 6 ||
274 ( modifiedOriginal.GetKineticEnergy()/GeV >= 1.0 &&
275 ( ( (originalIncident->GetDefinition() == G4KaonPlus::KaonPlus() ||
276 originalIncident->GetDefinition() == G4KaonMinus::KaonMinus() ||
277 originalIncident->GetDefinition() == G4KaonZeroLong::KaonZeroLong() ||
278 originalIncident->GetDefinition() == G4KaonZeroShort::KaonZeroShort() )
279 &&
280 rand1 < 0.5 )
281 || rand2 > twsup[vecLen] ) ) )
282
283 finishedGenXPt =
285 modifiedOriginal, originalIncident,
286 currentParticle, targetParticle,
287 originalTarget,
288 targetNucleus, incidentHasChanged,
289 targetHasChanged, leadFlag,
290 leadingStrangeParticle );
291 if( finishedGenXPt )
292 {
293 Rotate(vec, vecLen);
294 return;
295 }
296
297 G4bool finishedTwoClu = false;
298 if( modifiedOriginal.GetTotalMomentum()/MeV < 1.0 )
299 {
300 for(G4int i=0; i<vecLen; i++) delete vec[i];
301 vecLen = 0;
302 }
303 else
304 {
305 // Occaisionally, GenerateXandPt will fail in the annihilation channel.
306 // Restore current, target and secondaries to pre-GenerateXandPt state
307 // before trying annihilation in TwoCluster
308
309 if (!finishedGenXPt && annihilation) {
310 currentParticle = saveCurrent;
311 targetParticle = saveTarget;
312 for (G4int i = 0; i < vecLen; i++) delete vec[i];
313 vecLen = 0;
314 vec.Initialize( 0 );
315 for (G4int i = 0; i < G4int(savevec.size()); i++) {
317 *p = savevec[i];
318 vec.SetElement( vecLen++, p );
319 }
320 }
321
323 modifiedOriginal, currentParticle,
324 targetParticle, targetNucleus,
325 incidentHasChanged, targetHasChanged );
326 try
327 {
328 finishedTwoClu = theReactionDynamics.TwoCluster( vec, vecLen,
329 modifiedOriginal, originalIncident,
330 currentParticle, targetParticle,
331 originalTarget,
332 targetNucleus, incidentHasChanged,
333 targetHasChanged, leadFlag,
334 leadingStrangeParticle );
335 }
337 {
338 aC.Report(G4cout);
339 throw G4HadReentrentException(__FILE__, __LINE__, "Failing to calculate momenta");
340 }
341 }
342
343 if (finishedTwoClu) {
344 Rotate(vec, vecLen);
345 return;
346 }
347
348 theReactionDynamics.TwoBody(vec, vecLen, modifiedOriginal, originalTarget,
349 currentParticle, targetParticle,
350 targetNucleus, targetHasChanged);
351}
352
353
356{
357 G4double rotation = 2.*pi*G4UniformRand();
358 cache = rotation;
359 G4int i;
360 for (i = 0; i < vecLen; ++i) {
361 G4ThreeVector momentum = vec[i]->GetMomentum();
362 momentum = momentum.rotate(rotation, what);
363 vec[i]->SetMomentum(momentum);
364 }
365}
366
367
370 G4int& vecLen,
371 G4ReactionProduct& currentParticle,
372 G4ReactionProduct& targetParticle,
373 G4bool& incidentHasChanged)
374{
378 G4int i;
379 if (currentParticle.GetDefinition() == aKaonZL) {
380 if (G4UniformRand() <= 0.5) {
381 currentParticle.SetDefinition(aKaonZS);
382 incidentHasChanged = true;
383 }
384 } else if (currentParticle.GetDefinition() == aKaonZS) {
385 if (G4UniformRand() > 0.5) {
386 currentParticle.SetDefinition(aKaonZL);
387 incidentHasChanged = true;
388 }
389 }
390
391 if (targetParticle.GetDefinition() == aKaonZL) {
392 if (G4UniformRand() <= 0.5) targetParticle.SetDefinition(aKaonZS);
393 } else if (targetParticle.GetDefinition() == aKaonZS) {
394 if (G4UniformRand() > 0.5 ) targetParticle.SetDefinition(aKaonZL);
395 }
396
397 for (i = 0; i < vecLen; ++i) {
398 if (vec[i]->GetDefinition() == aKaonZL) {
399 if( G4UniformRand() <= 0.5) vec[i]->SetDefinition(aKaonZS);
400 } else if (vec[i]->GetDefinition() == aKaonZS) {
401 if (G4UniformRand() > 0.5) vec[i]->SetDefinition(aKaonZL);
402 }
403 }
404
405 if (incidentHasChanged) {
407 p0->SetDefinition( currentParticle.GetDefinition() );
408 p0->SetMomentum( currentParticle.GetMomentum() );
412 } else {
413 G4double p = currentParticle.GetMomentum().mag()/MeV;
414 G4ThreeVector pvec = currentParticle.GetMomentum();
415 if (p > DBL_MIN)
416 theParticleChange.SetMomentumChange(pvec.x()/p, pvec.y()/p, pvec.z()/p);
417 else
419
420 G4double aE = currentParticle.GetKineticEnergy();
421 if (std::fabs(aE) < .1*eV) aE = .1*eV;
423 }
424
425 if (targetParticle.GetMass() > 0.0) {
426 // targetParticle can be eliminated in TwoBody
428 p1->SetDefinition(targetParticle.GetDefinition() );
429 G4ThreeVector momentum = targetParticle.GetMomentum();
430 momentum = momentum.rotate(cache, what);
431 p1->SetMomentum(momentum);
433 }
434
436 for (i = 0; i < vecLen; ++i) {
437 p = new G4DynamicParticle(vec[i]->GetDefinition(), vec[i]->GetMomentum() );
439 delete vec[i];
440 }
441}
442
443
445{
446 G4IsoParticleChange* anIsoResult = theIsoResult;
447 if(theIsoResult) theOldIsoResult = theIsoResult;
448 theIsoResult = 0;
449 return anIsoResult;
450}
451
452
453void
455 const G4Nucleus& aNucleus)
456{
457 delete theOldIsoResult;
458 theOldIsoResult = 0;
459 delete theIsoResult;
460 theIsoResult = new G4IsoParticleChange;
461 G4IsoResult* anIsoResult = 0;
462 G4int nModels = theProductionModels.size();
463 if (nModels > 0) {
464 for (G4int i = 0; i < nModels; i++) {
465 anIsoResult = theProductionModels[i]->GetIsotope(theProjectile, aNucleus);
466 if (anIsoResult) break; // accept first result
467 }
468 if (!anIsoResult) anIsoResult = ExtractResidualNucleus(aNucleus);
469 } else {
470 // No production models active, use default iso production
471 anIsoResult = ExtractResidualNucleus(aNucleus);
472 }
473
474/*
475 G4cout << " contents of anIsoResult (from ExtractResidualNucleus) " << G4endl;
476 G4cout << " original projectile: "
477 << theProjectile->GetDefinition()->GetParticleName() << G4endl;
478 G4cout << " mother nucleus: "
479 << anIsoResult->GetMotherNucleus().GetA_asInt() << ","
480 << anIsoResult->GetMotherNucleus().GetZ_asInt() << G4endl;
481 G4cout << " extracted nucleus = " << anIsoResult->GetIsotope() << G4endl;
482 G4cout << " end contents of anIsoResult " << G4endl;
483*/
484 // Add all info explicitly and add typename from model called.
485 theIsoResult->SetIsotope(anIsoResult->GetIsotope());
486 theIsoResult->SetProductionTime(theProjectile->GetGlobalTime() );
487 theIsoResult->SetParentParticle(theProjectile->GetDefinition() );
488 theIsoResult->SetMotherNucleus(anIsoResult->GetMotherNucleus());
489 theIsoResult->SetProducer(this->GetModelName() );
490
491 delete anIsoResult;
492
493 // If isotope production is enabled the GetIsotopeProductionInfo()
494 // method must be called or else a memory leak will result
495 //
496 // The following code will fix the memory leak, but remove the
497 // isotope information:
498 //
499 // if(theIsoResult) {
500 // delete theIsoResult;
501 // theIsoResult = 0;
502 // }
503}
504
505
508{
509 G4double A = aNucleus.GetA_asInt();
510 G4double Z = aNucleus.GetZ_asInt();
511 G4double bufferA = 0;
512 G4double bufferZ = 0;
513
514 // Loop over theParticleChange, decrement A, Z accordingly, and
515 // cache the largest fragment
516 for (G4int i = 0; i < theParticleChange.GetNumberOfSecondaries(); ++i) {
518 const G4ParticleDefinition* part = aSecTrack->GetParticle()->GetParticleDefinition();
519 G4double Q = part->GetPDGCharge()/eplus;
520 G4double BN = part->GetBaryonNumber();
521 if (bufferA < BN) {
522 bufferA = BN;
523 bufferZ = Q;
524 }
525 Z -= Q;
526 A -= BN;
527 }
528
529 // if the fragment was part of the final state, it is
530 // assumed to be the heaviest secondary.
531 if (A < 0.1) {
532 A = bufferA;
533 Z = bufferZ;
534 }
535
536 // Prepare the IsoResult
537 std::ostringstream ost1;
538 ost1 << Z << "_" << A;
539 G4String biff = ost1.str();
540 G4IsoResult* theResult = new G4IsoResult(biff, aNucleus);
541
542 return theResult;
543}
544
545const std::pair<G4double, G4double> G4InelasticInteraction::GetFatalEnergyCheckLevels() const
546{
547 // max energy non-conservation is mass of heavy nucleus
548 return std::pair<G4double, G4double>(5*perCent,250*GeV);
549}
@ stopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
double z() const
double x() const
double y() const
double mag() const
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:28
Hep3Vector vect() const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
const G4ParticleDefinition * GetParticleDefinition() const
G4ParticleDefinition * GetDefinition() const
void SetMomentum(const G4ThreeVector &momentum)
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
void Initialize(G4int items)
Definition: G4FastVector.hh:63
void SetStatusChange(G4HadFinalStateStatus aS)
G4int GetNumberOfSecondaries() const
void AddSecondary(G4DynamicParticle *aP)
void SetEnergyChange(G4double anEnergy)
G4HadSecondary * GetSecondary(size_t i)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
void Report(std::ostream &aS)
G4DynamicParticle * GetParticle()
const G4String & GetModelName() const
G4IsoResult * ExtractResidualNucleus(const G4Nucleus &aNucleus)
static G4IsoParticleChange * GetIsotopeProductionInfo()
void Rotate(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen)
G4double Pmltpc(G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)
void GetNormalizationConstant(const G4double availableEnergy, G4double &n, G4double &anpn)
G4ReactionDynamics theReactionDynamics
G4InelasticInteraction(const G4String &name="LEInelastic")
void SetUpPions(const G4int np, const G4int nm, const G4int nz, G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen)
void CalculateMomenta(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4DynamicParticle *originalTarget, G4ReactionProduct &modifiedOriginal, G4Nucleus &targetNucleus, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool quasiElastic)
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
void DoIsotopeCounting(const G4HadProjectile *theProjectile, const G4Nucleus &aNucleus)
G4bool MarkLeadingStrangeParticle(const G4ReactionProduct &currentParticle, const G4ReactionProduct &targetParticle, G4ReactionProduct &leadParticle)
void SetUpChange(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
void SetProducer(const G4String &aProducer)
void SetIsotope(const G4String &anIsotope)
void SetProductionTime(const G4double &aProductionTime)
void SetMotherNucleus(const G4Nucleus &aTarget)
void SetParentParticle(const G4ParticleDefinition *aProjectile)
G4Nucleus GetMotherNucleus()
Definition: G4IsoResult.hh:51
G4String GetIsotope()
Definition: G4IsoResult.hh:50
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
static G4KaonZeroLong * KaonZeroLong()
static G4KaonZeroShort * KaonZeroShort()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double AnnihilationEvaporationEffects(G4double kineticEnergy, G4double ekOrg)
Definition: G4Nucleus.cc:323
G4double Cinema(G4double kineticEnergy)
Definition: G4Nucleus.cc:368
G4double GetPDGCharge() const
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4PionZero * PionZero()
Definition: G4PionZero.cc:104
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4bool TwoCluster(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &modifiedOriginal, const G4HadProjectile *originalIncident, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, const G4DynamicParticle *originalTarget, const G4Nucleus &targetNucleus, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool leadFlag, G4ReactionProduct &leadingStrangeParticle)
void SuppressChargedPions(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, const G4ReactionProduct &modifiedOriginal, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, const G4Nucleus &targetNucleus, G4bool &incidentHasChanged, G4bool &targetHasChanged)
void ProduceStrangeParticlePairs(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, const G4ReactionProduct &modifiedOriginal, const G4DynamicParticle *originalTarget, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged)
G4bool GenerateXandPt(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &modifiedOriginal, const G4HadProjectile *originalIncident, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, const G4DynamicParticle *originalTarget, const G4Nucleus &targetNucleus, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool leadFlag, G4ReactionProduct &leadingStrangeParticle)
void TwoBody(G4FastVector< G4ReactionProduct, GHADLISTSIZE > &vec, G4int &vecLen, G4ReactionProduct &modifiedOriginal, const G4DynamicParticle *originalTarget, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, const G4Nucleus &targetNucleus, G4bool &targetHasChanged)
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4double GetTotalMomentum() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void SetSide(const G4int sid)
void SetKineticEnergy(const G4double en)
G4ParticleDefinition * GetDefinition() const
void SetDefinition(G4ParticleDefinition *aParticleDefinition)
G4double GetMass() const
#define DBL_MIN
Definition: templates.hh:75