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
G4StatMF.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// Hadronic Process: Nuclear De-excitations
29// by V. Lara
30
31#include "G4StatMF.hh"
33#include "G4SystemOfUnits.hh"
34#include "G4Pow.hh"
36
38{
39 _secID = G4PhysicsModelCatalog::GetModelID("model_G4StatMF");
40}
41
43
45{
46 if (theFragment.GetExcitationEnergy() <= 0.0) {
47 return nullptr;
48 }
49
50 // Maximun average multiplicity: M_0 = 2.6 for A ~ 200
51 // and M_0 = 3.3 for A <= 110
52 G4double MaxAverageMultiplicity =
54
55
56 // We'll use two kinds of ensembles
57 G4StatMFMicroCanonical * theMicrocanonicalEnsemble = 0;
58 G4StatMFMacroCanonical * theMacrocanonicalEnsemble = 0;
59
60 //-------------------------------------------------------
61 // Direct simulation part (Microcanonical ensemble)
62 //-------------------------------------------------------
63
64 // Microcanonical ensemble initialization
65 theMicrocanonicalEnsemble = new G4StatMFMicroCanonical(theFragment);
66
67 G4int Iterations = 0;
68 G4int IterationsLimit = 100000;
69 G4double Temperature = 0.0;
70
71 G4bool FirstTime = true;
72 G4StatMFChannel * theChannel = 0;
73
74 G4bool ChannelOk;
75 do { // Try to de-excite as much as IterationLimit permits
76 do {
77
78 G4double theMeanMult = theMicrocanonicalEnsemble->GetMeanMultiplicity();
79 if (theMeanMult <= MaxAverageMultiplicity) {
80 // G4cout << "MICROCANONICAL" << G4endl;
81 // Choose fragments atomic numbers and charges from direct simulation
82 theChannel = theMicrocanonicalEnsemble->ChooseAandZ(theFragment);
83 _theEnsemble = theMicrocanonicalEnsemble;
84 } else {
85 //-----------------------------------------------------
86 // Non direct simulation part (Macrocanonical Ensemble)
87 //-----------------------------------------------------
88 if (FirstTime) {
89 // Macrocanonical ensemble initialization
90 theMacrocanonicalEnsemble = new G4StatMFMacroCanonical(theFragment);
91 _theEnsemble = theMacrocanonicalEnsemble;
92 FirstTime = false;
93 }
94 // G4cout << "MACROCANONICAL" << G4endl;
95 // Select calculated fragment total multiplicity,
96 // fragment atomic numbers and fragment charges.
97 theChannel = theMacrocanonicalEnsemble->ChooseAandZ(theFragment);
98 }
99
100 ChannelOk = theChannel->CheckFragments();
101 if (!ChannelOk) delete theChannel;
102
103 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
104 } while (!ChannelOk);
105
106
107 if (theChannel->GetMultiplicity() <= 1) {
108 G4FragmentVector * theResult = new G4FragmentVector;
109 theResult->push_back(new G4Fragment(theFragment));
110 delete theMicrocanonicalEnsemble;
111 if (theMacrocanonicalEnsemble != 0) delete theMacrocanonicalEnsemble;
112 delete theChannel;
113 return theResult;
114 }
115
116 //--------------------------------------
117 // Second part of simulation procedure.
118 //--------------------------------------
119
120 // Find temperature of breaking channel.
121 Temperature = _theEnsemble->GetMeanTemperature(); // Initial guess for Temperature
122
123 if (FindTemperatureOfBreakingChannel(theFragment,theChannel,Temperature)) break;
124
125 // Do not forget to delete this unusable channel, for which we failed to find the temperature,
126 // otherwise for very proton-reach nuclei it would lead to memory leak due to large
127 // number of iterations. N.B. "theChannel" is created in G4StatMFMacroCanonical::ChooseZ()
128
129 // G4cout << " Iteration # " << Iterations << " Mean Temperature = " << Temperature << G4endl;
130
131 delete theChannel;
132
133 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
134 } while (Iterations++ < IterationsLimit );
135
136 // If Iterations >= IterationsLimit means that we couldn't solve for temperature
137 if (Iterations >= IterationsLimit)
138 throw G4HadronicException(__FILE__, __LINE__, "G4StatMF::BreakItUp: Was not possible to solve for temperature of breaking channel");
139
140 G4FragmentVector * theResult = theChannel->
141 GetFragments(theFragment.GetA_asInt(),theFragment.GetZ_asInt(),Temperature);
142
143 // ~~~~~~ Energy conservation Patch !!!!!!!!!!!!!!!!!!!!!!
144 // Original nucleus 4-momentum in CM system
145 G4LorentzVector InitialMomentum(theFragment.GetMomentum());
146 InitialMomentum.boost(-InitialMomentum.boostVector());
147 G4double ScaleFactor = 0.0;
148 G4double SavedScaleFactor = 0.0;
149 do {
150 G4double FragmentsEnergy = 0.0;
151 G4FragmentVector::iterator j;
152 for (j = theResult->begin(); j != theResult->end(); j++)
153 FragmentsEnergy += (*j)->GetMomentum().e();
154 SavedScaleFactor = ScaleFactor;
155 ScaleFactor = InitialMomentum.e()/FragmentsEnergy;
156 G4ThreeVector ScaledMomentum(0.0,0.0,0.0);
157 for (j = theResult->begin(); j != theResult->end(); j++) {
158 ScaledMomentum = ScaleFactor * (*j)->GetMomentum().vect();
159 G4double Mass = (*j)->GetMomentum().m();
160 G4LorentzVector NewMomentum;
161 NewMomentum.setVect(ScaledMomentum);
162 NewMomentum.setE(std::sqrt(ScaledMomentum.mag2()+Mass*Mass));
163 (*j)->SetMomentum(NewMomentum);
164 }
165 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
166 } while (ScaleFactor > 1.0+1.e-5 && std::abs(ScaleFactor-SavedScaleFactor)/ScaleFactor > 1.e-10);
167 // ~~~~~~ End of patch !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
168
169 // Perform Lorentz boost
170 G4FragmentVector::iterator i;
171 for (i = theResult->begin(); i != theResult->end(); i++) {
172 G4LorentzVector FourMom = (*i)->GetMomentum();
173 FourMom.boost(theFragment.GetMomentum().boostVector());
174 (*i)->SetMomentum(FourMom);
175 (*i)->SetCreatorModelID(_secID);
176 }
177
178 // garbage collection
179 delete theMicrocanonicalEnsemble;
180 if (theMacrocanonicalEnsemble != 0) delete theMacrocanonicalEnsemble;
181 delete theChannel;
182
183 return theResult;
184}
185
186
187G4bool G4StatMF::FindTemperatureOfBreakingChannel(const G4Fragment & theFragment,
188 const G4StatMFChannel * aChannel,
189 G4double & Temperature)
190 // This finds temperature of breaking channel.
191{
192 G4int A = theFragment.GetA_asInt();
193 G4int Z = theFragment.GetZ_asInt();
194 G4double U = theFragment.GetExcitationEnergy();
195
196 G4double T = std::max(Temperature,0.0012*MeV);
197 G4double Ta = T;
198 G4double TotalEnergy = CalcEnergy(A,Z,aChannel,T);
199
200 G4double Da = (U - TotalEnergy)/U;
201 G4double Db = 0.0;
202
203 // bracketing the solution
204 if (Da == 0.0) {
205 Temperature = T;
206 return true;
207 } else if (Da < 0.0) {
208 do {
209 T *= 0.5;
210 if (T < 0.001*MeV) return false;
211
212 TotalEnergy = CalcEnergy(A,Z,aChannel,T);
213
214 Db = (U - TotalEnergy)/U;
215 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
216 } while (Db < 0.0);
217
218 } else {
219 do {
220 T *= 1.5;
221
222 TotalEnergy = CalcEnergy(A,Z,aChannel,T);
223
224 Db = (U - TotalEnergy)/U;
225 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
226 } while (Db > 0.0);
227 }
228
229 G4double eps = 1.0e-14 * std::abs(T-Ta);
230 //G4double eps = 1.0e-3 ;
231
232 // Start the bisection method
233 for (G4int j = 0; j < 1000; j++) {
234 G4double Tc = (Ta+T)*0.5;
235 if (std::abs(Ta-Tc) <= eps) {
236 Temperature = Tc;
237 return true;
238 }
239
240 T = Tc;
241 TotalEnergy = CalcEnergy(A,Z,aChannel,T);
242 G4double Dc = (U - TotalEnergy)/U;
243
244 if (Dc == 0.0) {
245 Temperature = Tc;
246 return true;
247 }
248 if (Da*Dc < 0.0) {
249 T = Tc;
250 Db = Dc;
251 } else {
252 Ta = Tc;
253 Da = Dc;
254 }
255 }
256
257 Temperature = (Ta+T)*0.5;
258 return false;
259}
260
261G4double G4StatMF::CalcEnergy(G4int A, G4int Z, const G4StatMFChannel * aChannel,
262 G4double T)
263{
265 G4double ChannelEnergy = aChannel->GetFragmentsEnergy(T);
266 return -MassExcess0 + G4StatMFParameters::GetCoulomb() + ChannelEnergy;
267}
268
269
270
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
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]
double mag2() const
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
void setVect(const Hep3Vector &)
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:312
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:322
G4int GetZ_asInt() const
Definition: G4Fragment.hh:289
G4int GetA_asInt() const
Definition: G4Fragment.hh:284
static G4double GetMassExcess(const G4int A, const G4int Z)
static G4int GetModelID(const G4int modelIndex)
G4bool CheckFragments(void)
G4double GetFragmentsEnergy(G4double T) const
size_t GetMultiplicity(void)
G4StatMFChannel * ChooseAandZ(const G4Fragment &theFragment)
G4StatMFChannel * ChooseAandZ(const G4Fragment &theFragment)
static G4double GetMaxAverageMultiplicity(G4int A)
static G4double GetCoulomb()
G4FragmentVector * BreakItUp(const G4Fragment &theNucleus) override
Definition: G4StatMF.cc:44
G4StatMF()
Definition: G4StatMF.cc:37
~G4StatMF()
Definition: G4StatMF.cc:42
G4double GetMeanTemperature(void) const
G4double GetMeanMultiplicity(void) const