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

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

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