Geant4  10.02.p03
G4NeutronRadCapture Class Reference

#include <G4NeutronRadCapture.hh>

Inheritance diagram for G4NeutronRadCapture:
Collaboration diagram for G4NeutronRadCapture:

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, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 

Private Member Functions

G4NeutronRadCaptureoperator= (const G4NeutronRadCapture &right)
 
 G4NeutronRadCapture (const G4NeutronRadCapture &)
 

Private Attributes

G4double lowestEnergyLimit
 
G4double minExcitation
 
G4VEvaporationChannelphotonEvaporation
 
G4IonTabletheTableOfIons
 
G4LorentzVector lab4mom
 

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 52 of file G4NeutronRadCapture.hh.

Constructor & Destructor Documentation

◆ G4NeutronRadCapture() [1/2]

G4NeutronRadCapture::G4NeutronRadCapture ( )

Definition at line 54 of file G4NeutronRadCapture.cc.

55  : G4HadronicInteraction("nRadCapture"),
56  lab4mom(0.,0.,0.,0.)
57 {
58  lowestEnergyLimit = 10*eV;
59  minExcitation = 1*keV;
60  SetMinEnergy( 0.0*GeV );
61  SetMaxEnergy( 100.*TeV );
62 
63  char* env = getenv("G4UsePhotonEvaporationOLD");
64  if(!env) { photonEvaporation = new G4PhotonEvaporation(); }
67 
69 }
virtual void SetICM(G4bool)
G4VEvaporationChannel * photonEvaporation
void SetMinEnergy(G4double anEnergy)
G4IonTable * GetIonTable() const
static const double GeV
Definition: G4SIunits.hh:214
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static const double eV
Definition: G4SIunits.hh:212
static G4ParticleTable * GetParticleTable()
void SetMaxEnergy(const G4double anEnergy)
static const double TeV
Definition: G4SIunits.hh:215
static const double keV
Definition: G4SIunits.hh:213
Here is the call graph for this function:

◆ ~G4NeutronRadCapture()

G4NeutronRadCapture::~G4NeutronRadCapture ( )
virtual

Definition at line 71 of file G4NeutronRadCapture.cc.

72 {
73  delete photonEvaporation;
74 }
G4VEvaporationChannel * photonEvaporation

◆ G4NeutronRadCapture() [2/2]

G4NeutronRadCapture::G4NeutronRadCapture ( const G4NeutronRadCapture )
private

Member Function Documentation

◆ ApplyYourself()

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

Implements G4HadronicInteraction.

Definition at line 76 of file G4NeutronRadCapture.cc.

78 {
81 
82  G4int A = theNucleus.GetA_asInt();
83  G4int Z = theNucleus.GetZ_asInt();
84 
85  G4double time = aTrack.GetGlobalTime();
86 
87  // Create initial state
89  lab4mom += aTrack.Get4Momentum();
90 
91  G4double M = lab4mom.mag();
92  ++A;
94  //G4cout << "Capture start: Z= " << Z << " A= " << A
95  // << " LabM= " << M << " Mcompound= " << mass << G4endl;
96 
97  // simplified method of 1 gamma emission
98  if(A <= 4) {
99 
101 
102  if(M - mass <= lowestEnergyLimit) {
103  return &theParticleChange;
104  }
105 
106  if (verboseLevel > 1) {
107  G4cout << "G4NeutronRadCapture::DoIt: Eini(MeV)="
108  << aTrack.GetKineticEnergy()/MeV << " Eexc(MeV)= "
109  << (M - mass)/MeV
110  << " Z= " << Z << " A= " << A << G4endl;
111  }
112  G4double e1 = (M - mass)*(M + mass)/(2*M);
113 
114  G4double cost = 2.0*G4UniformRand() - 1.0;
115  if(cost > 1.0) {cost = 1.0;}
116  else if(cost < -1.0) {cost = -1.0;}
117  G4double sint = std::sqrt((1. - cost)*(1.0 + cost));
119 
120  G4LorentzVector lv2(e1*sint*std::cos(phi),e1*sint*std::sin(phi),
121  e1*cost,e1);
122  lv2.boost(bst);
123  G4HadSecondary* news =
125  news->SetTime(time);
127  delete news;
128 
129  const G4ParticleDefinition* theDef = 0;
130 
131  lab4mom -= lv2;
132  if (Z == 1 && A == 2) {theDef = G4Deuteron::Deuteron();}
133  else if (Z == 1 && A == 3) {theDef = G4Triton::Triton();}
134  else if (Z == 2 && A == 3) {theDef = G4He3::He3();}
135  else if (Z == 2 && A == 4) {theDef = G4Alpha::Alpha();}
136  else { theDef = theTableOfIons->GetIon(Z,A,0.0); }
137 
138  if (verboseLevel > 1) {
139  G4cout << "Gamma 4-mom: " << lv2 << " "
140  << theDef->GetParticleName() << " " << lab4mom << G4endl;
141  }
142  if(theDef) {
143  news = new G4HadSecondary(new G4DynamicParticle(theDef, lab4mom));
144  news->SetTime(time);
146  delete news;
147  }
148 
149  // Use photon evaporation
150  } else {
151 
152  // protection against wrong kinematic
153  if(M < mass) {
154  G4double etot = std::max(mass, lab4mom.e());
155  G4double ptot = std::sqrt((etot - mass)*(etot + mass));
157  lab4mom.set(v.x()*ptot,v.y()*ptot,v.z()*ptot,etot);
158  }
159 
160  G4Fragment* aFragment = new G4Fragment(A, Z, lab4mom);
161 
162  if (verboseLevel > 1) {
163  G4cout << "G4NeutronRadCapture::ApplyYourself initial G4Fragmet:"
164  << G4endl;
165  G4cout << aFragment << G4endl;
166  }
167 
168  //
169  // Sample final state
170  //
172  if(!fv) { fv = new G4FragmentVector(); }
173  fv->push_back(aFragment);
174  size_t n = fv->size();
175 
176  if (verboseLevel > 1) {
177  G4cout << "G4NeutronRadCapture: " << n << " final particle" << G4endl;
178  }
179  for(size_t i=0; i<n; ++i) {
180 
181  G4Fragment* f = (*fv)[i];
182  G4double etot = f->GetMomentum().e();
183 
184  Z = f->GetZ_asInt();
185  A = f->GetA_asInt();
186 
187  const G4ParticleDefinition* theDef;
188  if(0 == Z && 0 == A) {theDef = f->GetParticleDefinition();}
189  else if (Z == 1 && A == 2) {theDef = G4Deuteron::Deuteron();}
190  else if (Z == 1 && A == 3) {theDef = G4Triton::Triton();}
191  else if (Z == 2 && A == 3) {theDef = G4He3::He3();}
192  else if (Z == 2 && A == 4) {theDef = G4Alpha::Alpha();}
193  else {
194  G4double eexc = f->GetExcitationEnergy();
195  if(eexc <= minExcitation) { eexc = 0.0; }
196  theDef = theTableOfIons->GetIon(Z, A, eexc);
197  /*
198  G4cout << "### Find ion Z= " << Z << " A= " << A
199  << " Eexc(MeV)= " << eexc/MeV << " "
200  << theDef << G4endl;
201  */
202  }
203  G4double ekin = std::max(0.0,etot - theDef->GetPDGMass());
204  if (verboseLevel > 1) {
205  G4cout << i << ". " << theDef->GetParticleName()
206  << " Ekin(MeV)= " << etot/MeV
207  << " p: " << f->GetMomentum().vect()
208  << G4endl;
209  }
210  G4HadSecondary* news = new G4HadSecondary(
211  new G4DynamicParticle(theDef,
212  f->GetMomentum().vect().unit(),
213  ekin));
214  G4double timeF = f->GetCreationTime();
215  if(timeF < 0.0) { timeF = 0.0; }
216  news->SetTime(time + timeF);
218  delete news;
219  delete f;
220  }
221  delete fv;
222  }
223  return &theParticleChange;
224 }
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:273
static const double MeV
Definition: G4SIunits.hh:211
static G4double GetNuclearMass(const G4double A, const G4double Z)
virtual G4FragmentVector * BreakUpFragment(G4Fragment *theNucleus)
const G4LorentzVector & Get4Momentum() const
G4int GetA_asInt() const
Definition: G4Fragment.hh:256
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
int G4int
Definition: G4Types.hh:78
G4VEvaporationChannel * photonEvaporation
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:410
void SetStatusChange(G4HadFinalStateStatus aS)
Hep3Vector vect() const
Char_t n[5]
const G4String & GetParticleName() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4int GetZ_asInt() const
Definition: G4Fragment.hh:261
Float_t Z
Hep3Vector unit() const
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:63
static G4Triton * Triton()
Definition: G4Triton.cc:95
G4double GetGlobalTime() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static const G4double e1
double x() const
void SetTime(G4double aT)
G4double GetKineticEnergy() const
void set(double x, double y, double z, double t)
double y() const
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:289
Hep3Vector boostVector() const
double z() const
#define G4endl
Definition: G4ios.hh:61
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
double G4double
Definition: G4Types.hh:76
static const double twopi
Definition: SystemOfUnits.h:54
static G4He3 * He3()
Definition: G4He3.cc:94
G4double GetCreationTime() const
Definition: G4Fragment.hh:420
Here is the call graph for this function:

◆ operator=()

G4NeutronRadCapture& G4NeutronRadCapture::operator= ( const G4NeutronRadCapture right)
private

Member Data Documentation

◆ lab4mom

G4LorentzVector G4NeutronRadCapture::lab4mom
private

Definition at line 72 of file G4NeutronRadCapture.hh.

◆ lowestEnergyLimit

G4double G4NeutronRadCapture::lowestEnergyLimit
private

Definition at line 68 of file G4NeutronRadCapture.hh.

◆ minExcitation

G4double G4NeutronRadCapture::minExcitation
private

Definition at line 69 of file G4NeutronRadCapture.hh.

◆ photonEvaporation

G4VEvaporationChannel* G4NeutronRadCapture::photonEvaporation
private

Definition at line 70 of file G4NeutronRadCapture.hh.

◆ theTableOfIons

G4IonTable* G4NeutronRadCapture::theTableOfIons
private

Definition at line 71 of file G4NeutronRadCapture.hh.


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