source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/processes/include/G4QCoherentChargeExchange.hh @ 1199

Last change on this file since 1199 was 1199, checked in by garnier, 15 years ago

nvx fichiers dans CVS

File size: 7.3 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// * 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: G4QCoherentChargeExchange.hh,v 1.1 2009/11/17 10:36:54 mkossov Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29//      ---------------- G4QCoherentChargeExchange header ----------------
30//                 by Mikhail Kossov, December 2003.
31//  Header of G4QCoherentChargeExchange class (hadron+A) of the CHIPS Simulation Branch in GEANT4
32// -------------------------------------------------------------------------------
33// This is a unique CHIPS class for the Hadron-Nuclear Elastic Scattering Prosesses
34// -------------------------------------------------------------------------------
35// At present (Jan-06) only proton-to-neutron & neutron-to-proton scattering on nuclei
36// are implemented. The scattering of mesons and nuclei on nuclei are possible...
37// The simulation is based on the CHIPS approximation of total elastic and differential
38// elastic cross sections from E=0 to the highest energyes.
39// -------------------------------------------------------------------------------
40// Short description: This class resolves an ambiguity in the definition of the
41// "inelastic" cross section. As it was shown in Ph.D.Thesis (M.Kosov,ITEP,1979)
42// it is more reasonable to subdivide the total cross-section in the coherent &
43// incoherent parts, but the measuring method for the "inelastic" cross-sections
44// consideres the lack of the projectile within the narrow forward solid angle
45// with the consequent extrapolation of these partial cross-sections, corresponding
46// to the particular solid angle, to the zero solid angle. The low angle region
47// is shadowed by the elastic (coherent) scattering. BUT the coherent charge
48// exchange (e.g. conversion p->n) is included by this procedure as a constant term
49// in the extrapolation, so the "inelastic" cross-section differes from the
50// incoherent cross-section by the value of the coherent charge exchange cross
51// section. Fortunately, this cross-sectoion drops ruther fast with energy increasing.
52// All Geant4 inelastic hadronic models (including CHIPS) simulate the incoherent
53// reactions. So the incoherent (including quasielastic) cross-section must be used
54// instead of the inelastic cross-section. For that the "inelastic" cross-section
55// must be reduced by the value of the coherent charge-exchange cross-section, which
56// is estimated (it must be tuned!) in this CHIPS class. The angular distribution
57// is made (at present) identical to the corresponding coherent-elastic scattering
58// -----------------------------------------------------------------------------------
59
60#ifndef G4QCoherentChargeExchange_hh
61#define G4QCoherentChargeExchange_hh
62
63// GEANT4 Headers
64#include "globals.hh"
65#include "G4ios.hh"
66#include "Randomize.hh"
67#include "G4VDiscreteProcess.hh"
68#include "G4Track.hh"
69#include "G4Step.hh"
70#include "G4ParticleTypes.hh"
71#include "G4VParticleChange.hh"
72#include "G4ParticleDefinition.hh"
73#include "G4DynamicParticle.hh"
74#include "G4ThreeVector.hh"
75#include "G4LorentzVector.hh"
76
77// CHIPS Headers
78#include "G4QuasiFreeRatios.hh"
79#include "G4QElasticCrossSection.hh"
80#include "G4QIsotope.hh"
81#include "G4QCHIPSWorld.hh"
82#include "G4QHadron.hh"
83#include <vector>
84
85class G4QCoherentChargeExchange : public G4VDiscreteProcess
86{
87public:
88
89  // Constructor
90  G4QCoherentChargeExchange(const G4String& processName ="CHIPS_CoherChargeExScattering");
91
92  // Destructor
93  ~G4QCoherentChargeExchange();
94
95  G4bool IsApplicable(const G4ParticleDefinition& particle);
96
97  G4double GetMeanFreePath(const G4Track& aTrack, G4double previousStepSize,
98                           G4ForceCondition* condition);
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
113  G4LorentzVector GetEnegryMomentumConservation();
114
115  G4int GetNumberOfNeutronsInTarget();
116
117private:
118
119  // Hide assignment operator as private
120  G4QCoherentChargeExchange& operator=(const G4QCoherentChargeExchange &right);
121
122  // Copy constructor
123  G4QCoherentChargeExchange(const G4QCoherentChargeExchange&);
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
Note: See TracBrowser for help on using the repository browser.