Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4EMDissociationCrossSection.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 // * *
21 // * Parts of this code which have been developed by QinetiQ Ltd *
22 // * under contract to the European Space Agency (ESA) are the *
23 // * intellectual property of ESA. Rights to use, copy, modify and *
24 // * redistribute this software for general public use are granted *
25 // * in compliance with any licensing, distribution and development *
26 // * policy adopted by the Geant4 Collaboration. This code has been *
27 // * written by QinetiQ Ltd for the European Space Agency, under ESA *
28 // * contract 17191/03/NL/LvH (Aurora Programme). *
29 // * *
30 // * By using, copying, modifying or distributing the software (or *
31 // * any work based on the software) you agree to acknowledge its *
32 // * use in resulting scientific publications, and indicate your *
33 // * acceptance of all terms of the Geant4 Software license. *
34 // ********************************************************************
35 //
36 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
37 //
38 // MODULE: G4EMDissociationCrossSection.cc
39 //
40 // Version: B.1
41 // Date: 15/04/04
42 // Author: P R Truscott
43 // Organisation: QinetiQ Ltd, UK
44 // Customer: ESA/ESTEC, NOORDWIJK
45 // Contract: 17191/03/NL/LvH
46 //
47 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
48 //
49 // CHANGE HISTORY
50 // --------------
51 //
52 // 17 October 2003, P R Truscott, QinetiQ Ltd, UK
53 // Created.
54 //
55 // 15 March 2004, P R Truscott, QinetiQ Ltd, UK
56 // Beta release
57 //
58 // 30 May 2005, J.P. Wellisch removed a compilation warning on gcc 3.4 for
59 // geant4 7.1.
60 // 09 November 2010, V.Ivanchenko make class applicable for Hydrogen but
61 // set cross section for Hydrogen to zero
62 //
63 // 17 August 2011, V.Ivanchenko, provide migration to new design of cross
64 // sections considering this cross section as element-wise
65 //
66 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
68 //
70 #include "G4PhysicalConstants.hh"
71 #include "G4SystemOfUnits.hh"
72 #include "G4ParticleTable.hh"
73 #include "G4IonTable.hh"
74 #include "G4HadTmpUtil.hh"
75 #include "globals.hh"
76 #include "G4NistManager.hh"
77 
78 
80  : G4VCrossSectionDataSet("Electromagnetic dissociation")
81 {
82  // This function makes use of the class which can sample the virtual photon
83  // spectrum, G4EMDissociationSpectrum.
84 
86 
87  // Define other constants.
88 
89  r0 = 1.18 * fermi;
90  J = 36.8 * MeV;
91  Qprime = 17.0 * MeV;
92  epsilon = 0.0768;
93  xd = 0.25;
94 }
95 
97 
99 {
100  delete thePhotonSpectrum;
101 }
103 //
104 G4bool
106  G4int /*ZZ*/, const G4Material*)
107 {
108 //
109 // The condition for the applicability of this class is that the projectile
110 // must be an ion and the target must have more than one nucleon. In reality
111 // the value of A for either the projectile or target could be much higher,
112 // since for cases where both he projectile and target are medium to small
113 // Z, the probability of the EMD process is, I think, VERY small.
114 //
115  if (G4ParticleTable::GetParticleTable()->GetIonTable()->IsIon(part->GetDefinition())) {
116  return true;
117  } else {
118  return false;
119  }
120 }
121 
123 //
125  (const G4DynamicParticle* theDynamicParticle, G4int Z,
126  const G4Material*)
127 {
128  // VI protection for Hydrogen
129  if(1 >= Z) { return 0.0; }
130 
131  //
132  // Get relevant information about the projectile and target (A, Z) and
133  // velocity of the projectile.
134  //
135  const G4ParticleDefinition *definitionP = theDynamicParticle->GetDefinition();
136  G4double AP = definitionP->GetBaryonNumber();
137  G4double ZP = definitionP->GetPDGCharge();
138  G4double b = theDynamicParticle->Get4Momentum().beta();
139 
141  G4double ZT = (G4double)Z;
142  G4double bmin = thePhotonSpectrum->GetClosestApproach(AP, ZP, AT, ZT, b);
143  //
144  //
145  // Calculate the cross-section for the projectile and then the target. The
146  // information is returned in a G4PhysicsFreeVector, which separates out the
147  // cross-sections for the E1 and E2 moments of the virtual photon field, and
148  // the energies (GDR and GQR).
149  //
150  G4PhysicsFreeVector *theProjectileCrossSections =
151  GetCrossSectionForProjectile (AP, ZP, AT, ZT, b, bmin);
152  G4double crossSection =
153  (*theProjectileCrossSections)[0]+(*theProjectileCrossSections)[1];
154  delete theProjectileCrossSections;
155  G4PhysicsFreeVector *theTargetCrossSections =
156  GetCrossSectionForTarget (AP, ZP, AT, ZT, b, bmin);
157  crossSection +=
158  (*theTargetCrossSections)[0]+(*theTargetCrossSections)[1];
159  delete theTargetCrossSections;
160  return crossSection;
161 }
163 //
166  G4double ZP, G4double /* AT */, G4double ZT, G4double b, G4double bmin)
167 {
168 //
169 //
170 // Use Wilson et al's approach to calculate the cross-sections due to the E1
171 // and E2 moments of the field at the giant dipole and quadrupole resonances
172 // respectively, Note that the algorithm is traditionally applied to the
173 // EMD break-up of the projectile in the field of the target, as is implemented
174 // here.
175 //
176 // Initialise variables and calculate the energies for the GDR and GQR.
177 //
178  G4double AProot3 = G4Pow::GetInstance()->powA(AP,1.0/3.0);
179  G4double u = 3.0 * J / Qprime / AProot3;
180  G4double R0 = r0 * AProot3;
181  G4double E_GDR = hbarc / std::sqrt(0.7*amu_c2*R0*R0/8.0/J*
182  (1.0 + u - (1.0 + epsilon + 3.0*u)/(1.0 + epsilon + u)*epsilon));
183  G4double E_GQR = 63.0 * MeV / AProot3;
184 //
185 //
186 // Determine the virtual photon spectra at these energies.
187 //
188  G4double ZTsq = ZT * ZT;
189  G4double nE1 = ZTsq *
190  thePhotonSpectrum->GetGeneralE1Spectrum(E_GDR, b, bmin);
191  G4double nE2 = ZTsq *
192  thePhotonSpectrum->GetGeneralE2Spectrum(E_GQR, b, bmin);
193 //
194 //
195 // Now calculate the cross-section of the projectile for interaction with the
196 // E1 and E2 fields.
197 //
198  G4double sE1 = 60.0 * millibarn * MeV * (AP-ZP)*ZP/AP;
199  G4double sE2 = 0.22 * microbarn / MeV * ZP * AProot3 * AProot3;
200  if (AP > 100.0) sE2 *= 0.9;
201  else if (AP > 40.0) sE2 *= 0.6;
202  else sE2 *= 0.3;
203 //
204 //
205 // ... and multiply with the intensity of the virtual photon spectra to get
206 // the probability of interaction.
207 //
208  G4PhysicsFreeVector *theCrossSectionVector = new G4PhysicsFreeVector(2);
209  theCrossSectionVector->PutValue(0, E_GDR, sE1*nE1);
210  theCrossSectionVector->PutValue(1, E_GQR, sE2*nE2*E_GQR*E_GQR);
211 
212  return theCrossSectionVector;
213 }
214 
216 //
219  G4double ZP, G4double AT, G4double ZT, G4double b, G4double bmin)
220 {
221 //
222 // This is a cheaky little member function to calculate the probability of
223 // EMD for the target in the field of the projectile ... just by reversing the
224 // A and Z's for the participants.
225 //
226  return GetCrossSectionForProjectile (AT, ZT, AP, ZP, b, bmin);
227 }
228 
230 //
231 G4double
233  G4double Z)
234 {
235 //
236 // This is a simple algorithm to choose whether a proton or neutron is ejected
237 // from the nucleus in the EMD interaction.
238 //
239  G4double p = 0.0;
240  if (Z < 6.0)
241  p = 0.5;
242  else if (Z < 8.0)
243  p = 0.6;
244  else if (Z < 14.0)
245  p = 0.7;
246  else
247  {
248  G4double p1 = (G4double) Z / (G4double) A;
249  G4double p2 = 1.95*G4Exp(-0.075*Z);
250  if (p1 < p2) p = p1;
251  else p = p2;
252  }
253 
254  return p;
255 }
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4ParticleTable * GetParticleTable()
static constexpr double MeV
Definition: G4SIunits.hh:214
const char * p
Definition: xmltok.h:285
static constexpr double hbarc
G4double GetPDGCharge() const
G4PhysicsFreeVector * GetCrossSectionForProjectile(G4double, G4double, G4double, G4double, G4double, G4double)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *)
G4double GetAtomicMassAmu(const G4String &symb) const
Float_t Z
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
static constexpr double amu_c2
TString part[npart]
Definition: Style.C:32
static constexpr double fermi
Definition: G4SIunits.hh:103
static constexpr double millibarn
Definition: G4SIunits.hh:106
static G4Pow * GetInstance()
Definition: G4Pow.cc:57
G4ParticleDefinition * GetDefinition() const
G4PhysicsFreeVector * GetCrossSectionForTarget(G4double, G4double, G4double, G4double, G4double, G4double)
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:242
G4double GetWilsonProbabilityForProtonDissociation(G4double, G4double)
G4LorentzVector Get4Momentum() const
double A(double temperature)
G4EMDissociationSpectrum * thePhotonSpectrum
static constexpr double microbarn
Definition: G4SIunits.hh:107
int G4int
Definition: G4Types.hh:78
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *)
G4double GetGeneralE1Spectrum(G4double, G4double, G4double)
G4double GetGeneralE2Spectrum(G4double, G4double, G4double)
void PutValue(size_t index, G4double energy, G4double dataValue)
static G4NistManager * Instance()