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

#include <G4NeutronRadCapture.hh>

+ Inheritance diagram for G4NeutronRadCapture:

Public Member Functions

 G4NeutronRadCapture ()
 
 ~G4NeutronRadCapture () override
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
 
void InitialiseModel () override
 
G4NeutronRadCaptureoperator= (const G4NeutronRadCapture &right)=delete
 
 G4NeutronRadCapture (const G4NeutronRadCapture &)=delete
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 51 of file G4NeutronRadCapture.hh.

Constructor & Destructor Documentation

◆ G4NeutronRadCapture() [1/2]

G4NeutronRadCapture::G4NeutronRadCapture ( )

Definition at line 56 of file G4NeutronRadCapture.cc.

57 : G4HadronicInteraction("nRadCapture"),
58 photonEvaporation(nullptr),lab4mom(0.,0.,0.,0.)
59{
60 lowestEnergyLimit = 10*CLHEP::eV;
61 minExcitation = 0.1*CLHEP::keV;
62 secID = -1;
64}
G4HadronicInteraction(const G4String &modelName="HadronicModel")
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()

Referenced by G4NeutronRadCapture(), and operator=().

◆ ~G4NeutronRadCapture()

G4NeutronRadCapture::~G4NeutronRadCapture ( )
override

Definition at line 66 of file G4NeutronRadCapture.cc.

67{
68 delete photonEvaporation;
69}

◆ G4NeutronRadCapture() [2/2]

G4NeutronRadCapture::G4NeutronRadCapture ( const G4NeutronRadCapture & )
delete

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4NeutronRadCapture::ApplyYourself ( const G4HadProjectile & aTrack,
G4Nucleus & targetNucleus )
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 83 of file G4NeutronRadCapture.cc.

85{
86 theParticleChange.Clear();
87 theParticleChange.SetStatusChange(stopAndKill);
88
89 G4int A = theNucleus.GetA_asInt();
90 G4int Z = theNucleus.GetZ_asInt();
91
92 G4double time = aTrack.GetGlobalTime();
93
94 // Create initial state
95 lab4mom.set(0.,0.,0.,G4NucleiProperties::GetNuclearMass(A, Z));
96 lab4mom += aTrack.Get4Momentum();
97
98 G4double M = lab4mom.mag();
99 ++A;
101 //G4cout << "Capture start: Z= " << Z << " A= " << A
102 // << " LabM= " << M << " Mcompound= " << mass << G4endl;
103
104 // simplified method of 1 gamma emission
105 if(A <= 4) {
106
107 G4ThreeVector bst = lab4mom.boostVector();
108
109 if(M - mass <= lowestEnergyLimit) {
110 return &theParticleChange;
111 }
112
113 if (verboseLevel > 1) {
114 G4cout << "G4NeutronRadCapture::DoIt: Eini(MeV)="
115 << aTrack.GetKineticEnergy()/MeV << " Eexc(MeV)= "
116 << (M - mass)/MeV
117 << " Z= " << Z << " A= " << A << G4endl;
118 }
119 G4double e1 = (M - mass)*(M + mass)/(2*M);
121 lv2.boost(bst);
122 G4HadSecondary* news =
123 new G4HadSecondary(new G4DynamicParticle(G4Gamma::Gamma(), lv2));
124 news->SetTime(time);
125 news->SetCreatorModelID(secID);
126 theParticleChange.AddSecondary(*news);
127 delete news;
128
129 const G4ParticleDefinition* theDef = 0;
130
131 lab4mom -= lv2;
132 if (Z == 1 && A == 2) {theDef = G4Deuteron::Deuteron();}
133 else if (Z == 1 && A == 3) {theDef = G4Triton::Triton();}
134 else if (Z == 2 && A == 3) {theDef = G4He3::He3();}
135 else if (Z == 2 && A == 4) {theDef = G4Alpha::Alpha();}
136 else { theDef = theTableOfIons->GetIon(Z,A,0.0,noFloat,0); }
137
138 if (verboseLevel > 1) {
139 G4cout << "Gamma 4-mom: " << lv2 << " "
140 << theDef->GetParticleName() << " " << lab4mom << G4endl;
141 }
142 if(theDef) {
143 news = new G4HadSecondary(new G4DynamicParticle(theDef, lab4mom));
144 news->SetTime(time);
145 news->SetCreatorModelID(secID);
146 theParticleChange.AddSecondary(*news);
147 delete news;
148 }
149
150 // Use photon evaporation
151 } else {
152
153 // protection against wrong kinematic
154 if(M < mass) {
155 G4double etot = std::max(mass, lab4mom.e());
156 G4double ptot = std::sqrt((etot - mass)*(etot + mass));
157 G4ThreeVector v = lab4mom.vect().unit();
158 lab4mom.set(v.x()*ptot,v.y()*ptot,v.z()*ptot,etot);
159 }
160
161 G4Fragment* aFragment = new G4Fragment(A, Z, lab4mom);
162
163 if (verboseLevel > 1) {
164 G4cout << "G4NeutronRadCapture::ApplyYourself initial G4Fragmet:"
165 << G4endl;
166 G4cout << aFragment << G4endl;
167 }
168
169 //
170 // Sample final state
171 //
172 G4FragmentVector* fv = photonEvaporation->BreakUpFragment(aFragment);
173 if (nullptr == fv) { fv = new G4FragmentVector(); }
174 fv->push_back(aFragment);
175 std::size_t n = fv->size();
176
177 if (verboseLevel > 1) {
178 G4cout << "G4NeutronRadCapture: " << n << " final particles" << G4endl;
179 }
180 for(std::size_t i=0; i<n; ++i) {
181
182 G4Fragment* f = (*fv)[i];
183 G4double etot = f->GetMomentum().e();
184
185 Z = f->GetZ_asInt();
186 A = f->GetA_asInt();
187
188 const G4ParticleDefinition* theDef;
189 if(0 == Z && 0 == A) {theDef = f->GetParticleDefinition();}
190 else if (Z == 1 && A == 2) {theDef = G4Deuteron::Deuteron();}
191 else if (Z == 1 && A == 3) {theDef = G4Triton::Triton();}
192 else if (Z == 2 && A == 3) {theDef = G4He3::He3();}
193 else if (Z == 2 && A == 4) {theDef = G4Alpha::Alpha();}
194 else {
195 G4double eexc = f->GetExcitationEnergy();
196 if(eexc <= minExcitation) { eexc = 0.0; }
197 theDef = theTableOfIons->GetIon(Z, A, eexc, noFloat, 0);
198 /*
199 G4cout << "### NC Find ion Z= " << Z << " A= " << A
200 << " Eexc(MeV)= " << eexc/MeV << " "
201 << theDef << G4endl;
202 */
203 }
204 G4double ekin = std::max(0.0,etot - theDef->GetPDGMass());
205 if (verboseLevel > 1) {
206 G4cout << i << ". " << theDef->GetParticleName()
207 << " Ekin(MeV)= " << etot/MeV
208 << " p: " << f->GetMomentum().vect()
209 << G4endl;
210 }
211 G4HadSecondary* news = new G4HadSecondary(
212 new G4DynamicParticle(theDef,
213 f->GetMomentum().vect().unit(),
214 ekin));
215 G4double timeF = f->GetCreationTime();
216 if(timeF < 0.0) { timeF = 0.0; }
217 news->SetTime(time + timeF);
218 news->SetCreatorModelID(secID);
219 theParticleChange.AddSecondary(*news);
220 delete news;
221 delete f;
222 }
223 delete fv;
224 }
225 //G4cout << "Capture done" << G4endl;
226 return &theParticleChange;
227}
std::vector< G4Fragment * > G4FragmentVector
Definition G4Fragment.hh:65
@ stopAndKill
#define noFloat
Definition G4Ions.hh:119
CLHEP::HepLorentzVector G4LorentzVector
#define M(row, col)
G4ThreeVector G4RandomDirection()
CLHEP::Hep3Vector G4ThreeVector
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
void set(double x, double y, double z)
Hep3Vector vect() const
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
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
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
static G4He3 * He3()
Definition G4He3.cc:90
static G4double GetNuclearMass(const G4double A, const G4double Z)
const G4String & GetParticleName() const
static G4Triton * Triton()
Definition G4Triton.cc:90

◆ InitialiseModel()

void G4NeutronRadCapture::InitialiseModel ( )
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 71 of file G4NeutronRadCapture.cc.

72{
73 if(photonEvaporation != nullptr) { return; }
74 G4DeexPrecoParameters* param =
76 minExcitation = param->GetMinExcitation();
78 photonEvaporation = new G4PhotonEvaporation();
79 photonEvaporation->Initialise();
80 photonEvaporation->SetICM(true);
81}
const G4String & GetModelName() const
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
static G4int GetModelID(const G4int modelIndex)

◆ operator=()

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

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