Geant4 10.7.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 ()
 
virtual ~G4NeutronRadCapture ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) final
 
virtual void InitialiseModel () final
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
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 &)
 
virtual void InitialiseModel ()
 
 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()

G4NeutronRadCapture::G4NeutronRadCapture ( )
explicit

Definition at line 55 of file G4NeutronRadCapture.cc.

56 : G4HadronicInteraction("nRadCapture"),
57 photonEvaporation(nullptr),lab4mom(0.,0.,0.,0.)
58{
59 lowestEnergyLimit = 10*CLHEP::eV;
60 minExcitation = 0.1*CLHEP::keV;
61 SetMinEnergy( 0.0*CLHEP::GeV );
63
64 electron = G4Electron::Electron();
65 icID = -1;
66
68}
static G4Electron * Electron()
Definition: G4Electron.cc:93
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()

◆ ~G4NeutronRadCapture()

G4NeutronRadCapture::~G4NeutronRadCapture ( )
virtual

Definition at line 70 of file G4NeutronRadCapture.cc.

71{
72 delete photonEvaporation;
73}

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 88 of file G4NeutronRadCapture.cc.

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

◆ InitialiseModel()

void G4NeutronRadCapture::InitialiseModel ( )
finalvirtual

Reimplemented from G4HadronicInteraction.

Definition at line 75 of file G4NeutronRadCapture.cc.

76{
77 if(photonEvaporation != nullptr) { return; }
78 G4DeexPrecoParameters* param =
80 minExcitation = param->GetMinExcitation();
81 icID = param->GetInternalConversionID();
82
83 photonEvaporation = new G4PhotonEvaporation();
84 photonEvaporation->Initialise();
85 photonEvaporation->SetICM(true);
86}
G4double GetMinExcitation() const
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
virtual void SetICM(G4bool)

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