Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
G4LENDModel Class Reference

#include <G4LENDModel.hh>

Inheritance diagram for G4LENDModel:
G4HadronicInteraction G4LENDCapture G4LENDElastic G4LENDFission G4LENDInelastic

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 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
 
virtual G4int GetVerboseLevel () const
 
virtual 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, G4double
GetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () 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 ( G4String  name = "LENDModel")

Definition at line 45 of file G4LENDModel.cc.

References python.hepunit::eV, G4LENDManager::GetInstance(), lend_manager, python.hepunit::MeV, proj, G4HadronicInteraction::SetMaxEnergy(), and G4HadronicInteraction::SetMinEnergy().

46 :G4HadronicInteraction( name )
47 {
48 
49  proj = NULL; //will be set in an inherited class
50 
51  SetMinEnergy( 0.*eV );
52  SetMaxEnergy( 20.*MeV );
53 
54  //default_evaluation = "endl99";
55  default_evaluation = "ENDF.B-VII.0";
56 
57  allow_nat = false;
58  allow_any = false;
59 
61 
62 }
G4ParticleDefinition * proj
Definition: G4LENDModel.hh:77
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4LENDManager * GetInstance()
void SetMaxEnergy(const G4double anEnergy)
G4LENDManager * lend_manager
Definition: G4LENDModel.hh:78
G4LENDModel::~G4LENDModel ( )

Definition at line 64 of file G4LENDModel.cc.

References usedTarget_map.

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

Member Function Documentation

void G4LENDModel::AllowAnyCandidateTarget ( )
inline

Definition at line 62 of file G4LENDModel.hh.

References recreate_used_target_map().

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

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

Definition at line 61 of file G4LENDModel.hh.

References recreate_used_target_map().

61 { allow_nat = true; recreate_used_target_map(); };
void recreate_used_target_map()
Definition: G4LENDModel.cc:74
G4HadFinalState * G4LENDModel::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus aTargetNucleus 
)
virtual

Implements G4HadronicInteraction.

Reimplemented in G4LENDFission, G4LENDCapture, G4LENDElastic, and G4LENDInelastic.

Definition at line 172 of file G4LENDModel.cc.

References G4HadFinalState::AddSecondary(), G4ParticleTable::FindIon(), G4UniformRand, G4HadProjectile::Get4Momentum(), G4Nucleus::GetA_asInt(), G4HadProjectile::GetDefinition(), G4GIDI_target::getElasticFinalState(), G4Nucleus::GetIsotope(), G4HadProjectile::GetKineticEnergy(), G4ReactionProduct::GetKineticEnergy(), G4Isotope::Getm(), G4ReactionProduct::GetMass(), G4HadProjectile::GetMaterial(), G4ReactionProduct::GetMomentum(), G4LENDManager::GetNucleusEncoding(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGMass(), G4Material::GetTemperature(), G4ReactionProduct::GetTotalEnergy(), G4ReactionProduct::GetTotalMomentum(), G4Nucleus::GetZ_asInt(), python.hepunit::k_Boltzmann, lend_manager, G4ReactionProduct::Lorentz(), python.hepunit::MeV, CLHEP::Hep3Vector::phi(), python.hepunit::second, G4DynamicParticle::SetDefinition(), G4HadFinalState::SetEnergyChange(), G4ReactionProduct::SetKineticEnergy(), G4ReactionProduct::SetMass(), G4ReactionProduct::SetMomentum(), G4DynamicParticle::SetMomentum(), G4HadFinalState::SetMomentumChange(), G4ReactionProduct::SetTotalEnergy(), CLHEP::Hep3Vector::setX(), CLHEP::Hep3Vector::setY(), CLHEP::Hep3Vector::setZ(), G4INCL::DeJongSpin::shoot(), CLHEP::Hep3Vector::theta(), theTarget, python.hepunit::twopi, CLHEP::Hep3Vector::unit(), usedTarget_map, test::v, and CLHEP::HepLorentzVector::vect().

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

Definition at line 64 of file G4LENDModel.hh.

References recreate_used_target_map().

void recreate_used_target_map()
Definition: G4LENDModel.cc:74
void G4LENDModel::ChangeDefaultEvaluation ( G4String  name)
inline

Definition at line 60 of file G4LENDModel.hh.

References recreate_used_target_map().

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

60 { default_evaluation = name; recreate_used_target_map(); };
void recreate_used_target_map()
Definition: G4LENDModel.cc:74
const XML_Char * name
void G4LENDModel::create_used_target_map ( )
protected

Definition at line 90 of file G4LENDModel.cc.

References G4LENDUsedTarget::AllowAny(), G4LENDUsedTarget::AllowNat(), G4cout, G4endl, G4Element::GetElementTable(), G4Element::GetIsotope(), G4NistElementBuilder::GetIsotopeAbundance(), G4Isotope::Getm(), G4HadronicInteraction::GetModelName(), G4Isotope::GetN(), G4LENDManager::GetNistElementBuilder(), G4NistElementBuilder::GetNistFirstIsotopeN(), G4LENDManager::GetNucleusEncoding(), G4Element::GetNumberOfElements(), G4Element::GetNumberOfIsotopes(), G4NistElementBuilder::GetNumberOfNistIsotopes(), G4Isotope::GetZ(), G4Element::GetZ(), int(), lend_manager, proj, G4LENDManager::RequestChangeOfVerboseLevel(), usedTarget_map, and G4HadronicInteraction::verboseLevel.

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

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

Definition at line 74 of file G4LENDModel.cc.

References create_used_target_map(), and usedTarget_map.

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

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

Field Documentation

G4LENDManager* G4LENDModel::lend_manager
protected
G4ParticleDefinition* G4LENDModel::proj
protected
std::map< G4int , G4LENDUsedTarget* > G4LENDModel::usedTarget_map
protected

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