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

#include <G4GEMChannelVI.hh>

+ Inheritance diagram for G4GEMChannelVI:

Public Member Functions

 G4GEMChannelVI (G4int theA, G4int theZ)
 
 ~G4GEMChannelVI () final
 
G4double GetEmissionProbability (G4Fragment *theNucleus) final
 
G4FragmentEmittedFragment (G4Fragment *theNucleus) final
 
void Dump () const final
 
- 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 40 of file G4GEMChannelVI.hh.

Constructor & Destructor Documentation

◆ G4GEMChannelVI()

G4GEMChannelVI::G4GEMChannelVI ( G4int  theA,
G4int  theZ 
)
explicit

Definition at line 42 of file G4GEMChannelVI.cc.

43 : A(theA), Z(theZ), secID(-1)
44{
46 pairingCorrection = nData->GetPairingCorrection();
47 const G4LevelManager* lManager = nullptr;
48 if(A > 4) { lManager = nData->GetLevelManager(Z, A); }
50 evapMass2 = evapMass*evapMass;
51
52 cBarrier = new G4CoulombBarrier(A, Z);
53 fProbability = new G4GEMProbabilityVI(A, Z, lManager);
54
55 resA = resZ = fragZ = fragA = 0;
56 mass = resMass = 0.0;
57 secID = G4PhysicsModelCatalog::GetModelID("model_G4GEMChannelVI");
58}
const G4LevelManager * GetLevelManager(G4int Z, G4int A)
G4PairingCorrection * GetPairingCorrection()
static G4NuclearLevelData * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
static G4int GetModelID(const G4int modelIndex)

◆ ~G4GEMChannelVI()

G4GEMChannelVI::~G4GEMChannelVI ( )
final

Definition at line 60 of file G4GEMChannelVI.cc.

61{
62 delete cBarrier;
63 delete fProbability;
64}

Member Function Documentation

◆ Dump()

void G4GEMChannelVI::Dump ( ) const
finalvirtual

Reimplemented from G4VEvaporationChannel.

Definition at line 126 of file G4GEMChannelVI.cc.

127{}

◆ EmittedFragment()

G4Fragment * G4GEMChannelVI::EmittedFragment ( G4Fragment theNucleus)
finalvirtual

Reimplemented from G4VEvaporationChannel.

Definition at line 96 of file G4GEMChannelVI.cc.

97{
98 // assumed, that TotalProbability(...) was already called
99 // if value iz zero no possiblity to sample final state
100 G4Fragment* evFragment = nullptr;
101 G4LorentzVector lv0 = theNucleus->GetMomentum();
102 if(resA <= 4 || fProbability->GetProbability() == 0.0) {
103 G4double ekin =
104 std::max(0.5*(mass*mass - resMass*resMass + evapMass2)/mass
105 - evapMass, 0.0);
106 G4LorentzVector lv(std::sqrt(ekin*(ekin + 2.0*evapMass))
107 *G4RandomDirection(), ekin + evapMass);
108 lv.boost(lv0.boostVector());
109 evFragment = new G4Fragment(A, Z, lv);
110 lv0 -= lv;
111 } else {
112 evFragment = fProbability->SampleEvaporationFragment();
113 G4LorentzVector lv = evFragment->GetMomentum();
114 lv.boost(lv0.boostVector());
115 evFragment->SetMomentum(lv);
116 lv0 -= lv;
117 }
118 if(evFragment != nullptr) { evFragment->SetCreatorModelID(secID); }
119 theNucleus->SetZandA_asInt(resZ, resA);
120 theNucleus->SetMomentum(lv0);
121 theNucleus->SetCreatorModelID(secID);
122
123 return evFragment;
124}
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
G4Fragment * SampleEvaporationFragment()

◆ GetEmissionProbability()

G4double G4GEMChannelVI::GetEmissionProbability ( G4Fragment theNucleus)
finalvirtual

Implements G4VEvaporationChannel.

Definition at line 66 of file G4GEMChannelVI.cc.

67{
68 fProbability->ResetProbability();
69 fragZ = fragment->GetZ_asInt();
70 fragA = fragment->GetA_asInt();
71 resZ = fragZ - Z;
72 resA = fragA - A;
73 if(resA < A || resA < resZ || resZ < 0 || (resA == A && resZ < Z)) {
74 return 0.0;
75 }
76
77 const G4double exc = fragment->GetExcitationEnergy();
78 const G4double delta0 =
79 std::max(pairingCorrection->GetPairingCorrection(fragA, fragZ),0.0);
80 if(exc < delta0) { return 0.0; }
81
82 resMass = G4NucleiProperties::GetNuclearMass(resA, resZ);
83 const G4double fragM = fragment->GetGroundStateMass() + exc;
84 const G4double CB = cBarrier->GetCoulombBarrier(resA, resZ, exc);
85
86 const G4double delta1 =
87 std::max(0.0,pairingCorrection->GetPairingCorrection(resA,resZ));
88 if(fragM <= resMass + CB + delta1) { return 0.0; }
89
90 fProbability->SetDecayKinematics(resZ, resA, resMass, fragM);
91 G4double prob = fProbability->ComputeTotalProbability(*fragment, CB);
92 //G4cout<<"G4EvaporationChannel: probability= "<< prob <<G4endl;
93 return prob;
94}
G4double ComputeTotalProbability(const G4Fragment &, G4double CB)
G4double GetPairingCorrection(G4int A, G4int Z) const
virtual G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U=0.0) const =0
void SetDecayKinematics(G4int rZ, G4int rA, G4double rmass, G4double fmass)

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