Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UnstableFragmentBreakUp.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// -------------------------------------------------------------------
28// GEANT 4 class file
29//
30// CERN, Geneva, Switzerland
31//
32// File name: G4UnstableFragmentBreakUp
33//
34// Author: V.Ivanchenko
35//
36// Creation date: 11 May 2010
37//
38//Modifications:
39//
40// -------------------------------------------------------------------
41//
42
44#include "Randomize.hh"
45#include "G4RandomDirection.hh"
47#include "G4LorentzVector.hh"
48#include "G4Fragment.hh"
49#include "G4FragmentVector.hh"
50#include "G4NucleiProperties.hh"
51#include "G4NuclearLevelData.hh"
52#include "G4LevelManager.hh"
53#include "Randomize.hh"
54
55const G4int G4UnstableFragmentBreakUp::Zfr[] = {0, 1, 1, 1, 2, 2};
56const G4int G4UnstableFragmentBreakUp::Afr[] = {1, 1, 2, 3, 3, 4};
57
59{
61 for(G4int i=0; i<6; ++i) {
62 masses[i] = G4NucleiProperties::GetNuclearMass(Afr[i], Zfr[i]);
63 }
64}
65
67{}
68
70 G4Fragment* nucleus)
71{
72 //G4cout << "G4UnstableFragmentBreakUp::EmittedFragment" << G4endl;
73 G4Fragment* frag = nullptr;
74
75 G4int Z = nucleus->GetZ_asInt();
76 G4int A = nucleus->GetA_asInt();
77
78 G4LorentzVector lv = nucleus->GetMomentum();
79 G4double time = nucleus->GetCreationTime();
80
81 G4double mass1(0.0), mass2(0.0);
82
83 // look for the decay channel with normal masses
84 // without Coulomb barrier and paring corrections
85 // 1 - recoil, 2 - emitted light ion
86 if(fVerbose > 1) {
87 G4cout << "#Unstable decay " << " Z= " << Z << " A= " << A
88 << " Eex(MeV)= " << nucleus->GetExcitationEnergy() << G4endl;
89 }
90 const G4double tolerance = 10*CLHEP::eV;
91 const G4double dmlimit = 0.2*CLHEP::MeV;
92 G4double mass = lv.mag();
93 G4double exca = -1000.0;
94 G4bool isChannel = false;
95 G4int idx = -1;
96 for(G4int i=0; i<6; ++i) {
97 G4int Zres = Z - Zfr[i];
98 G4int Ares = A - Afr[i];
99 if(Zres >= 0 && Ares >= Zres && Ares >= Afr[i]) {
100 if(Ares <= 4) {
101 for(G4int j=0; j<6; ++j) {
102 if(Zres == Zfr[j] && Ares == Afr[j]) {
103 /*
104 G4cout << "i= " << i << " j= " << j << " Zres= " << Zres
105 << " Ares= " << Ares << " dm= " << mass - masses[i] - masses[j]
106 << G4endl;
107 */
108 G4double delm = mass - masses[i] - masses[j];
109 if(delm > exca) {
110 mass2 = masses[i]; // emitted
111 mass1 = masses[j]; // recoil
112 exca = delm;
113 idx = i;
114 if(delm > 0.0) { isChannel = true; }
115 break;
116 }
117 }
118 }
119 }
120 if(isChannel) { break; }
121 // no simple channel
123 G4double e = mass - mres - masses[i];
124 // select excited state
125 const G4LevelManager* lman = fLevelData->GetLevelManager(Zres, Ares);
126 if(lman && e >= 0.0) {
127 mass2 = masses[i];
128 mass1 = mres + e*G4UniformRand();
129 idx = i;
130 isChannel = true;
131 break;
132 }
133 // if physical channel is not identified
134 // check excitation energy
135 if(e > exca) {
136 mass2 = masses[i];
137 mass1 = mres;
138 if(e > 0.0) { mass1 += e; }
139 exca = e;
140 idx = i;
141 }
142 }
143 }
144 G4double massmin = mass1 + mass2;
145 if(mass < massmin) {
146 if(mass + dmlimit < massmin) { return false; }
147 if(fVerbose > 1) {
148 G4cout << "#Unstable decay correction: Z= " << Z << " A= " << A
149 << " idx= " << idx
150 << " deltaM(MeV)= " << mass - massmin
151 << G4endl;
152 }
153 mass = massmin;
154 G4double e = std::max(lv.e(), mass + tolerance);
155 G4double mom = std::sqrt((e - mass)*(e + mass));
156 G4ThreeVector dir = lv.vect().unit();
157 lv.set(dir*mom, e);
158 }
159
160 // compute energy of light fragment
161 G4double e2 = 0.5*((mass - mass1)*(mass + mass1) + mass2*mass2)/mass;
162 e2 = std::max(e2, mass2);
163 G4double mom = std::sqrt((e2 - mass2)*(e2 + mass2));
164
165 // sample decay
166 G4ThreeVector bst = lv.boostVector();
168 G4LorentzVector mom2 = G4LorentzVector(v*mom, e2);
169 mom2.boost(bst);
170 frag = new G4Fragment(Afr[idx], Zfr[idx], mom2);
171 frag->SetCreationTime(time);
172 results->push_back(frag);
173
174 // residual
175 lv -= mom2;
176 Z -= Zfr[idx];
177 A -= Afr[idx];
178
179 nucleus->SetZandA_asInt(Z, A);
180 nucleus->SetMomentum(lv);
181 return true;
182}
183
185{
186 return 0.0;
187}
double A(double temperature)
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:63
CLHEP::HepLorentzVector G4LorentzVector
G4ThreeVector G4RandomDirection()
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
void set(double x, double y, double z, double t)
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:275
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:299
G4double GetCreationTime() const
Definition: G4Fragment.hh:440
G4int GetZ_asInt() const
Definition: G4Fragment.hh:263
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:445
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:304
void SetZandA_asInt(G4int Znew, G4int Anew)
Definition: G4Fragment.hh:268
G4int GetA_asInt() const
Definition: G4Fragment.hh:258
const G4LevelManager * GetLevelManager(G4int Z, G4int A)
static G4NuclearLevelData * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
virtual G4bool BreakUpChain(G4FragmentVector *, G4Fragment *) final
virtual G4double GetEmissionProbability(G4Fragment *fragment) final