Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4LENDInelastic Class Reference

#include <G4LENDInelastic.hh>

Inheritance diagram for G4LENDInelastic:
Collaboration diagram for G4LENDInelastic:

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 &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, G4double
GetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void InitialiseModel ()
 

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 47 of file G4LENDInelastic.hh.

Constructor & Destructor Documentation

G4LENDInelastic::G4LENDInelastic ( G4ParticleDefinition pd)
inline

Definition at line 52 of file G4LENDInelastic.hh.

53  :G4LENDModel( "LENDInelastic" )
54  {
55  proj = pd;
56 
57 // theModelName = "LENDInelastic for ";
58 // theModelName += proj->GetParticleName();
60 
63  G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(p);
64  if(!pre) { pre = new G4PreCompoundModel(); }
65  preco = pre;
66  };
const char * p
Definition: xmltok.h:285
G4ParticleDefinition * proj
Definition: G4LENDModel.hh:79
G4LENDModel(G4String name="LENDModel")
Definition: G4LENDModel.cc:47
G4HadronicInteraction * FindModel(const G4String &name)
void create_used_target_map()
Definition: G4LENDModel.cc:93
static G4HadronicInteractionRegistry * Instance()

Here is the call graph for this function:

G4LENDInelastic::~G4LENDInelastic ( )
inline

Definition at line 68 of file G4LENDInelastic.hh.

68 {;};

Member Function Documentation

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

Reimplemented from G4LENDModel.

Definition at line 32 of file G4LENDInelastic.cc.

33 {
34  //return preco->ApplyYourself( aTrack, aTarg );
35 
36  G4ThreeVector proj_p = aTrack.Get4Momentum().vect();
37 
38  G4double temp = aTrack.GetMaterial()->GetTemperature();
39 
40  //G4int iZ = int ( aTarg.GetZ() );
41  //G4int iA = int ( aTarg.GetN() );
42  //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
43  G4int iZ = aTarg.GetZ_asInt();
44  G4int iA = aTarg.GetA_asInt();
45  //G4int iM = aTarg.GetM_asInt();
46  G4int iM = 0;
47  if ( aTarg.GetIsotope() != NULL ) {
48  iM = aTarg.GetIsotope()->Getm();
49  }
50  //G4cout << "target: Z = " << iZ << " N = " << iA << G4endl;
51 
52  G4double ke = aTrack.GetKineticEnergy();
53  //G4cout << "projectile: KE = " << ke/MeV << " [MeV]" << G4endl;
54 
55  G4HadFinalState* theResult = &theParticleChange;
56  theResult->Clear();
57 
58  G4GIDI_target* aTarget = usedTarget_map.find( lend_manager->GetNucleusEncoding( iZ , iA , iM ) )->second->GetTarget();
59  //std::vector<G4GIDI_Product>* products = aTarget->getOthersFinalState( ke*MeV, temp, MyRNG, NULL );
60 
61  std::vector<G4GIDI_Product>* products;
62  for ( G4int i = 0 ; i != 1024 ; i++ ) {
63  products = aTarget->getOthersFinalState( ke*MeV, temp, MyRNG, NULL );
64  if ( products != NULL ) break;
65  }
66  //return preco->ApplyYourself( aTrack, aTarg );
67 
68  G4int iTotZ = iZ + aTrack.GetDefinition()->GetAtomicNumber();
69  G4int iTotA = iA + aTrack.GetDefinition()->GetAtomicMass();
70 
71  if ( products != NULL )
72  {
73  //G4cout << "Using LENDModel" << G4endl;
74 
75  G4ThreeVector psum(0);
76  G4bool needResidual = true;
77  int totN = 0;
78  int totZ = 0;
79  for ( G4int j = 0; j < int( products->size() ); j++ )
80  {
81 
82  G4int jZ = (*products)[j].Z;
83  G4int jA = (*products)[j].A;
84  G4int jm = (*products)[j].m;
85  //TK
86  //We need coordination LEND *products)[j].m and G4IonTable(Z,A,m)
87  //Excitation energy of isomer level is might (probably) different each other.
88  //
89 
90  //G4cout << "ZA = " << 1000 * (*products)[j].Z + (*products)[j].A << " EK = "
91  // << (*products)[j].kineticEnergy
92  // << " px " << (*products)[j].px
93  // << " py " << (*products)[j].py
94  // << " pz " << (*products)[j].pz
95  // << G4endl;
96 
97  iTotZ -= jZ;
98  iTotA -= jA;
99 
100  G4DynamicParticle* theSec = new G4DynamicParticle;
101 
102  if ( jA == 1 && jZ == 1 )
103  {
104  theSec->SetDefinition( G4Proton::Proton() );
105  totN += 1;
106  totZ += 1;
107  }
108  else if ( jA == 1 && jZ == 0 )
109  {
110  theSec->SetDefinition( G4Neutron::Neutron() );
111  totN += 1;
112  }
113  else if ( jZ > 0 )
114  {
115  if ( jA != 0 )
116  {
117  theSec->SetDefinition( G4IonTable::GetIonTable()->GetIon( jZ , jA , jm ) );
118  totN += jA;
119  totZ += jZ;
120  }
121  else
122  {
123  theSec->SetDefinition( G4IonTable::GetIonTable()->GetIon( jZ , iA+aTrack.GetDefinition()->GetAtomicMass()-totN , jm ) );
124  iTotZ -= jZ;
125  iTotA -= iA+aTrack.GetDefinition()->GetAtomicMass()-totN;
126  needResidual=false;
127  }
128  }
129  else
130  {
131  theSec->SetDefinition( G4Gamma::Gamma() );
132  }
133 
134  G4ThreeVector p( (*products)[j].px*MeV , (*products)[j].py*MeV , (*products)[j].pz*MeV );
135  psum += p;
136  if ( p.mag() == 0 ) p = proj_p - psum;
137 
138  theSec->SetMomentum( p );
139 
140  theResult->AddSecondary( theSec );
141  }
142 
143  if ( !( iTotZ == 0 && iTotA == 0 ) ) {
144 
145  if ( iTotZ >= 0 && iTotA > 0 ) {
146  if ( needResidual ) {
147  G4DynamicParticle* residual = new G4DynamicParticle;
148  if ( iTotZ > 0 ) {
149  residual->SetDefinition( G4IonTable::GetIonTable()->GetIon( iTotZ , iTotA ) );
150  } else if ( iTotA == 1 ) {
151  residual->SetDefinition( G4Neutron::Neutron() );
152  } else {
153  //G4cout << "Charge or Baryon Number Error #3 iTotZ = " << iTotZ << ", iTotA = " << iTotA << G4endl;
154  ;
155  }
156  residual->SetMomentum( proj_p - psum );
157  theResult->AddSecondary( residual );
158  } else {
159  //G4cout << "Charge or Baryon Number Error #1 iTotZ = " << iTotZ << ", iTotA = " << iTotA << G4endl;
160  ;
161  }
162  } else {
163 
164  if ( needResidual ) {
165  //G4cout << "Charge or Baryon Number Error #2 iTotZ = " << iTotZ << ", iTotA = " << iTotA << G4endl;
166  ;
167  }
168  }
169 
170  }
171 
172  }
173  else {
174  //G4cout << "Using PreCompoundModel" << G4endl;
175  if ( aTrack.GetDefinition() == G4Proton::Proton() ||
176  aTrack.GetDefinition() == G4Neutron::Neutron() ) {
177  theResult = preco->ApplyYourself( aTrack, aTarg );
178  } else {
179  return theResult;
180  }
181  }
182  delete products;
183 
184  theResult->SetStatusChange( stopAndKill );
185 
186  return theResult;
187 
188 }
void SetMomentum(const G4ThreeVector &momentum)
double MyRNG(void *)
Definition: G4LENDModel.cc:45
const char * p
Definition: xmltok.h:285
static constexpr double second
Definition: G4SIunits.hh:157
std::vector< G4GIDI_Product > * getOthersFinalState(double e_in, double temperature, double(*rng)(void *), void *rngState)
int G4int
Definition: G4Types.hh:78
G4int GetAtomicNumber() const
void SetStatusChange(G4HadFinalStateStatus aS)
G4int GetNucleusEncoding(G4int iZ, G4int iA, G4int iM)
Hep3Vector vect() const
const G4ParticleDefinition * GetDefinition() const
std::map< G4int, G4LENDUsedTarget * > usedTarget_map
Definition: G4LENDModel.hh:81
bool G4bool
Definition: G4Types.hh:79
G4double GetKineticEnergy() const
typedef int(XMLCALL *XML_NotStandaloneHandler)(void *userData)
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
G4int GetAtomicMass() const
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
G4double GetTemperature() const
Definition: G4Material.hh:182
static constexpr double MeV
Definition: G4SIunits.hh:214
const G4Material * GetMaterial() const
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4LENDManager * lend_manager
Definition: G4LENDModel.hh:80
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0

Here is the call graph for this function:


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