Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4VPreCompoundFragment.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id$
27//
28// J. M. Quesada (August 2008). Based on previous work by V. Lara
29//
30// Modified:
31// 20.08.2010 V.Ivanchenko added G4Pow and G4PreCompoundParameters pointers
32// use int Z and A and cleanup
33
35#include "G4SystemOfUnits.hh"
37#include "G4NucleiProperties.hh"
38
40 const G4ParticleDefinition* part, G4VCoulombBarrier* aCoulombBarrier)
41 : particle(part), theCoulombBarrierPtr(aCoulombBarrier),
42 theRestNucleusA(0),theRestNucleusZ(0),theBindingEnergy(0.0),
43 theMaximalKineticEnergy(-MeV),theRestNucleusMass(0.0),
44 theReducedMass(0.0),theMomentum(0.,0.,0.,0.),
45 theEmissionProbability(0.0),theCoulombBarrier(0.0),
46 OPTxs(3),useSICB(false)
47{
48 theA = particle->GetBaryonNumber();
49 theZ = G4int(particle->GetPDGCharge()/eplus + 0.1);
50 theMass = particle->GetPDGMass();
53 theRestNucleusA13 = 0;
54}
55
57{}
58
59std::ostream&
60operator << (std::ostream &out, const G4VPreCompoundFragment &theFragment)
61{
62 out << &theFragment;
63 return out;
64}
65
66std::ostream&
67operator << (std::ostream &out, const G4VPreCompoundFragment *theFragment)
68{
69 out
70 << "PreCompoundModel Emitted Fragment: Z= " << theFragment->GetZ()
71 << " A= " << theFragment->GetA()
72 << " Mass(GeV)= " << theFragment->GetNuclearMass()/CLHEP::GeV;
73 return out;
74}
75
76void
78{
79 theRestNucleusA = aFragment.GetA_asInt() - theA;
80 theRestNucleusZ = aFragment.GetZ_asInt() - theZ;
81
82 if ((theRestNucleusA < theRestNucleusZ) ||
83 (theRestNucleusA < theA) ||
84 (theRestNucleusZ < theZ))
85 {
86 // In order to be sure that emission probability will be 0.
87 theMaximalKineticEnergy = 0.0;
88 return;
89 }
90
91 theRestNucleusA13 = g4pow->Z13(theRestNucleusA);
92
93 // Calculate Coulomb barrier
94 theCoulombBarrier = theCoulombBarrierPtr->
95 GetCoulombBarrier(theRestNucleusA,theRestNucleusZ,
96 aFragment.GetExcitationEnergy());
97
98 // Calculate masses
99 theRestNucleusMass =
100 G4NucleiProperties::GetNuclearMass(theRestNucleusA, theRestNucleusZ);
101 theReducedMass = theRestNucleusMass*theMass/(theRestNucleusMass + theMass);
102
103 // Compute Binding Energies for fragments
104 // needed to separate a fragment from the nucleus
105 theBindingEnergy =
106 theRestNucleusMass + theMass - aFragment.GetGroundStateMass();
107
108 // Compute Maximal Kinetic Energy which can be carried by fragments
109 // after separation - the true assimptotic value
110 G4double Ecm = aFragment.GetMomentum().m();
111 theMaximalKineticEnergy =
112 ((Ecm-theRestNucleusMass)*(Ecm+theRestNucleusMass) + theMass*theMass)/(2.0*Ecm)-theMass;
113}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
std::ostream & operator<<(std::ostream &out, const G4VPreCompoundFragment &theFragment)
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:240
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:235
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:251
G4int GetZ_asInt() const
Definition: G4Fragment.hh:223
G4int GetA_asInt() const
Definition: G4Fragment.hh:218
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetPDGCharge() const
static G4Pow * GetInstance()
Definition: G4Pow.cc:50
G4double Z13(G4int Z)
Definition: G4Pow.hh:110
static G4PreCompoundParameters * GetAddress()
G4VPreCompoundFragment(const G4ParticleDefinition *, G4VCoulombBarrier *aCoulombBarrier)
G4double GetNuclearMass() const
G4double GetCoulombBarrier() const
void Initialize(const G4Fragment &aFragment)
G4PreCompoundParameters * theParameters