source: trunk/source/processes/hadronic/cross_sections/src/G4EMDissociationCrossSection.cc@ 1007

Last change on this file since 1007 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 10.2 KB
RevLine 
[819]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 gant4 7.1.
59//
60// %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
61////////////////////////////////////////////////////////////////////////////////
62//
63#include "G4EMDissociationCrossSection.hh"
64#include "G4PhysicsFreeVector.hh"
65#include "G4ParticleTable.hh"
66#include "G4IonTable.hh"
67#include "globals.hh"
68////////////////////////////////////////////////////////////////////////////////
69//
70G4EMDissociationCrossSection::G4EMDissociationCrossSection ()
71{
72//
73//
74// This function makes use of the class which can sample the virtual photon
75// spectrum, G4EMDissociationSpectrum.
76//
77 thePhotonSpectrum = new G4EMDissociationSpectrum();
78//
79//
80// Define other constants.
81//
82 r0 = 1.18 * fermi;
83 J = 36.8 * MeV;
84 Qprime = 17.0 * MeV;
85 epsilon = 0.0768;
86 xd = 0.25;
87}
88////////////////////////////////////////////////////////////////////////////////
89//
90G4EMDissociationCrossSection::~G4EMDissociationCrossSection()
91{
92 delete thePhotonSpectrum;
93}
94///////////////////////////////////////////////////////////////////////////////
95//
96G4bool G4EMDissociationCrossSection::IsZAApplicable
97(const G4DynamicParticle* theDynamicParticle, G4double /*ZZ*/, G4double AA)
98{
99//
100// The condition for the applicability of this class is that the projectile
101// must be an ion and the target must have more than one nucleon. In reality
102// the value of A for either the projectile or target could be much higher,
103// since for cases where both he projectile and target are medium to small
104// Z, the probability of the EMD process is, I think, VERY small.
105//
106 if (G4ParticleTable::GetParticleTable()->GetIonTable()->
107 IsIon(theDynamicParticle->GetDefinition()) && AA > 1.0)
108 return true;
109 else
110 return false;
111}
112
113G4bool G4EMDissociationCrossSection::IsApplicable
114 (const G4DynamicParticle* theDynamicParticle, const G4Element* theElement)
115{
116 return IsZAApplicable(theDynamicParticle, 0., theElement->GetN());
117}
118
119//////////////////////////////////////////////////////////////////////////////
120//
121G4double G4EMDissociationCrossSection::GetCrossSection
122 (const G4DynamicParticle* theDynamicParticle, const G4Element* theElement,
123 G4double temperature)
124{
125 G4int nIso = theElement->GetNumberOfIsotopes();
126 G4double crossSection = 0;
127
128 if (nIso) {
129 G4double sig;
130 G4IsotopeVector* isoVector = theElement->GetIsotopeVector();
131 G4double* abundVector = theElement->GetRelativeAbundanceVector();
132 G4double ZZ;
133 G4double AA;
134
135 for (G4int i = 0; i < nIso; i++) {
136 ZZ = G4double( (*isoVector)[i]->GetZ() );
137 AA = G4double( (*isoVector)[i]->GetN() );
138 sig = GetIsoZACrossSection(theDynamicParticle, ZZ, AA, temperature);
139 crossSection += sig*abundVector[i];
140 }
141
142 } else {
143 crossSection =
144 GetIsoZACrossSection(theDynamicParticle, theElement->GetZ(),
145 theElement->GetN(), temperature);
146 }
147
148 return crossSection;
149}
150
151
152G4double G4EMDissociationCrossSection::GetIsoZACrossSection
153 (const G4DynamicParticle *theDynamicParticle, G4double ZZ, G4double AA,
154 G4double /*temperature*/)
155{
156//
157//
158// Get relevant information about the projectile and target (A, Z) and
159// velocity of the projectile.
160//
161 G4ParticleDefinition *definitionP = theDynamicParticle->GetDefinition();
162 G4double AP = definitionP->GetBaryonNumber();
163 G4double ZP = definitionP->GetPDGCharge();
164 G4double b = theDynamicParticle->Get4Momentum().beta();
165// G4double bsq = b * b;
166
167 G4double AT = AA;
168 G4double ZT = ZZ;
169 G4double bmin = thePhotonSpectrum->GetClosestApproach(AP, ZP, AT, ZT, b);
170//
171//
172// Calculate the cross-section for the projectile and then the target. The
173// information is returned in a G4PhysicsFreeVector, which separates out the
174// cross-sections for the E1 and E2 moments of the virtual photon field, and
175// the energies (GDR and GQR).
176//
177 G4PhysicsFreeVector *theProjectileCrossSections =
178 GetCrossSectionForProjectile (AP, ZP, AT, ZT, b, bmin);
179 G4double crossSection =
180 (*theProjectileCrossSections)[0]+(*theProjectileCrossSections)[1];
181 delete theProjectileCrossSections;
182 G4PhysicsFreeVector *theTargetCrossSections =
183 GetCrossSectionForTarget (AP, ZP, AT, ZT, b, bmin);
184 crossSection +=
185 (*theTargetCrossSections)[0]+(*theTargetCrossSections)[1];
186 delete theTargetCrossSections;
187
188 return crossSection;
189}
190////////////////////////////////////////////////////////////////////////////////
191//
192G4PhysicsFreeVector *
193 G4EMDissociationCrossSection::GetCrossSectionForProjectile (G4double AP,
194 G4double ZP, G4double /* AT */, G4double ZT, G4double b, G4double bmin)
195{
196//
197//
198// Use Wilson et al's approach to calculate the cross-sections due to the E1
199// and E2 moments of the field at the giant dipole and quadrupole resonances
200// respectively, Note that the algorithm is traditionally applied to the
201// EMD break-up of the projectile in the field of the target, as is implemented
202// here.
203//
204// Initialise variables and calculate the energies for the GDR and GQR.
205//
206 G4double AProot3 = std::pow(AP,1.0/3.0);
207 G4double u = 3.0 * J / Qprime / AProot3;
208 G4double R0 = r0 * AProot3;
209 G4double E_GDR = hbarc / std::sqrt(0.7*amu_c2*R0*R0/8.0/J*
210 (1.0 + u - (1.0 + epsilon + 3.0*u)/(1.0 + epsilon + u)*epsilon));
211 G4double E_GQR = 63.0 * MeV / AProot3;
212//
213//
214// Determine the virtual photon spectra at these energies.
215//
216 G4double ZTsq = ZT * ZT;
217 G4double nE1 = ZTsq *
218 thePhotonSpectrum->GetGeneralE1Spectrum(E_GDR, b, bmin);
219 G4double nE2 = ZTsq *
220 thePhotonSpectrum->GetGeneralE2Spectrum(E_GQR, b, bmin);
221//
222//
223// Now calculate the cross-section of the projectile for interaction with the
224// E1 and E2 fields.
225//
226 G4double sE1 = 60.0 * millibarn * MeV * (AP-ZP)*ZP/AP;
227 G4double sE2 = 0.22 * microbarn / MeV * ZP * AProot3 * AProot3;
228 if (AP > 100.0) sE2 *= 0.9;
229 else if (AP > 40.0) sE2 *= 0.6;
230 else sE2 *= 0.3;
231//
232//
233// ... and multiply with the intensity of the virtual photon spectra to get
234// the probability of interaction.
235//
236 G4PhysicsFreeVector *theCrossSectionVector = new G4PhysicsFreeVector(2);
237 theCrossSectionVector->PutValue(0, E_GDR, sE1*nE1);
238 theCrossSectionVector->PutValue(1, E_GQR, sE2*nE2*E_GQR*E_GQR);
239
240 return theCrossSectionVector;
241}
242////////////////////////////////////////////////////////////////////////////////
243//
244G4PhysicsFreeVector *
245 G4EMDissociationCrossSection::GetCrossSectionForTarget (G4double AP,
246 G4double ZP, G4double AT, G4double ZT, G4double b, G4double bmin)
247{
248//
249//
250// This is a cheaky little member function to calculate the probability of
251// EMD for the target in the field of the projectile ... just by reversing the
252// A and Z's for the participants.
253//
254 return GetCrossSectionForProjectile (AT, ZT, AP, ZP, b, bmin);
255}
256////////////////////////////////////////////////////////////////////////////////
257//
258G4double
259 G4EMDissociationCrossSection::GetWilsonProbabilityForProtonDissociation
260 (G4double A, G4double Z)
261{
262//
263//
264// This is a simple algorithm to choose whether a proton or neutron is ejected
265// from the nucleus in the EMD interaction.
266//
267 G4double p = 0.0;
268 if (Z < 6.0)
269 p = 0.5;
270 else if (Z < 8.0)
271 p = 0.6;
272 else if (Z < 14.0)
273 p = 0.7;
274 else
275 {
276 G4double p1 = (G4double) Z / (G4double) A;
277 G4double p2 = 1.95*std::exp(-0.075*Z);
278 if (p1 < p2) p = p1;
279 else p = p2;
280 }
281
282 return p;
283}
284////////////////////////////////////////////////////////////////////////////////
285//
Note: See TracBrowser for help on using the repository browser.