Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4LENDFission Class Reference

#include <G4LENDFission.hh>

Inheritance diagram for G4LENDFission:
G4LENDModel G4HadronicInteraction

Public Member Functions

 G4LENDFission (G4ParticleDefinition *pd)
 
 ~G4LENDFission ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
- Public Member Functions inherited from G4LENDModel
 G4LENDModel (G4String name="LENDModel")
 
 ~G4LENDModel ()
 
void ChangeDefaultEvaluation (G4String name)
 
void AllowNaturalAbundanceTarget ()
 
void AllowAnyCandidateTarget ()
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
- 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 &, 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
 
virtual G4int GetVerboseLevel () const
 
virtual 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, G4double
GetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 

Additional Inherited Members

- Protected Member Functions inherited from G4LENDModel
void create_used_target_map ()
 
void recreate_used_target_map ()
 
- 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 ( G4ParticleDefinition pd)
inline

Definition at line 51 of file G4LENDFission.hh.

References G4LENDModel::create_used_target_map(), and G4LENDModel::proj.

52  :G4LENDModel( "LENDFission" )
53  {
54  proj = pd;
55 
56 // theModelName = "LENDFission for ";
57 // theModelName += proj->GetParticleName();
59  };
G4ParticleDefinition * proj
Definition: G4LENDModel.hh:77
G4LENDModel(G4String name="LENDModel")
Definition: G4LENDModel.cc:45
void create_used_target_map()
Definition: G4LENDModel.cc:90
G4LENDFission::~G4LENDFission ( )
inline

Definition at line 61 of file G4LENDFission.hh.

61 {;};

Member Function Documentation

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

Reimplemented from G4LENDModel.

Definition at line 31 of file G4LENDFission.cc.

References G4HadFinalState::AddSecondary(), G4HadFinalState::Clear(), G4Gamma::Gamma(), G4Nucleus::GetA_asInt(), G4GIDI_target::getFissionFinalState(), G4Nucleus::GetIsotope(), G4HadProjectile::GetKineticEnergy(), G4Isotope::Getm(), G4HadProjectile::GetMaterial(), G4LENDManager::GetNucleusEncoding(), G4ParticleTable::GetParticleTable(), G4Material::GetTemperature(), G4Nucleus::GetZ_asInt(), int(), G4LENDModel::lend_manager, python.hepunit::MeV, G4Neutron::Neutron(), python.hepunit::second, G4DynamicParticle::SetDefinition(), G4DynamicParticle::SetMomentum(), G4HadFinalState::SetStatusChange(), stopAndKill, G4HadronicInteraction::theParticleChange, and G4LENDModel::usedTarget_map.

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 
47  G4HadFinalState* theResult = &theParticleChange;
48  theResult->Clear();
49 
50  G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA , iM ) )->second->GetTarget();
51  std::vector<G4GIDI_Product>* products = aTarget->getFissionFinalState( ke*MeV, temp, NULL, NULL );
52  if ( products != NULL )
53  {
54  for ( G4int j = 0; j < int( products->size() ); j++ )
55  {
56  G4int jZ = (*products)[j].Z;
57  G4int jA = (*products)[j].A;
58 
59  //G4cout << "Z = " << (*products)[j].Z
60  // << ", A = " << (*products)[j].A
61  // << ", EK = " << (*products)[j].kineticEnergy << " [MeV]"
62  // << ", px = " << (*products)[j].px
63  // << ", py = " << (*products)[j].py
64  // << ", pz = " << (*products)[j].pz
65  // << G4endl;
66 
68 
69  if ( jZ > 0 )
70  {
71  //Ex j?
72  theSec->SetDefinition( G4ParticleTable::GetParticleTable()->FindIon( jZ, jA , 0, 0 ) );
73  }
74  else if ( jA == 1 && jZ == 0 )
75  {
76  theSec->SetDefinition( G4Neutron::Neutron() );
77  }
78  else
79  {
80  theSec->SetDefinition( G4Gamma::Gamma() );
81  }
82 
83  theSec->SetMomentum( G4ThreeVector( (*products)[j].px*MeV , (*products)[j].py*MeV , (*products)[j].pz*MeV ) );
84  //G4cout << theSec->GetDefinition()->GetParticleName() << G4endl;
85  theResult->AddSecondary( theSec );
86  }
87  }
88  delete products;
89 
90  theResult->SetStatusChange( stopAndKill );
91 
92  return theResult;
93 
94 }
void SetMomentum(const G4ThreeVector &momentum)
CLHEP::Hep3Vector G4ThreeVector
typedef int(XMLCALL *XML_NotStandaloneHandler)(void *userData)
int G4int
Definition: G4Types.hh:78
void SetStatusChange(G4HadFinalStateStatus aS)
std::vector< G4GIDI_Product > * getFissionFinalState(double e_in, double temperature, double(*rng)(void *), void *rngState)
G4int GetNucleusEncoding(G4int iZ, G4int iA, G4int iM)
std::map< G4int, G4LENDUsedTarget * > usedTarget_map
Definition: G4LENDModel.hh:79
G4double GetKineticEnergy() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4ParticleTable * GetParticleTable()
G4double GetTemperature() const
Definition: G4Material.hh:180
const G4Material * GetMaterial() const
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4LENDManager * lend_manager
Definition: G4LENDModel.hh:78
void AddSecondary(G4DynamicParticle *aP)

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