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

#include <G4LENDElastic.hh>

+ Inheritance diagram for G4LENDElastic:

Public Member Functions

 G4LENDElastic (G4ParticleDefinition *pd)
 
 ~G4LENDElastic ()
 
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 &)
 
void DumpLENDTargetInfo (G4bool force=false)
 
- 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 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
 
G4int secID
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 44 of file G4LENDElastic.hh.

Constructor & Destructor Documentation

◆ G4LENDElastic()

G4LENDElastic::G4LENDElastic ( G4ParticleDefinition * pd)
inline

Definition at line 49 of file G4LENDElastic.hh.

50 :G4LENDModel( "LENDElastic" )
51 {
52 proj = pd;
53
54 //theModelName = "LEND Elastic Model for ";
55 //theModelName += proj->GetParticleName();
57 };
G4LENDModel(G4String name="LENDModel")
void create_used_target_map()
G4ParticleDefinition * proj

◆ ~G4LENDElastic()

G4LENDElastic::~G4LENDElastic ( )
inline

Definition at line 59 of file G4LENDElastic.hh.

59{;};

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4LENDModel.

Definition at line 37 of file G4LENDElastic.cc.

38{
39
40 G4double temp = aTrack.GetMaterial()->GetTemperature();
41
42 //G4int iZ = int ( aTarg.GetZ() );
43 //G4int iA = int ( aTarg.GetN() );
44 //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
45 G4int iZ = aTarg.GetZ_asInt();
46 G4int iA = aTarg.GetA_asInt();
47 G4int iM = 0;
48 if ( aTarg.GetIsotope() != NULL ) {
49 iM = aTarg.GetIsotope()->Getm();
50 }
51
52 G4double ke = aTrack.GetKineticEnergy();
53
55 theResult->Clear();
56
58 if ( aTarget == NULL ) return returnUnchanged( aTrack , theResult );
59
60 G4double aMu = aTarget->getElasticFinalState( ke*MeV, temp, MyRNG , NULL );
61
62 G4double phi = twopi*G4UniformRand();
63 G4double theta = std::acos( aMu );
64 //G4double sinth = std::sin( theta );
65
66 G4ReactionProduct theNeutron( aTrack.GetDefinition() );
67 theNeutron.SetMomentum( aTrack.Get4Momentum().vect() );
68 theNeutron.SetKineticEnergy( ke );
69
70 //G4ParticleDefinition* pd = G4IonTable::GetIonTable()->GetIon( iZ , iA , iM );
71 //TK 170509 Fix for the case of excited isomer target
72 G4double EE = 0.0;
73 if ( iM != 0 ) {
75 }
76 G4ParticleDefinition* target_pd = G4IonTable::GetIonTable()->GetIon( iZ , iA , EE );
77 G4ReactionProduct theTarget( target_pd );
78
79 G4double mass = target_pd->GetPDGMass();
80
81// add Thermal motion
82 G4double kT = k_Boltzmann*temp;
83 G4ThreeVector v ( G4RandGauss::shoot() * std::sqrt( kT*mass )
84 , G4RandGauss::shoot() * std::sqrt( kT*mass )
85 , G4RandGauss::shoot() * std::sqrt( kT*mass ) );
86 theTarget.SetMomentum( v );
87
88 G4ThreeVector the3Neutron = theNeutron.GetMomentum();
89 G4double nEnergy = theNeutron.GetTotalEnergy();
90 G4ThreeVector the3Target = theTarget.GetMomentum();
91 G4double tEnergy = theTarget.GetTotalEnergy();
92 G4ReactionProduct theCMS;
93 G4double totE = nEnergy+tEnergy;
94 G4ThreeVector the3CMS = the3Target+the3Neutron;
95 theCMS.SetMomentum(the3CMS);
96 G4double cmsMom = std::sqrt(the3CMS*the3CMS);
97 G4double sqrts = std::sqrt((totE-cmsMom)*(totE+cmsMom));
98 theCMS.SetMass(sqrts);
99 theCMS.SetTotalEnergy(totE);
100
101 theNeutron.Lorentz(theNeutron, theCMS);
102 theTarget.Lorentz(theTarget, theCMS);
103 G4double en = theNeutron.GetTotalMomentum(); // already in CMS.
104 G4ThreeVector cms3Mom=theNeutron.GetMomentum(); // for neutron direction in CMS
105 G4double cms_theta=cms3Mom.theta();
106 G4double cms_phi=cms3Mom.phi();
107 G4ThreeVector tempVector;
108 tempVector.setX( std::cos(theta)*std::sin(cms_theta)*std::cos(cms_phi)
109 +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::cos(cms_phi)
110 -std::sin(theta)*std::sin(phi)*std::sin(cms_phi) );
111 tempVector.setY( std::cos(theta)*std::sin(cms_theta)*std::sin(cms_phi)
112 +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::sin(cms_phi)
113 +std::sin(theta)*std::sin(phi)*std::cos(cms_phi) );
114 tempVector.setZ( std::cos(theta)*std::cos(cms_theta)
115 -std::sin(theta)*std::cos(phi)*std::sin(cms_theta) );
116 tempVector *= en;
117 theNeutron.SetMomentum(tempVector);
118 theTarget.SetMomentum(-tempVector);
119 G4double tP = theTarget.GetTotalMomentum();
120 G4double tM = theTarget.GetMass();
121 theTarget.SetTotalEnergy(std::sqrt((tP+tM)*(tP+tM)-2.*tP*tM));
122
123
124 theNeutron.Lorentz(theNeutron, -1.*theCMS);
125 theTarget.Lorentz(theTarget, -1.*theCMS);
126
127//110913 Add Protection for very low energy (1e-6eV) scattering
128 if ( theNeutron.GetKineticEnergy() <= 0 )
129 {
130 theNeutron.SetTotalEnergy ( theNeutron.GetMass() * ( 1 + G4Pow::GetInstance()->powA( 10 , -15.65 ) ) );
131 }
132
133 if ( theTarget.GetKineticEnergy() < 0 )
134 {
135 theTarget.SetTotalEnergy ( theTarget.GetMass() * ( 1 + G4Pow::GetInstance()->powA( 10 , -15.65 ) ) );
136 }
137//110913 END
138
139 theResult->SetEnergyChange(theNeutron.GetKineticEnergy());
140 theResult->SetMomentumChange(theNeutron.GetMomentum().unit());
141 G4DynamicParticle* theRecoil = new G4DynamicParticle;
142
143 theRecoil->SetDefinition( target_pd );
144 theRecoil->SetMomentum( theTarget.GetMomentum() );
145 theResult->AddSecondary( theRecoil, secID );
146
147 return theResult;
148
149}
double MyRNG(void *)
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
double phi() const
double theta() const
void setY(double)
void setZ(double)
void setX(double)
Hep3Vector vect() const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4ThreeVector &momentum)
double getElasticFinalState(double e_in, double temperature, double(*rng)(void *), void *rngState)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
static G4IonTable * GetIonTable()
G4double GetExcitationEnergyOfExcitedIsomer(G4int, G4int, G4int)
G4int GetNucleusEncoding(G4int iZ, G4int iA, G4int iM)
static G4LENDManager * GetInstance()
G4LENDManager * lend_manager
G4HadFinalState * returnUnchanged(const G4HadProjectile &aTrack, G4HadFinalState *theResult)
G4GIDI_target * get_target_from_map(G4int nuclear_code)
G4double GetTemperature() const
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition G4Pow.hh:230
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
void SetMass(const G4double mas)

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