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

#include <G4GEMChannel.hh>

+ Inheritance diagram for G4GEMChannel:

Public Member Functions

 G4GEMChannel (G4int theA, G4int theZ, const G4String &aName, G4GEMProbability *aEmissionStrategy)
 
virtual ~G4GEMChannel ()
 
virtual G4double GetEmissionProbability (G4Fragment *theNucleus)
 
virtual G4FragmentEmittedFragment (G4Fragment *theNucleus)
 
virtual void Dump () const
 
void SetLevelDensityParameter (G4VLevelDensityParameter *aLevelDensity)
 
- Public Member Functions inherited from G4VEvaporationChannel
 G4VEvaporationChannel (const G4String &aName="")
 
virtual ~G4VEvaporationChannel ()
 
virtual G4double GetEmissionProbability (G4Fragment *theNucleus)=0
 
virtual void Initialise ()
 
virtual G4double GetLifeTime (G4Fragment *theNucleus)
 
virtual G4FragmentEmittedFragment (G4Fragment *theNucleus)
 
virtual G4bool BreakUpChain (G4FragmentVector *theResult, G4Fragment *theNucleus)
 
G4FragmentVectorBreakUpFragment (G4Fragment *theNucleus)
 
virtual void Dump () const
 
virtual void SetICM (G4bool)
 
virtual void RDMForced (G4bool)
 
void SetOPTxs (G4int opt)
 
void UseSICB (G4bool use)
 

Additional Inherited Members

- Protected Attributes inherited from G4VEvaporationChannel
G4int OPTxs
 
G4bool useSICB
 

Detailed Description

Definition at line 47 of file G4GEMChannel.hh.

Constructor & Destructor Documentation

◆ G4GEMChannel()

G4GEMChannel::G4GEMChannel ( G4int  theA,
G4int  theZ,
const G4String aName,
G4GEMProbability aEmissionStrategy 
)
explicit

Definition at line 48 of file G4GEMChannel.cc.

49 :
51 A(theA),
52 Z(theZ),
53 EmissionProbability(0.0),
54 MaximalKineticEnergy(-CLHEP::GeV),
55 theEvaporationProbabilityPtr(aEmissionStrategy)
56{
57 theCoulombBarrierPtr = new G4GEMCoulombBarrier(theA, theZ);
58 theEvaporationProbabilityPtr->SetCoulomBarrier(theCoulombBarrierPtr);
59 theLevelDensityPtr = new G4EvaporationLevelDensityParameter;
60 MyOwnLevelDensity = true;
61 EvaporatedMass = G4NucleiProperties::GetNuclearMass(A, Z);
62 ResidualMass = CoulombBarrier = 0.0;
63 fG4pow = G4Pow::GetInstance();
64 ResidualZ = ResidualA = 0;
66}
void SetCoulomBarrier(const G4VCoulombBarrier *aCoulombBarrierStrategy)
static G4NuclearLevelData * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
static G4Pow * GetInstance()
Definition: G4Pow.cc:41

◆ ~G4GEMChannel()

G4GEMChannel::~G4GEMChannel ( )
virtual

Definition at line 68 of file G4GEMChannel.cc.

69{
70 if (MyOwnLevelDensity) { delete theLevelDensityPtr; }
71 delete theCoulombBarrierPtr;
72}

Member Function Documentation

◆ Dump()

void G4GEMChannel::Dump ( ) const
virtual

Reimplemented from G4VEvaporationChannel.

Definition at line 249 of file G4GEMChannel.cc.

250{
251 theEvaporationProbabilityPtr->Dump();
252}

◆ EmittedFragment()

G4Fragment * G4GEMChannel::EmittedFragment ( G4Fragment theNucleus)
virtual

Reimplemented from G4VEvaporationChannel.

Definition at line 127 of file G4GEMChannel.cc.

128{
129 G4Fragment* evFragment = 0;
130 G4double evEnergy = SampleKineticEnergy(*theNucleus) + EvaporatedMass;
131
132 G4ThreeVector momentum = G4RandomDirection()*
133 std::sqrt((evEnergy - EvaporatedMass)*(evEnergy + EvaporatedMass));
134
135 G4LorentzVector EvaporatedMomentum(momentum, evEnergy);
136 G4LorentzVector ResidualMomentum = theNucleus->GetMomentum();
137 EvaporatedMomentum.boost(ResidualMomentum.boostVector());
138
139 evFragment = new G4Fragment(A, Z, EvaporatedMomentum);
140 ResidualMomentum -= EvaporatedMomentum;
141 theNucleus->SetZandA_asInt(ResidualZ, ResidualA);
142 theNucleus->SetMomentum(ResidualMomentum);
143
144 return evFragment;
145}
G4ThreeVector G4RandomDirection()
double G4double
Definition: G4Types.hh:83
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:299
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:304
void SetZandA_asInt(G4int Znew, G4int Anew)
Definition: G4Fragment.hh:268

◆ GetEmissionProbability()

G4double G4GEMChannel::GetEmissionProbability ( G4Fragment theNucleus)
virtual

Implements G4VEvaporationChannel.

Definition at line 74 of file G4GEMChannel.cc.

75{
76 G4int anA = fragment->GetA_asInt();
77 G4int aZ = fragment->GetZ_asInt();
78 ResidualA = anA - A;
79 ResidualZ = aZ - Z;
80 /*
81 G4cout << "G4GEMChannel: Z= " << Z << " A= " << A
82 << " FragmentZ= " << aZ << " FragmentA= " << anA
83 << " Zres= " << ResidualZ << " Ares= " << ResidualA
84 << G4endl;
85 */
86 // We only take into account channels which are physically allowed
87 EmissionProbability = 0.0;
88
89 // Only channels which are physically allowed are taken into account
90 if (ResidualA >= ResidualZ && ResidualZ >= 0 && ResidualA >= A) {
91
92 //Effective excitation energy
93 G4double ExEnergy = fragment->GetExcitationEnergy()
94 - fNucData->GetPairingCorrection(aZ, anA);
95 if(ExEnergy > 0.0) {
96 ResidualMass = G4NucleiProperties::GetNuclearMass(ResidualA, ResidualZ);
97 G4double FragmentMass = fragment->GetGroundStateMass();
98 G4double Etot = FragmentMass + ExEnergy;
99 // Coulomb Barrier calculation
100 CoulombBarrier =
101 theCoulombBarrierPtr->GetCoulombBarrier(ResidualA,ResidualZ,ExEnergy);
102 /*
103 G4cout << "Eexc(MeV)= " << ExEnergy/MeV
104 << " CoulBarrier(MeV)= " << CoulombBarrier/MeV << G4endl;
105 */
106 if(Etot > ResidualMass + EvaporatedMass + CoulombBarrier) {
107
108 // Maximal Kinetic Energy
109 MaximalKineticEnergy = ((Etot-ResidualMass)*(Etot+ResidualMass)
110 + EvaporatedMass*EvaporatedMass)/(2.0*Etot)
111 - EvaporatedMass - CoulombBarrier;
112
113 //G4cout << "CBarrier(MeV)= " << CoulombBarrier/MeV << G4endl;
114
115 if (MaximalKineticEnergy > 0.0) {
116 // Total emission probability for this channel
117 EmissionProbability = theEvaporationProbabilityPtr->
118 EmissionProbability(*fragment, MaximalKineticEnergy);
119 }
120 }
121 }
122 }
123 //G4cout << "Prob= " << EmissionProbability << G4endl;
124 return EmissionProbability;
125}
int G4int
Definition: G4Types.hh:85
G4PairingCorrection * GetPairingCorrection()
virtual G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U) const =0

◆ SetLevelDensityParameter()

void G4GEMChannel::SetLevelDensityParameter ( G4VLevelDensityParameter aLevelDensity)
inline

Definition at line 63 of file G4GEMChannel.hh.

64 {
65 if (MyOwnLevelDensity) { delete theLevelDensityPtr; }
66 theLevelDensityPtr = aLevelDensity;
67 MyOwnLevelDensity = false;
68 }

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