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

#include <G4LENDFission.hh>

+ Inheritance diagram for G4LENDFission:

Public Member Functions

 G4LENDFission (G4ParticleDefinition *pd)
 
 ~G4LENDFission ()
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
- Public Member Functions inherited from G4LENDModel
 G4LENDModel (G4String name="LENDModel")
 
 ~G4LENDModel ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
void ChangeDefaultEvaluation (G4String name)
 
void AllowNaturalAbundanceTarget ()
 
void AllowAnyCandidateTarget ()
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
void DumpLENDTargetInfo (G4bool force=false)
 
- 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 G4LENDModel
void create_used_target_map ()
 
void recreate_used_target_map ()
 
G4GIDI_targetget_target_from_map (G4int nuclear_code)
 
G4HadFinalStatereturnUnchanged (const G4HadProjectile &aTrack, G4HadFinalState *theResult)
 
- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4LENDModel
G4ParticleDefinitionproj
 
G4LENDManagerlend_manager
 
std::map< G4int, G4LENDUsedTarget * > usedTarget_map
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 46 of file G4LENDFission.hh.

Constructor & Destructor Documentation

◆ G4LENDFission()

G4LENDFission::G4LENDFission ( G4ParticleDefinition pd)
inline

Definition at line 51 of file G4LENDFission.hh.

52 :G4LENDModel( "LENDFission" )
53 {
54 proj = pd;
55
56// theModelName = "LENDFission for ";
57// theModelName += proj->GetParticleName();
59 };
void create_used_target_map()
Definition: G4LENDModel.cc:93
G4ParticleDefinition * proj
Definition: G4LENDModel.hh:83

◆ ~G4LENDFission()

G4LENDFission::~G4LENDFission ( )
inline

Definition at line 61 of file G4LENDFission.hh.

61{;};

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4LENDFission::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus aTargetNucleus 
)
virtual

Reimplemented from G4LENDModel.

Definition at line 31 of file G4LENDFission.cc.

32{
33
34 G4double temp = aTrack.GetMaterial()->GetTemperature();
35
36 //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
37 G4int iZ = aTarg.GetZ_asInt();
38 G4int iA = aTarg.GetA_asInt();
39 //G4int iM = aTarg.GetM_asInt();
40 G4int iM = 0;
41 if ( aTarg.GetIsotope() != NULL ) {
42 iM = aTarg.GetIsotope()->Getm();
43 }
44
45 G4double ke = aTrack.GetKineticEnergy();
46
48 theResult->Clear();
49
51 if ( aTarget == NULL ) return returnUnchanged( aTrack , theResult );
52 std::vector<G4GIDI_Product>* products = aTarget->getFissionFinalState( ke*MeV, temp, MyRNG, NULL );
53 if ( products != NULL )
54 {
55 for ( G4int j = 0; j < int( products->size() ); j++ )
56 {
57 G4int jZ = (*products)[j].Z;
58 G4int jA = (*products)[j].A;
59 G4int jM = (*products)[j].m;
60
61 //G4cout << "Z = " << (*products)[j].Z
62 // << ", A = " << (*products)[j].A
63 // << ", EK = " << (*products)[j].kineticEnergy << " [MeV]"
64 // << ", px = " << (*products)[j].px
65 // << ", py = " << (*products)[j].py
66 // << ", pz = " << (*products)[j].pz
67 // << ", birthTimeSec = " << (*products)[j].birthTimeSec << " [second]"
68 // << G4endl;
69
71
72 if ( jZ > 0 )
73 {
74 theSec->SetDefinition( G4IonTable::GetIonTable()->GetIon( jZ, jA , jM ) );
75 }
76 else if ( jA == 1 && jZ == 0 )
77 {
79 }
80 else
81 {
82 theSec->SetDefinition( G4Gamma::Gamma() );
83 }
84
85 theSec->SetMomentum( G4ThreeVector( (*products)[j].px*MeV , (*products)[j].py*MeV , (*products)[j].pz*MeV ) );
86 //G4cout << theSec->GetDefinition()->GetParticleName() << G4endl;
87 theResult->AddSecondary( theSec );
88 //Set time for delayed neutrons
89 //Current implementation is a little tricky,
90 if ( (*products)[j].birthTimeSec != 0 ) {
91 G4double time = (*products)[j].birthTimeSec*second + aTrack.GetGlobalTime();
92 theResult->GetSecondary(theResult->GetNumberOfSecondaries()-1)->SetTime(time);
93 }
94 }
95 }
96 delete products;
97
98 theResult->SetStatusChange( stopAndKill );
99
100 return theResult;
101
102}
@ stopAndKill
double MyRNG(void *)
Definition: G4LENDModel.cc:45
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4ThreeVector &momentum)
std::vector< G4GIDI_Product > * getFissionFinalState(double e_in, double temperature, double(*rng)(void *), void *rngState)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
const G4Material * GetMaterial() const
G4double GetKineticEnergy() const
G4double GetGlobalTime() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
static G4IonTable * GetIonTable()
Definition: G4IonTable.cc:170
G4int GetNucleusEncoding(G4int iZ, G4int iA, G4int iM)
G4LENDManager * lend_manager
Definition: G4LENDModel.hh:84
G4HadFinalState * returnUnchanged(const G4HadProjectile &aTrack, G4HadFinalState *theResult)
Definition: G4LENDModel.cc:253
G4GIDI_target * get_target_from_map(G4int nuclear_code)
Definition: G4LENDModel.cc:267
G4double GetTemperature() const
Definition: G4Material.hh:180
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103

◆ GetFatalEnergyCheckLevels()

const std::pair< G4double, G4double > G4LENDFission::GetFatalEnergyCheckLevels ( ) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 103 of file G4LENDFission.cc.

104{
105 // max energy non-conservation is mass of heavy nucleus
106 //return std::pair<G4double, G4double>(5*perCent,250*GeV);
107 return std::pair<G4double, G4double>(5*perCent,DBL_MAX);
108}
#define DBL_MAX
Definition: templates.hh:62

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