Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronRadCapture.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//
27//
28// Physics model class G4NeutronRadCapture
29// Created: 31 August 2009
30// Author V.Ivanchenko
31//
32// Modified:
33// 09.09.2010 V.Ivanchenko added usage of G4PhotonEvaporation
34//
35
37#include "G4SystemOfUnits.hh"
39#include "G4Fragment.hh"
40#include "G4FragmentVector.hh"
41#include "G4NucleiProperties.hh"
44#include "G4DynamicParticle.hh"
45#include "G4ParticleTable.hh"
46#include "G4IonTable.hh"
47#include "G4Electron.hh"
48#include "G4Deuteron.hh"
49#include "G4Triton.hh"
50#include "G4He3.hh"
51#include "G4Alpha.hh"
52#include "G4RandomDirection.hh"
55
57 : G4HadronicInteraction("nRadCapture"),
58 photonEvaporation(nullptr),lab4mom(0.,0.,0.,0.)
59{
60 lowestEnergyLimit = 10*CLHEP::eV;
61 minExcitation = 0.1*CLHEP::keV;
62 SetMinEnergy( 0.0*CLHEP::GeV );
64
65 electron = G4Electron::Electron();
66 icID = -1;
67 secID = -1;
69}
70
72{
73 delete photonEvaporation;
74}
75
77{
78 if(photonEvaporation != nullptr) { return; }
79 G4DeexPrecoParameters* param =
81 minExcitation = param->GetMinExcitation();
82 icID = G4PhysicsModelCatalog::GetModelID("model_e-InternalConversion");
84 photonEvaporation = new G4PhotonEvaporation();
85 photonEvaporation->Initialise();
86 photonEvaporation->SetICM(true);
87}
88
90 const G4HadProjectile& aTrack, G4Nucleus& theNucleus)
91{
94
95 G4int A = theNucleus.GetA_asInt();
96 G4int Z = theNucleus.GetZ_asInt();
97
98 G4double time = aTrack.GetGlobalTime();
99
100 // Create initial state
101 lab4mom.set(0.,0.,0.,G4NucleiProperties::GetNuclearMass(A, Z));
102 lab4mom += aTrack.Get4Momentum();
103
104 G4double M = lab4mom.mag();
105 ++A;
107 //G4cout << "Capture start: Z= " << Z << " A= " << A
108 // << " LabM= " << M << " Mcompound= " << mass << G4endl;
109
110 // simplified method of 1 gamma emission
111 if(A <= 4) {
112
113 G4ThreeVector bst = lab4mom.boostVector();
114
115 if(M - mass <= lowestEnergyLimit) {
116 return &theParticleChange;
117 }
118
119 if (verboseLevel > 1) {
120 G4cout << "G4NeutronRadCapture::DoIt: Eini(MeV)="
121 << aTrack.GetKineticEnergy()/MeV << " Eexc(MeV)= "
122 << (M - mass)/MeV
123 << " Z= " << Z << " A= " << A << G4endl;
124 }
125 G4double e1 = (M - mass)*(M + mass)/(2*M);
127 lv2.boost(bst);
128 G4HadSecondary* news =
130 news->SetTime(time);
131 news->SetCreatorModelID(secID);
133 delete news;
134
135 const G4ParticleDefinition* theDef = 0;
136
137 lab4mom -= lv2;
138 if (Z == 1 && A == 2) {theDef = G4Deuteron::Deuteron();}
139 else if (Z == 1 && A == 3) {theDef = G4Triton::Triton();}
140 else if (Z == 2 && A == 3) {theDef = G4He3::He3();}
141 else if (Z == 2 && A == 4) {theDef = G4Alpha::Alpha();}
142 else { theDef = theTableOfIons->GetIon(Z,A,0.0,noFloat,0); }
143
144 if (verboseLevel > 1) {
145 G4cout << "Gamma 4-mom: " << lv2 << " "
146 << theDef->GetParticleName() << " " << lab4mom << G4endl;
147 }
148 if(theDef) {
149 news = new G4HadSecondary(new G4DynamicParticle(theDef, lab4mom));
150 news->SetTime(time);
151 news->SetCreatorModelID(secID);
153 delete news;
154 }
155
156 // Use photon evaporation
157 } else {
158
159 // protection against wrong kinematic
160 if(M < mass) {
161 G4double etot = std::max(mass, lab4mom.e());
162 G4double ptot = std::sqrt((etot - mass)*(etot + mass));
163 G4ThreeVector v = lab4mom.vect().unit();
164 lab4mom.set(v.x()*ptot,v.y()*ptot,v.z()*ptot,etot);
165 }
166
167 G4Fragment* aFragment = new G4Fragment(A, Z, lab4mom);
168
169 if (verboseLevel > 1) {
170 G4cout << "G4NeutronRadCapture::ApplyYourself initial G4Fragmet:"
171 << G4endl;
172 G4cout << aFragment << G4endl;
173 }
174
175 //
176 // Sample final state
177 //
178 G4FragmentVector* fv = photonEvaporation->BreakUpFragment(aFragment);
179 if(!fv) { fv = new G4FragmentVector(); }
180 fv->push_back(aFragment);
181 size_t n = fv->size();
182
183 if (verboseLevel > 1) {
184 G4cout << "G4NeutronRadCapture: " << n << " final particle icID= " << icID << G4endl;
185 }
186 for(size_t i=0; i<n; ++i) {
187
188 G4Fragment* f = (*fv)[i];
189 G4double etot = f->GetMomentum().e();
190
191 Z = f->GetZ_asInt();
192 A = f->GetA_asInt();
193
194 const G4ParticleDefinition* theDef;
195 if(0 == Z && 0 == A) {theDef = f->GetParticleDefinition();}
196 else if (Z == 1 && A == 2) {theDef = G4Deuteron::Deuteron();}
197 else if (Z == 1 && A == 3) {theDef = G4Triton::Triton();}
198 else if (Z == 2 && A == 3) {theDef = G4He3::He3();}
199 else if (Z == 2 && A == 4) {theDef = G4Alpha::Alpha();}
200 else {
201 G4double eexc = f->GetExcitationEnergy();
202 if(eexc <= minExcitation) { eexc = 0.0; }
203 theDef = theTableOfIons->GetIon(Z, A, eexc, noFloat, 0);
204 /*
205 G4cout << "### NC Find ion Z= " << Z << " A= " << A
206 << " Eexc(MeV)= " << eexc/MeV << " "
207 << theDef << G4endl;
208 */
209 }
210 G4double ekin = std::max(0.0,etot - theDef->GetPDGMass());
211 if (verboseLevel > 1) {
212 G4cout << i << ". " << theDef->GetParticleName()
213 << " Ekin(MeV)= " << etot/MeV
214 << " p: " << f->GetMomentum().vect()
215 << G4endl;
216 }
217 G4HadSecondary* news = new G4HadSecondary(
218 new G4DynamicParticle(theDef,
219 f->GetMomentum().vect().unit(),
220 ekin));
221 G4double timeF = f->GetCreationTime();
222 if(timeF < 0.0) { timeF = 0.0; }
223 news->SetTime(time + timeF);
224 if(theDef == electron) {
225 news->SetCreatorModelID(icID);
226 } else {
227 news->SetCreatorModelID(secID);
228 }
230 delete news;
231 delete f;
232 }
233 delete fv;
234 }
235 //G4cout << "Capture done" << G4endl;
236 return &theParticleChange;
237}
238
std::vector< G4Fragment * > G4FragmentVector
Definition G4Fragment.hh:65
@ stopAndKill
#define noFloat
Definition G4Ions.hh:119
#define M(row, col)
G4ThreeVector G4RandomDirection()
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double z() const
Hep3Vector unit() const
double x() const
double y() const
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
void set(double x, double y, double z, double t)
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
static G4Electron * Electron()
Definition G4Electron.cc:91
G4double GetExcitationEnergy() const
const G4LorentzVector & GetMomentum() const
G4double GetCreationTime() const
G4int GetZ_asInt() const
const G4ParticleDefinition * GetParticleDefinition() const
G4int GetA_asInt() const
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
void SetMinEnergy(G4double anEnergy)
const G4String & GetModelName() const
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
static G4He3 * He3()
Definition G4He3.cc:90
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) final
virtual void InitialiseModel() final
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4int GetModelID(const G4int modelIndex)
static G4Triton * Triton()
Definition G4Triton.cc:90
G4FragmentVector * BreakUpFragment(G4Fragment *theNucleus)
virtual void SetICM(G4bool)