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

#include <G4LENDModel.hh>

+ Inheritance diagram for G4LENDModel:

Public Member Functions

 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
 

Protected Member Functions

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

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 49 of file G4LENDModel.hh.

Constructor & Destructor Documentation

◆ G4LENDModel()

G4LENDModel::G4LENDModel ( G4String  name = "LENDModel")

Definition at line 45 of file G4LENDModel.cc.

47{
48
49 SetMinEnergy( 0.*eV );
50 SetMaxEnergy( 20.*MeV );
51
52 default_evaluation = "endl99";
53 allow_nat = false;
54 allow_any = false;
55
57
58}
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
static G4LENDManager * GetInstance()
G4LENDManager * lend_manager
Definition: G4LENDModel.hh:78

◆ ~G4LENDModel()

G4LENDModel::~G4LENDModel ( )

Definition at line 60 of file G4LENDModel.cc.

61{
62 for ( std::map< G4int , G4LENDUsedTarget* >::iterator
63 it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ )
64 {
65 delete it->second;
66 }
67}
std::map< G4int, G4LENDUsedTarget * > usedTarget_map
Definition: G4LENDModel.hh:79

Member Function Documentation

◆ AllowAnyCandidateTarget()

void G4LENDModel::AllowAnyCandidateTarget ( )
inline

Definition at line 62 of file G4LENDModel.hh.

62{ allow_any = true; recreate_used_target_map(); };
void recreate_used_target_map()
Definition: G4LENDModel.cc:70

Referenced by G4NeutronLENDBuilder::Build().

◆ AllowNaturalAbundanceTarget()

void G4LENDModel::AllowNaturalAbundanceTarget ( )
inline

Definition at line 61 of file G4LENDModel.hh.

61{ allow_nat = true; recreate_used_target_map(); };

◆ ApplyYourself()

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

Implements G4HadronicInteraction.

Reimplemented in G4LENDCapture, G4LENDElastic, G4LENDFission, and G4LENDInelastic.

Definition at line 166 of file G4LENDModel.cc.

167{
168
169 G4double temp = aTrack.GetMaterial()->GetTemperature();
170
171 //G4int iZ = int ( aTarg.GetZ() );
172 //G4int iA = int ( aTarg.GetN() );
173 //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
174 G4int iZ = aTarg.GetZ_asInt();
175 G4int iA = aTarg.GetA_asInt();
176
177 G4double ke = aTrack.GetKineticEnergy();
178
179 G4HadFinalState* theResult = new G4HadFinalState();
180
181 G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA ) )->second->GetTarget();
182
183 G4double aMu = aTarget->getElasticFinalState( ke*MeV, temp, NULL, NULL );
184
185 G4double phi = twopi*G4UniformRand();
186 G4double theta = std::acos( aMu );
187 //G4double sinth = std::sin( theta );
188
189 G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>( aTrack.GetDefinition() ) );
190 theNeutron.SetMomentum( aTrack.Get4Momentum().vect() );
191 theNeutron.SetKineticEnergy( ke );
192
193 G4ReactionProduct theTarget( G4ParticleTable::GetParticleTable()->FindIon( iZ , iA , 0 , iZ ) );
194
195 G4double mass = G4ParticleTable::GetParticleTable()->FindIon( iZ , iA , 0 , iZ )->GetPDGMass();
196
197// add Thermal motion
198 G4double kT = k_Boltzmann*temp;
199 G4ThreeVector v ( G4RandGauss::shoot() * std::sqrt( kT*mass )
200 , G4RandGauss::shoot() * std::sqrt( kT*mass )
201 , G4RandGauss::shoot() * std::sqrt( kT*mass ) );
202
203 theTarget.SetMomentum( v );
204
205
206 G4ThreeVector the3Neutron = theNeutron.GetMomentum();
207 G4double nEnergy = theNeutron.GetTotalEnergy();
208 G4ThreeVector the3Target = theTarget.GetMomentum();
209 G4double tEnergy = theTarget.GetTotalEnergy();
210 G4ReactionProduct theCMS;
211 G4double totE = nEnergy+tEnergy;
212 G4ThreeVector the3CMS = the3Target+the3Neutron;
213 theCMS.SetMomentum(the3CMS);
214 G4double cmsMom = std::sqrt(the3CMS*the3CMS);
215 G4double sqrts = std::sqrt((totE-cmsMom)*(totE+cmsMom));
216 theCMS.SetMass(sqrts);
217 theCMS.SetTotalEnergy(totE);
218
219 theNeutron.Lorentz(theNeutron, theCMS);
220 theTarget.Lorentz(theTarget, theCMS);
221 G4double en = theNeutron.GetTotalMomentum(); // already in CMS.
222 G4ThreeVector cms3Mom=theNeutron.GetMomentum(); // for neutron direction in CMS
223 G4double cms_theta=cms3Mom.theta();
224 G4double cms_phi=cms3Mom.phi();
225 G4ThreeVector tempVector;
226 tempVector.setX(std::cos(theta)*std::sin(cms_theta)*std::cos(cms_phi)
227 +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::cos(cms_phi)
228 -std::sin(theta)*std::sin(phi)*std::sin(cms_phi) );
229 tempVector.setY(std::cos(theta)*std::sin(cms_theta)*std::sin(cms_phi)
230 +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::sin(cms_phi)
231 +std::sin(theta)*std::sin(phi)*std::cos(cms_phi) );
232 tempVector.setZ(std::cos(theta)*std::cos(cms_theta)
233 -std::sin(theta)*std::cos(phi)*std::sin(cms_theta) );
234 tempVector *= en;
235 theNeutron.SetMomentum(tempVector);
236 theTarget.SetMomentum(-tempVector);
237 G4double tP = theTarget.GetTotalMomentum();
238 G4double tM = theTarget.GetMass();
239 theTarget.SetTotalEnergy(std::sqrt((tP+tM)*(tP+tM)-2.*tP*tM));
240 theNeutron.Lorentz(theNeutron, -1.*theCMS);
241 theTarget.Lorentz(theTarget, -1.*theCMS);
242
243 theResult->SetEnergyChange(theNeutron.GetKineticEnergy());
244 theResult->SetMomentumChange(theNeutron.GetMomentum().unit());
245 G4DynamicParticle* theRecoil = new G4DynamicParticle;
246
247 theRecoil->SetDefinition( G4ParticleTable::GetParticleTable()->FindIon( iZ, iA , 0, iZ ) );
248 theRecoil->SetMomentum( theTarget.GetMomentum() );
249
250 theResult->AddSecondary( theRecoil );
251
252 return theResult;
253
254}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
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)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4int GetNucleusEncoding(G4int iZ, G4int iA)
G4double GetTemperature() const
Definition: G4Material.hh:181
G4ParticleDefinition * FindIon(G4int atomicNumber, G4int atomicMass, G4double excitationEnergy)
static G4ParticleTable * GetParticleTable()
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
void SetMass(const G4double mas)

◆ BuildPhysicsTable()

void G4LENDModel::BuildPhysicsTable ( const G4ParticleDefinition )
inline

Definition at line 64 of file G4LENDModel.hh.

◆ ChangeDefaultEvaluation()

void G4LENDModel::ChangeDefaultEvaluation ( G4String  name)
inline

Definition at line 60 of file G4LENDModel.hh.

60{ default_evaluation = name; recreate_used_target_map(); };

Referenced by G4NeutronLENDBuilder::Build().

◆ create_used_target_map()

void G4LENDModel::create_used_target_map ( )
protected

Definition at line 86 of file G4LENDModel.cc.

87{
88
90
91 size_t numberOfElements = G4Element::GetNumberOfElements();
92 static const G4ElementTable* theElementTable = G4Element::GetElementTable();
93
94 for ( size_t i = 0 ; i < numberOfElements ; ++i )
95 {
96
97 const G4Element* anElement = (*theElementTable)[i];
98 G4int numberOfIsotope = anElement->GetNumberOfIsotopes();
99
100 if ( numberOfIsotope > 0 )
101 {
102 // User Defined Abundances
103 for ( G4int i_iso = 0 ; i_iso < numberOfIsotope ; i_iso++ )
104 {
105 G4int iZ = anElement->GetIsotope( i_iso )->GetZ();
106 G4int iA = anElement->GetIsotope( i_iso )->GetN();
107
108 G4LENDUsedTarget* aTarget = new G4LENDUsedTarget ( proj , default_evaluation , iZ , iA );
109 if ( allow_nat == true ) aTarget->AllowNat();
110 if ( allow_any == true ) aTarget->AllowAny();
111 usedTarget_map.insert( std::pair< G4int , G4LENDUsedTarget* > ( lend_manager->GetNucleusEncoding( iZ , iA ) , aTarget ) );
112 }
113 }
114 else
115 {
116 // Natural Abundances
118 G4int iZ = int ( anElement->GetZ() );
119 //G4cout << nistElementBuild->GetNumberOfNistIsotopes( int ( anElement->GetZ() ) ) << G4endl;
120 G4int numberOfNistIso = nistElementBuild->GetNumberOfNistIsotopes( int ( anElement->GetZ() ) );
121
122 for ( G4int ii = 0 ; ii < numberOfNistIso ; ii++ )
123 {
124 //G4cout << nistElementBuild->GetIsotopeAbundance( iZ , nistElementBuild->GetNistFirstIsotopeN( iZ ) + i ) << G4endl;
125 if ( nistElementBuild->GetIsotopeAbundance( iZ , nistElementBuild->GetNistFirstIsotopeN( iZ ) + ii ) > 0 )
126 {
127 G4int iMass = nistElementBuild->GetNistFirstIsotopeN( iZ ) + ii;
128 //G4cout << iZ << " " << nistElementBuild->GetNistFirstIsotopeN( iZ ) + i << " " << nistElementBuild->GetIsotopeAbundance ( iZ , iMass ) << G4endl;
129
130 G4LENDUsedTarget* aTarget = new G4LENDUsedTarget ( proj , default_evaluation , iZ , iMass );
131 if ( allow_nat == true ) aTarget->AllowNat();
132 if ( allow_any == true ) aTarget->AllowAny();
133 usedTarget_map.insert( std::pair< G4int , G4LENDUsedTarget* > ( lend_manager->GetNucleusEncoding( iZ , iMass ) , aTarget ) );
134
135 }
136
137 }
138
139 }
140 }
141
142
143
144 G4cout << "Dump UsedTarget for " << GetModelName() << G4endl;
145 G4cout << "Requested Evaluation, Z , A -> Actual Evaluation, Z , A(0=Nat) , Pointer of Target" << G4endl;
146 for ( std::map< G4int , G4LENDUsedTarget* >::iterator
147 it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ )
148 {
149 G4cout
150 << " " << it->second->GetWantedEvaluation()
151 << ", " << it->second->GetWantedZ()
152 << ", " << it->second->GetWantedA()
153 << " -> " << it->second->GetActualEvaluation()
154 << ", " << it->second->GetActualZ()
155 << ", " << it->second->GetActualA()
156 << ", " << it->second->GetTarget()
157 << G4endl;
158 }
159
160}
std::vector< G4Element * > G4ElementTable
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double GetZ() const
Definition: G4Element.hh:131
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
const G4String & GetModelName() const
G4int GetZ() const
Definition: G4Isotope.hh:91
G4int GetN() const
Definition: G4Isotope.hh:94
G4bool RequestChangeOfVerboseLevel(G4int)
G4NistElementBuilder * GetNistElementBuilder()
G4ParticleDefinition * proj
Definition: G4LENDModel.hh:77
G4int GetNumberOfNistIsotopes(G4int Z)
G4double GetIsotopeAbundance(G4int Z, G4int N)
G4int GetNistFirstIsotopeN(G4int Z)

Referenced by G4LENDCapture::G4LENDCapture(), G4LENDElastic::G4LENDElastic(), G4LENDFission::G4LENDFission(), G4LENDInelastic::G4LENDInelastic(), and recreate_used_target_map().

◆ recreate_used_target_map()

void G4LENDModel::recreate_used_target_map ( )
protected

Definition at line 70 of file G4LENDModel.cc.

71{
72
73 for ( std::map< G4int , G4LENDUsedTarget* >::iterator
74 it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ )
75 {
76 delete it->second;
77 }
78 usedTarget_map.clear();
79
81
82}
void create_used_target_map()
Definition: G4LENDModel.cc:86

Referenced by AllowAnyCandidateTarget(), AllowNaturalAbundanceTarget(), BuildPhysicsTable(), and ChangeDefaultEvaluation().

Member Data Documentation

◆ lend_manager

◆ proj

◆ usedTarget_map


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