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

#include <G4EvaporationChannel.hh>

+ Inheritance diagram for G4EvaporationChannel:

Public Member Functions

 G4EvaporationChannel (G4int theA, G4int theZ, const G4String &aName, G4EvaporationProbability *aEmissionStrategy, G4VCoulombBarrier *aCoulombBarrier)
 
virtual ~G4EvaporationChannel ()
 
void SetEmissionStrategy (G4EvaporationProbability *aEmissionStrategy)
 
void SetCoulombBarrierStrategy (G4VCoulombBarrier *aCoulombBarrier)
 
virtual G4double GetEmissionProbability (G4Fragment *fragment)
 
G4FragmentVectorBreakUp (const G4Fragment &theNucleus)
 
G4double GetMaximalKineticEnergy (void) const
 
- Public Member Functions inherited from G4VEvaporationChannel
 G4VEvaporationChannel (const G4String &aName="Anonymous")
 
virtual ~G4VEvaporationChannel ()
 
virtual G4FragmentEmittedFragment (G4Fragment *theNucleus)
 
virtual G4FragmentVectorBreakUpFragment (G4Fragment *theNucleus)
 
virtual G4FragmentVectorBreakUp (const G4Fragment &theNucleus)=0
 
virtual G4double GetEmissionProbability (G4Fragment *theNucleus)=0
 
G4String GetName () const
 
void SetName (const G4String &aName)
 
void SetOPTxs (G4int opt)
 
void UseSICB (G4bool use)
 

Protected Member Functions

 G4EvaporationChannel ()
 

Additional Inherited Members

- Protected Attributes inherited from G4VEvaporationChannel
G4int OPTxs
 
G4bool useSICB
 

Detailed Description

Definition at line 47 of file G4EvaporationChannel.hh.

Constructor & Destructor Documentation

◆ G4EvaporationChannel() [1/2]

G4EvaporationChannel::G4EvaporationChannel ( G4int  theA,
G4int  theZ,
const G4String aName,
G4EvaporationProbability aEmissionStrategy,
G4VCoulombBarrier aCoulombBarrier 
)

Definition at line 51 of file G4EvaporationChannel.cc.

54 :
56 theA(anA),
57 theZ(aZ),
58 theEvaporationProbabilityPtr(aEmissionStrategy),
59 theCoulombBarrierPtr(aCoulombBarrier),
60 EmissionProbability(0.0),
61 MaximalKineticEnergy(-1000.0)
62{
63 ResidualA = 0;
64 ResidualZ = 0;
65 ResidualMass = CoulombBarrier=0.0;
66 EvaporatedMass = G4NucleiProperties::GetNuclearMass(theA, theZ);
67 theLevelDensityPtr = new G4EvaporationLevelDensityParameter;
68}
static G4double GetNuclearMass(const G4double A, const G4double Z)

◆ ~G4EvaporationChannel()

G4EvaporationChannel::~G4EvaporationChannel ( )
virtual

Definition at line 85 of file G4EvaporationChannel.cc.

86{
87 delete theLevelDensityPtr;
88}

◆ G4EvaporationChannel() [2/2]

G4EvaporationChannel::G4EvaporationChannel ( )
protected

Definition at line 70 of file G4EvaporationChannel.cc.

70 :
72 theA(0),
73 theZ(0),
74 theEvaporationProbabilityPtr(0),
75 theCoulombBarrierPtr(0),
76 EmissionProbability(0.0),
77 MaximalKineticEnergy(-1000.0)
78{
79 ResidualA = 0;
80 ResidualZ = 0;
81 EvaporatedMass = ResidualMass = CoulombBarrier = 0.0;
82 theLevelDensityPtr = new G4EvaporationLevelDensityParameter;
83}

Member Function Documentation

◆ BreakUp()

G4FragmentVector * G4EvaporationChannel::BreakUp ( const G4Fragment theNucleus)
virtual

Implements G4VEvaporationChannel.

Definition at line 151 of file G4EvaporationChannel.cc.

152{
153 /*
154 G4double Ecm = GetKineticEnergy(theNucleus) + ResidualMass + EvaporatedMass;
155
156 G4double EvaporatedEnergy =
157 ((Ecm-ResidualMass)*(Ecm+ResidualMass) + EvaporatedMass*EvaporatedMass)/(2*Ecm);
158 */
159 G4double EvaporatedEnergy = GetKineticEnergy(theNucleus) + EvaporatedMass;
160
161 G4ThreeVector momentum(IsotropicVector
162 (std::sqrt((EvaporatedEnergy - EvaporatedMass)*
163 (EvaporatedEnergy + EvaporatedMass))));
164
165 G4LorentzVector EvaporatedMomentum(momentum,EvaporatedEnergy);
166 G4LorentzVector ResidualMomentum = theNucleus.GetMomentum();
167 EvaporatedMomentum.boost(ResidualMomentum.boostVector());
168
169 G4Fragment * EvaporatedFragment = new G4Fragment(theA,theZ,EvaporatedMomentum);
170 ResidualMomentum -= EvaporatedMomentum;
171
172 G4Fragment * ResidualFragment = new G4Fragment(ResidualA, ResidualZ, ResidualMomentum);
173
174 G4FragmentVector * theResult = new G4FragmentVector;
175
176#ifdef debug
177 G4double Efinal = ResidualMomentum.e() + EvaporatedMomentum.e();
178 G4ThreeVector Pfinal = ResidualMomentum.vect() + EvaporatedMomentum.vect();
179 if (std::abs(Efinal-theNucleus.GetMomentum().e()) > 1.0*keV) {
180 G4cout << "@@@@@@@@@@@@@@@@@@@@@ G4Evaporation Chanel: ENERGY @@@@@@@@@@@@@@@@" << G4endl;
181 G4cout << "Initial : " << theNucleus.GetMomentum().e()/MeV << " MeV Final :"
182 <<Efinal/MeV << " MeV Delta: " << (Efinal-theNucleus.GetMomentum().e())/MeV
183 << " MeV" << G4endl;
184 }
185 if (std::abs(Pfinal.x()-theNucleus.GetMomentum().x()) > 1.0*keV ||
186 std::abs(Pfinal.y()-theNucleus.GetMomentum().y()) > 1.0*keV ||
187 std::abs(Pfinal.z()-theNucleus.GetMomentum().z()) > 1.0*keV ) {
188 G4cout << "@@@@@@@@@@@@@@@@@@@@@ G4Evaporation Chanel: MOMENTUM @@@@@@@@@@@@@@@@" << G4endl;
189 G4cout << "Initial : " << theNucleus.GetMomentum().vect() << " MeV Final :"
190 <<Pfinal/MeV << " MeV Delta: " << Pfinal-theNucleus.GetMomentum().vect()
191 << " MeV" << G4endl;
192
193 }
194#endif
195 theResult->push_back(EvaporatedFragment);
196 theResult->push_back(ResidualFragment);
197 return theResult;
198}
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
double z() const
double x() const
double y() const
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:251

◆ GetEmissionProbability()

G4double G4EvaporationChannel::GetEmissionProbability ( G4Fragment fragment)
virtual

Implements G4VEvaporationChannel.

Definition at line 93 of file G4EvaporationChannel.cc.

94{
95 //for inverse cross section choice
96 theEvaporationProbabilityPtr->SetOPTxs(OPTxs);
97 // for superimposed Coulomb Barrier for inverse cross sections
98 theEvaporationProbabilityPtr->UseSICB(useSICB);
99
100 G4int FragmentA = fragment->GetA_asInt();
101 G4int FragmentZ = fragment->GetZ_asInt();
102 ResidualA = FragmentA - theA;
103 ResidualZ = FragmentZ - theZ;
104 //G4cout << "G4EvaporationChannel::Initialize Z= " << theZ << " A= " << theA
105 // << " FragZ= " << FragmentZ << " FragA= " << FragmentA << G4endl;
106
107 //Effective excitation energy
108 G4double ExEnergy = fragment->GetExcitationEnergy() -
110
111 // Only channels which are physically allowed are taken into account
112 if (ResidualA <= 0 || ResidualZ <= 0 || ResidualA < ResidualZ ||
113 (ResidualA == ResidualZ && ResidualA > 1) || ExEnergy <= 0.0) {
114 CoulombBarrier = ResidualMass = 0.0;
115 MaximalKineticEnergy = -1000.0*MeV;
116 EmissionProbability = 0.0;
117 } else {
118 ResidualMass = G4NucleiProperties::GetNuclearMass(ResidualA, ResidualZ);
119 G4double FragmentMass = fragment->GetGroundStateMass();
120 CoulombBarrier = theCoulombBarrierPtr->GetCoulombBarrier(ResidualA,ResidualZ,ExEnergy);
121 // Maximal Kinetic Energy
122 MaximalKineticEnergy = CalcMaximalKineticEnergy(FragmentMass + ExEnergy);
123 //MaximalKineticEnergy = ExEnergy + fragment.GetGroundStateMass()
124 // - EvaporatedMass - ResidualMass;
125
126 // Emission probability
127 // Protection for the case Tmax<V. If not set in this way we could end up in an
128 // infinite loop in the method GetKineticEnergy if OPTxs!=0 && useSICB=true.
129 // Of course for OPTxs=0 we have the Coulomb barrier
130
131 G4double limit;
132 if (OPTxs==0 || (OPTxs!=0 && useSICB))
133 limit= CoulombBarrier;
134 else limit=0.;
135 limit = std::max(limit, FragmentMass - ResidualMass - EvaporatedMass);
136
137 // The threshold for charged particle emission must be set to 0 if Coulomb
138 //cutoff is included in the cross sections
139 if (MaximalKineticEnergy <= limit) EmissionProbability = 0.0;
140 else {
141 // Total emission probability for this channel
142 EmissionProbability = theEvaporationProbabilityPtr->
143 EmissionProbability(*fragment, MaximalKineticEnergy);
144 }
145 }
146 //G4cout << "G4EvaporationChannel:: probability= " << EmissionProbability << G4endl;
147
148 return EmissionProbability;
149}
int G4int
Definition: G4Types.hh:66
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:240
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:235
G4int GetZ_asInt() const
Definition: G4Fragment.hh:223
G4int GetA_asInt() const
Definition: G4Fragment.hh:218
static G4PairingCorrection * GetInstance()
G4double GetPairingCorrection(G4int A, G4int Z) const
virtual G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U) const =0

◆ GetMaximalKineticEnergy()

G4double G4EvaporationChannel::GetMaximalKineticEnergy ( void  ) const
inline

Definition at line 76 of file G4EvaporationChannel.hh.

77 { return MaximalKineticEnergy; }

◆ SetCoulombBarrierStrategy()

void G4EvaporationChannel::SetCoulombBarrierStrategy ( G4VCoulombBarrier aCoulombBarrier)
inline

Definition at line 61 of file G4EvaporationChannel.hh.

62 {theCoulombBarrierPtr = aCoulombBarrier;}

◆ SetEmissionStrategy()

void G4EvaporationChannel::SetEmissionStrategy ( G4EvaporationProbability aEmissionStrategy)
inline

Definition at line 58 of file G4EvaporationChannel.hh.

59 {theEvaporationProbabilityPtr = aEmissionStrategy;}

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