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

#include <G4NeutronRadCapture.hh>

Inheritance diagram for G4NeutronRadCapture:
G4HadronicInteraction

Public Member Functions

 G4NeutronRadCapture ()
 
virtual ~G4NeutronRadCapture ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
- 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 G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 51 of file G4NeutronRadCapture.hh.

Constructor & Destructor Documentation

G4NeutronRadCapture::G4NeutronRadCapture ( )

Definition at line 53 of file G4NeutronRadCapture.cc.

References python.hepunit::eV, G4ParticleTable::GetIonTable(), G4ParticleTable::GetParticleTable(), python.hepunit::GeV, python.hepunit::keV, G4HadronicInteraction::SetMaxEnergy(), G4HadronicInteraction::SetMinEnergy(), and python.hepunit::TeV.

54  : G4HadronicInteraction("nRadCapture")
55 {
56  lowestEnergyLimit = 0.1*eV;
57  minExcitation = 1*keV;
58  SetMinEnergy( 0.0*GeV );
59  SetMaxEnergy( 100.*TeV );
60  photonEvaporation = new G4PhotonEvaporation();
61  //photonEvaporation = 0;
62  theTableOfIons = G4ParticleTable::GetParticleTable()->GetIonTable();
63 }
void SetMinEnergy(G4double anEnergy)
G4IonTable * GetIonTable() const
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4ParticleTable * GetParticleTable()
void SetMaxEnergy(const G4double anEnergy)
G4NeutronRadCapture::~G4NeutronRadCapture ( )
virtual

Definition at line 65 of file G4NeutronRadCapture.cc.

66 {
67  delete photonEvaporation;
68 }

Member Function Documentation

G4HadFinalState * G4NeutronRadCapture::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
virtual

Implements G4HadronicInteraction.

Definition at line 70 of file G4NeutronRadCapture.cc.

References G4HadFinalState::AddSecondary(), G4Alpha::Alpha(), CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), G4PhotonEvaporation::BreakUpFragment(), G4HadFinalState::Clear(), G4Deuteron::Deuteron(), CLHEP::HepLorentzVector::e(), G4cout, G4endl, G4UniformRand, G4Gamma::Gamma(), G4HadProjectile::Get4Momentum(), G4Fragment::GetA_asInt(), G4Nucleus::GetA_asInt(), G4Fragment::GetCreationTime(), G4Ions::GetExcitationEnergy(), G4Fragment::GetExcitationEnergy(), G4HadProjectile::GetGlobalTime(), G4IonTable::GetIon(), G4HadProjectile::GetKineticEnergy(), G4Fragment::GetMomentum(), G4NucleiProperties::GetNuclearMass(), G4Fragment::GetParticleDefinition(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMass(), G4Fragment::GetZ_asInt(), G4Nucleus::GetZ_asInt(), G4He3::He3(), CLHEP::HepLorentzVector::mag(), python.hepunit::MeV, n, G4HadFinalState::SetStatusChange(), G4HadSecondary::SetTime(), stopAndKill, G4HadronicInteraction::theParticleChange, G4Triton::Triton(), CLHEP::Hep3Vector::unit(), CLHEP::HepLorentzVector::vect(), and G4HadronicInteraction::verboseLevel.

72 {
75 
76  G4int A = theNucleus.GetA_asInt();
77  G4int Z = theNucleus.GetZ_asInt();
78 
79  G4double time = aTrack.GetGlobalTime();
80 
81  // Create initial state
83  G4LorentzVector lv0(0.0,0.0,0.0,m1);
84  G4LorentzVector lv1 = aTrack.Get4Momentum() + lv0;
85 
86  // simplified method of 1 gamma emission
87  if(A <= 3) {
88 
89  G4ThreeVector bst = lv1.boostVector();
90  G4double M = lv1.mag();
91 
92  ++A;
94  if(M - mass <= lowestEnergyLimit) {
95  return &theParticleChange;
96  }
97 
98  if (verboseLevel > 1) {
99  G4cout << "G4NeutronRadCapture::DoIt: Eini(MeV)="
100  << aTrack.GetKineticEnergy()/MeV << " Eexc(MeV)= "
101  << (M - mass)/MeV
102  << " Z= " << Z << " A= " << A << G4endl;
103  }
104  G4double e1 = (M - mass)*(M + mass)/(2*M);
105  G4double cost = 2.0*G4UniformRand() - 1.0;
106  if(cost > 1.0) {cost = 1.0;}
107  else if(cost < -1.0) {cost = -1.0;}
108  G4double sint = std::sqrt((1. - cost)*(1.0 + cost));
109  G4double phi = G4UniformRand()*CLHEP::twopi;
110  G4LorentzVector lv2(e1*sint*std::cos(phi),e1*sint*std::sin(phi),e1*cost,e1);
111  lv2.boost(bst);
112  G4HadSecondary* news =
114  news->SetTime(time);
116  delete news;
117 
118  G4ParticleDefinition* theDef = 0;
119 
120  lv1 -= lv2;
121  if (Z == 1 && A == 2) {theDef = G4Deuteron::Deuteron();}
122  else if (Z == 1 && A == 3) {theDef = G4Triton::Triton();}
123  else if (Z == 2 && A == 3) {theDef = G4He3::He3();}
124  else if (Z == 2 && A == 4) {theDef = G4Alpha::Alpha();}
125  else { theDef = theTableOfIons->GetIon(Z,A,0); }
126 
127  if (verboseLevel > 1) {
128  G4cout << "Gamma 4-mom: " << lv2 << " "
129  << theDef->GetParticleName() << " " << lv1 << G4endl;
130  }
131  if(theDef) {
132  news = new G4HadSecondary(new G4DynamicParticle(theDef, lv1));
133  news->SetTime(time);
135  delete news;
136  }
137 
138  // Use photon evaporation
139  } else {
140 
141  G4Fragment* aFragment = new G4Fragment(A+1, Z, lv1);
142 
143  if (verboseLevel > 1) {
144  G4cout << "G4NeutronRadCapture::ApplyYourself initial G4Fragmet:" << G4endl;
145  G4cout << aFragment << G4endl;
146  }
147 
148  //
149  // Sample final state
150  //
151  G4FragmentVector* fv = photonEvaporation->BreakUpFragment(aFragment);
152  if(!fv) { fv = new G4FragmentVector(); }
153  fv->push_back(aFragment);
154  size_t n = fv->size();
155 
156  if (verboseLevel > 1) {
157  G4cout << "G4NeutronRadCapture: " << n << " final particle" << G4endl;
158  }
159  for(size_t i=0; i<n; ++i) {
160 
161  G4Fragment* f = (*fv)[i];
162  G4double etot = f->GetMomentum().e();
163 
164  Z = f->GetZ_asInt();
165  A = f->GetA_asInt();
166 
167  G4ParticleDefinition* theDef = 0;
168  if(0 == Z && 0 == A) {theDef = f->GetParticleDefinition();}
169  else if (Z == 1 && A == 2) {theDef = G4Deuteron::Deuteron();}
170  else if (Z == 1 && A == 3) {theDef = G4Triton::Triton();}
171  else if (Z == 2 && A == 3) {theDef = G4He3::He3();}
172  else if (Z == 2 && A == 4) {theDef = G4Alpha::Alpha();}
173  else {
174  G4double eexc = f->GetExcitationEnergy();
175  G4double excitation = eexc;
176  G4int level = 0;
177  theDef = theTableOfIons->GetIon(Z, A, level);
178  /*
179  G4cout << "### Find ion Z= " << theFragmentZ << " A= " << theFragmentA
180  << " Eexc(MeV)= " << excitation/MeV << " "
181  << theKindOfFragment << G4endl;
182  */
183  // production of an isomer
184  if(eexc > minExcitation) {
185  G4double elevel1 = 0.0;
186  G4double elevel2 = 0.0;
188  for(level=1; level<9; ++level) {
189  ion = theTableOfIons->GetIon(Z, A, level);
190  //G4cout << level << " " << ion << G4endl;
191  if(ion) {
192  G4Ions* ip = dynamic_cast<G4Ions*>(ion);
193  if(ip) {
194  elevel2 = ip->GetExcitationEnergy();
195  //G4cout<<" Level "<<level<<" E(MeV)= "<<elevel2/MeV<<G4endl;
196  // close level
197  if(std::fabs(eexc - elevel2) < minExcitation) {
198  excitation = eexc - elevel2;
199  theDef = ion;
200  break;
201  // previous level was closer
202  } else if(elevel2 - eexc >= eexc - elevel1) {
203  excitation = eexc - elevel1;
204  break;
205  // will check next level and save current
206  } else {
207  theDef = ion;
208  excitation = eexc - elevel2;
209  elevel1 = elevel2;
210  }
211  }
212  } else {
213  break;
214  }
215  }
216  }
217  // correction of total energy for ground state isotopes
218  etot += excitation;
219  etot -= theDef->GetPDGMass();
220  if(etot < 0.0) { etot = 0.0; }
221  }
222  if (verboseLevel > 1) {
223  G4cout << i << ". " << theDef->GetParticleName()
224  << " Ekin(MeV)= " << etot/MeV
225  << " p: " << f->GetMomentum().vect()
226  << G4endl;
227  }
228  if(theDef) {
229  G4HadSecondary* news =
230  new G4HadSecondary(new G4DynamicParticle(theDef,
231  f->GetMomentum().vect().unit(),
232  etot));
233  G4double timeF = f->GetCreationTime();
234  if(timeF < 0.0) { timeF = 0.0; }
235  news->SetTime(time + timeF);
237  delete news;
238  }
239  delete f;
240  }
241  delete fv;
242  }
243  return &theParticleChange;
244 }
Hep3Vector boostVector() const
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetExcitationEnergy() const
Definition: G4Ions.hh:113
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
void SetStatusChange(G4HadFinalStateStatus aS)
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4int GetA_asInt() const
Definition: G4Fragment.hh:238
Definition: G4Ions.hh:51
G4double GetCreationTime() const
Definition: G4Fragment.hh:398
double mag() const
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:271
G4double GetKineticEnergy() const
G4double GetGlobalTime() const
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
static G4Triton * Triton()
Definition: G4Triton.cc:95
const G4int n
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4LorentzVector & Get4Momentum() const
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
void SetTime(G4double aT)
G4double GetPDGMass() const
Hep3Vector unit() const
G4int GetZ_asInt() const
Definition: G4Fragment.hh:243
G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:388
#define G4endl
Definition: G4ios.hh:61
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
virtual G4FragmentVector * BreakUpFragment(G4Fragment *theNucleus)
double G4double
Definition: G4Types.hh:76
static G4He3 * He3()
Definition: G4He3.cc:94
void AddSecondary(G4DynamicParticle *aP)
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:255

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