Geant4 11.4.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 &)
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

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 //default_evaluation = "ENDF/B-8.0";
61 //default_evaluation = "ENDF/B-7.1";
62 default_evaluation = "";
63
64 allow_nat = false;
65 allow_any = false;
66
68
70}
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
const G4String & GetModelName() const
void SetMaxEnergy(const G4double anEnergy)
static G4LENDManager * GetInstance()
G4LENDManager * lend_manager
G4ParticleDefinition * proj
static G4int GetModelID(const G4int modelIndex)

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

◆ ~G4LENDModel()

G4LENDModel::~G4LENDModel ( )

Definition at line 72 of file G4LENDModel.cc.

73{
74 for ( std::map< G4int , G4LENDUsedTarget* >::iterator
75 it = usedTarget_map.begin() ; it != usedTarget_map.end() ; it ++ )
76 {
77 delete it->second;
78 }
79}
std::map< G4int, G4LENDUsedTarget * > usedTarget_map

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()

◆ AllowNaturalAbundanceTarget()

void G4LENDModel::AllowNaturalAbundanceTarget ( )
inline

Definition at line 63 of file G4LENDModel.hh.

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

Referenced by G4HadronElasticPhysicsLEND::ConstructProcess().

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

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

Definition at line 163 of file G4LENDModel.cc.

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

Reimplemented in G4LENDorBERTModel.

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 G4HadronElasticPhysicsLEND::ConstructProcess().

◆ create_used_target_map()

void G4LENDModel::create_used_target_map ( )
protected

Definition at line 98 of file G4LENDModel.cc.

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

Referenced by G4LENDCombinedModel::BuildPhysicsTable(), G4LENDGammaModel::BuildPhysicsTable(), G4LENDorBERTModel::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 280 of file G4LENDModel.cc.

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

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 272 of file G4LENDModel.cc.

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

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 82 of file G4LENDModel.cc.

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

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

◆ returnUnchanged()

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

Definition at line 258 of file G4LENDModel.cc.

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