Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4StatMFFragment.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: G4StatMFFragment.cc 91834 2015-08-07 07:24:22Z gcosmo $
28 //
29 // Hadronic Process: Nuclear De-excitations
30 // by V. Lara
31 
32 #include "G4StatMFFragment.hh"
33 #include "G4PhysicalConstants.hh"
34 #include "G4HadronicException.hh"
35 #include "G4Pow.hh"
36 
37 // Copy constructor
39 {
40  throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::copy_constructor meant to not be accessable");
41 }
42 
43 // Operators
44 
45 G4StatMFFragment & G4StatMFFragment::
46 operator=(const G4StatMFFragment & )
47 {
48  throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::operator= meant to not be accessable");
49  return *this;
50 }
51 
53 {
54 // throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::operator== meant to not be accessable");
55  return false;
56 }
57 
59 {
60 // throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::operator!= meant to not be accessable");
61  return true;
62 }
63 
65 {
66  G4double res = 0.0;
67  if (theZ >= 1) {
69  }
70  return res;
71 }
72 
74 {
75  if (theA < 1 || theZ < 0 || theZ > theA) {
76  G4cout << "G4StatMFFragment::GetEnergy: A = " << theA
77  << ", Z = " << theZ << G4endl;
78  throw G4HadronicException(__FILE__, __LINE__,
79  "G4StatMFFragment::GetEnergy: Wrong values for A and Z!");
80  }
81  G4double BulkEnergy = G4NucleiProperties::GetMassExcess(theA,theZ);
82 
83  if (theA < 4) return BulkEnergy - GetCoulombEnergy();
84 
85  G4double SurfaceEnergy;
86  if (G4StatMFParameters::DBetaDT(T) == 0.0) SurfaceEnergy = 0.0;
87  else SurfaceEnergy = 2.5*G4Pow::GetInstance()->Z23(theA)*T*T*
91 
92  G4double ExchangeEnergy = theA*T*T/GetInvLevelDensity();
93  if (theA != 4) ExchangeEnergy += SurfaceEnergy;
94  return BulkEnergy + ExchangeEnergy - GetCoulombEnergy();
95 }
96 
98 {
99  G4double res = 0.0;
100  if (theA > 1) {
101  res = G4StatMFParameters::GetEpsilon0()*(1.0+3.0/(theA - 1.0));
102  }
103  return res;
104 }
105 
107 {
108  G4double U = CalcExcitationEnergy(T);
109  G4double M = GetNuclearMass();
110  G4LorentzVector FourMomentum(_momentum,std::sqrt(_momentum.mag2()+(M+U)*(M+U)));
111  G4Fragment * theFragment = new G4Fragment(theA, theZ, FourMomentum);
112  return theFragment;
113 }
114 
115 G4double G4StatMFFragment::CalcExcitationEnergy(const G4double T)
116 {
117  if (theA <= 3) return 0.0;
118 
119  G4double BulkEnergy = theA*T*T/GetInvLevelDensity();
120 
121  // if it is an alpha particle: done
122  if (theA == 4) return BulkEnergy;
123 
124  // Term connected with surface energy
125  G4double SurfaceEnergy = 0.0;
127  if (std::abs(q) > 1.0e-20) {
128  SurfaceEnergy = 2.5*G4Pow::GetInstance()->Z23(theA)
130  }
131  return BulkEnergy + SurfaceEnergy;
132 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4Fragment * GetFragment(const G4double T)
G4double GetCoulombEnergy(void) const
G4double GetNuclearMass(void)
G4bool operator!=(const G4StatMFFragment &right) const
static G4double GetMassExcess(const G4int A, const G4int Z)
G4double GetEnergy(const G4double T) const
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
static G4double DBetaDT(G4double T)
double mag2() const
G4bool operator==(const G4StatMFFragment &right) const
static G4double GetCoulomb()
G4double Z23(G4int Z) const
Definition: G4Pow.hh:154
#define G4endl
Definition: G4ios.hh:61
static G4double GetEpsilon0()
static G4double GetBeta0()
double G4double
Definition: G4Types.hh:76
G4StatMFFragment(G4int anA, G4int aZ)
static G4double Beta(G4double T)
G4double GetInvLevelDensity(void) const
static G4double GetCriticalTemp()