Geant4 11.1.1
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
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}
59
61{
62 delete cBarrier;
63 delete fProbability;
64}
65
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}
95
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}
125
127{}
128
129
130
G4ThreeVector G4RandomDirection()
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:317
void SetZandA_asInt(G4int Znew, G4int Anew, G4int Lnew=0)
Definition: G4Fragment.hh:295
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:312
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:322
void SetCreatorModelID(G4int value)
Definition: G4Fragment.hh:433
G4int GetZ_asInt() const
Definition: G4Fragment.hh:289
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:327
G4int GetA_asInt() const
Definition: G4Fragment.hh:284
void Dump() const final
G4GEMChannelVI(G4int theA, G4int theZ)
G4double GetEmissionProbability(G4Fragment *theNucleus) final
~G4GEMChannelVI() final
G4Fragment * EmittedFragment(G4Fragment *theNucleus) final
G4Fragment * SampleEvaporationFragment()
G4double ComputeTotalProbability(const G4Fragment &, G4double CB)
const G4LevelManager * GetLevelManager(G4int Z, G4int A)
G4PairingCorrection * GetPairingCorrection()
static G4NuclearLevelData * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double GetPairingCorrection(G4int A, G4int Z) const
static G4int GetModelID(const G4int modelIndex)
virtual G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U=0.0) const =0
void SetDecayKinematics(G4int rZ, G4int rA, G4double rmass, G4double fmass)