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

#include <G4NeutronFissionVI.hh>

+ Inheritance diagram for G4NeutronFissionVI:

Public Member Functions

 G4NeutronFissionVI ()
 
 ~G4NeutronFissionVI () override
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
 
void InitialiseModel () override
 
G4NeutronFissionVIoperator= (const G4NeutronFissionVI &right)=delete
 
 G4NeutronFissionVI (const G4NeutronFissionVI &)=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 53 of file G4NeutronFissionVI.hh.

Constructor & Destructor Documentation

◆ G4NeutronFissionVI() [1/2]

G4NeutronFissionVI::G4NeutronFissionVI ( )

Definition at line 59 of file G4NeutronFissionVI.cc.

60 : G4HadronicInteraction("nFissionVI"),
62 minExcitation(0.1*CLHEP::keV),
63 emaxT(fManagerHP->GetMaxEnergyDoppler()),
64 lab4mom(0.,0.,0.,0.)
65{
66 SetMinEnergy( 0.0*CLHEP::GeV );
69}
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
G4double GetMaxEnergyDoppler() const
static G4ParticleHPManager * GetInstance()
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()

◆ ~G4NeutronFissionVI()

G4NeutronFissionVI::~G4NeutronFissionVI ( )
override

Definition at line 71 of file G4NeutronFissionVI.cc.

72{
73 if (fLocalHandler) { delete fHandler; }
74}

◆ G4NeutronFissionVI() [2/2]

G4NeutronFissionVI::G4NeutronFissionVI ( const G4NeutronFissionVI & )
delete

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 97 of file G4NeutronFissionVI.cc.

99{
102 G4double T = aTrack.GetMaterial()->GetTemperature();
103
104 G4int A = theNucleus.GetA_asInt();
105 G4int Z = theNucleus.GetZ_asInt();
106
107 G4double time = aTrack.GetGlobalTime();
108 G4double ekin = aTrack.GetKineticEnergy();
109
110 // Create initial state
112
113 // no Doppler broading
114 G4double factT = T/CLHEP::STP_Temperature;
115
116 if (ekin >= emaxT*factT || fManagerHP->GetNeglectDoppler()) {
117 lab4mom.set(0.,0.,0.,mass);
118
119 } else {
120 G4double lambda = 1.0/(CLHEP::k_Boltzmann*T);
121 G4double erand = G4RandGamma::shoot(2.0, lambda);
122 auto mom = G4RandomDirection()*std::sqrt(2*mass*erand);
123 lab4mom.set(mom.x(), mom.y(), mom.z(), mass + erand);
124 }
125
126 lab4mom += aTrack.Get4Momentum();
127
128 G4double M = lab4mom.mag();
129 ++A;
131 //G4cout << "Fission start: Z= " << Z << " A= " << A
132 // << " LabM= " << M << " Mcompound= " << mass << G4endl;
133
134 // protection against wrong kinematic
135 if (M < mass) {
136 G4double etot = std::max(mass, lab4mom.e());
137 G4double ptot = std::sqrt((etot - mass)*(etot + mass));
138 G4ThreeVector v = lab4mom.vect().unit();
139 lab4mom.set(v.x()*ptot,v.y()*ptot,v.z()*ptot,etot);
140 }
141
142 G4Fragment* aFragment = new G4Fragment(A, Z, lab4mom);
143
144 if (verboseLevel > 1) {
145 G4cout << "G4NeutronFissionVI::ApplyYourself initial G4Fragmet:"
146 << G4endl;
147 G4cout << aFragment << G4endl;
148 }
149
150 //
151 // Sample final state
152 //
153 fFission->GetEmissionProbability(aFragment);
154 G4Fragment* frag = fFission->EmittedFragment(aFragment);
155 G4ReactionProductVector* final = fHandler->BreakItUp(*aFragment);
156 if (nullptr != frag) {
157 G4ReactionProductVector* v = fHandler->BreakItUp(*frag);
158 for (auto & p : *v) {
159 final->push_back(p);
160 }
161 delete v;
162 delete frag;
163 }
164
165 if (verboseLevel > 1) {
166 G4cout << "G4NeutronFissionVI: " << final->size()
167 << " final particle secID= " << secID << G4endl;
168 }
169 for (auto const & ptr : *final) {
170 G4double etot = ptr->GetTotalEnergy();
171 const G4ParticleDefinition* theDef = ptr->GetDefinition();
172 ekin = std::max(0.0, etot - theDef->GetPDGMass());
173 if (verboseLevel > 1) {
174 G4cout << theDef->GetParticleName()
175 << " Ekin(MeV)= " << ekin/MeV
176 << " p: " << ptr->GetMomentum()
177 << G4endl;
178 }
179 G4HadSecondary* news = new G4HadSecondary(
180 new G4DynamicParticle(theDef, ptr->GetMomentum().unit(), ekin));
181 G4double timeF = std::max(ptr->GetFormationTime(), 0.0);
182 news->SetTime(time + timeF);
183 news->SetCreatorModelID(secID);
185 delete news;
186 delete ptr;
187 }
188 delete final;
189 delete aFragment;
190
191 //G4cout << "Fission done" << G4endl;
192 return &theParticleChange;
193}
@ stopAndKill
#define M(row, col)
G4ThreeVector G4RandomDirection()
std::vector< G4ReactionProduct * > G4ReactionProductVector
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 vect() const
void set(double x, double y, double z, double t)
G4ReactionProductVector * BreakItUp(const G4Fragment &theInitialState)
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
const G4Material * GetMaterial() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetGlobalTime() const
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
G4double GetTemperature() const
static G4double GetNuclearMass(const G4double A, const G4double Z)
const G4String & GetParticleName() const
G4bool GetNeglectDoppler() const
virtual G4Fragment * EmittedFragment(G4Fragment *theNucleus)
virtual G4double GetEmissionProbability(G4Fragment *theNucleus)=0

◆ InitialiseModel()

void G4NeutronFissionVI::InitialiseModel ( )
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 76 of file G4NeutronFissionVI.cc.

77{
78 if (fFission != nullptr && fHandler != nullptr) { return; }
81 if (nullptr != p) {
82 fHandler = (static_cast<G4VPreCompoundModel*>(p))->GetExcitationHandler();
83 }
84 if (nullptr == fHandler) {
85 fHandler = new G4ExcitationHandler();
86 fLocalHandler = true;
87 }
88 fHandler->Initialise();
89 fFission = fHandler->GetEvaporation()->GetFissionChannel();
90
91 G4DeexPrecoParameters* param =
93 minExcitation = param->GetMinExcitation();
95}
G4VEvaporation * GetEvaporation()
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
const G4String & GetModelName() const
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
static G4int GetModelID(const G4int modelIndex)
G4VEvaporationChannel * GetFissionChannel() const

◆ operator=()

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

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