Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GEMChannelVI.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// GEM de-excitation model
28// by V. Ivanchenko (July 2019)
29//
30
31#include "G4GEMChannelVI.hh"
32#include "G4GEMProbabilityVI.hh"
33#include "G4VCoulombBarrier.hh"
34#include "G4CoulombBarrier.hh"
36#include "G4NuclearLevelData.hh"
37#include "G4LevelManager.hh"
38#include "G4NucleiProperties.hh"
39#include "G4RandomDirection.hh"
41#include "Randomize.hh"
42
43namespace
44{
45 const G4double minExc = 1.0*CLHEP::MeV;
46 const G4int nProbMax = 10;
47}
48
50 : A(theA), Z(theZ)
51{
53 pairingCorrection = nData->GetPairingCorrection();
54 const G4LevelManager* lManager = nullptr;
55 if (A > 4) { lManager = nData->GetLevelManager(Z, A); }
56 fEvapMass = G4NucleiProperties::GetNuclearMass(A, Z);
57 fEvapMass2 = fEvapMass*fEvapMass;
58
59 cBarrier = new G4CoulombBarrier(A, Z);
60 fProbability = new G4GEMProbabilityVI(A, Z, lManager);
61
62 fCoeff = CLHEP::millibarn/((CLHEP::pi*CLHEP::hbarc)*(CLHEP::pi*CLHEP::hbarc));
63
64 secID = G4PhysicsModelCatalog::GetModelID("model_G4GEMChannelVI");
65 if (Z == 0 && A == 1) {
66 indexC = 0;
67 fCoeff *= 2.0;
68 } else if (Z == 1 && A == 1) {
69 indexC = 1;
70 fCoeff *= 2.0;
71 } else if (Z == 1 && A == 2) {
72 indexC = 2;
73 fCoeff *= 3.0;
74 } else if (Z == 1 && A == 3) {
75 indexC = 3;
76 fCoeff *= 2.0;
77 } else if (Z == 2 && A == 3) {
78 indexC = 4;
79 fCoeff *= 2.0;
80 } else if (Z == 2 && A == 4) {
81 indexC = 5;
82 } else {
83 indexC = 6;
84 }
85}
86
88{
89 delete cBarrier;
90 delete fProbability;
91}
92
94{
95 fProbability->Initialise();
97}
98
100{
101 fProbability->ResetProbability();
102 fragZ = fragment->GetZ_asInt();
103 fragA = fragment->GetA_asInt();
104 resZ = fragZ - Z;
105 resA = fragA - A;
106 if(resA < A || resA < resZ || resZ < 0 || (resA == A && resZ < Z)) {
107 return 0.0;
108 }
109
110 fExc = fragment->GetExcitationEnergy();
111 fMass = fragment->GetGroundStateMass() + fExc;
112 fResMass = G4NucleiProperties::GetNuclearMass(resA, resZ);
113
114 // limit for the case when both evaporation and residual
115 // fragments are in ground states
116 if (fMass <= fEvapMass + fResMass) { return 0.0; }
117
118 if (Z > 0) {
119 bCoulomb = cBarrier->GetCoulombBarrier(resA, resZ, 0.0);
120 }
121 G4double de = fMass - fEvapMass - fResMass - bCoulomb;
122 nProb = (G4int)(de/minExc);
123 if (nProb <= 1 || indexC < 6 || resA <= 4) {
124 nProb = 1;
125 } else {
126 nProb = std::min(nProb, nProbMax);
127 }
128 if (2 < fVerbose) {
129 G4cout << "## G4GEMChannelVI::GetEmissionProbability fragZ="
130 << fragZ << " fragA=" << fragA << " Z=" << Z << " A=" << A
131 << " Eex(MeV)=" << fExc << " nProb=" << nProb
132 << G4endl;
133 }
134 fProbability->SetDecayKinematics(resZ, resA, fResMass, fMass);
135 G4double sump = 0.0;
136 for (G4int i=0; i<nProb; ++i) {
137 G4double exc = std::min(minExc*i, de);
138 G4double m1 = fEvapMass + exc;
139 G4double e2 = 0.5*((fMass-fResMass)*(fMass+fResMass) + m1*m1)/fMass - m1;
140 G4double m2 = fMass - m1 - 0.5*bCoulomb;
141 if (m2 < fResMass) {
142 nProb = i;
143 break;
144 }
145 G4double e1 = std::max(0.5*((fMass-m2)*(fMass+m2) + m1*m1)/fMass - m1, 0.0);
146 if (e1 >= e2) {
147 nProb = i;
148 break;
149 }
150 sump += fProbability->TotalProbability(*fragment, e1, e2, bCoulomb, fExc, exc);
151 fEData[i].exc = exc;
152 fEData[i].ekin1 = e1;
153 fEData[i].ekin2 = e2;
154 fEData[i].prob = sump;
155 }
156 return sump;
157}
158
160{
161 // assumed, that TotalProbability(...) was already called
162 // if value iz zero no possiblity to sample final state
163 G4Fragment* evFragment = nullptr;
164 G4LorentzVector lv0 = theNucleus->GetMomentum();
165 G4double ekin;
166 G4double exc = 0.0;
167 G4double probMax = std::max(fEData[nProb - 1].prob, 0.0);
168 if (0.0 >= probMax) {
169 ekin = std::max(0.5*(fMass*fMass - fResMass*fResMass + fEvapMass2)
170 /fMass - fEvapMass, 0.0);
171 } else if (1 == nProb) {
172 ekin = fProbability->SampleEnergy(fEData[0].ekin1, fEData[0].ekin2,
173 bCoulomb, fExc, 0.0);
174 } else {
175 G4double p = G4UniformRand()*probMax;
176 G4int i{1};
177 for (; i<nProb; ++i) {
178 if (p <= fEData[i].prob) { break; }
179 }
180 G4double e1 = fEData[i - 1].exc;
181 G4double e2 = fEData[i].exc;
182 G4double p1 = fEData[i - 1].prob;
183 G4double p2 = fEData[i].prob;
184 exc = e1 + (e2 - e1)*(p - p1)/(p2 - p1);
185 ekin = fProbability->SampleEnergy(fEData[i].ekin1, fEData[i].ekin2,
186 bCoulomb, fExc, exc);
187 }
188 G4double m1 = fEvapMass + exc;
189 G4LorentzVector lv(std::sqrt(ekin*(ekin + 2.0*m1))
190 *G4RandomDirection(), ekin + m1);
191 lv.boost(lv0.boostVector());
192 evFragment = new G4Fragment(A, Z, lv);
193 lv0 -= lv;
194 evFragment->SetCreatorModelID(secID);
195 theNucleus->SetZandA_asInt(resZ, resA);
196 theNucleus->SetMomentum(lv0);
197 theNucleus->SetCreatorModelID(secID);
198
199 return evFragment;
200}
201
203{}
204
205
206
CLHEP::HepLorentzVector G4LorentzVector
G4ThreeVector G4RandomDirection()
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
G4double GetGroundStateMass() const
void SetZandA_asInt(G4int Znew, G4int Anew, G4int Lnew=0)
G4double GetExcitationEnergy() const
const G4LorentzVector & GetMomentum() const
void SetCreatorModelID(G4int value)
G4int GetZ_asInt() const
void SetMomentum(const G4LorentzVector &value)
G4int GetA_asInt() const
G4GEMChannelVI(G4int theA, G4int theZ)
void Dump() const override
~G4GEMChannelVI() override
G4double GetEmissionProbability(G4Fragment *theNucleus) override
void Initialise() override
G4Fragment * EmittedFragment(G4Fragment *theNucleus) override
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)