Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4QCoherentChargeExchange.hh
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// * 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$
27//
28// ---------------- G4QCoherentChargeExchange header ----------------
29// by Mikhail Kossov, December 2003.
30// Header of G4QCoherentChargeExchange class (hA) of the CHIPS Simulation Branch
31// -------------------------------------------------------------------------------
32// This is a unique CHIPS class for the Hadron-Nuclear Elastic Scattering Prosesses
33// -------------------------------------------------------------------------------
34// At present (Jan-06) only proton-to-neutron & neutron-to-proton scattering on nuclei
35// are implemented. The scattering of mesons and nuclei on nuclei are possible...
36// The simulation is based on the CHIPS approximation of total elastic and differential
37// elastic cross sections from E=0 to the highest energyes.
38// -------------------------------------------------------------------------------
39// Short description: This class resolves an ambiguity in the definition of the
40// "inelastic" cross section. As it was shown in Ph.D.Thesis (M.Kosov,ITEP,1979)
41// it is more reasonable to subdivide the total cross-section in the coherent &
42// incoherent parts, but the measuring method for the "inelastic" cross-sections
43// consideres the lack of the projectile within the narrow forward solid angle
44// with the consequent extrapolation of these partial cross-sections, corresponding
45// to the particular solid angle, to the zero solid angle. The low angle region
46// is shadowed by the elastic (coherent) scattering. BUT the coherent charge
47// exchange (e.g. conversion p->n) is included by this procedure as a constant term
48// in the extrapolation, so the "inelastic" cross-section differes from the
49// incoherent cross-section by the value of the coherent charge exchange cross
50// section. Fortunately, this cross-sectoion drops ruther fast with energy increasing.
51// All Geant4 inelastic hadronic models (including CHIPS) simulate the incoherent
52// reactions. So the incoherent (including quasielastic) cross-section must be used
53// instead of the inelastic cross-section. For that the "inelastic" cross-section
54// must be reduced by the value of the coherent charge-exchange cross-section, which
55// is estimated (it must be tuned!) in this CHIPS class. The angular distribution
56// is made (at present) identical to the corresponding coherent-elastic scattering
57// -----------------------------------------------------------------------------------
58
59#ifndef G4QCoherentChargeExchange_hh
60#define G4QCoherentChargeExchange_hh
61
62// GEANT4 Headers
63#include "globals.hh"
64#include "G4ios.hh"
65#include "Randomize.hh"
66#include "G4VDiscreteProcess.hh"
67#include "G4Track.hh"
68#include "G4Step.hh"
69#include "G4ParticleTypes.hh"
70#include "G4VParticleChange.hh"
72#include "G4DynamicParticle.hh"
73#include "G4ThreeVector.hh"
74#include "G4LorentzVector.hh"
75
76// CHIPS Headers
77#include "G4QuasiFreeRatios.hh"
80#include "G4QIsotope.hh"
81#include "G4QCHIPSWorld.hh"
82#include "G4QHadron.hh"
83#include <vector>
84
86{
87public:
88
89 // Constructor
90 G4QCoherentChargeExchange(const G4String& processName ="CHIPS_CoherChargeExScattering");
91
92 // Destructor
94
96
97 G4double GetMeanFreePath(const G4Track& aTrack, G4double previousStepSize,
99 // It returns the MeanFreePath of the process for the current track :
100 // (energy, material)
101 // The previousStepSize and G4ForceCondition* are not used.
102 // This function overloads a virtual function of the base class.
103 // It is invoked by the ProcessManager of the Particle.
104
105
106 G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep);
107 // It computes the final state of the process (at end of step),
108 // returned as a ParticleChange object.
109 // This function overloads a virtual function of the base class.
110 // It is invoked by the ProcessManager of the Particle.
111
112
114
116
117private:
118
119 // Hide assignment operator as private
121
122 // Copy constructor
124
125 // Calculate XS/t: oxs=true - only CS; xst=true - calculate XS, xst=false(oxs=f/t) - t/tm
126 G4double CalculateXSt(G4bool oxs, G4bool xst, G4double p, G4int Z, G4int N, G4int pPDG);
127
128 // BODY
129 // Static Parameters --------------------------------------------------------------------
130 static G4int nPartCWorld; // The#of particles for hadronization (limit of A of fragm.)
131 //--------------------------------- End of static parameters ---------------------------
132 // Working parameters
133 G4VQCrossSection* theCS;
134 G4LorentzVector EnMomConservation; // Residual of Energy/Momentum Cons.
135 G4int nOfNeutrons; // #of neutrons in the target nucleus
136
137 // Modifires for the reaction
138 G4double Time; // Time shift of the capture reaction
139 G4double EnergyDeposition; // Energy deposited in the reaction
140 static std::vector <G4int> ElementZ; // Z of the element(i) in theLastCalc
141 static std::vector <G4double> ElProbInMat; // SumProbabilityElements in Material
142 static std::vector <std::vector<G4int>*> ElIsoN; // N of isotope(j) of Element(i)
143 static std::vector <std::vector<G4double>*> IsoProbInEl;// SumProbabIsotopes in Element i
144};
145#endif
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
G4bool IsApplicable(const G4ParticleDefinition &particle)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4LorentzVector GetEnegryMomentumConservation()
Definition: G4Step.hh:78