Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4FermiBreakUpVI.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// FermiBreakUp de-excitation model
28// by V. Ivanchenko (July 2016)
29//
30
31#include "G4FermiBreakUpVI.hh"
34#include "G4FermiChannels.hh"
35#include "G4FermiPair.hh"
36#include "G4RandomDirection.hh"
38#include "Randomize.hh"
40
41G4FermiFragmentsPoolVI* G4FermiBreakUpVI::thePool = nullptr;
42
43#ifdef G4MULTITHREADED
44G4Mutex G4FermiBreakUpVI::FermiBreakUpVIMutex = G4MUTEX_INITIALIZER;
45#endif
46
48 : theDecay(nullptr), rndmEngine(nullptr), maxZ(9), maxA(17), secID(-1)
49{
50 frag.reserve(10);
51 lvect.reserve(10);
52 Z = A = spin = 0;
53 secID = G4PhysicsModelCatalog::GetModelID("model_G4FermiBreakUpVI");
54 mass = elim = excitation = 0.0;
55 tolerance = CLHEP::MeV;
56 frag1 = frag2 = nullptr;
57 prob.resize(12,0.0);
58 Initialise();
59}
60
62{
64 delete thePool;
65 thePool = nullptr;
66 }
67}
68
70{
71 if(verbose > 1) {
72 G4cout << "### G4FermiBreakUpVI::Initialise(): " << thePool << G4endl;
73 }
74 if(thePool == nullptr) { InitialisePool(); }
75 theDecay = thePool->FermiDecayProbability();
76 elim = thePool->GetEnergyLimit();
77}
78
79void G4FermiBreakUpVI::InitialisePool()
80{
81#ifdef G4MULTITHREADED
82 G4MUTEXLOCK(&G4FermiBreakUpVI::FermiBreakUpVIMutex);
83#endif
84 if(thePool == nullptr) {
85 thePool = new G4FermiFragmentsPoolVI();
86 }
87#ifdef G4MULTITHREADED
88 G4MUTEXUNLOCK(&G4FermiBreakUpVI::FermiBreakUpVIMutex);
89#endif
90}
91
93{
94 return (ZZ < maxZ && AA < maxA && AA > 0 && eexc <= elim
95 && thePool->HasChannels(ZZ, AA, eexc));
96}
97
99 G4Fragment* theNucleus)
100{
101 if(verbose > 1) {
102 G4cout << "### G4FermiBreakUpVI::BreakFragment start new fragment "
103 << G4endl;
104 G4cout << *theNucleus << G4endl;
105 }
106
107 // initial fragment
108 Z = theNucleus->GetZ_asInt();
109 A = theNucleus->GetA_asInt();
110 excitation = theNucleus->GetExcitationEnergy();
111 mass = theNucleus->GetGroundStateMass() + excitation;
112 spin = -1;
113
114 lv0 = theNucleus->GetMomentum();
115 rndmEngine = G4Random::getTheEngine();
116
117 // sample first decay of an initial state
118 // if not possible to decay - exit
119 if(!SampleDecay()) {
120 return;
121 }
122
123 G4double time = theNucleus->GetCreationTime();
124 delete theNucleus;
125
126 static const G4int imax = 100;
127
128 // loop over vector of Fermi fragments
129 // vector may grow at each iteraction
130 for(size_t i=0; i<frag.size(); ++i) {
131 Z = frag[i]->GetZ();
132 A = frag[i]->GetA();
133 spin = frag[i]->GetSpin();
134 mass = frag[i]->GetTotalEnergy();
135 lv0 = lvect[i];
136 if(verbose > 1) {
137 G4cout << "# FermiFrag " << i << ". Z= " << Z << " A= " << A
138 << " mass= " << mass << " exc= "
139 << frag[i]->GetExcitationEnergy() << G4endl;
140 }
141 // stable fragment
142 if(!SampleDecay()) {
143 if(verbose > 1) { G4cout << " New G4Fragment" << G4endl; }
144 G4Fragment* f = new G4Fragment(A, Z, lv0);
145 f->SetSpin(0.5*spin);
146 f->SetCreationTime(time);
147 f->SetCreatorModelID(secID);
148 theResult->push_back(f);
149 }
150 // limit the loop
151 if(i == imax) {
152 break;
153 }
154 }
155 frag.clear();
156 lvect.clear();
157}
158
159G4bool G4FermiBreakUpVI::SampleDecay()
160{
161 const G4FermiChannels* chan = thePool->ClosestChannels(Z, A, mass);
162 if(nullptr == chan) { return false; }
163 size_t nn = chan->GetNumberOfChannels();
164 if(verbose > 1) {
165 G4cout << "== SampleDecay " << nn << " channels Eex= "
166 << chan->GetExcitation() << G4endl;
167 }
168 if(0 == nn) { return false; }
169
170 const G4FermiPair* fpair = nullptr;
171
172 // one unstable fragment
173 if(1 == nn) {
174 fpair = chan->GetPair(0);
175
176 // more pairs
177 } else {
178
179 // in static probabilities may be used
180 if(std::abs(excitation - chan->GetExcitation()) < tolerance) {
181 fpair = chan->SamplePair(rndmEngine->flat());
182
183 } else {
184
185 // recompute probabilities
186 const std::vector<const G4FermiPair*>& pvect = chan->GetChannels();
187 if(nn > 12) { prob.resize(nn, 0.0); }
188 G4double ptot = 0.0;
189 if(verbose > 2) {
190 G4cout << "Start recompute probabilities" << G4endl;
191 }
192 for(size_t i=0; i<nn; ++i) {
193 ptot += theDecay->ComputeProbability(Z, A, -1, mass,
194 pvect[i]->GetFragment1(),
195 pvect[i]->GetFragment2());
196 prob[i] = ptot;
197 if(verbose > 2) {
198 G4cout << i << ". " << prob[i]
199 << " Z1= " << pvect[i]->GetFragment1()->GetZ()
200 << " A1= " << pvect[i]->GetFragment1()->GetA()
201 << " Z2= " << pvect[i]->GetFragment2()->GetZ()
202 << " A2= " << pvect[i]->GetFragment2()->GetA()
203 << G4endl;
204 }
205 }
206 ptot *= rndmEngine->flat();
207 for(size_t i=0; i<nn; ++i) {
208 if(ptot <= prob[i] || i+1 == nn) {
209 fpair = pvect[i];
210 break;
211 }
212 }
213 }
214 }
215 if(!fpair) { return false; }
216
217 frag1 = fpair->GetFragment1();
218 frag2 = fpair->GetFragment2();
219
220 G4double mass1 = frag1->GetTotalEnergy();
221 G4double mass2 = frag2->GetTotalEnergy();
222 if(verbose > 2) {
223 G4cout << " M= " << mass << " M1= " << mass1 << " M2= " << mass2
224 << " Exc1= " << frag1->GetExcitationEnergy()
225 << " Exc2= " << frag2->GetExcitationEnergy() << G4endl;
226 }
227 // sample fragment1
228 G4double e1 = 0.5*(mass*mass - mass2*mass2 + mass1*mass1)/mass;
229 //G4cout << "ekin1= " << e1 - mass1 << G4endl;
230 G4double p1(0.0);
231 if(e1 > mass1) {
232 p1 = std::sqrt((e1 - mass1)*(e1 + mass1));
233 } else {
234 e1 = mass1;
235 }
237 G4LorentzVector lv1 = G4LorentzVector(v*p1, e1);
238
239 // compute kinematics
240 boostVector = lv0.boostVector();
241 lv1.boost(boostVector);
242
243 lv0 -= lv1;
244
245 G4double e2 = lv0.e();
246 if(e2 < mass2) {
247 lv0.set(0.,0.,0.,mass2);
248 }
249
250 frag.push_back(frag1);
251 frag.push_back(frag2);
252 lvect.push_back(lv1);
253 lvect.push_back(lv0);
254
255 return true;
256}
257
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
CLHEP::HepLorentzVector G4LorentzVector
G4ThreeVector G4RandomDirection()
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
std::mutex G4Mutex
Definition: G4Threading.hh:81
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
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
void set(double x, double y, double z, double t)
virtual double flat()=0
void Initialise() final
void BreakFragment(G4FragmentVector *, G4Fragment *theNucleus) final
G4bool IsApplicable(G4int ZZ, G4int AA, G4double etot) const final
G4double GetExcitation() const
const G4FermiPair * SamplePair(G4double rand) const
const G4FermiPair * GetPair(size_t idx) const
size_t GetNumberOfChannels() const
const std::vector< const G4FermiPair * > & GetChannels() const
G4double ComputeProbability(G4int Z, G4int A, G4int spin, G4double TotalE, const G4FermiFragment *f1, const G4FermiFragment *f2) const
G4double GetExcitationEnergy() const
G4double GetTotalEnergy(void) const
const G4FermiDecayProbability * FermiDecayProbability() const
G4bool HasChannels(G4int Z, G4int A, G4double exc) const
const G4FermiChannels * ClosestChannels(G4int Z, G4int A, G4double mass) const
const G4FermiFragment * GetFragment2() const
Definition: G4FermiPair.hh:101
const G4FermiFragment * GetFragment1() const
Definition: G4FermiPair.hh:96
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:317
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:312
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:322
void SetCreatorModelID(G4int value)
Definition: G4Fragment.hh:433
G4double GetCreationTime() const
Definition: G4Fragment.hh:479
G4int GetZ_asInt() const
Definition: G4Fragment.hh:289
void SetCreationTime(G4double time)
Definition: G4Fragment.hh:484
void SetSpin(G4double value)
Definition: G4Fragment.hh:443
G4int GetA_asInt() const
Definition: G4Fragment.hh:284
static G4int GetModelID(const G4int modelIndex)
G4bool IsMasterThread()
Definition: G4Threading.cc:124