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

#include <G4LENDInelastic.hh>

Inheritance diagram for G4LENDInelastic:
G4LENDModel G4HadronicInteraction

Public Member Functions

 G4LENDInelastic (G4ParticleDefinition *pd)
 
 ~G4LENDInelastic ()
 
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 &)
 
- 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
 

Additional Inherited Members

- Protected Member Functions inherited from G4LENDModel
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 inherited from G4LENDModel
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 44 of file G4LENDInelastic.hh.

Constructor & Destructor Documentation

G4LENDInelastic::G4LENDInelastic ( G4ParticleDefinition pd)
inline

Definition at line 49 of file G4LENDInelastic.hh.

References G4LENDModel::create_used_target_map(), and G4LENDModel::proj.

50  :G4LENDModel( "LENDInelastic" )
51  {
52  proj = pd;
53 
54 // theModelName = "LENDInelastic for ";
55 // theModelName += proj->GetParticleName();
57  };
G4ParticleDefinition * proj
Definition: G4LENDModel.hh:77
G4LENDModel(G4String name="LENDModel")
Definition: G4LENDModel.cc:45
void create_used_target_map()
Definition: G4LENDModel.cc:90
G4LENDInelastic::~G4LENDInelastic ( )
inline

Definition at line 59 of file G4LENDInelastic.hh.

59 {;};

Member Function Documentation

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

Reimplemented from G4LENDModel.

Definition at line 31 of file G4LENDInelastic.cc.

References G4HadFinalState::AddSecondary(), G4HadFinalState::Clear(), G4Gamma::Gamma(), G4HadProjectile::Get4Momentum(), G4Nucleus::GetA_asInt(), G4Nucleus::GetIsotope(), G4HadProjectile::GetKineticEnergy(), G4Isotope::Getm(), G4HadProjectile::GetMaterial(), G4LENDManager::GetNucleusEncoding(), G4GIDI_target::getOthersFinalState(), G4ParticleTable::GetParticleTable(), G4Material::GetTemperature(), G4Nucleus::GetZ_asInt(), int(), G4LENDModel::lend_manager, CLHEP::Hep3Vector::mag(), python.hepunit::MeV, G4Neutron::Neutron(), G4Proton::Proton(), python.hepunit::second, G4DynamicParticle::SetDefinition(), G4DynamicParticle::SetMomentum(), G4HadFinalState::SetStatusChange(), stopAndKill, G4HadronicInteraction::theParticleChange, G4LENDModel::usedTarget_map, and CLHEP::HepLorentzVector::vect().

32 {
33 
34  G4ThreeVector proj_p = aTrack.Get4Momentum().vect();
35 
36  G4double temp = aTrack.GetMaterial()->GetTemperature();
37 
38  //G4int iZ = int ( aTarg.GetZ() );
39  //G4int iA = int ( aTarg.GetN() );
40  //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
41  G4int iZ = aTarg.GetZ_asInt();
42  G4int iA = aTarg.GetA_asInt();
43  //G4int iM = aTarg.GetM_asInt();
44  G4int iM = 0;
45  if ( aTarg.GetIsotope() != NULL ) {
46  iM = aTarg.GetIsotope()->Getm();
47  }
48  //G4cout << "target: Z = " << iZ << " N = " << iA << G4endl;
49 
50  G4double ke = aTrack.GetKineticEnergy();
51  //G4cout << "projectile: KE = " << ke/MeV << " [MeV]" << G4endl;
52 
53  G4HadFinalState* theResult = &theParticleChange;
54  theResult->Clear();
55 
56  G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA , iM ) )->second->GetTarget();
57  std::vector<G4GIDI_Product>* products = aTarget->getOthersFinalState( ke*MeV, temp, NULL, NULL );
58  if ( products != NULL )
59  {
60 
61  G4ThreeVector psum(0);
62 
63  int totN = 0;
64  for ( G4int j = 0; j < int( products->size() ); j++ )
65  {
66 
67  G4int jZ = (*products)[j].Z;
68  G4int jA = (*products)[j].A;
69 
70  //G4cout << "ZA = " << 1000 * (*products)[j].Z + (*products)[j].A << " EK = "
71  // << (*products)[j].kineticEnergy
72  // << " px " << (*products)[j].px
73  // << " py " << (*products)[j].py
74  // << " pz " << (*products)[j].pz
75  // << G4endl;
76 
78 
79  if ( jA == 1 && jZ == 1 )
80  {
81  theSec->SetDefinition( G4Proton::Proton() );
82  totN += 1;
83  }
84  else if ( jA == 1 && jZ == 0 )
85  {
86  theSec->SetDefinition( G4Neutron::Neutron() );
87  totN += 1;
88  }
89  else if ( jZ > 0 )
90  {
91  if ( jA != 0 )
92  {
93  theSec->SetDefinition( G4ParticleTable::GetParticleTable()->FindIon( jZ , jA , 0 , 0 ) );
94  totN += jA;
95  }
96  else
97  {
98  theSec->SetDefinition( G4ParticleTable::GetParticleTable()->FindIon( jZ , iA+1-totN , 0 , 0 ) );
99  }
100  }
101  else
102  {
103  theSec->SetDefinition( G4Gamma::Gamma() );
104  }
105 
106  G4ThreeVector p( (*products)[j].px*MeV , (*products)[j].py*MeV , (*products)[j].pz*MeV );
107  psum += p;
108  if ( p.mag() == 0 ) p = proj_p - psum;
109 
110  theSec->SetMomentum( p );
111 
112  theResult->AddSecondary( theSec );
113  }
114  }
115  delete products;
116 
117  theResult->SetStatusChange( stopAndKill );
118 
119  return theResult;
120 
121 }
void SetMomentum(const G4ThreeVector &momentum)
typedef int(XMLCALL *XML_NotStandaloneHandler)(void *userData)
const char * p
Definition: xmltok.h:285
std::vector< G4GIDI_Product > * getOthersFinalState(double e_in, double temperature, double(*rng)(void *), void *rngState)
int G4int
Definition: G4Types.hh:78
void SetStatusChange(G4HadFinalStateStatus aS)
G4int GetNucleusEncoding(G4int iZ, G4int iA, G4int iM)
Hep3Vector vect() const
std::map< G4int, G4LENDUsedTarget * > usedTarget_map
Definition: G4LENDModel.hh:79
G4double GetKineticEnergy() const
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4LorentzVector & Get4Momentum() 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 AddSecondary(G4DynamicParticle *aP)

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