source: trunk/source/processes/hadronic/models/de_excitation/multifragmentation/src/G4StatMF.cc @ 846

Last change on this file since 846 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 9.8 KB
Line 
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,v 1.5 2006/06/29 20:24:43 gunter Exp $
28// GEANT4 tag $Name:  $
29//
30// Hadronic Process: Nuclear De-excitations
31// by V. Lara
32
33#include "G4StatMF.hh"
34
35
36
37// Default constructor
38G4StatMF::G4StatMF() : _theEnsemble(0) {}
39
40
41// Destructor
42G4StatMF::~G4StatMF() {} //{if (_theEnsemble != 0) delete _theEnsemble;}
43
44
45// Copy constructor
46G4StatMF::G4StatMF(const G4StatMF & ) : G4VMultiFragmentation()
47{
48    throw G4HadronicException(__FILE__, __LINE__, "G4StatMF::copy_constructor meant to not be accessable");
49}
50
51
52// Operators
53
54G4StatMF & G4StatMF::operator=(const G4StatMF & )
55{
56    throw G4HadronicException(__FILE__, __LINE__, "G4StatMF::operator= meant to not be accessable");
57    return *this;
58}
59
60
61G4bool G4StatMF::operator==(const G4StatMF & )
62{
63    throw G4HadronicException(__FILE__, __LINE__, "G4StatMF::operator== meant to not be accessable");
64    return false;
65}
66 
67
68G4bool G4StatMF::operator!=(const G4StatMF & )
69{
70    throw G4HadronicException(__FILE__, __LINE__, "G4StatMF::operator!= meant to not be accessable");
71    return true;
72}
73
74
75
76
77
78G4FragmentVector * G4StatMF::BreakItUp(const G4Fragment &theFragment)
79{
80  //    G4FragmentVector * theResult = new G4FragmentVector;
81
82  if (theFragment.GetExcitationEnergy() <= 0.0) {
83    G4FragmentVector * theResult = new G4FragmentVector;
84    theResult->push_back(new G4Fragment(theFragment));
85    return 0;
86  }
87
88
89  // Maximun average multiplicity: M_0 = 2.6 for A ~ 200
90  // and M_0 = 3.3 for A <= 110
91  G4double MaxAverageMultiplicity = 
92    G4StatMFParameters::GetMaxAverageMultiplicity(static_cast<G4int>(theFragment.GetA()));
93
94       
95    // We'll use two kinds of ensembles
96  G4StatMFMicroCanonical * theMicrocanonicalEnsemble = 0;
97  G4StatMFMacroCanonical * theMacrocanonicalEnsemble = 0;
98
99       
100        //-------------------------------------------------------
101        // Direct simulation part (Microcanonical ensemble)
102        //-------------------------------------------------------
103 
104        // Microcanonical ensemble initialization
105  theMicrocanonicalEnsemble = new G4StatMFMicroCanonical(theFragment);
106
107  G4int Iterations = 0;
108  G4double Temperature = 0.0;
109 
110  G4bool FirstTime = true;
111  G4StatMFChannel * theChannel = 0;
112
113  G4bool ChannelOk;
114  do {  // Try to de-excite as much as 10 times
115    do {
116     
117      G4double theMeanMult = theMicrocanonicalEnsemble->GetMeanMultiplicity();
118      if (theMeanMult <= MaxAverageMultiplicity) {
119        // G4cout << "MICROCANONICAL" << G4endl;
120        // Choose fragments atomic numbers and charges from direct simulation
121        theChannel = theMicrocanonicalEnsemble->ChooseAandZ(theFragment);
122        _theEnsemble = theMicrocanonicalEnsemble;
123      } else {
124        //-----------------------------------------------------
125        // Non direct simulation part (Macrocanonical Ensemble)
126        //-----------------------------------------------------
127        if (FirstTime) {
128          // Macrocanonical ensemble initialization
129          theMacrocanonicalEnsemble = new G4StatMFMacroCanonical(theFragment);
130          _theEnsemble = theMacrocanonicalEnsemble;
131          FirstTime = false;
132        }
133        // G4cout << "MACROCANONICAL" << G4endl;
134        // Select calculated fragment total multiplicity,
135        // fragment atomic numbers and fragment charges.
136        theChannel = theMacrocanonicalEnsemble->ChooseAandZ(theFragment);
137      }
138     
139      if (!(ChannelOk = theChannel->CheckFragments())) delete theChannel; 
140     
141    } while (!ChannelOk);
142   
143   
144    if (theChannel->GetMultiplicity() <= 1) {
145      G4FragmentVector * theResult = new G4FragmentVector;
146      theResult->push_back(new G4Fragment(theFragment));
147      delete theMicrocanonicalEnsemble;
148      if (theMacrocanonicalEnsemble != 0) delete theMacrocanonicalEnsemble;
149      delete theChannel;
150      return theResult;
151    }
152   
153    //--------------------------------------
154    // Second part of simulation procedure.
155    //--------------------------------------
156   
157    // Find temperature of breaking channel.
158    Temperature = _theEnsemble->GetMeanTemperature(); // Initial value for Temperature
159   
160    if (FindTemperatureOfBreakingChannel(theFragment,theChannel,Temperature)) break;
161   
162  } while (Iterations++ < 10);
163 
164 
165  // If Iterations >= 10 means that we couldn't solve for temperature
166  if (Iterations >= 10) 
167    throw G4HadronicException(__FILE__, __LINE__, "G4StatMF::BreakItUp: Was not possible to solve for temperature of breaking channel");
168 
169 
170  G4FragmentVector * theResult = theChannel->
171    GetFragments(theFragment.GetA(),theFragment.GetZ(),Temperature);
172 
173 
174       
175  // ~~~~~~ Energy conservation Patch !!!!!!!!!!!!!!!!!!!!!!
176  // Original nucleus 4-momentum in CM system
177  G4LorentzVector InitialMomentum(theFragment.GetMomentum());
178  InitialMomentum.boost(-InitialMomentum.boostVector());
179  G4double ScaleFactor = 0.0;
180  G4double SavedScaleFactor = 0.0;
181  do {
182    G4double FragmentsEnergy = 0.0;
183    G4FragmentVector::iterator j;
184    for (j = theResult->begin(); j != theResult->end(); j++) 
185      FragmentsEnergy += (*j)->GetMomentum().e();
186    SavedScaleFactor = ScaleFactor;
187    ScaleFactor = InitialMomentum.e()/FragmentsEnergy;
188    G4ThreeVector ScaledMomentum(0.0,0.0,0.0);
189    for (j = theResult->begin(); j != theResult->end(); j++) {
190      ScaledMomentum = ScaleFactor * (*j)->GetMomentum().vect();
191      G4double Mass = (*j)->GetMomentum().m();
192      G4LorentzVector NewMomentum;
193      NewMomentum.setVect(ScaledMomentum);
194      NewMomentum.setE(std::sqrt(ScaledMomentum.mag2()+Mass*Mass));
195      (*j)->SetMomentum(NewMomentum);           
196    }
197  } while (ScaleFactor > 1.0+1.e-5 && std::abs(ScaleFactor-SavedScaleFactor)/ScaleFactor > 1.e-10);
198  // ~~~~~~ End of patch !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
199 
200  // Perform Lorentz boost
201  G4FragmentVector::iterator i;
202  for (i = theResult->begin(); i != theResult->end(); i++) {
203    G4LorentzVector FourMom = (*i)->GetMomentum();
204    FourMom.boost(theFragment.GetMomentum().boostVector());
205    (*i)->SetMomentum(FourMom);
206#ifdef PRECOMPOUND_TEST
207    (*i)->SetCreatorModel(G4String("G4StatMF"));
208#endif
209  }
210 
211  // garbage collection
212  delete theMicrocanonicalEnsemble;
213  if (theMacrocanonicalEnsemble != 0) delete theMacrocanonicalEnsemble;
214  delete theChannel;
215 
216  return theResult;
217}
218
219
220G4bool G4StatMF::FindTemperatureOfBreakingChannel(const G4Fragment & theFragment,
221                                                  const G4StatMFChannel * aChannel,
222                                                  G4double & Temperature)
223  // This finds temperature of breaking channel.
224{
225  G4double A = theFragment.GetA();
226  G4double Z = theFragment.GetZ();
227  G4double U = theFragment.GetExcitationEnergy();
228 
229  G4double T = std::max(Temperature,0.0012*MeV);
230 
231  G4double Ta = T;
232  G4double Tb = T;
233 
234 
235  G4double TotalEnergy = CalcEnergy(A,Z,aChannel,T);
236 
237  G4double Da = (U - TotalEnergy)/U;
238  G4double Db = 0.0;
239 
240  // bracketing the solution
241  if (Da == 0.0) {
242    Temperature = T;
243    return true;
244  } else if (Da < 0.0) {
245    do {
246      Tb -= 0.5 * std::abs(Tb);
247      T = Tb;
248      if (Tb < 0.001*MeV) return false;
249     
250      TotalEnergy = CalcEnergy(A,Z,aChannel,T);
251     
252      Db = (U - TotalEnergy)/U;
253    } while (Db < 0.0);
254   
255  } else {
256    do {
257      Tb += 0.5 * std::abs(Tb);
258      T = Tb;
259     
260      TotalEnergy = CalcEnergy(A,Z,aChannel,T);
261     
262      Db = (U - TotalEnergy)/U;
263    } while (Db > 0.0); 
264  }
265 
266  G4double eps = 1.0e-14 * std::abs(Tb-Ta);
267  //G4double eps = 1.0e-3 ;
268 
269  // Start the bisection method
270  for (G4int j = 0; j < 1000; j++) {
271    G4double Tc =  (Ta+Tb)/2.0;
272    if (std::abs(Ta-Tc) <= eps) {
273      Temperature = Tc;
274      return true;
275    }
276   
277    T = Tc;
278   
279    TotalEnergy = CalcEnergy(A,Z,aChannel,T);
280   
281    G4double Dc = (U - TotalEnergy)/U; 
282   
283    if (Dc == 0.0) {
284      Temperature  = Tc;
285      return true;
286    }
287   
288    if (Da*Dc < 0.0) {
289      Tb = Tc;
290      Db = Dc;
291    } else {
292      Ta = Tc;
293      Da = Dc;
294    }
295  }
296 
297  Temperature  = (Ta+Tb)/2.0;
298  return false;
299}
300
301
302
303G4double G4StatMF::CalcEnergy(const G4double A, const G4double Z, const G4StatMFChannel * aChannel,
304                              const G4double T)
305{
306    G4double MassExcess0 = G4NucleiProperties::GetMassExcess(static_cast<G4int>(A),static_cast<G4int>(Z));
307       
308    G4double Coulomb = (3./5.)*(elm_coupling*Z*Z)*std::pow(1.0+G4StatMFParameters::GetKappaCoulomb(),1./3.)/
309        (G4StatMFParameters::Getr0()*std::pow(A,1./3.));
310
311    G4double ChannelEnergy = aChannel->GetFragmentsEnergy(T);
312       
313    return -MassExcess0 + Coulomb + ChannelEnergy;
314
315}
316
317
318
Note: See TracBrowser for help on using the repository browser.