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
G4LEAlphaInelastic.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: Alpha 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
31#include <iostream>
32
33#include "G4LEAlphaInelastic.hh"
34#include "G4SystemOfUnits.hh"
35#include "Randomize.hh"
36#include "G4Electron.hh"
37
40{
41 SetMinEnergy(0.0*GeV);
42 SetMaxEnergy(10.*TeV);
43 G4cout << "WARNING: model G4LEAlphaInelastic is being deprecated and will\n"
44 << "disappear in Geant4 version 10.0" << G4endl;
45}
46
47
48void G4LEAlphaInelastic::ModelDescription(std::ostream& outFile) const
49{
50 outFile << "G4LEAlphaInelastic is one of the Low Energy Parameterized\n"
51 << "(LEP) models used to implement inelastic alpha scattering\n"
52 << "from nuclei. It is a re-engineered version of the GHEISHA\n"
53 << "code of H. Fesefeldt. It divides the initial collision\n"
54 << "products into backward- and forward-going clusters which are\n"
55 << "then decayed into final state hadrons. The model does not\n"
56 << "conserve energy on an event-by-event basis. It may be\n"
57 << "applied to alphas with initial energies between 0 and 10\n"
58 << "TeV.\n";
59}
60
61
64 G4Nucleus& targetNucleus)
65{
67 const G4HadProjectile* originalIncident = &aTrack;
68
69 G4double A = targetNucleus.GetA_asInt();
70 G4double Z = targetNucleus.GetZ_asInt();
71
72 G4double kineticEnergy = aTrack.Get4Momentum().e()-aTrack.GetDefinition()->GetPDGMass();
73 if (verboseLevel > 1) {
74 const G4Material *targetMaterial = aTrack.GetMaterial();
75 G4cout << "G4LEAlphaInelastic::ApplyYourself called" << G4endl;
76 G4cout << "kinetc energy = " <<kineticEnergy/MeV << "MeV, ";
77 G4cout << "target material = " << targetMaterial->GetName() << G4endl;
78 }
79
80 // Work-around for lack of model above 100 MeV
81 if (kineticEnergy/MeV > 100. || kineticEnergy <= 0.1*MeV) {
85 return &theParticleChange;
86 }
87 G4double theAtomicMass = targetNucleus.AtomicMass( A, Z );
88 G4double massVec[9];
89 massVec[0] = targetNucleus.AtomicMass( A+4.0, Z+2.0 );
90 massVec[1] = targetNucleus.AtomicMass( A+3.0, Z+2.0 );
91 massVec[2] = targetNucleus.AtomicMass( A+3.0, Z+1.0 );
92 massVec[3] = targetNucleus.AtomicMass( A+2.0, Z+1.0 );
93 massVec[4] = targetNucleus.AtomicMass( A+1.0, Z+1.0 );
94 massVec[5] = theAtomicMass;
95 massVec[6] = targetNucleus.AtomicMass( A+2.0, Z+2.0 );
96 massVec[7] = massVec[3];
97 massVec[8] = targetNucleus.AtomicMass( A+2.0, Z );
98
99 G4FastVector<G4ReactionProduct,4> vec; // vec will contain the secondary particles
100 G4int vecLen = 0;
101 vec.Initialize( 0 );
102
103 theReactionDynamics.NuclearReaction(vec, vecLen, &aTrack,
104 targetNucleus, theAtomicMass, massVec);
105
106 G4double p = vec[0]->GetMomentum().mag();
107 theParticleChange.SetMomentumChange( vec[0]->GetMomentum() *(1./p));
108 theParticleChange.SetEnergyChange( vec[0]->GetKineticEnergy() );
109 delete vec[0];
110
111 if (vecLen <= 1)
112 {
116 return &theParticleChange;
117 }
118
120 for (G4int i = 1; i < vecLen; ++i) {
121 pd = new G4DynamicParticle();
122 pd->SetDefinition( vec[i]->GetDefinition() );
123 pd->SetMomentum( vec[i]->GetMomentum() );
125 delete vec[i];
126 }
127
128 if (isotopeProduction) DoIsotopeCounting(originalIncident, targetNucleus);
129 return &theParticleChange;
130}
@ 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
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
G4ReactionDynamics theReactionDynamics
void DoIsotopeCounting(const G4HadProjectile *theProjectile, const G4Nucleus &aNucleus)
G4LEAlphaInelastic(const G4String &name="G4LEAlphaInelastic")
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)