Geant4 9.6.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 (const G4int theA, const G4int theZ, const G4String &aName, G4GEMProbability *aEmissionStrategy, G4VCoulombBarrier *aCoulombBarrier)
 
virtual ~G4GEMChannel ()
 
virtual G4double GetEmissionProbability (G4Fragment *theNucleus)
 
virtual G4FragmentVectorBreakUp (const G4Fragment &theNucleus)
 
void SetLevelDensityParameter (G4VLevelDensityParameter *aLevelDensity)
 
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)
 

Additional Inherited Members

- Protected Attributes inherited from G4VEvaporationChannel
G4int OPTxs
 
G4bool useSICB
 

Detailed Description

Definition at line 49 of file G4GEMChannel.hh.

Constructor & Destructor Documentation

◆ G4GEMChannel()

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

Definition at line 44 of file G4GEMChannel.cc.

46 :
48 A(theA),
49 Z(theZ),
50 theEvaporationProbabilityPtr(aEmissionStrategy),
51 theCoulombBarrierPtr(aCoulombBarrier),
52 EmissionProbability(0.0),
53 MaximalKineticEnergy(-CLHEP::GeV)
54{
55 theLevelDensityPtr = new G4EvaporationLevelDensityParameter;
56 MyOwnLevelDensity = true;
57 EvaporatedMass = G4NucleiProperties::GetNuclearMass(A, Z);
58 ResidualMass = CoulombBarrier = 0.0;
59 fG4pow = G4Pow::GetInstance();
60 ResidualZ = ResidualA = 0;
61}
static G4double GetNuclearMass(const G4double A, const G4double Z)
static G4Pow * GetInstance()
Definition: G4Pow.cc:50

◆ ~G4GEMChannel()

G4GEMChannel::~G4GEMChannel ( )
virtual

Definition at line 63 of file G4GEMChannel.cc.

64{
65 if (MyOwnLevelDensity) { delete theLevelDensityPtr; }
66}

Member Function Documentation

◆ BreakUp()

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

Implements G4VEvaporationChannel.

Definition at line 135 of file G4GEMChannel.cc.

136{
137 G4double EvaporatedKineticEnergy = CalcKineticEnergy(theNucleus);
138 G4double EvaporatedEnergy = EvaporatedKineticEnergy + EvaporatedMass;
139
140 G4ThreeVector momentum(IsotropicVector(std::sqrt(EvaporatedKineticEnergy*
141 (EvaporatedKineticEnergy+2.0*EvaporatedMass))));
142
143 momentum.rotateUz(theNucleus.GetMomentum().vect().unit());
144
145 G4LorentzVector EvaporatedMomentum(momentum,EvaporatedEnergy);
146 EvaporatedMomentum.boost(theNucleus.GetMomentum().boostVector());
147 G4Fragment * EvaporatedFragment = new G4Fragment(A,Z,EvaporatedMomentum);
148 // ** And now the residual nucleus **
149 G4double theExEnergy = theNucleus.GetExcitationEnergy();
150 G4double theMass = theNucleus.GetGroundStateMass();
151 G4double ResidualEnergy =
152 theMass + (theExEnergy - EvaporatedKineticEnergy) - EvaporatedMass;
153
154 G4LorentzVector ResidualMomentum(-momentum,ResidualEnergy);
155 ResidualMomentum.boost(theNucleus.GetMomentum().boostVector());
156
157 G4Fragment * ResidualFragment = new G4Fragment( ResidualA, ResidualZ, ResidualMomentum );
158
159 G4FragmentVector * theResult = new G4FragmentVector;
160
161 theResult->push_back(EvaporatedFragment);
162 theResult->push_back(ResidualFragment);
163 return theResult;
164}
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
double G4double
Definition: G4Types.hh:64
Hep3Vector unit() const
Hep3Vector boostVector() const
Hep3Vector vect() const
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:240
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:235
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:251

◆ GetEmissionProbability()

G4double G4GEMChannel::GetEmissionProbability ( G4Fragment theNucleus)
virtual

Implements G4VEvaporationChannel.

Definition at line 68 of file G4GEMChannel.cc.

69{
70 G4int anA = fragment->GetA_asInt();
71 G4int aZ = fragment->GetZ_asInt();
72 ResidualA = anA - A;
73 ResidualZ = aZ - Z;
74 //G4cout << "G4GEMChannel::Initialize: Z= " << aZ << " A= " << anA
75 // << " Zres= " << ResidualZ << " Ares= " << ResidualA << G4endl;
76
77 // We only take into account channels which are physically allowed
78 if (ResidualA <= 0 || ResidualZ <= 0 || ResidualA < ResidualZ ||
79 (ResidualA == ResidualZ && ResidualA > 1))
80 {
81 CoulombBarrier = 0.0;
82 MaximalKineticEnergy = -CLHEP::GeV;
83 EmissionProbability = 0.0;
84 }
85 else
86 {
87 // Effective excitation energy
88 // JMQ 071009: pairing in ExEnergy should be the one of parent compound nucleus
89 // FIXED the bug causing reported crash by VI (negative Probabilities
90 // due to inconsistency in Coulomb barrier calculation (CoulombBarrier and -Beta
91 // param for protons must be the same)
92 // G4double ExEnergy = fragment.GetExcitationEnergy() -
93 // G4PairingCorrection::GetInstance()->GetPairingCorrection(ResidualA,ResidualZ);
94 G4double ExEnergy = fragment->GetExcitationEnergy() -
96
97 //G4cout << "Eexc(MeV)= " << ExEnergy/MeV << G4endl;
98
99 if( ExEnergy <= 0.0) {
100 CoulombBarrier = 0.0;
101 MaximalKineticEnergy = -1000.0*MeV;
102 EmissionProbability = 0.0;
103
104 } else {
105
106 ResidualMass = G4NucleiProperties::GetNuclearMass(ResidualA, ResidualZ);
107
108 // Coulomb Barrier calculation
109 CoulombBarrier = theCoulombBarrierPtr->GetCoulombBarrier(ResidualA,ResidualZ,ExEnergy);
110 //G4cout << "CBarrier(MeV)= " << CoulombBarrier/MeV << G4endl;
111
112 //Maximal kinetic energy (JMQ : at the Coulomb barrier)
113 MaximalKineticEnergy =
114 CalcMaximalKineticEnergy(fragment->GetGroundStateMass()+ExEnergy);
115 //G4cout << "MaxE(MeV)= " << MaximalKineticEnergy/MeV << G4endl;
116
117 // Emission probability
118 if (MaximalKineticEnergy <= 0.0)
119 {
120 EmissionProbability = 0.0;
121 }
122 else
123 {
124 // Total emission probability for this channel
125 EmissionProbability =
126 theEvaporationProbabilityPtr->EmissionProbability(*fragment,
127 MaximalKineticEnergy);
128 }
129 }
130 }
131 //G4cout << "Prob= " << EmissionProbability << G4endl;
132 return EmissionProbability;
133}
int G4int
Definition: G4Types.hh:66
G4double EmissionProbability(const G4Fragment &fragment, G4double anEnergy)
static G4PairingCorrection * GetInstance()
G4double GetPairingCorrection(G4int A, G4int Z) const
virtual G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U) const =0

◆ GetMaximalKineticEnergy()

G4double G4GEMChannel::GetMaximalKineticEnergy ( void  ) const
inline

Definition at line 71 of file G4GEMChannel.hh.

72 { return MaximalKineticEnergy; }

◆ SetLevelDensityParameter()

void G4GEMChannel::SetLevelDensityParameter ( G4VLevelDensityParameter aLevelDensity)
inline

Definition at line 64 of file G4GEMChannel.hh.

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

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