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

#include <G4InuclEvaporation.hh>

+ Inheritance diagram for G4InuclEvaporation:

Public Member Functions

 G4InuclEvaporation ()
 
 ~G4InuclEvaporation ()
 
G4FragmentVectorBreakItUp (const G4Fragment &theNucleus)
 
void setVerboseLevel (const G4int verbose)
 
- Public Member Functions inherited from G4VEvaporation
 G4VEvaporation ()
 
virtual ~G4VEvaporation ()
 
virtual void BreakFragment (G4FragmentVector *, G4Fragment *theNucleus)
 
virtual void InitialiseChannels ()
 
virtual void SetPhotonEvaporation (G4VEvaporationChannel *ptr)
 
void SetFermiBreakUp (G4VFermiBreakUp *ptr)
 
G4VFermiBreakUpGetFermiBreakUp () const
 
G4VEvaporationChannelGetPhotonEvaporation ()
 
G4VEvaporationChannelGetFissionChannel ()
 
G4VEvaporationChannelGetChannel (size_t idx)
 
void SetOPTxs (G4int opt)
 
void UseSICB (G4bool use)
 
size_t GetNumberOfChannels () const
 
 G4VEvaporation (const G4VEvaporation &right)=delete
 
const G4VEvaporationoperator= (const G4VEvaporation &right)=delete
 
G4bool operator== (const G4VEvaporation &right) const =delete
 
G4bool operator!= (const G4VEvaporation &right) const =delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEvaporation
void CleanChannels ()
 
- Protected Attributes inherited from G4VEvaporation
G4VEvaporationChannelthePhotonEvaporation
 
G4VFermiBreakUptheFBU
 
G4int OPTxs
 
G4bool useSICB
 
std::vector< G4VEvaporationChannel * > * theChannels
 
G4VEvaporationFactorytheChannelFactory
 

Detailed Description

Definition at line 43 of file G4InuclEvaporation.hh.

Constructor & Destructor Documentation

◆ G4InuclEvaporation()

G4InuclEvaporation::G4InuclEvaporation ( )

Definition at line 72 of file G4InuclEvaporation.cc.

73 : verboseLevel(0), evaporator(new G4EvaporationInuclCollider) {}

◆ ~G4InuclEvaporation()

G4InuclEvaporation::~G4InuclEvaporation ( )

Definition at line 79 of file G4InuclEvaporation.cc.

79 {
80 delete evaporator;
81}

Member Function Documentation

◆ BreakItUp()

G4FragmentVector * G4InuclEvaporation::BreakItUp ( const G4Fragment theNucleus)

Definition at line 100 of file G4InuclEvaporation.cc.

100 {
101 G4FragmentVector* theResult = new G4FragmentVector;
102
103 if (theNucleus.GetExcitationEnergy() <= 0.0) { // Check that Excitation Energy > 0
104 theResult->push_back(new G4Fragment(theNucleus));
105 return theResult;
106 }
107
108 G4int A = theNucleus.GetA_asInt();
109 G4int Z = theNucleus.GetZ_asInt();
110 G4double mTar = G4NucleiProperties::GetNuclearMass(A, Z); // Mass of the target nucleus
111
112 G4ThreeVector momentum = theNucleus.GetMomentum().vect();
113 G4double exitationE = theNucleus.GetExcitationEnergy();
114
115 G4double mass = mTar;
116 G4ThreeVector boostToLab( momentum/mass );
117
118 if ( verboseLevel > 2 )
119 G4cout << " G4InuclEvaporation : initial kinematics : boostToLab vector = " << boostToLab << G4endl
120 << " excitation energy : " << exitationE << G4endl;
121
122 if (verboseLevel > 2) {
123 G4cout << "G4InuclEvaporation::BreakItUp >>> A: " << A << " Z: " << Z
124 << " exitation E: " << exitationE << " mass: " << mTar/GeV << " GeV"
125 << G4endl;
126 };
127
128 G4InuclNuclei* nucleus = new G4InuclNuclei(A, Z);
129 nucleus->setExitationEnergy(exitationE);
130
131 G4CollisionOutput output;
132 evaporator->collide(0, nucleus, output);
133
134 const std::vector<G4InuclNuclei>& outgoingNuclei = output.getOutgoingNuclei();
135 const std::vector<G4InuclElementaryParticle>& particles = output.getOutgoingParticles();
136
137 G4int i=1;
138
139 if (!particles.empty()) {
140 G4int outgoingType;
141 particleIterator ipart = particles.cbegin();
142 for (; ipart != particles.cend(); ++ipart) {
143 outgoingType = ipart->type();
144
145 if (verboseLevel > 2) {
146 G4cout << "Evaporated particle: " << i << " of type: "
147 << outgoingType << G4endl;
148 ++i;
149 }
150
151 G4LorentzVector vlab = ipart->getMomentum().boost(boostToLab);
152
153 // TEMPORARY: Remove constness on PD until G4Fragment is fixed
154 theResult->push_back( new G4Fragment(vlab, ipart->getDefinition()) );
155 }
156 }
157
158 if (!outgoingNuclei.empty()) {
159 nucleiIterator ifrag = outgoingNuclei.cbegin();
160 for (i=1; ifrag != outgoingNuclei.cend(); ++ifrag) {
161 if (verboseLevel > 2) {
162 G4cout << " Nuclei fragment: " << i << G4endl; ++i;
163 }
164
165 G4LorentzVector vlab = ifrag->getMomentum().boost(boostToLab);
166
167 G4int fragA = ifrag->getA();
168 G4int fragZ = ifrag->getZ();
169 if (verboseLevel > 2) {
170 G4cout << "boosted v" << vlab << G4endl;
171 }
172 theResult->push_back( new G4Fragment(fragA, fragZ, vlab) );
173 }
174 }
175
176 return theResult;
177}
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:64
std::vector< G4InuclNuclei >::const_iterator nucleiIterator
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
const std::vector< G4InuclNuclei > & getOutgoingNuclei() const
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:312
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:322
G4int GetZ_asInt() const
Definition: G4Fragment.hh:289
G4int GetA_asInt() const
Definition: G4Fragment.hh:284
void setExitationEnergy(G4double e)
static G4double GetNuclearMass(const G4double A, const G4double Z)
virtual void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &globalOutput)

◆ setVerboseLevel()

void G4InuclEvaporation::setVerboseLevel ( const G4int  verbose)

Definition at line 96 of file G4InuclEvaporation.cc.

96 {
97 verboseLevel = verbose;
98}

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