Geant4 9.6.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 ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
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)
 
const G4HadronicInteractionGetMyPointer () const
 
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
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) 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
 

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 50 of file G4NeutronRadCapture.hh.

Constructor & Destructor Documentation

◆ G4NeutronRadCapture()

G4NeutronRadCapture::G4NeutronRadCapture ( )

Definition at line 53 of file G4NeutronRadCapture.cc.

54 : G4HadronicInteraction("nRadCapture")
55{
56 lowestEnergyLimit = 0.1*eV;
57 SetMinEnergy( 0.0*GeV );
58 SetMaxEnergy( 100.*TeV );
59 photonEvaporation = new G4PhotonEvaporation();
60 //photonEvaporation = 0;
61}
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)

◆ ~G4NeutronRadCapture()

G4NeutronRadCapture::~G4NeutronRadCapture ( )
virtual

Definition at line 63 of file G4NeutronRadCapture.cc.

64{
65 delete photonEvaporation;
66}

Member Function Documentation

◆ ApplyYourself()

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

Implements G4HadronicInteraction.

Definition at line 68 of file G4NeutronRadCapture.cc.

70{
73
74 G4int A = theNucleus.GetA_asInt();
75 G4int Z = theNucleus.GetZ_asInt();
76
77 // Create initial state
79 G4LorentzVector lv0(0.0,0.0,0.0,m1);
80 G4LorentzVector lv1 = aTrack.Get4Momentum() + lv0;
81
82 // simplified method of 1 gamma emission
83 if(A <= 3) {
84
85 G4ThreeVector bst = lv1.boostVector();
86 G4double M = lv1.mag();
87
88 ++A;
90 if(M - mass <= lowestEnergyLimit) {
91 return &theParticleChange;
92 }
93
94 if (verboseLevel > 1) {
95 G4cout << "G4NeutronRadCapture::DoIt: Eini(MeV)="
96 << aTrack.GetKineticEnergy()/MeV << " Eexc(MeV)= "
97 << (M - mass)/MeV
98 << " Z= " << Z << " A= " << A << G4endl;
99 }
100 G4double e1 = (M - mass)*(M + mass)/(2*M);
101 G4double cost = 2.0*G4UniformRand() - 1.0;
102 if(cost > 1.0) {cost = 1.0;}
103 else if(cost < -1.0) {cost = -1.0;}
104 G4double sint = std::sqrt((1. - cost)*(1.0 + cost));
105 G4double phi = G4UniformRand()*CLHEP::twopi;
106 G4LorentzVector lv2(e1*sint*std::cos(phi),e1*sint*std::sin(phi),e1*cost,e1);
107 lv2.boost(bst);
109 G4ParticleDefinition* theDef = 0;
110
111 lv1 -= lv2;
112 if (Z == 1 && A == 2) {theDef = G4Deuteron::Deuteron();}
113 else if (Z == 1 && A == 3) {theDef = G4Triton::Triton();}
114 else if (Z == 2 && A == 3) {theDef = G4He3::He3();}
115 else if (Z == 2 && A == 4) {theDef = G4Alpha::Alpha();}
116 else
117 {
118 theDef =
120 }
121
122 if (verboseLevel > 1) {
123 G4cout << "Gamma 4-mom: " << lv2 << " "
124 << theDef->GetParticleName() << " " << lv1 << G4endl;
125 }
126 if(theDef) {
128 }
129
130 // Use photon evaporation
131 } else {
132
133 G4Fragment* aFragment = new G4Fragment(A+1, Z, lv1);
134
135 if (verboseLevel > 1) {
136 G4cout << "G4NeutronRadCapture::ApplyYourself initial G4Fragmet:" << G4endl;
137 G4cout << aFragment << G4endl;
138 }
139
140 //
141 // Sample final state
142 //
143 G4FragmentVector* fv = photonEvaporation->BreakUpFragment(aFragment);
144 if(!fv) { fv = new G4FragmentVector(); }
145 fv->push_back(aFragment);
146 size_t n = fv->size();
147
148 if (verboseLevel > 1) {
149 G4cout << "G4NeutronRadCapture: " << n << " final particle" << G4endl;
150 }
151 for(size_t i=0; i<n; ++i) {
152 G4Fragment* f = (*fv)[i];
153 G4LorentzVector lvres = f->GetMomentum();
154 Z = f->GetZ_asInt();
155 A = f->GetA_asInt();
156
157 G4ParticleDefinition* theDef = 0;
158 if(0 == Z && 0 == A) {theDef = f->GetParticleDefinition();}
159 else
160 {
162 }
163
164 if (verboseLevel > 1) {
165 G4cout << i << ". " << theDef->GetParticleName()
166 << " " << lvres << G4endl;
167 }
168 if(theDef) {
170 }
171 delete f;
172 }
173 delete fv;
174 }
175 return &theParticleChange;
176}
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
@ stopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector boostVector() const
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:251
G4int GetZ_asInt() const
Definition: G4Fragment.hh:223
G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:368
G4int GetA_asInt() const
Definition: G4Fragment.hh:218
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP)
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
static G4He3 * He3()
Definition: G4He3.cc:94
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int J=0)
Definition: G4IonTable.cc:267
static G4double GetNuclearMass(const G4double A, const G4double Z)
const G4String & GetParticleName() const
static G4ParticleTable * GetParticleTable()
G4IonTable * GetIonTable()
virtual G4FragmentVector * BreakUpFragment(G4Fragment *theNucleus)
static G4Triton * Triton()
Definition: G4Triton.cc:95

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