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

#include <G4LENDInelastic.hh>

+ Inheritance diagram for G4LENDInelastic:

Public Member Functions

 G4LENDInelastic (G4ParticleDefinition *pd)
 
 ~G4LENDInelastic ()
 
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 &)
 
- 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 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 44 of file G4LENDInelastic.hh.

Constructor & Destructor Documentation

◆ G4LENDInelastic()

G4LENDInelastic::G4LENDInelastic ( G4ParticleDefinition pd)
inline

Definition at line 49 of file G4LENDInelastic.hh.

50 :G4LENDModel( "LENDInelastic" )
51 {
52 proj = pd;
53
54// theModelName = "LENDInelastic for ";
55// theModelName += proj->GetParticleName();
57 };
void create_used_target_map()
Definition: G4LENDModel.cc:86
G4ParticleDefinition * proj
Definition: G4LENDModel.hh:77

◆ ~G4LENDInelastic()

G4LENDInelastic::~G4LENDInelastic ( )
inline

Definition at line 59 of file G4LENDInelastic.hh.

59{;};

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4LENDModel.

Definition at line 31 of file G4LENDInelastic.cc.

32{
33
34 G4ThreeVector proj_p = aTrack.Get4Momentum().vect();
35
36 G4double temp = aTrack.GetMaterial()->GetTemperature();
37
38 //G4int iZ = int ( aTarg.GetZ() );
39 //G4int iA = int ( aTarg.GetN() );
40 //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
41 G4int iZ = aTarg.GetZ_asInt();
42 G4int iA = aTarg.GetA_asInt();
43 //G4cout << "target: Z = " << iZ << " N = " << iA << G4endl;
44
45 G4double ke = aTrack.GetKineticEnergy();
46 //G4cout << "projectile: KE = " << ke/MeV << " [MeV]" << G4endl;
47
49 theResult->Clear();
50
51 G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA ) )->second->GetTarget();
52 std::vector<G4GIDI_Product>* products = aTarget->getOthersFinalState( ke*MeV, temp, NULL, NULL );
53 if ( products != NULL )
54 {
55
56 G4ThreeVector psum(0);
57
58 int totN = 0;
59 for ( G4int j = 0; j < int( products->size() ); j++ )
60 {
61
62 G4int jZ = (*products)[j].Z;
63 G4int jA = (*products)[j].A;
64
65 //G4cout << "ZA = " << 1000 * (*products)[j].Z + (*products)[j].A << " EK = "
66 // << (*products)[j].kineticEnergy
67 // << " px " << (*products)[j].px
68 // << " py " << (*products)[j].py
69 // << " pz " << (*products)[j].pz
70 // << G4endl;
71
73
74 if ( jA == 1 && jZ == 1 )
75 {
77 totN += 1;
78 }
79 else if ( jA == 1 && jZ == 0 )
80 {
82 totN += 1;
83 }
84 else if ( jZ > 0 )
85 {
86 if ( jA != 0 )
87 {
88 theSec->SetDefinition( G4ParticleTable::GetParticleTable()->FindIon( jZ , jA , 0 , 0 ) );
89 totN += jA;
90 }
91 else
92 {
93 theSec->SetDefinition( G4ParticleTable::GetParticleTable()->FindIon( jZ , iA+1-totN , 0 , 0 ) );
94 }
95 }
96 else
97 {
98 theSec->SetDefinition( G4Gamma::Gamma() );
99 }
100
101 G4ThreeVector p( (*products)[j].px*MeV , (*products)[j].py*MeV , (*products)[j].pz*MeV );
102 psum += p;
103 if ( p.mag() == 0 ) p = proj_p - psum;
104
105 theSec->SetMomentum( p );
106
107 theResult->AddSecondary( theSec );
108 }
109 }
110 delete products;
111
112 theResult->SetStatusChange( stopAndKill );
113
114 return theResult;
115
116}
@ stopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
Hep3Vector vect() const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4ThreeVector &momentum)
std::vector< G4GIDI_Product > * getOthersFinalState(double e_in, double temperature, double(*rng)(void *), void *rngState)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP)
const G4Material * GetMaterial() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4int GetNucleusEncoding(G4int iZ, G4int iA)
std::map< G4int, G4LENDUsedTarget * > usedTarget_map
Definition: G4LENDModel.hh:79
G4LENDManager * lend_manager
Definition: G4LENDModel.hh:78
G4double GetTemperature() const
Definition: G4Material.hh:181
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4ParticleDefinition * FindIon(G4int atomicNumber, G4int atomicMass, G4double excitationEnergy)
static G4ParticleTable * GetParticleTable()
static G4Proton * Proton()
Definition: G4Proton.cc:93

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