Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4EvaporationChannel.cc
이 파일의 문서화 페이지로 가기
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 // $Id: G4EvaporationChannel.cc 105799 2017-08-21 07:35:55Z gcosmo $
27 //
28 //J.M. Quesada (August2008). Based on:
29 //
30 // Hadronic Process: Nuclear De-excitations
31 // by V. Lara (Oct 1998)
32 //
33 // Modified:
34 // 03-09-2008 J.M. Quesada for external choice of inverse cross section option
35 // 06-09-2008 J.M. Quesada Also external choices have been added for superimposed
36 // Coulomb barrier (if useSICB is set true, by default is false)
37 // 17-11-2010 V.Ivanchenko in constructor replace G4VEmissionProbability by
38 // G4EvaporationProbability and do not new and delete probability
39 // object at each call; use G4Pow
40 
41 #include "G4EvaporationChannel.hh"
42 #include "G4PairingCorrection.hh"
43 #include "G4NuclearLevelData.hh"
44 #include "G4NucleiProperties.hh"
45 #include "G4Pow.hh"
46 #include "G4Log.hh"
47 #include "G4Exp.hh"
48 #include "G4PhysicalConstants.hh"
49 #include "G4SystemOfUnits.hh"
50 #include "Randomize.hh"
51 #include "G4RandomDirection.hh"
52 #include "G4Alpha.hh"
53 
55  const G4String & aName,
57  G4VCoulombBarrier* barrier):
58  G4VEvaporationChannel(aName),
59  theA(anA),
60  theZ(aZ),
61  theProbability(aprob),
62  theCoulombBarrier(barrier)
63 {
64  ResA = ResZ = 0;
68 }
69 
71 {}
72 
74 {
77 }
78 
80 {
81  G4int FragA = fragment->GetA_asInt();
82  G4int FragZ = fragment->GetZ_asInt();
83  ResA = FragA - theA;
84  ResZ = FragZ - theZ;
85 
86  G4double FragmentMass = fragment->GetGroundStateMass();
87  G4double ExEnergy = fragment->GetExcitationEnergy();
88  Mass = FragmentMass + ExEnergy;
89  //G4cout << "G4EvaporationChannel::Initialize Z= " << theZ << " A= " << theA
90  // << " FragZ= " << FragZ << " FragA= " << FragA << G4endl;
91  EmissionProbability = 0.0;
92 
93  // Only channels which are physically allowed are taken into account
94  if (ResA >= ResZ && ResZ > 0 && ResA >= theA) {
95 
96  //Effective excitation energy
98 
99  CoulombBarrier = (0 == theZ) ? 0.0 :
101 
102  G4double delta0 =
104  G4double delta1 =
106  ResMass += delta1;
107  /*
108  G4cout << "ExEnergy= " << ExEnergy << " Ec= " << CoulombBarrier
109  << " delta0= " << delta0 << " delta1= " << delta1
110  << " Free= " << Mass - ResMass - EvapMass
111  << G4endl;
112  */
113  // for OPTxs >0 penetration under the barrier is taken into account
114  // G4double elim = (0 == OPTxs) ? CoulombBarrier : CoulombBarrier*0.5;
115  static const G4double dCB = 3.5*CLHEP::MeV;
116  G4double elim = (0 == OPTxs) ? CoulombBarrier : CoulombBarrier - dCB*theZ;
117  if(ExEnergy >= delta0 && Mass >= ResMass + EvapMass + elim) {
118  G4double xm2 = (Mass - EvapMass)*(Mass - EvapMass);
119  G4double xm = Mass - EvapMass - elim;
120  MinKinEnergy = (0.0 >= elim) ? 0.0 : std::max(0.5*(xm2 - xm*xm)/Mass, 0.0);
121  MaxKinEnergy = std::max(0.5*(xm2 - ResMass*ResMass)/Mass, 0.0);
122  //G4cout << "Emin= " << MinKinEnergy << " Emax= " << MaxKinEnergy
123  // << " xm= " << xm << G4endl;
125  TotalProbability(*fragment, MinKinEnergy, MaxKinEnergy, CoulombBarrier);
126  }
127  }
128  //G4cout << "G4EvaporationChannel:: probability= "
129  // << EmissionProbability << G4endl;
130  return EmissionProbability;
131 }
132 
134 {
135  G4Fragment* evFragment = nullptr;
136  G4double ekin = 0.0;
137  if(ResA <= 4 &&
138  ((ResA == 4 && ResZ == 2) || (ResA == 3 && ResZ == 2) ||
139  (ResA == 3 && ResZ == 1) || (ResA == 2 && ResZ == 1) ||
140  (ResA == 1 && ResZ == 1) || (ResA == 1 && ResZ == 0) )) {
142  ekin = 0.5*(Mass*Mass - mres*mres + EvapMass*EvapMass)/Mass - EvapMass;
143  } else {
146  }
147  G4LorentzVector lv0 = theNucleus->GetMomentum();
148  G4LorentzVector lv(std::sqrt(ekin*(ekin + 2.0*EvapMass))*G4RandomDirection(),
149  ekin + EvapMass);
150  lv.boost(lv0.boostVector());
151 
152  evFragment = new G4Fragment(theA, theZ, lv);
153  lv0 -= lv;
154  theNucleus->SetZandA_asInt(ResZ, ResA);
155  theNucleus->SetMomentum(lv0);
156 
157  return evFragment;
158 }
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:300
G4EvaporationChannel(G4int A, G4int Z, const G4String &aName, G4EvaporationProbability *, G4VCoulombBarrier *)
void SetZandA_asInt(G4int Znew, G4int Anew)
Definition: G4Fragment.hh:269
static G4NuclearLevelData * GetInstance()
virtual G4Fragment * EmittedFragment(G4Fragment *theNucleus)
G4int GetA_asInt() const
Definition: G4Fragment.hh:259
G4ThreeVector G4RandomDirection()
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:281
G4PairingCorrection * pairingCorrection
double G4double
Definition: G4Types.hh:76
static constexpr double MeV
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:305
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetZ_asInt() const
Definition: G4Fragment.hh:264
G4EvaporationProbability * theProbability
virtual G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U) const =0
int G4int
Definition: G4Types.hh:78
G4double GetPairingCorrection(G4int A, G4int Z) const
Hep3Vector boostVector() const
G4VCoulombBarrier * theCoulombBarrier
G4PairingCorrection * GetPairingCorrection()
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:276
G4double SampleKineticEnergy(G4double minKineticEnergy, G4double maxKineticEnergy, G4double CoulombBarrier=0.0)
HepLorentzVector & boost(double, double, double)
virtual G4double GetEmissionProbability(G4Fragment *fragment)