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

#include <G4NeutronRadCaptureHP.hh>

+ Inheritance diagram for G4NeutronRadCaptureHP:

Public Member Functions

 G4NeutronRadCaptureHP ()
 
 ~G4NeutronRadCaptureHP () override
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4NeutronRadCaptureHPoperator= (const G4NeutronRadCaptureHP &right)=delete
 
 G4NeutronRadCaptureHP (const G4NeutronRadCaptureHP &)=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 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 54 of file G4NeutronRadCaptureHP.hh.

Constructor & Destructor Documentation

◆ G4NeutronRadCaptureHP() [1/2]

G4NeutronRadCaptureHP::G4NeutronRadCaptureHP ( )

Definition at line 60 of file G4NeutronRadCaptureHP.cc.

61 : G4HadronicInteraction("nRadCaptureHP"),
62 electron(G4Electron::Electron()),
64 lowestEnergyLimit(1.0e-11*CLHEP::eV),
65 minExcitation(0.1*CLHEP::keV),
66 emax(20*CLHEP::MeV),
67 emaxT(fManagerHP->GetMaxEnergyDoppler()),
68 lab4mom(0.,0.,0.,0.)
69{
72}
static G4Electron * Electron()
Definition G4Electron.cc:91
G4HadronicInteraction(const G4String &modelName="HadronicModel")
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
static G4ParticleHPManager * GetInstance()
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()

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

◆ ~G4NeutronRadCaptureHP()

G4NeutronRadCaptureHP::~G4NeutronRadCaptureHP ( )
override

Definition at line 74 of file G4NeutronRadCaptureHP.cc.

75{
76 if (fLocalPE) { delete photonEvaporation; }
77}

◆ G4NeutronRadCaptureHP() [2/2]

G4NeutronRadCaptureHP::G4NeutronRadCaptureHP ( const G4NeutronRadCaptureHP & )
delete

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 103 of file G4NeutronRadCaptureHP.cc.

105{
106 theParticleChange.Clear();
107 G4double ekin = aTrack.GetKineticEnergy();
108 if (ekin > emax) {
109 return &theParticleChange;
110 }
111
112 theParticleChange.SetStatusChange(stopAndKill);
113 G4double T = aTrack.GetMaterial()->GetTemperature();
114
115 G4int A = theNucleus.GetA_asInt();
116 G4int Z = theNucleus.GetZ_asInt();
117
118 G4double time = aTrack.GetGlobalTime();
119
120 // Create initial state
122
123 // no Doppler broading
124 G4double factT = T/CLHEP::STP_Temperature;
125
126 if (ekin >= emaxT*factT || fManagerHP->GetNeglectDoppler()) {
127 lab4mom.set(0.,0.,0.,mass);
128 } else {
129 G4double lambda = 1.0/(CLHEP::k_Boltzmann*T);
130 G4double erand = G4RandGamma::shoot(2.0, lambda);
131 auto mom = G4RandomDirection()*std::sqrt(2*mass*erand);
132 lab4mom.set(mom.x(), mom.y(), mom.z(), mass + erand);
133 }
134
135 lab4mom += aTrack.Get4Momentum();
136
137 G4double M = lab4mom.mag();
138 ++A;
140 //G4cout << "Capture start: Z= " << Z << " A= " << A
141 // << " LabM= " << M << " Mcompound= " << mass << G4endl;
142
143 // simplified method of 1 gamma emission
144 if (A <= 4) {
145
146 if (verboseLevel > 1) {
147 G4cout << "G4NeutronRadCaptureHP::DoIt: Eini(MeV)="
148 << ekin/MeV << " Eexc(MeV)= "
149 << (M - mass)/MeV
150 << " Z= " << Z << " A= " << A << G4endl;
151 }
152 if (M - mass > lowestEnergyLimit) {
153 G4ThreeVector bst = lab4mom.boostVector();
154 G4double e1 = (M - mass)*(M + mass)/(2*M);
156 lv2.boost(bst);
157 if (verboseLevel > 1) {
158 G4cout << "Gamma 4-mom: " << lv2 << " Escm(MeV)=" << e1/CLHEP::MeV << G4endl;
159 }
160 lab4mom -= lv2;
161 G4HadSecondary* news =
162 new G4HadSecondary(new G4DynamicParticle(G4Gamma::Gamma(), lv2));
163 news->SetTime(time);
164 news->SetCreatorModelID(secID);
165 theParticleChange.AddSecondary(*news);
166 delete news;
167 }
168
169 const G4ParticleDefinition* theDef = nullptr;
170
171 if (Z == 1 && A == 2) {theDef = G4Deuteron::Deuteron(); }
172 else if (Z == 1 && A == 3) {theDef = G4Triton::Triton(); }
173 else if (Z == 2 && A == 3) {theDef = G4He3::He3(); }
174 else if (Z == 2 && A == 4) {theDef = G4Alpha::Alpha(); }
175 else { theDef = theTableOfIons->GetIon(Z, A, 0.0, noFloat, 0); }
176
177 if (nullptr != theDef) {
178 G4HadSecondary* news =
179 new G4HadSecondary(new G4DynamicParticle(theDef, lab4mom));
180 news->SetTime(time);
181 news->SetCreatorModelID(secID);
182 theParticleChange.AddSecondary(*news);
183 delete news;
184 }
185
186 // Use photon evaporation
187 } else {
188
189 // protection against wrong kinematic
190 if (M < mass) {
191 G4double etot = std::max(mass, lab4mom.e());
192 G4double ptot = std::sqrt((etot - mass)*(etot + mass));
193 G4ThreeVector v = lab4mom.vect().unit();
194 lab4mom.set(v.x()*ptot,v.y()*ptot,v.z()*ptot,etot);
195 }
196
197 G4Fragment* aFragment = new G4Fragment(A, Z, lab4mom);
198
199 if (verboseLevel > 1) {
200 G4cout << "G4NeutronRadCaptureHP::ApplyYourself initial G4Fragmet:"
201 << G4endl;
202 G4cout << aFragment << G4endl;
203 }
204
205 //
206 // Sample final state
207 //
208 G4FragmentVector* fv = photonEvaporation->BreakUpFragment(aFragment);
209 if (nullptr == fv) { fv = new G4FragmentVector(); }
210 fv->push_back(aFragment);
211
212 if (verboseLevel > 1) {
213 G4cout << "G4NeutronRadCaptureHP: " << fv->size() << " final particle icID= "
214 << icID << G4endl;
215 }
216 for (auto const & f : *fv) {
217 G4double etot = f->GetMomentum().e();
218
219 Z = f->GetZ_asInt();
220 A = f->GetA_asInt();
221
222 const G4ParticleDefinition* theDef;
223 if (0 == Z && 0 == A) { theDef = f->GetParticleDefinition(); }
224 else if (Z == 1 && A == 2) { theDef = G4Deuteron::Deuteron(); }
225 else if (Z == 1 && A == 3) { theDef = G4Triton::Triton(); }
226 else if (Z == 2 && A == 3) { theDef = G4He3::He3(); }
227 else if (Z == 2 && A == 4) { theDef = G4Alpha::Alpha(); }
228 else {
229 G4double eexc = f->GetExcitationEnergy();
230 if (eexc <= minExcitation) { eexc = 0.0; }
231 theDef = theTableOfIons->GetIon(Z, A, eexc, noFloat, 0);
232 /*
233 G4cout << "### NC Find ion Z= " << Z << " A= " << A
234 << " Eexc(MeV)= " << eexc/MeV << " "
235 << theDef << G4endl;
236 */
237 }
238 ekin = std::max(0.0, etot - theDef->GetPDGMass());
239 if (verboseLevel > 1) {
240 G4cout << theDef->GetParticleName()
241 << " Ekin(MeV)= " << ekin/MeV
242 << " p: " << f->GetMomentum().vect()
243 << G4endl;
244 }
245 G4HadSecondary* news =
246 new G4HadSecondary(new G4DynamicParticle(theDef,
247 f->GetMomentum().vect().unit(),
248 ekin));
249 G4double timeF = std::max(f->GetCreationTime(), 0.0);
250 news->SetTime(time + timeF);
251 if (theDef == electron) {
252 news->SetCreatorModelID(icID);
253 } else {
254 news->SetCreatorModelID(secID);
255 }
256 theParticleChange.AddSecondary(*news);
257 delete news;
258 delete f;
259 }
260 delete fv;
261 }
262 //G4cout << "Capture done" << G4endl;
263 return &theParticleChange;
264}
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)
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
const G4Material * GetMaterial() const
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
G4double GetTemperature() const
static G4double GetNuclearMass(const G4double A, const G4double Z)
const G4String & GetParticleName() const
static G4Triton * Triton()
Definition G4Triton.cc:90

◆ BuildPhysicsTable()

void G4NeutronRadCaptureHP::BuildPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 79 of file G4NeutronRadCaptureHP.cc.

80{
81 if (photonEvaporation != nullptr) { return; }
84 if (nullptr != p) {
85 auto handler =
86 (static_cast<G4VPreCompoundModel*>(p))->GetExcitationHandler();
87 if (nullptr != handler)
88 photonEvaporation = handler->GetPhotonEvaporation();
89 }
90 G4DeexPrecoParameters* param =
92 minExcitation = param->GetMinExcitation();
93 icID = G4PhysicsModelCatalog::GetModelID("model_e-InternalConversion");
95 if (nullptr == photonEvaporation) {
96 photonEvaporation = new G4PhotonEvaporation();
97 fLocalPE = true;
98 }
99 photonEvaporation->Initialise();
100 photonEvaporation->SetICM(true);
101}
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
const G4String & GetModelName() const
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
static G4int GetModelID(const G4int modelIndex)

◆ operator=()

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

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