Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Nucleus Class Reference

#include <G4Nucleus.hh>

+ Inheritance diagram for G4Nucleus:

Public Member Functions

 G4Nucleus ()
 
 G4Nucleus (const G4double A, const G4double Z)
 
 G4Nucleus (const G4int A, const G4int Z)
 
 G4Nucleus (const G4Material *aMaterial)
 
 ~G4Nucleus ()
 
 G4Nucleus (const G4Nucleus &right)
 
G4Nucleusoperator= (const G4Nucleus &right)
 
G4bool operator== (const G4Nucleus &right) const
 
G4bool operator!= (const G4Nucleus &right) const
 
void ChooseParameters (const G4Material *aMaterial)
 
void SetParameters (const G4double A, const G4double Z)
 
void SetParameters (const G4int A, const G4int Z)
 
G4int GetA_asInt () const
 
G4int GetN_asInt () const
 
G4int GetZ_asInt () const
 
const G4IsotopeGetIsotope ()
 
void SetIsotope (const G4Isotope *iso)
 
G4DynamicParticleReturnTargetParticle () const
 
G4double AtomicMass (const G4double A, const G4double Z) const
 
G4double AtomicMass (const G4int A, const G4int Z) const
 
G4double GetThermalPz (const G4double mass, const G4double temp) const
 
G4ReactionProduct GetThermalNucleus (G4double aMass, G4double temp=-1) const
 
G4ReactionProduct GetBiasedThermalNucleus (G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
 
G4double Cinema (G4double kineticEnergy)
 
G4double EvaporationEffects (G4double kineticEnergy)
 
G4double AnnihilationEvaporationEffects (G4double kineticEnergy, G4double ekOrg)
 
G4double GetPNBlackTrackEnergy () const
 
G4double GetDTABlackTrackEnergy () const
 
G4double GetAnnihilationPNBlackTrackEnergy () const
 
G4double GetAnnihilationDTABlackTrackEnergy () const
 
G4ThreeVector GetFermiMomentum ()
 
G4ReactionProductVectorFragmentate ()
 
void AddExcitationEnergy (G4double anEnergy)
 
void AddMomentum (const G4ThreeVector aMomentum)
 
G4double GetEnergyDeposit ()
 

Detailed Description

Definition at line 50 of file G4Nucleus.hh.

Constructor & Destructor Documentation

◆ G4Nucleus() [1/5]

G4Nucleus::G4Nucleus ( )

Definition at line 49 of file G4Nucleus.cc.

50 : theA(0), theZ(0), aEff(0.0), zEff(0)
51{
52 pnBlackTrackEnergy = 0.0;
53 dtaBlackTrackEnergy = 0.0;
54 pnBlackTrackEnergyfromAnnihilation = 0.0;
55 dtaBlackTrackEnergyfromAnnihilation = 0.0;
56 excitationEnergy = 0.0;
57 momentum = G4ThreeVector(0.,0.,0.);
58 fermiMomentum = 1.52*hbarc/fermi;
59 theTemp = 293.16*kelvin;
60 fIsotope = 0;
61}
CLHEP::Hep3Vector G4ThreeVector

◆ G4Nucleus() [2/5]

G4Nucleus::G4Nucleus ( const G4double  A,
const G4double  Z 
)

Definition at line 63 of file G4Nucleus.cc.

64{
65 SetParameters( A, Z );
66 pnBlackTrackEnergy = 0.0;
67 dtaBlackTrackEnergy = 0.0;
68 pnBlackTrackEnergyfromAnnihilation = 0.0;
69 dtaBlackTrackEnergyfromAnnihilation = 0.0;
70 excitationEnergy = 0.0;
71 momentum = G4ThreeVector(0.,0.,0.);
72 fermiMomentum = 1.52*hbarc/fermi;
73 theTemp = 293.16*kelvin;
74 fIsotope = 0;
75}
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:198

◆ G4Nucleus() [3/5]

G4Nucleus::G4Nucleus ( const G4int  A,
const G4int  Z 
)

Definition at line 77 of file G4Nucleus.cc.

78{
79 SetParameters( A, Z );
80 pnBlackTrackEnergy = 0.0;
81 dtaBlackTrackEnergy = 0.0;
82 pnBlackTrackEnergyfromAnnihilation = 0.0;
83 dtaBlackTrackEnergyfromAnnihilation = 0.0;
84 excitationEnergy = 0.0;
85 momentum = G4ThreeVector(0.,0.,0.);
86 fermiMomentum = 1.52*hbarc/fermi;
87 theTemp = 293.16*kelvin;
88 fIsotope = 0;
89}

◆ G4Nucleus() [4/5]

G4Nucleus::G4Nucleus ( const G4Material aMaterial)

Definition at line 91 of file G4Nucleus.cc.

92{
93 ChooseParameters( aMaterial );
94 pnBlackTrackEnergy = 0.0;
95 dtaBlackTrackEnergy = 0.0;
96 pnBlackTrackEnergyfromAnnihilation = 0.0;
97 dtaBlackTrackEnergyfromAnnihilation = 0.0;
98 excitationEnergy = 0.0;
99 momentum = G4ThreeVector(0.,0.,0.);
100 fermiMomentum = 1.52*hbarc/fermi;
101 theTemp = aMaterial->GetTemperature();
102 fIsotope = 0;
103}
G4double GetTemperature() const
Definition: G4Material.hh:181
void ChooseParameters(const G4Material *aMaterial)
Definition: G4Nucleus.cc:158

◆ ~G4Nucleus()

G4Nucleus::~G4Nucleus ( )

Definition at line 105 of file G4Nucleus.cc.

105{}

◆ G4Nucleus() [5/5]

G4Nucleus::G4Nucleus ( const G4Nucleus right)
inline

Definition at line 61 of file G4Nucleus.hh.

62 { *this = right; }

Member Function Documentation

◆ AddExcitationEnergy()

void G4Nucleus::AddExcitationEnergy ( G4double  anEnergy)

Definition at line 435 of file G4Nucleus.cc.

436 {
437 excitationEnergy+=anEnergy;
438 }

Referenced by G4KaonMinusAbsorptionAtRest::AtRestDoIt(), and G4BertiniEvaporation::BreakItUp().

◆ AddMomentum()

void G4Nucleus::AddMomentum ( const G4ThreeVector  aMomentum)

Definition at line 430 of file G4Nucleus.cc.

431 {
432 momentum+=(aMomentum);
433 }

◆ AnnihilationEvaporationEffects()

G4double G4Nucleus::AnnihilationEvaporationEffects ( G4double  kineticEnergy,
G4double  ekOrg 
)

Definition at line 323 of file G4Nucleus.cc.

324 {
325 // Nuclear evaporation as a function of atomic number and kinetic
326 // energy (MeV) of primary particle. Modified for annihilation effects.
327 //
328 if( aEff < 1.5 || ekOrg < 0.)
329 {
330 pnBlackTrackEnergyfromAnnihilation = 0.0;
331 dtaBlackTrackEnergyfromAnnihilation = 0.0;
332 return 0.0;
333 }
334 G4double ek = kineticEnergy/GeV;
335 G4float ekin = std::min( 4.0, std::max( 0.1, ek ) );
336 const G4float atno = std::min( 120., aEff );
337 const G4float gfa = 2.0*((aEff-1.0)/70.)*std::exp(-(aEff-1.0)/70.);
338
339 G4float cfa = std::max( 0.15, 0.35 + ((0.35-0.05)/2.3)*std::log(ekin) );
340 G4float exnu = 7.716 * cfa * std::exp(-cfa)
341 * ((atno-1.0)/120.)*std::exp(-(atno-1.0)/120.);
342 G4float fpdiv = std::max( 0.5, 1.0-0.25*ekin*ekin );
343
344 pnBlackTrackEnergyfromAnnihilation = exnu*fpdiv;
345 dtaBlackTrackEnergyfromAnnihilation = exnu*(1.0-fpdiv);
346
347 G4double ran1 = -6.0;
348 G4double ran2 = -6.0;
349 for( G4int i=0; i<12; ++i ) {
350 ran1 += G4UniformRand();
351 ran2 += G4UniformRand();
352 }
353 pnBlackTrackEnergyfromAnnihilation *= 1.0 + ran1*gfa;
354 dtaBlackTrackEnergyfromAnnihilation *= 1.0 + ran2*gfa;
355
356 pnBlackTrackEnergyfromAnnihilation = std::max( 0.0, pnBlackTrackEnergyfromAnnihilation);
357 dtaBlackTrackEnergyfromAnnihilation = std::max( 0.0, dtaBlackTrackEnergyfromAnnihilation);
358 G4double blackSum = pnBlackTrackEnergyfromAnnihilation+dtaBlackTrackEnergyfromAnnihilation;
359 if (blackSum >= ekOrg/GeV) {
360 pnBlackTrackEnergyfromAnnihilation *= ekOrg/GeV/blackSum;
361 dtaBlackTrackEnergyfromAnnihilation *= ekOrg/GeV/blackSum;
362 }
363
364 return (pnBlackTrackEnergyfromAnnihilation+dtaBlackTrackEnergyfromAnnihilation)*GeV;
365 }
double G4double
Definition: G4Types.hh:64
float G4float
Definition: G4Types.hh:65
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53

Referenced by G4RPGInelastic::CalculateMomenta(), and G4InelasticInteraction::CalculateMomenta().

◆ AtomicMass() [1/2]

G4double G4Nucleus::AtomicMass ( const G4double  A,
const G4double  Z 
) const

Definition at line 240 of file G4Nucleus.cc.

241 {
242 // Now returns (atomic mass - electron masses)
244 }
static G4double GetNuclearMass(const G4double A, const G4double Z)

Referenced by G4WilsonAbrasionModel::ApplyYourself(), G4LEAlphaInelastic::ApplyYourself(), G4LEDeuteronInelastic::ApplyYourself(), and G4LETritonInelastic::ApplyYourself().

◆ AtomicMass() [2/2]

G4double G4Nucleus::AtomicMass ( const G4int  A,
const G4int  Z 
) const

Definition at line 247 of file G4Nucleus.cc.

248 {
249 // Now returns (atomic mass - electron masses)
251 }

◆ ChooseParameters()

void G4Nucleus::ChooseParameters ( const G4Material aMaterial)

Definition at line 158 of file G4Nucleus.cc.

159{
160 G4double random = G4UniformRand();
161 G4double sum = aMaterial->GetTotNbOfAtomsPerVolume();
162 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
163 G4double running(0);
164 // G4Element* element(0);
165 G4Element* element = (*theElementVector)[aMaterial->GetNumberOfElements()-1];
166
167 for (unsigned int i = 0; i < aMaterial->GetNumberOfElements(); ++i) {
168 running += aMaterial->GetVecNbOfAtomsPerVolume()[i];
169 if (running > random*sum) {
170 element = (*theElementVector)[i];
171 break;
172 }
173 }
174
175 if (element->GetNumberOfIsotopes() > 0) {
176 G4double randomAbundance = G4UniformRand();
177 G4double sumAbundance = element->GetRelativeAbundanceVector()[0];
178 unsigned int iso=0;
179 while (iso < element->GetNumberOfIsotopes() &&
180 sumAbundance < randomAbundance) {
181 ++iso;
182 sumAbundance += element->GetRelativeAbundanceVector()[iso];
183 }
184 theA=element->GetIsotope(iso)->GetN();
185 theZ=element->GetIsotope(iso)->GetZ();
186 aEff=theA;
187 zEff=theZ;
188 } else {
189 aEff = element->GetN();
190 zEff = element->GetZ();
191 theZ = G4int(zEff + 0.5);
192 theA = G4int(aEff + 0.5);
193 }
194}
std::vector< G4Element * > G4ElementVector
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
G4double GetZ() const
Definition: G4Element.hh:131
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4double GetN() const
Definition: G4Element.hh:134
G4int GetZ() const
Definition: G4Isotope.hh:91
G4int GetN() const
Definition: G4Isotope.hh:94
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:208
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:205

Referenced by G4Nucleus().

◆ Cinema()

G4double G4Nucleus::Cinema ( G4double  kineticEnergy)

Definition at line 368 of file G4Nucleus.cc.

369 {
370 // derived from original FORTRAN code CINEMA by H. Fesefeldt (14-Oct-1987)
371 //
372 // input: kineticEnergy (MeV)
373 // returns modified kinetic energy (MeV)
374 //
375 static const G4double expxu = 82.; // upper bound for arg. of exp
376 static const G4double expxl = -expxu; // lower bound for arg. of exp
377
378 G4double ek = kineticEnergy/GeV;
379 G4double ekLog = std::log( ek );
380 G4double aLog = std::log( aEff );
381 G4double em = std::min( 1.0, 0.2390 + 0.0408*aLog*aLog );
382 G4double temp1 = -ek * std::min( 0.15, 0.0019*aLog*aLog*aLog );
383 G4double temp2 = std::exp( std::max( expxl, std::min( expxu, -(ekLog-em)*(ekLog-em)*2.0 ) ) );
384 G4double result = 0.0;
385 if( std::abs( temp1 ) < 1.0 )
386 {
387 if( temp2 > 1.0e-10 )result = temp1*temp2;
388 }
389 else result = temp1*temp2;
390 if( result < -ek )result = -ek;
391 return result*GeV;
392 }

Referenced by G4LEAntiKaonZeroInelastic::ApplyYourself(), G4LEAntiLambdaInelastic::ApplyYourself(), G4LEAntiNeutronInelastic::ApplyYourself(), G4LEAntiOmegaMinusInelastic::ApplyYourself(), G4LEAntiProtonInelastic::ApplyYourself(), G4LEAntiSigmaMinusInelastic::ApplyYourself(), G4LEAntiSigmaPlusInelastic::ApplyYourself(), G4LEAntiXiMinusInelastic::ApplyYourself(), G4LEAntiXiZeroInelastic::ApplyYourself(), G4LEKaonMinusInelastic::ApplyYourself(), G4LEKaonPlusInelastic::ApplyYourself(), G4LEKaonZeroInelastic::ApplyYourself(), G4LELambdaInelastic::ApplyYourself(), G4LENeutronInelastic::ApplyYourself(), G4LEOmegaMinusInelastic::ApplyYourself(), G4LEPionMinusInelastic::ApplyYourself(), G4LEPionPlusInelastic::ApplyYourself(), G4LEProtonInelastic::ApplyYourself(), G4LESigmaMinusInelastic::ApplyYourself(), G4LESigmaPlusInelastic::ApplyYourself(), G4LEXiMinusInelastic::ApplyYourself(), G4LEXiZeroInelastic::ApplyYourself(), G4RPGAntiKZeroInelastic::ApplyYourself(), G4RPGAntiLambdaInelastic::ApplyYourself(), G4RPGAntiNeutronInelastic::ApplyYourself(), G4RPGAntiOmegaMinusInelastic::ApplyYourself(), G4RPGAntiProtonInelastic::ApplyYourself(), G4RPGAntiSigmaMinusInelastic::ApplyYourself(), G4RPGAntiSigmaPlusInelastic::ApplyYourself(), G4RPGAntiXiMinusInelastic::ApplyYourself(), G4RPGAntiXiZeroInelastic::ApplyYourself(), G4RPGKMinusInelastic::ApplyYourself(), G4RPGKPlusInelastic::ApplyYourself(), G4RPGKZeroInelastic::ApplyYourself(), G4RPGLambdaInelastic::ApplyYourself(), G4RPGNeutronInelastic::ApplyYourself(), G4RPGOmegaMinusInelastic::ApplyYourself(), G4RPGPiMinusInelastic::ApplyYourself(), G4RPGPiPlusInelastic::ApplyYourself(), G4RPGProtonInelastic::ApplyYourself(), G4RPGSigmaMinusInelastic::ApplyYourself(), G4RPGSigmaPlusInelastic::ApplyYourself(), G4RPGXiMinusInelastic::ApplyYourself(), G4RPGXiZeroInelastic::ApplyYourself(), G4RPGInelastic::CalculateMomenta(), and G4InelasticInteraction::CalculateMomenta().

◆ EvaporationEffects()

G4double G4Nucleus::EvaporationEffects ( G4double  kineticEnergy)

Definition at line 264 of file G4Nucleus.cc.

265 {
266 // derived from original FORTRAN code EXNU by H. Fesefeldt (10-Dec-1986)
267 //
268 // Nuclear evaporation as function of atomic number
269 // and kinetic energy (MeV) of primary particle
270 //
271 // returns kinetic energy (MeV)
272 //
273 if( aEff < 1.5 )
274 {
275 pnBlackTrackEnergy = dtaBlackTrackEnergy = 0.0;
276 return 0.0;
277 }
278 G4double ek = kineticEnergy/GeV;
279 G4float ekin = std::min( 4.0, std::max( 0.1, ek ) );
280 const G4float atno = std::min( 120., aEff );
281 const G4float gfa = 2.0*((aEff-1.0)/70.)*std::exp(-(aEff-1.0)/70.);
282 //
283 // 0.35 value at 1 GeV
284 // 0.05 value at 0.1 GeV
285 //
286 G4float cfa = std::max( 0.15, 0.35 + ((0.35-0.05)/2.3)*std::log(ekin) );
287 G4float exnu = 7.716 * cfa * std::exp(-cfa)
288 * ((atno-1.0)/120.)*std::exp(-(atno-1.0)/120.);
289 G4float fpdiv = std::max( 0.5, 1.0-0.25*ekin*ekin );
290 //
291 // pnBlackTrackEnergy is the kinetic energy (in GeV) available for
292 // proton/neutron black track particles
293 // dtaBlackTrackEnergy is the kinetic energy (in GeV) available for
294 // deuteron/triton/alpha black track particles
295 //
296 pnBlackTrackEnergy = exnu*fpdiv;
297 dtaBlackTrackEnergy = exnu*(1.0-fpdiv);
298
299 if( G4int(zEff+0.1) != 82 )
300 {
301 G4double ran1 = -6.0;
302 G4double ran2 = -6.0;
303 for( G4int i=0; i<12; ++i )
304 {
305 ran1 += G4UniformRand();
306 ran2 += G4UniformRand();
307 }
308 pnBlackTrackEnergy *= 1.0 + ran1*gfa;
309 dtaBlackTrackEnergy *= 1.0 + ran2*gfa;
310 }
311 pnBlackTrackEnergy = std::max( 0.0, pnBlackTrackEnergy );
312 dtaBlackTrackEnergy = std::max( 0.0, dtaBlackTrackEnergy );
313 while( pnBlackTrackEnergy+dtaBlackTrackEnergy >= ek )
314 {
315 pnBlackTrackEnergy *= 1.0 - 0.5*G4UniformRand();
316 dtaBlackTrackEnergy *= 1.0 - 0.5*G4UniformRand();
317 }
318// G4cout << "EvaporationEffects "<<kineticEnergy<<" "
319// <<pnBlackTrackEnergy+dtaBlackTrackEnergy<<endl;
320 return (pnBlackTrackEnergy+dtaBlackTrackEnergy)*GeV;
321 }

Referenced by G4LEAntiKaonZeroInelastic::ApplyYourself(), G4LEAntiLambdaInelastic::ApplyYourself(), G4LEAntiNeutronInelastic::ApplyYourself(), G4LEAntiOmegaMinusInelastic::ApplyYourself(), G4LEAntiProtonInelastic::ApplyYourself(), G4LEAntiSigmaMinusInelastic::ApplyYourself(), G4LEAntiSigmaPlusInelastic::ApplyYourself(), G4LEAntiXiMinusInelastic::ApplyYourself(), G4LEAntiXiZeroInelastic::ApplyYourself(), G4LEKaonMinusInelastic::ApplyYourself(), G4LEKaonPlusInelastic::ApplyYourself(), G4LEKaonZeroInelastic::ApplyYourself(), G4LELambdaInelastic::ApplyYourself(), G4LENeutronInelastic::ApplyYourself(), G4LEOmegaMinusInelastic::ApplyYourself(), G4LEPionMinusInelastic::ApplyYourself(), G4LEPionPlusInelastic::ApplyYourself(), G4LEProtonInelastic::ApplyYourself(), G4LESigmaMinusInelastic::ApplyYourself(), G4LESigmaPlusInelastic::ApplyYourself(), G4LEXiMinusInelastic::ApplyYourself(), G4LEXiZeroInelastic::ApplyYourself(), G4RPGAntiKZeroInelastic::ApplyYourself(), G4RPGAntiLambdaInelastic::ApplyYourself(), G4RPGAntiNeutronInelastic::ApplyYourself(), G4RPGAntiOmegaMinusInelastic::ApplyYourself(), G4RPGAntiProtonInelastic::ApplyYourself(), G4RPGAntiSigmaMinusInelastic::ApplyYourself(), G4RPGAntiSigmaPlusInelastic::ApplyYourself(), G4RPGAntiXiMinusInelastic::ApplyYourself(), G4RPGAntiXiZeroInelastic::ApplyYourself(), G4RPGKMinusInelastic::ApplyYourself(), G4RPGKPlusInelastic::ApplyYourself(), G4RPGKZeroInelastic::ApplyYourself(), G4RPGLambdaInelastic::ApplyYourself(), G4RPGNeutronInelastic::ApplyYourself(), G4RPGOmegaMinusInelastic::ApplyYourself(), G4RPGPiMinusInelastic::ApplyYourself(), G4RPGPiPlusInelastic::ApplyYourself(), G4RPGProtonInelastic::ApplyYourself(), G4RPGSigmaMinusInelastic::ApplyYourself(), G4RPGSigmaPlusInelastic::ApplyYourself(), G4RPGXiMinusInelastic::ApplyYourself(), and G4RPGXiZeroInelastic::ApplyYourself().

◆ Fragmentate()

G4ReactionProductVector * G4Nucleus::Fragmentate ( )

Definition at line 424 of file G4Nucleus.cc.

425 {
426 // needs implementation!
427 return NULL;
428 }

◆ GetA_asInt()

G4int G4Nucleus::GetA_asInt ( ) const
inline

Definition at line 109 of file G4Nucleus.hh.

110 { return theA; }

Referenced by G4RPGReaction::AddBlackTrackParticles(), G4WilsonAbrasionModel::ApplyYourself(), G4EMDissociation::ApplyYourself(), G4LENDCapture::ApplyYourself(), G4LENDElastic::ApplyYourself(), G4LENDFission::ApplyYourself(), G4LENDInelastic::ApplyYourself(), G4LENDModel::ApplyYourself(), G4ElectroVDNuclearModel::ApplyYourself(), G4ChiralInvariantPhaseSpace::ApplyYourself(), G4ChargeExchange::ApplyYourself(), G4HadronElastic::ApplyYourself(), G4LEnp::ApplyYourself(), G4LEpp::ApplyYourself(), G4NeutronRadCapture::ApplyYourself(), G4HEAntiKaonZeroInelastic::ApplyYourself(), G4HEAntiLambdaInelastic::ApplyYourself(), G4HEAntiNeutronInelastic::ApplyYourself(), G4HEAntiOmegaMinusInelastic::ApplyYourself(), G4HEAntiProtonInelastic::ApplyYourself(), G4HEAntiSigmaMinusInelastic::ApplyYourself(), G4HEAntiSigmaPlusInelastic::ApplyYourself(), G4HEAntiXiMinusInelastic::ApplyYourself(), G4HEAntiXiZeroInelastic::ApplyYourself(), G4HEKaonMinusInelastic::ApplyYourself(), G4HEKaonPlusInelastic::ApplyYourself(), G4HEKaonZeroInelastic::ApplyYourself(), G4HEKaonZeroLongInelastic::ApplyYourself(), G4HEKaonZeroShortInelastic::ApplyYourself(), G4HELambdaInelastic::ApplyYourself(), G4HENeutronInelastic::ApplyYourself(), G4HEOmegaMinusInelastic::ApplyYourself(), G4HEPionMinusInelastic::ApplyYourself(), G4HEPionPlusInelastic::ApplyYourself(), G4HEProtonInelastic::ApplyYourself(), G4HESigmaMinusInelastic::ApplyYourself(), G4HESigmaPlusInelastic::ApplyYourself(), G4HEXiMinusInelastic::ApplyYourself(), G4HEXiZeroInelastic::ApplyYourself(), G4LCapture::ApplyYourself(), G4LEAlphaInelastic::ApplyYourself(), G4LEDeuteronInelastic::ApplyYourself(), G4LElastic::ApplyYourself(), G4LETritonInelastic::ApplyYourself(), G4LFission::ApplyYourself(), G4QMDReaction::ApplyYourself(), G4EmCaptureCascade::ApplyYourself(), G4MuMinusCapturePrecompound::ApplyYourself(), G4MuonMinusBoundDecay::ApplyYourself(), G4BinaryCascade::ApplyYourself(), G4BinaryLightIonReaction::ApplyYourself(), G4CascadeInterface::ApplyYourself(), G4ParaFissionModel::ApplyYourself(), G4INCLXXInterface::ApplyYourself(), G4LowEIonFragmentation::ApplyYourself(), G4PreCompoundModel::ApplyYourself(), G4TheoFSGenerator::ApplyYourself(), G4HadronStoppingProcess::AtRestDoIt(), G4KaonMinusAbsorptionAtRest::AtRestDoIt(), G4BertiniEvaporation::BreakItUp(), G4RPGInelastic::CalculateMomenta(), G4InelasticInteraction::CalculateMomenta(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4HadronicProcess::CheckResult(), G4CascadeInterface::createTarget(), G4InelasticInteraction::ExtractResidualNucleus(), G4ReactionDynamics::GenerateXandPt(), G4QuasiElasticChannel::GetFraction(), G4FTFModel::Init(), G4RPGReaction::NuclearReaction(), G4ReactionDynamics::NuclearReaction(), G4HadronicProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4WHadronElasticProcess::PostStepDoIt(), G4RPGFragmentation::ReactionStage(), G4RPGPionSuppression::ReactionStage(), G4RPGTwoBody::ReactionStage(), G4RPGTwoCluster::ReactionStage(), G4RPGReaction::Rotate(), G4ProjectileDiffractiveChannel::Scatter(), G4QuasiElasticChannel::Scatter(), G4HadronicWhiteBoard::SetTargetNucleus(), G4ReactionDynamics::SuppressChargedPions(), G4ReactionDynamics::TwoBody(), and G4ReactionDynamics::TwoCluster().

◆ GetAnnihilationDTABlackTrackEnergy()

G4double G4Nucleus::GetAnnihilationDTABlackTrackEnergy ( ) const
inline

Definition at line 159 of file G4Nucleus.hh.

160 { return dtaBlackTrackEnergyfromAnnihilation; }

Referenced by G4ReactionDynamics::GenerateXandPt(), G4RPGFragmentation::ReactionStage(), G4RPGTwoCluster::ReactionStage(), and G4ReactionDynamics::TwoCluster().

◆ GetAnnihilationPNBlackTrackEnergy()

G4double G4Nucleus::GetAnnihilationPNBlackTrackEnergy ( ) const
inline

Definition at line 156 of file G4Nucleus.hh.

157 { return pnBlackTrackEnergyfromAnnihilation; }

Referenced by G4ReactionDynamics::GenerateXandPt(), G4RPGFragmentation::ReactionStage(), G4RPGTwoCluster::ReactionStage(), and G4ReactionDynamics::TwoCluster().

◆ GetBiasedThermalNucleus()

G4ReactionProduct G4Nucleus::GetBiasedThermalNucleus ( G4double  aMass,
G4ThreeVector  aVelocity,
G4double  temp = -1 
) const

Definition at line 107 of file G4Nucleus.cc.

109{
110 G4double velMag = aVelocity.mag();
111 G4ReactionProduct result;
112 G4double value = 0;
113 G4double random = 1;
114 G4double norm = 3.*std::sqrt(k_Boltzmann*temp*aMass*G4Neutron::Neutron()->GetPDGMass());
115 norm /= G4Neutron::Neutron()->GetPDGMass();
116 norm *= 5.;
117 norm += velMag;
118 norm /= velMag;
119 while(value/norm<random)
120 {
121 result = GetThermalNucleus(aMass, temp);
122 G4ThreeVector targetVelocity = 1./result.GetMass()*result.GetMomentum();
123 value = (targetVelocity+aVelocity).mag()/velMag;
124 random = G4UniformRand();
125 }
126 return result;
127}
double mag() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:130
G4ThreeVector GetMomentum() const
G4double GetMass() const

Referenced by G4FissionLibrary::ApplyYourself(), G4NeutronHPCaptureFS::ApplyYourself(), G4NeutronHPElasticFS::ApplyYourself(), G4NeutronHPFissionFS::ApplyYourself(), G4NeutronHPInelasticBaseFS::BaseApply(), G4NeutronHPInelasticCompFS::CompositeApply(), and G4NeutronHPThermalBoost::GetThermalEnergy().

◆ GetDTABlackTrackEnergy()

G4double G4Nucleus::GetDTABlackTrackEnergy ( ) const
inline

◆ GetEnergyDeposit()

G4double G4Nucleus::GetEnergyDeposit ( )
inline

Definition at line 184 of file G4Nucleus.hh.

184{return excitationEnergy; }

Referenced by G4WilsonAbrasionModel::ApplyYourself(), G4KaonMinusAbsorptionAtRest::AtRestDoIt(), and G4BertiniEvaporation::BreakItUp().

◆ GetFermiMomentum()

G4ThreeVector G4Nucleus::GetFermiMomentum ( )

Definition at line 398 of file G4Nucleus.cc.

399 {
400 // chv: .. we assume zero temperature!
401
402 // momentum is equally distributed in each phasespace volume dpx, dpy, dpz.
403 G4double ranflat1=
404 CLHEP::RandFlat::shoot((G4double)0.,(G4double)fermiMomentum);
405 G4double ranflat2=
406 CLHEP::RandFlat::shoot((G4double)0.,(G4double)fermiMomentum);
407 G4double ranflat3=
408 CLHEP::RandFlat::shoot((G4double)0.,(G4double)fermiMomentum);
409 G4double ranmax = (ranflat1>ranflat2? ranflat1: ranflat2);
410 ranmax = (ranmax>ranflat3? ranmax : ranflat3);
411
412 // Isotropic momentum distribution
413 G4double costheta = 2.*G4UniformRand() - 1.0;
414 G4double sintheta = std::sqrt(1.0 - costheta*costheta);
415 G4double phi = 2.0*pi*G4UniformRand();
416
417 G4double pz=costheta*ranmax;
418 G4double px=sintheta*std::cos(phi)*ranmax;
419 G4double py=sintheta*std::sin(phi)*ranmax;
420 G4ThreeVector p(px,py,pz);
421 return p;
422 }
static double shoot()
Definition: RandFlat.cc:59
const G4double pi

◆ GetIsotope()

const G4Isotope * G4Nucleus::GetIsotope ( )
inline

Definition at line 119 of file G4Nucleus.hh.

120 { return fIsotope; }

Referenced by G4HadronicProcess::GetTargetIsotope().

◆ GetN_asInt()

◆ GetPNBlackTrackEnergy()

G4double G4Nucleus::GetPNBlackTrackEnergy ( ) const
inline

◆ GetThermalNucleus()

G4ReactionProduct G4Nucleus::GetThermalNucleus ( G4double  aMass,
G4double  temp = -1 
) const

Definition at line 130 of file G4Nucleus.cc.

131{
132 G4double currentTemp = temp;
133 if(currentTemp < 0) currentTemp = theTemp;
134 G4ReactionProduct theTarget;
135 theTarget.SetMass(targetMass*G4Neutron::Neutron()->GetPDGMass());
136 G4double px, py, pz;
137 px = GetThermalPz(theTarget.GetMass(), currentTemp);
138 py = GetThermalPz(theTarget.GetMass(), currentTemp);
139 pz = GetThermalPz(theTarget.GetMass(), currentTemp);
140 theTarget.SetMomentum(px, py, pz);
141 G4double tMom = std::sqrt(px*px+py*py+pz*pz);
142 G4double tEtot = std::sqrt((tMom+theTarget.GetMass())*
143 (tMom+theTarget.GetMass())-
144 2.*tMom*theTarget.GetMass());
145 if(1-tEtot/theTarget.GetMass()>0.001)
146 {
147 theTarget.SetTotalEnergy(tEtot);
148 }
149 else
150 {
151 theTarget.SetKineticEnergy(tMom*tMom/(2.*theTarget.GetMass()));
152 }
153 return theTarget;
154}
G4double GetThermalPz(const G4double mass, const G4double temp) const
Definition: G4Nucleus.cc:254
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
void SetKineticEnergy(const G4double en)
void SetMass(const G4double mas)

Referenced by GetBiasedThermalNucleus(), G4NeutronHPCaptureData::GetCrossSection(), G4NeutronHPElasticData::GetCrossSection(), G4NeutronHPFissionData::GetCrossSection(), G4NeutronHPInelasticData::GetCrossSection(), G4NeutronHPorLCaptureData::GetCrossSection(), G4NeutronHPorLEInelasticData::GetCrossSection(), G4NeutronHPorLElasticData::GetCrossSection(), and G4NeutronHPorLFissionData::GetCrossSection().

◆ GetThermalPz()

G4double G4Nucleus::GetThermalPz ( const G4double  mass,
const G4double  temp 
) const

Definition at line 254 of file G4Nucleus.cc.

255 {
256 G4double result = G4RandGauss::shoot();
257 result *= std::sqrt(k_Boltzmann*temp*mass); // Das ist impuls (Pz),
258 // nichtrelativistische rechnung
259 // Maxwell verteilung angenommen
260 return result;
261 }

Referenced by GetThermalNucleus().

◆ GetZ_asInt()

G4int G4Nucleus::GetZ_asInt ( ) const
inline

Definition at line 115 of file G4Nucleus.hh.

116 { return theZ; }

Referenced by G4RPGReaction::AddBlackTrackParticles(), G4LightMedia::AntiLambdaExchange(), G4LightMedia::AntiNeutronExchange(), G4LightMedia::AntiOmegaMinusExchange(), G4LightMedia::AntiProtonExchange(), G4LightMedia::AntiSigmaMinusExchange(), G4LightMedia::AntiSigmaPlusExchange(), G4LightMedia::AntiXiMinusExchange(), G4LightMedia::AntiXiZeroExchange(), G4WilsonAbrasionModel::ApplyYourself(), G4EMDissociation::ApplyYourself(), G4ElectroNuclearReaction::ApplyYourself(), G4ProtonAntiProtonReaction::ApplyYourself(), G4LENDCapture::ApplyYourself(), G4LENDElastic::ApplyYourself(), G4LENDFission::ApplyYourself(), G4LENDInelastic::ApplyYourself(), G4LENDModel::ApplyYourself(), G4ElectroVDNuclearModel::ApplyYourself(), G4NeutronHPorLCaptureModel::ApplyYourself(), G4NeutronHPorLEInelasticModel::ApplyYourself(), G4NeutronHPorLElasticModel::ApplyYourself(), G4NeutronHPorLFissionModel::ApplyYourself(), G4NeutronHPThermalScattering::ApplyYourself(), G4ChiralInvariantPhaseSpace::ApplyYourself(), G4ChargeExchange::ApplyYourself(), G4HadronElastic::ApplyYourself(), G4LEnp::ApplyYourself(), G4LEpp::ApplyYourself(), G4NeutronRadCapture::ApplyYourself(), G4HEAntiKaonZeroInelastic::ApplyYourself(), G4HEAntiLambdaInelastic::ApplyYourself(), G4HEAntiNeutronInelastic::ApplyYourself(), G4HEAntiOmegaMinusInelastic::ApplyYourself(), G4HEAntiProtonInelastic::ApplyYourself(), G4HEAntiSigmaMinusInelastic::ApplyYourself(), G4HEAntiSigmaPlusInelastic::ApplyYourself(), G4HEAntiXiMinusInelastic::ApplyYourself(), G4HEAntiXiZeroInelastic::ApplyYourself(), G4HEKaonMinusInelastic::ApplyYourself(), G4HEKaonPlusInelastic::ApplyYourself(), G4HEKaonZeroInelastic::ApplyYourself(), G4HEKaonZeroLongInelastic::ApplyYourself(), G4HEKaonZeroShortInelastic::ApplyYourself(), G4HELambdaInelastic::ApplyYourself(), G4HENeutronInelastic::ApplyYourself(), G4HEOmegaMinusInelastic::ApplyYourself(), G4HEPionMinusInelastic::ApplyYourself(), G4HEPionPlusInelastic::ApplyYourself(), G4HEProtonInelastic::ApplyYourself(), G4HESigmaMinusInelastic::ApplyYourself(), G4HESigmaPlusInelastic::ApplyYourself(), G4HEXiMinusInelastic::ApplyYourself(), G4HEXiZeroInelastic::ApplyYourself(), G4LCapture::ApplyYourself(), G4LEAlphaInelastic::ApplyYourself(), G4LEDeuteronInelastic::ApplyYourself(), G4LElastic::ApplyYourself(), G4LETritonInelastic::ApplyYourself(), G4LFission::ApplyYourself(), G4QMDReaction::ApplyYourself(), G4EmCaptureCascade::ApplyYourself(), G4MuMinusCapturePrecompound::ApplyYourself(), G4MuonMinusBoundDecay::ApplyYourself(), G4BinaryCascade::ApplyYourself(), G4BinaryLightIonReaction::ApplyYourself(), G4ParaFissionModel::ApplyYourself(), G4INCLXXInterface::ApplyYourself(), G4LowEIonFragmentation::ApplyYourself(), G4PreCompoundModel::ApplyYourself(), G4TheoFSGenerator::ApplyYourself(), G4ProtonAntiProtonAtRestChips::AtRestDoIt(), G4HadronStoppingProcess::AtRestDoIt(), G4KaonMinusAbsorptionAtRest::AtRestDoIt(), G4BertiniEvaporation::BreakItUp(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4HadronicProcess::CheckResult(), G4CascadeInterface::createTarget(), G4InelasticInteraction::ExtractResidualNucleus(), G4ReactionDynamics::GenerateXandPt(), G4ProjectileDiffractiveChannel::GetFraction(), G4QuasiElasticChannel::GetFraction(), G4FTFModel::Init(), G4LightMedia::KaonPlusExchange(), G4LightMedia::KaonZeroShortExchange(), G4LightMedia::LambdaExchange(), G4LightMedia::NeutronExchange(), G4LightMedia::OmegaMinusExchange(), G4LightMedia::PionPlusExchange(), G4HadronicProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4WHadronElasticProcess::PostStepDoIt(), G4LightMedia::ProtonExchange(), G4RPGFragmentation::ReactionStage(), G4RPGPionSuppression::ReactionStage(), G4RPGTwoBody::ReactionStage(), G4RPGTwoCluster::ReactionStage(), G4ProjectileDiffractiveChannel::Scatter(), G4QuasiElasticChannel::Scatter(), G4HadronicWhiteBoard::SetTargetNucleus(), G4LightMedia::SigmaMinusExchange(), G4LightMedia::SigmaPlusExchange(), G4ReactionDynamics::SuppressChargedPions(), G4ReactionDynamics::TwoBody(), G4ReactionDynamics::TwoCluster(), G4LightMedia::XiMinusExchange(), and G4LightMedia::XiZeroExchange().

◆ operator!=()

G4bool G4Nucleus::operator!= ( const G4Nucleus right) const
inline

Definition at line 89 of file G4Nucleus.hh.

90 { return ( this != (G4Nucleus *) &right ); }

◆ operator=()

G4Nucleus & G4Nucleus::operator= ( const G4Nucleus right)
inline

Definition at line 64 of file G4Nucleus.hh.

65 {
66 if (this != &right) {
67 theA=right.theA;
68 theZ=right.theZ;
69 aEff=right.aEff;
70 zEff=right.zEff;
71 fIsotope = right.fIsotope;
72 pnBlackTrackEnergy=right.pnBlackTrackEnergy;
73 dtaBlackTrackEnergy=right.dtaBlackTrackEnergy;
74 pnBlackTrackEnergyfromAnnihilation =
75 right.pnBlackTrackEnergyfromAnnihilation;
76 dtaBlackTrackEnergyfromAnnihilation =
77 right.dtaBlackTrackEnergyfromAnnihilation;
78 theTemp = right.theTemp;
79 excitationEnergy = right.excitationEnergy;
80 momentum = right.momentum;
81 fermiMomentum = right.fermiMomentum;
82 }
83 return *this;
84 }

◆ operator==()

G4bool G4Nucleus::operator== ( const G4Nucleus right) const
inline

Definition at line 86 of file G4Nucleus.hh.

87 { return ( this == (G4Nucleus *) &right ); }

◆ ReturnTargetParticle()

G4DynamicParticle * G4Nucleus::ReturnTargetParticle ( ) const

Definition at line 227 of file G4Nucleus.cc.

228 {
229 // choose a proton or a neutron as the target particle
230
231 G4DynamicParticle *targetParticle = new G4DynamicParticle;
232 if( G4UniformRand() < zEff/aEff )
233 targetParticle->SetDefinition( G4Proton::Proton() );
234 else
235 targetParticle->SetDefinition( G4Neutron::Neutron() );
236 return targetParticle;
237 }
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
static G4Proton * Proton()
Definition: G4Proton.cc:93

Referenced by G4LightMedia::AntiLambdaExchange(), G4LightMedia::AntiNeutronExchange(), G4LightMedia::AntiOmegaMinusExchange(), G4LightMedia::AntiProtonExchange(), G4LightMedia::AntiSigmaMinusExchange(), G4LightMedia::AntiSigmaPlusExchange(), G4LightMedia::AntiXiMinusExchange(), G4LightMedia::AntiXiZeroExchange(), G4LEnp::ApplyYourself(), G4LEpp::ApplyYourself(), G4LEAntiKaonZeroInelastic::ApplyYourself(), G4LEAntiLambdaInelastic::ApplyYourself(), G4LEAntiNeutronInelastic::ApplyYourself(), G4LEAntiOmegaMinusInelastic::ApplyYourself(), G4LEAntiProtonInelastic::ApplyYourself(), G4LEAntiSigmaMinusInelastic::ApplyYourself(), G4LEAntiSigmaPlusInelastic::ApplyYourself(), G4LEAntiXiMinusInelastic::ApplyYourself(), G4LEAntiXiZeroInelastic::ApplyYourself(), G4LEKaonMinusInelastic::ApplyYourself(), G4LEKaonPlusInelastic::ApplyYourself(), G4LEKaonZeroInelastic::ApplyYourself(), G4LELambdaInelastic::ApplyYourself(), G4LENeutronInelastic::ApplyYourself(), G4LEOmegaMinusInelastic::ApplyYourself(), G4LEPionMinusInelastic::ApplyYourself(), G4LEPionPlusInelastic::ApplyYourself(), G4LEProtonInelastic::ApplyYourself(), G4LESigmaMinusInelastic::ApplyYourself(), G4LESigmaPlusInelastic::ApplyYourself(), G4LEXiMinusInelastic::ApplyYourself(), G4LEXiZeroInelastic::ApplyYourself(), G4RPGAntiKZeroInelastic::ApplyYourself(), G4RPGAntiLambdaInelastic::ApplyYourself(), G4RPGAntiNeutronInelastic::ApplyYourself(), G4RPGAntiOmegaMinusInelastic::ApplyYourself(), G4RPGAntiProtonInelastic::ApplyYourself(), G4RPGAntiSigmaMinusInelastic::ApplyYourself(), G4RPGAntiSigmaPlusInelastic::ApplyYourself(), G4RPGAntiXiMinusInelastic::ApplyYourself(), G4RPGAntiXiZeroInelastic::ApplyYourself(), G4RPGKMinusInelastic::ApplyYourself(), G4RPGKPlusInelastic::ApplyYourself(), G4RPGKZeroInelastic::ApplyYourself(), G4RPGLambdaInelastic::ApplyYourself(), G4RPGNeutronInelastic::ApplyYourself(), G4RPGOmegaMinusInelastic::ApplyYourself(), G4RPGPiMinusInelastic::ApplyYourself(), G4RPGPiPlusInelastic::ApplyYourself(), G4RPGProtonInelastic::ApplyYourself(), G4RPGSigmaMinusInelastic::ApplyYourself(), G4RPGSigmaPlusInelastic::ApplyYourself(), G4RPGXiMinusInelastic::ApplyYourself(), G4RPGXiZeroInelastic::ApplyYourself(), G4LightMedia::KaonPlusExchange(), G4LightMedia::KaonZeroShortExchange(), G4LightMedia::LambdaExchange(), G4LightMedia::NeutronExchange(), G4LightMedia::OmegaMinusExchange(), G4LightMedia::PionPlusExchange(), G4LightMedia::ProtonExchange(), G4LightMedia::SigmaMinusExchange(), G4LightMedia::SigmaPlusExchange(), G4LightMedia::XiMinusExchange(), and G4LightMedia::XiZeroExchange().

◆ SetIsotope()

void G4Nucleus::SetIsotope ( const G4Isotope iso)
inline

Definition at line 122 of file G4Nucleus.hh.

123 {
124 fIsotope = iso;
125 if(iso) {
126 theZ = iso->GetZ();
127 theA = iso->GetN();
128 aEff = theA;
129 zEff = theZ;
130 }
131 }

Referenced by G4CrossSectionDataStore::SampleZandA().

◆ SetParameters() [1/2]

void G4Nucleus::SetParameters ( const G4double  A,
const G4double  Z 
)

Definition at line 198 of file G4Nucleus.cc.

199{
200 theZ = G4lrint(Z);
201 theA = G4lrint(A);
202 if (theA<1 || theZ<0 || theZ>theA) {
203 throw G4HadronicException(__FILE__, __LINE__,
204 "G4Nucleus::SetParameters called with non-physical parameters");
205 }
206 aEff = A; // atomic weight
207 zEff = Z; // atomic number
208 fIsotope = 0;
209}
int G4lrint(double ad)
Definition: templates.hh:163

Referenced by G4NeutronHPCapture::ApplyYourself(), G4NeutronHPElastic::ApplyYourself(), G4NeutronHPFission::ApplyYourself(), G4NeutronHPInelastic::ApplyYourself(), G4Nucleus(), and G4ElementSelector::SelectZandA().

◆ SetParameters() [2/2]

void G4Nucleus::SetParameters ( const G4int  A,
const G4int  Z 
)

Definition at line 212 of file G4Nucleus.cc.

213{
214 theZ = Z;
215 theA = A;
216 if( theA<1 || theZ<0 || theZ>theA )
217 {
218 throw G4HadronicException(__FILE__, __LINE__,
219 "G4Nucleus::SetParameters called with non-physical parameters");
220 }
221 aEff = A; // atomic weight
222 zEff = Z; // atomic number
223 fIsotope = 0;
224}

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