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