Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EMDissociationCrossSection.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// * *
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// %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
67//////////////////////////////////////////////////////////////////////////////
68//
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
85 thePhotonSpectrum = new G4EMDissociationSpectrum();
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
96//////////////////////////////////////////////////////////////////////////////
97
99{
100 delete thePhotonSpectrum;
101}
102/////////////////////////////////////////////////////////////////////////////
103//
104G4bool
106 G4int /*ZZ*/, const G4Material*)
107{
108 return true;
109}
110
111//////////////////////////////////////////////////////////////////////////////
112//
114 (const G4DynamicParticle* theDynamicParticle, G4int Z, const G4Material*)
115{
116 // VI protection for Hydrogen
117 if(1 >= Z) { return 0.0; }
118
119 // Zero cross-section for particles with kinetic energy less than 2 MeV to prevent
120 // possible abort signal from bad arithmetic in GetCrossSectionForProjectile
121 if ( theDynamicParticle->GetKineticEnergy() < 2.0*CLHEP::MeV ) { return 0.0; }
122
123 //
124 // Get relevant information about the projectile and target (A, Z) and
125 // velocity of the projectile.
126 //
127 const G4ParticleDefinition *definitionP = theDynamicParticle->GetDefinition();
128 G4double AP = definitionP->GetBaryonNumber();
129 G4double ZP = definitionP->GetPDGCharge();
130 G4double b = theDynamicParticle->GetBeta();
131 if (b <= 0.0) { return 0.0; }
132
134 G4double ZT = (G4double)Z;
135 G4double bmin = thePhotonSpectrum->GetClosestApproach(AP, ZP, AT, ZT, b);
136 //
137 //
138 // Calculate the cross-section for the projectile and then the target. The
139 // information is returned in a G4PhysicsFreeVector, which separates out the
140 // cross-sections for the E1 and E2 moments of the virtual photon field, and
141 // the energies (GDR and GQR).
142 //
143 G4PhysicsFreeVector *theProjectileCrossSections =
144 GetCrossSectionForProjectile (AP, ZP, AT, ZT, b, bmin);
145 G4double crossSection =
146 (*theProjectileCrossSections)[0]+(*theProjectileCrossSections)[1];
147 delete theProjectileCrossSections;
148 G4PhysicsFreeVector *theTargetCrossSections =
149 GetCrossSectionForTarget (AP, ZP, AT, ZT, b, bmin);
150 crossSection +=
151 (*theTargetCrossSections)[0]+(*theTargetCrossSections)[1];
152 delete theTargetCrossSections;
153 return crossSection;
154}
155////////////////////////////////////////////////////////////////////////////////
156//
159 G4double ZP, G4double /* AT */, G4double ZT, G4double b, G4double bmin)
160{
161//
162//
163// Use Wilson et al's approach to calculate the cross-sections due to the E1
164// and E2 moments of the field at the giant dipole and quadrupole resonances
165// respectively, Note that the algorithm is traditionally applied to the
166// EMD break-up of the projectile in the field of the target, as is implemented
167// here.
168//
169// Initialise variables and calculate the energies for the GDR and GQR.
170//
171 G4double AProot3 = G4Pow::GetInstance()->A13(AP);
172 G4double u = 3.0 * J / Qprime / AProot3;
173 G4double R0 = r0 * AProot3;
174 G4double E_GDR = hbarc / std::sqrt(0.7*amu_c2*R0*R0/8.0/J*
175 (1.0 + u - (1.0 + epsilon + 3.0*u)/(1.0 + epsilon + u)*epsilon));
176 G4double E_GQR = 63.0 * MeV / AProot3;
177//
178//
179// Determine the virtual photon spectra at these energies.
180//
181 G4double ZTsq = ZT * ZT;
182 G4double nE1 = ZTsq *
183 thePhotonSpectrum->GetGeneralE1Spectrum(E_GDR, b, bmin);
184 G4double nE2 = ZTsq *
185 thePhotonSpectrum->GetGeneralE2Spectrum(E_GQR, b, bmin);
186//
187//
188// Now calculate the cross-section of the projectile for interaction with the
189// E1 and E2 fields.
190//
191 G4double sE1 = 60.0 * millibarn * MeV * (AP-ZP)*ZP/AP;
192 G4double sE2 = 0.22 * microbarn / MeV * ZP * AProot3 * AProot3;
193 if (AP > 100.0) sE2 *= 0.9;
194 else if (AP > 40.0) sE2 *= 0.6;
195 else sE2 *= 0.3;
196//
197//
198// ... and multiply with the intensity of the virtual photon spectra to get
199// the probability of interaction.
200//
201 G4PhysicsFreeVector *theCrossSectionVector = new G4PhysicsFreeVector(2);
202 theCrossSectionVector->PutValue(0, E_GDR, sE1*nE1);
203 theCrossSectionVector->PutValue(1, E_GQR, sE2*nE2*E_GQR*E_GQR);
204
205 return theCrossSectionVector;
206}
207
208////////////////////////////////////////////////////////////////////////////////
209//
212 G4double ZP, G4double AT, G4double ZT, G4double b, G4double bmin)
213{
214//
215// This is a cheaky little member function to calculate the probability of
216// EMD for the target in the field of the projectile ... just by reversing the
217// A and Z's for the participants.
218//
219 return GetCrossSectionForProjectile (AT, ZT, AP, ZP, b, bmin);
220}
221
222////////////////////////////////////////////////////////////////////////////////
223//
226 G4double Z)
227{
228//
229// This is a simple algorithm to choose whether a proton or neutron is ejected
230// from the nucleus in the EMD interaction.
231//
232 G4double p = 0.0;
233 if (Z < 2.0)
234 p = 0.0; // To avoid to remove one proton from hydrogen isotopes
235 else if (Z < 6.0)
236 p = 0.5;
237 else if (Z < 8.0)
238 p = 0.6;
239 else if (Z < 14.0)
240 p = 0.7;
241 else
242 {
243 G4double p1 = (G4double) Z / (G4double) A;
244 G4double p2 = 1.95*G4Exp(-0.075*Z);
245 if (p1 < p2) p = p1;
246 else p = p2;
247 }
248
249 return p;
250}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:132
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetBeta() const
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *) override
G4PhysicsFreeVector * GetCrossSectionForProjectile(G4double, G4double, G4double, G4double, G4double, G4double)
G4double GetWilsonProbabilityForProtonDissociation(G4double, G4double)
G4PhysicsFreeVector * GetCrossSectionForTarget(G4double, G4double, G4double, G4double, G4double, G4double)
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *) override
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
void PutValue(const std::size_t index, const G4double e, const G4double value)
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double A13(G4double A) const
Definition G4Pow.cc:116
G4VCrossSectionDataSet(const G4String &nam="")