Geant4 11.1.1
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 &)
 
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
 

Protected Member Functions

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

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

Constructor & Destructor Documentation

◆ G4LENDModel()

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

Definition at line 48 of file G4LENDModel.cc.

49 :G4HadronicInteraction( name ), secID(-1)
50{
51
52 proj = NULL; //will be set in an inherited class
53
54 SetMinEnergy( 0.*eV );
55 SetMaxEnergy( 20.*MeV );
56
57 //default_evaluation = "endl99";
58 //default_evaluation = "ENDF.B-VII.0";
59 default_evaluation = "ENDF/BVII.1";
60
61 allow_nat = false;
62 allow_any = false;
63
65
67}
void SetMinEnergy(G4double anEnergy)
const G4String & GetModelName() const
void SetMaxEnergy(const G4double anEnergy)
static G4LENDManager * GetInstance()
G4LENDManager * lend_manager
Definition: G4LENDModel.hh:84
G4ParticleDefinition * proj
Definition: G4LENDModel.hh:83
static G4int GetModelID(const G4int modelIndex)

◆ ~G4LENDModel()

G4LENDModel::~G4LENDModel ( )

Definition at line 69 of file G4LENDModel.cc.

70{
71 for ( std::map< G4int , G4LENDUsedTarget* >::iterator
72 it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ )
73 {
74 delete it->second;
75 }
76}
std::map< G4int, G4LENDUsedTarget * > usedTarget_map
Definition: G4LENDModel.hh:85

Member Function Documentation

◆ AllowAnyCandidateTarget()

void G4LENDModel::AllowAnyCandidateTarget ( )
inline

Definition at line 64 of file G4LENDModel.hh.

64{ allow_any = true; recreate_used_target_map(); };
void recreate_used_target_map()
Definition: G4LENDModel.cc:79

◆ AllowNaturalAbundanceTarget()

void G4LENDModel::AllowNaturalAbundanceTarget ( )
inline

Definition at line 63 of file G4LENDModel.hh.

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

Referenced by G4NeutronLENDBuilder::Build(), and G4HadronElasticPhysicsLEND::ConstructProcess().

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Reimplemented in G4LENDorBERTModel, G4LENDCapture, G4LENDCombinedModel, G4LENDElastic, G4LENDFission, G4LENDGammaModel, and G4LENDInelastic.

Definition at line 160 of file G4LENDModel.cc.

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

Referenced by G4LENDCombinedModel::ApplyYourself(), and G4LENDGammaModel::ApplyYourself().

◆ BuildPhysicsTable()

void G4LENDModel::BuildPhysicsTable ( const G4ParticleDefinition )
inlinevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 66 of file G4LENDModel.hh.

◆ ChangeDefaultEvaluation()

void G4LENDModel::ChangeDefaultEvaluation ( G4String  name)
inline

Definition at line 62 of file G4LENDModel.hh.

62{ default_evaluation = name; recreate_used_target_map(); };
const char * name(G4int ptype)

Referenced by G4NeutronLENDBuilder::Build(), and G4HadronElasticPhysicsLEND::ConstructProcess().

◆ create_used_target_map()

void G4LENDModel::create_used_target_map ( )
protected

Definition at line 95 of file G4LENDModel.cc.

96{
97
99
100 std::size_t numberOfElements = G4Element::GetNumberOfElements();
101 static const G4ElementTable* theElementTable = G4Element::GetElementTable();
102
103 for ( std::size_t i = 0 ; i < numberOfElements ; ++i )
104 {
105
106 const G4Element* anElement = (*theElementTable)[i];
107 G4int numberOfIsotope = (G4int)anElement->GetNumberOfIsotopes();
108
109 if ( numberOfIsotope > 0 )
110 {
111 // User Defined Abundances
112 for ( G4int i_iso = 0 ; i_iso < numberOfIsotope ; ++i_iso )
113 {
114 G4int iZ = anElement->GetIsotope( i_iso )->GetZ();
115 G4int iA = anElement->GetIsotope( i_iso )->GetN();
116 G4int iIsomer = anElement->GetIsotope( i_iso )->Getm();
117
118 G4LENDUsedTarget* aTarget = new G4LENDUsedTarget ( proj , default_evaluation , iZ , iA , iIsomer );
119 if ( allow_nat == true ) aTarget->AllowNat();
120 if ( allow_any == true ) aTarget->AllowAny();
121 usedTarget_map.insert( std::pair< G4int , G4LENDUsedTarget* > ( lend_manager->GetNucleusEncoding( iZ , iA , iIsomer ) , aTarget ) );
122 }
123 }
124 else
125 {
126 // Natural Abundances
128 G4int iZ = int ( anElement->GetZ() );
129 //G4cout << nistElementBuild->GetNumberOfNistIsotopes( int ( anElement->GetZ() ) ) << G4endl;
130 G4int numberOfNistIso = nistElementBuild->GetNumberOfNistIsotopes( int ( anElement->GetZ() ) );
131
132 for ( G4int ii = 0 ; ii < numberOfNistIso ; ii++ )
133 {
134 //G4cout << nistElementBuild->GetIsotopeAbundance( iZ , nistElementBuild->GetNistFirstIsotopeN( iZ ) + i ) << G4endl;
135 if ( nistElementBuild->GetIsotopeAbundance( iZ , nistElementBuild->GetNistFirstIsotopeN( iZ ) + ii ) > 0 )
136 {
137 G4int iMass = nistElementBuild->GetNistFirstIsotopeN( iZ ) + ii;
138 //G4cout << iZ << " " << nistElementBuild->GetNistFirstIsotopeN( iZ ) + i << " " << nistElementBuild->GetIsotopeAbundance ( iZ , iMass ) << G4endl;
139 G4int iIsomer = 0;
140
141 G4LENDUsedTarget* aTarget = new G4LENDUsedTarget ( proj , default_evaluation , iZ , iMass );
142 if ( allow_nat == true ) aTarget->AllowNat();
143 if ( allow_any == true ) aTarget->AllowAny();
144 usedTarget_map.insert( std::pair< G4int , G4LENDUsedTarget* > ( lend_manager->GetNucleusEncoding( iZ , iMass , iIsomer ) , aTarget ) );
145
146 }
147
148 }
149
150 }
151 }
152
154}
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
G4double GetZ() const
Definition: G4Element.hh:131
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4int GetZ() const
Definition: G4Isotope.hh:90
G4int Getm() const
Definition: G4Isotope.hh:99
G4int GetN() const
Definition: G4Isotope.hh:93
G4bool RequestChangeOfVerboseLevel(G4int)
G4NistElementBuilder * GetNistElementBuilder()
void DumpLENDTargetInfo(G4bool force=false)
Definition: G4LENDModel.cc:277
G4int GetNumberOfNistIsotopes(G4int Z) const
G4double GetIsotopeAbundance(G4int Z, G4int N) const
G4int GetNistFirstIsotopeN(G4int Z) const

Referenced by G4LENDorBERTModel::BuildPhysicsTable(), G4LENDCombinedModel::BuildPhysicsTable(), G4LENDGammaModel::BuildPhysicsTable(), DumpLENDTargetInfo(), G4LENDCapture::G4LENDCapture(), G4LENDElastic::G4LENDElastic(), G4LENDFission::G4LENDFission(), G4LENDInelastic::G4LENDInelastic(), and recreate_used_target_map().

◆ DumpLENDTargetInfo()

void G4LENDModel::DumpLENDTargetInfo ( G4bool  force = false)

Definition at line 277 of file G4LENDModel.cc.

277 {
278
279 if ( lend_manager->GetVerboseLevel() >= 1 || force ) {
280 if ( usedTarget_map.size() == 0 ) create_used_target_map();
281 G4cout << "Dumping UsedTarget of " << GetModelName() << " for " << proj->GetParticleName() << G4endl;
282 G4cout << "Requested Evaluation, Z , A -> Actual Evaluation, Z , A(0=Nat) " << G4endl;
283 for ( std::map< G4int , G4LENDUsedTarget* >::iterator
284 it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ ) {
285 G4cout
286 << " " << it->second->GetWantedEvaluation()
287 << ", " << it->second->GetWantedZ()
288 << ", " << it->second->GetWantedA()
289 << " -> " << it->second->GetActualEvaluation()
290 << ", " << it->second->GetActualZ()
291 << ", " << it->second->GetActualA()
292 << G4endl;
293 }
294 }
295}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4int GetVerboseLevel()
void create_used_target_map()
Definition: G4LENDModel.cc:95
const G4String & GetParticleName() const

Referenced by G4LENDBertiniGammaElectroNuclearBuilder::Build(), G4HadronElasticPhysicsLEND::ConstructProcess(), and create_used_target_map().

◆ get_target_from_map()

G4GIDI_target * G4LENDModel::get_target_from_map ( G4int  nuclear_code)
protected

Definition at line 269 of file G4LENDModel.cc.

269 {
270 G4GIDI_target* target = NULL;
271 if ( usedTarget_map.find( nuclear_code ) != usedTarget_map.end() ) {
272 target = usedTarget_map.find( nuclear_code )->second->GetTarget();
273 }
274 return target;
275}

Referenced by G4LENDCapture::ApplyYourself(), G4LENDElastic::ApplyYourself(), G4LENDFission::ApplyYourself(), G4LENDInelastic::ApplyYourself(), G4LENDCombinedModel::HasData(), and G4LENDGammaModel::HasData().

◆ recreate_used_target_map()

void G4LENDModel::recreate_used_target_map ( )
protected

Definition at line 79 of file G4LENDModel.cc.

80{
81
82 for ( std::map< G4int , G4LENDUsedTarget* >::iterator
83 it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ )
84 {
85 delete it->second;
86 }
87 usedTarget_map.clear();
88
90
91}

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

◆ returnUnchanged()

G4HadFinalState * G4LENDModel::returnUnchanged ( const G4HadProjectile aTrack,
G4HadFinalState theResult 
)
protected

Definition at line 255 of file G4LENDModel.cc.

255 {
256 if ( lend_manager->GetVerboseLevel() >= 1 ) {
257 G4String message;
258 message = "Produce unchanged final state is requested in ";
259 message += this->GetModelName();
260 message += ". Cross section and model likely have an inconsistency.";
261 G4Exception( "G4LENDModel::returnUnchanged(,)" , "LENDModel-01" , JustWarning ,
262 message );
263 }
264 theResult->SetEnergyChange( aTrack.GetKineticEnergy() );
265 theResult->SetMomentumChange( aTrack.Get4Momentum().getV().unit() );
266 return theResult;
267}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
Hep3Vector unit() const
Hep3Vector getV() const

Referenced by G4LENDCapture::ApplyYourself(), G4LENDElastic::ApplyYourself(), and G4LENDFission::ApplyYourself().

Member Data Documentation

◆ lend_manager

◆ proj

◆ secID

◆ usedTarget_map

std::map< G4int , G4LENDUsedTarget* > G4LENDModel::usedTarget_map
protected

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