source: trunk/source/processes/hadronic/models/de_excitation/multifragmentation/src/G4StatMFFragment.cc@ 1199

Last change on this file since 1199 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 5.3 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: G4StatMFFragment.cc,v 1.7 2008/07/25 11:20:47 vnivanch Exp $
28// GEANT4 tag $Name: geant4-09-03-cand-01 $
29//
30// Hadronic Process: Nuclear De-excitations
31// by V. Lara
32
33#include "G4StatMFFragment.hh"
34#include "G4HadronicException.hh"
35
36
37// Copy constructor
38G4StatMFFragment::G4StatMFFragment(const G4StatMFFragment & )
39{
40 throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::copy_constructor meant to not be accessable");
41}
42
43// Operators
44
45G4StatMFFragment & G4StatMFFragment::
46operator=(const G4StatMFFragment & )
47{
48 throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::operator= meant to not be accessable");
49 return *this;
50}
51
52
53G4bool G4StatMFFragment::operator==(const G4StatMFFragment & ) const
54{
55// throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::operator== meant to not be accessable");
56 return false;
57}
58
59
60G4bool G4StatMFFragment::operator!=(const G4StatMFFragment & ) const
61{
62// throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::operator!= meant to not be accessable");
63 return true;
64}
65
66
67
68G4double G4StatMFFragment::GetCoulombEnergy(void) const
69{
70 if (theZ <= 0.1) return 0.0;
71 G4double Coulomb = (3./5.)*(elm_coupling*theZ*theZ)*
72 std::pow(1.0+G4StatMFParameters::GetKappaCoulomb(),1./3.)/
73 (G4StatMFParameters::Getr0()*std::pow(theA,1./3.));
74
75 return Coulomb;
76}
77
78
79G4double G4StatMFFragment::GetEnergy(const G4double T) const
80{
81 if (theA < 1 || theZ < 0 || theZ > theA) {
82 G4cerr << "G4StatMFFragment::GetEnergy: A = " << theA
83 << ", Z = " << theZ << G4endl;
84 throw G4HadronicException(__FILE__, __LINE__,
85 "G4StatMFFragment::GetEnergy: Wrong values for A and Z!");
86 }
87 G4double BulkEnergy = G4NucleiProperties::GetMassExcess(static_cast<G4int>(theA),
88 static_cast<G4int>(theZ));
89
90 if (theA < 4) return BulkEnergy - GetCoulombEnergy();
91
92 G4double SurfaceEnergy;
93 if (G4StatMFParameters::DBetaDT(T) == 0.0) SurfaceEnergy = 0.0;
94 else SurfaceEnergy = (5./2.)*std::pow(theA,2.0/3.0)*T*T*
95 G4StatMFParameters::GetBeta0()/
96 (G4StatMFParameters::GetCriticalTemp()*
97 G4StatMFParameters::GetCriticalTemp());
98
99
100 G4double ExchangeEnergy = theA*T*T/GetInvLevelDensity();
101 if (theA != 4) ExchangeEnergy += SurfaceEnergy;
102
103 return BulkEnergy + ExchangeEnergy - GetCoulombEnergy();
104
105}
106
107
108G4double G4StatMFFragment::GetInvLevelDensity(void) const
109{
110 // Calculate Inverse Density Level
111 // Epsilon0*(1 + 3 /(Af - 1))
112 if (theA == 1) return 0.0;
113 else return
114 G4StatMFParameters::GetEpsilon0()*(1.0+3.0/(theA - 1.0));
115}
116
117
118
119G4Fragment * G4StatMFFragment::GetFragment(const G4double T)
120{
121 G4double U = CalcExcitationEnergy(T);
122
123 G4double M = GetNuclearMass();
124
125 G4LorentzVector FourMomentum(_momentum,std::sqrt(_momentum.mag2()+(M+U)*(M+U)));
126
127 G4Fragment * theFragment = new G4Fragment(static_cast<G4int>(theA),static_cast<G4int>(theZ),FourMomentum);
128
129 return theFragment;
130}
131
132
133G4double G4StatMFFragment::CalcExcitationEnergy(const G4double T)
134{
135 if (theA <= 3) return 0.0;
136
137 G4double BulkEnergy = theA*T*T/GetInvLevelDensity();
138
139 // if it is an alpha particle: done
140 if (theA == 4) return BulkEnergy;
141
142 // Term connected with surface energy
143 G4double SurfaceEnergy = 0.0;
144 if (std::abs(G4StatMFParameters::DBetaDT(T)) > 1.0e-20)
145// SurfaceEnergy = (5./2.)*std::pow(theA,2.0/3.0)*T*T*G4StatMFParameters::GetBeta0()/
146// (G4StatMFParameters::GetCriticalTemp()*G4StatMFParameters::GetCriticalTemp());
147 SurfaceEnergy = (5./2.)*std::pow(theA,2.0/3.0)*(G4StatMFParameters::Beta(T) -
148 T*G4StatMFParameters::DBetaDT(T) - G4StatMFParameters::GetBeta0());
149
150 return BulkEnergy + SurfaceEnergy;
151}
152
153
Note: See TracBrowser for help on using the repository browser.