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