Geant4 11.1.1
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 ()=default
 
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 G4double ComputeInverseXSection (G4Fragment *theNucleus, G4double kinEnergy)
 
virtual G4double ComputeProbability (G4Fragment *theNucleus, G4double kinEnergy)
 
virtual void Dump () const
 
virtual void SetICM (G4bool)
 
virtual void RDMForced (G4bool)
 
void SetOPTxs (G4int opt)
 
void UseSICB (G4bool use)
 
 G4VEvaporationChannel (const G4VEvaporationChannel &right)=delete
 
const G4VEvaporationChanneloperator= (const G4VEvaporationChannel &right)=delete
 
G4bool operator== (const G4VEvaporationChannel &right) const =delete
 
G4bool operator!= (const G4VEvaporationChannel &right) const =delete
 

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 49 of file G4GEMChannel.cc.

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

◆ ~G4GEMChannel()

G4GEMChannel::~G4GEMChannel ( )
virtual

Definition at line 71 of file G4GEMChannel.cc.

72{
73 if (MyOwnLevelDensity) { delete theLevelDensityPtr; }
74 delete theCoulombBarrierPtr;
75}

Member Function Documentation

◆ Dump()

void G4GEMChannel::Dump ( ) const
virtual

Reimplemented from G4VEvaporationChannel.

Definition at line 254 of file G4GEMChannel.cc.

255{
256 theEvaporationProbabilityPtr->Dump();
257}

◆ EmittedFragment()

G4Fragment * G4GEMChannel::EmittedFragment ( G4Fragment theNucleus)
virtual

Reimplemented from G4VEvaporationChannel.

Definition at line 130 of file G4GEMChannel.cc.

131{
132 G4Fragment* evFragment = 0;
133 G4double evEnergy = SampleKineticEnergy(*theNucleus) + EvaporatedMass;
134
135 G4ThreeVector momentum = G4RandomDirection()*
136 std::sqrt((evEnergy - EvaporatedMass)*(evEnergy + EvaporatedMass));
137
138 G4LorentzVector EvaporatedMomentum(momentum, evEnergy);
139 G4LorentzVector ResidualMomentum = theNucleus->GetMomentum();
140 EvaporatedMomentum.boost(ResidualMomentum.boostVector());
141
142 evFragment = new G4Fragment(A, Z, EvaporatedMomentum);
143 if ( evFragment != nullptr ) { evFragment->SetCreatorModelID(secID); }
144 ResidualMomentum -= EvaporatedMomentum;
145 theNucleus->SetZandA_asInt(ResidualZ, ResidualA);
146 theNucleus->SetMomentum(ResidualMomentum);
147 theNucleus->SetCreatorModelID(secID);
148
149 return evFragment;
150}
G4ThreeVector G4RandomDirection()
double G4double
Definition: G4Types.hh:83
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
void SetZandA_asInt(G4int Znew, G4int Anew, G4int Lnew=0)
Definition: G4Fragment.hh:295
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:322
void SetCreatorModelID(G4int value)
Definition: G4Fragment.hh:433
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:327

◆ GetEmissionProbability()

G4double G4GEMChannel::GetEmissionProbability ( G4Fragment theNucleus)
virtual

Implements G4VEvaporationChannel.

Definition at line 77 of file G4GEMChannel.cc.

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