Geant4-11
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"
55
56const G4int G4UnstableFragmentBreakUp::Zfr[] = {0, 1, 1, 1, 2, 2};
57const G4int G4UnstableFragmentBreakUp::Afr[] = {1, 1, 2, 3, 3, 4};
58
60{
62 for(G4int i=0; i<6; ++i) {
64 }
65 fSecID = G4PhysicsModelCatalog::GetModelID("model_G4UnstableFragmentBreakUp");
66}
67
69{}
70
72 G4Fragment* nucleus)
73{
74 //G4cout << "G4UnstableFragmentBreakUp::EmittedFragment" << G4endl;
75 G4Fragment* frag = nullptr;
76
77 G4int Z = nucleus->GetZ_asInt();
78 G4int A = nucleus->GetA_asInt();
79
80 G4LorentzVector lv = nucleus->GetMomentum();
81 G4double time = nucleus->GetCreationTime();
82
83 G4double mass1(0.0), mass2(0.0);
84
85 // look for the decay channel with normal masses
86 // without Coulomb barrier and paring corrections
87 // 1 - recoil, 2 - emitted light ion
88 if(fVerbose > 1) {
89 G4cout << "#Unstable decay " << " Z= " << Z << " A= " << A
90 << " Eex(MeV)= " << nucleus->GetExcitationEnergy() << G4endl;
91 }
92 const G4double tolerance = 10*CLHEP::eV;
93 const G4double dmlimit = 0.2*CLHEP::MeV;
94 G4double mass = lv.mag();
95 G4double exca = -1000.0;
96 G4bool isChannel = false;
97 G4int idx = -1;
98 for(G4int i=0; i<6; ++i) {
99 G4int Zres = Z - Zfr[i];
100 G4int Ares = A - Afr[i];
101 if(Zres >= 0 && Ares >= Zres && Ares >= Afr[i]) {
102 if(Ares <= 4) {
103 for(G4int j=0; j<6; ++j) {
104 if(Zres == Zfr[j] && Ares == Afr[j]) {
105 /*
106 G4cout << "i= " << i << " j= " << j << " Zres= " << Zres
107 << " Ares= " << Ares << " dm= " << mass - masses[i] - masses[j]
108 << G4endl;
109 */
110 G4double delm = mass - masses[i] - masses[j];
111 if(delm > exca) {
112 mass2 = masses[i]; // emitted
113 mass1 = masses[j]; // recoil
114 exca = delm;
115 idx = i;
116 if(delm > 0.0) { isChannel = true; }
117 break;
118 }
119 }
120 }
121 }
122 if(isChannel) { break; }
123 // no simple channel
125 G4double e = mass - mres - masses[i];
126 // select excited state
127 const G4LevelManager* lman = fLevelData->GetLevelManager(Zres, Ares);
128 if(lman && e >= 0.0) {
129 mass2 = masses[i];
130 mass1 = mres + e*G4UniformRand();
131 idx = i;
132 isChannel = true;
133 break;
134 }
135 // if physical channel is not identified
136 // check excitation energy
137 if(e > exca) {
138 mass2 = masses[i];
139 mass1 = mres;
140 if(e > 0.0) { mass1 += e; }
141 exca = e;
142 idx = i;
143 }
144 }
145 }
146 G4double massmin = mass1 + mass2;
147 if(mass < massmin) {
148 if(mass + dmlimit < massmin) { return false; }
149 if(fVerbose > 1) {
150 G4cout << "#Unstable decay correction: Z= " << Z << " A= " << A
151 << " idx= " << idx
152 << " deltaM(MeV)= " << mass - massmin
153 << G4endl;
154 }
155 mass = massmin;
156 G4double e = std::max(lv.e(), mass + tolerance);
157 G4double mom = std::sqrt((e - mass)*(e + mass));
158 G4ThreeVector dir = lv.vect().unit();
159 lv.set(dir*mom, e);
160 }
161
162 // compute energy of light fragment
163 G4double e2 = 0.5*((mass - mass1)*(mass + mass1) + mass2*mass2)/mass;
164 e2 = std::max(e2, mass2);
165 G4double mom = std::sqrt((e2 - mass2)*(e2 + mass2));
166
167 // sample decay
168 G4ThreeVector bst = lv.boostVector();
170 G4LorentzVector mom2 = G4LorentzVector(v*mom, e2);
171 mom2.boost(bst);
172 frag = new G4Fragment(Afr[idx], Zfr[idx], mom2);
173 frag->SetCreationTime(time);
175 results->push_back(frag);
176
177 // residual
178 lv -= mom2;
179 Z -= Zfr[idx];
180 A -= Afr[idx];
181
182 nucleus->SetZandA_asInt(Z, A);
183 nucleus->SetMomentum(lv);
184 nucleus->SetCreatorModelID(fSecID);
185 return true;
186}
187
189{
190 return 0.0;
191}
static const G4double e2[44]
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:64
CLHEP::HepLorentzVector G4LorentzVector
G4ThreeVector G4RandomDirection()
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#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:299
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:323
void SetCreatorModelID(G4int value)
Definition: G4Fragment.hh:428
G4double GetCreationTime() const
Definition: G4Fragment.hh:464
G4int GetZ_asInt() const
Definition: G4Fragment.hh:276
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:469
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:328
void SetZandA_asInt(G4int Znew, G4int Anew)
Definition: G4Fragment.hh:281
G4int GetA_asInt() const
Definition: G4Fragment.hh:271
const G4LevelManager * GetLevelManager(G4int Z, G4int A)
static G4NuclearLevelData * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
static G4int GetModelID(const G4int modelIndex)
virtual G4bool BreakUpChain(G4FragmentVector *, G4Fragment *) final
virtual G4double GetEmissionProbability(G4Fragment *fragment) final
static constexpr double MeV
static constexpr double eV
T max(const T t1, const T t2)
brief Return the largest of the two arguments