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

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

nvx fichiers dans CVS

File size: 5.8 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: G4QElastic.hh,v 1.1 2009/11/17 10:36:54 mkossov Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29//      ---------------- G4QElastic header ----------------
30//                 by Mikhail Kossov, December 2003.
31//  Header of G4QElastic 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-proton scattering is implemented The interaction with
36// nuclei are planned only. The scattering of 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: At present this is a process for nucleon-nucleus
41// elastic scattering. Mesons and hyperons exist only for the Hydrogen
42// target (see G4QuasiFreeRatios).
43// ---------------------------------------------------------------
44
45#ifndef G4QElastic_hh
46#define G4QElastic_hh
47
48// GEANT4 Headers
49#include "globals.hh"
50#include "G4ios.hh"
51#include "Randomize.hh"
52#include "G4VDiscreteProcess.hh"
53#include "G4Track.hh"
54#include "G4Step.hh"
55#include "G4ParticleTypes.hh"
56#include "G4VParticleChange.hh"
57#include "G4ParticleDefinition.hh"
58#include "G4DynamicParticle.hh"
59#include "G4ThreeVector.hh"
60#include "G4LorentzVector.hh"
61
62// CHIPS Headers
63#include "G4QElasticCrossSection.hh"
64#include "G4QIsotope.hh"
65#include "G4QCHIPSWorld.hh"
66#include "G4QHadron.hh"
67#include <vector>
68
69class G4QElastic : public G4VDiscreteProcess
70{
71public:
72
73  // Constructor
74  G4QElastic(const G4String& processName ="CHIPSElasticScattering");
75
76  // Destructor
77  ~G4QElastic();
78
79  G4bool IsApplicable(const G4ParticleDefinition& particle);
80
81  G4double GetMeanFreePath(const G4Track& aTrack, G4double previousStepSize,
82                           G4ForceCondition* condition);
83  // It returns the MeanFreePath of the process for the current track :
84  // (energy, material)
85  // The previousStepSize and G4ForceCondition* are not used.
86  // This function overloads a virtual function of the base class.       
87  // It is invoked by the ProcessManager of the Particle.
88 
89
90  G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep); 
91  // It computes the final state of the process (at end of step),
92  // returned as a ParticleChange object.       
93  // This function overloads a virtual function of the base class.
94  // It is invoked by the ProcessManager of the Particle.
95
96
97  G4LorentzVector GetEnegryMomentumConservation();
98
99  G4int GetNumberOfNeutronsInTarget();
100
101private:
102
103  // Hide assignment operator as private
104  G4QElastic& operator=(const G4QElastic &right);
105
106  // Copy constructor
107  G4QElastic(const G4QElastic&);
108
109  // Calculate XS/t: oxs=true - only CS; xst=true - calculate XS, xst=false(oxs=f/t) - t/tm
110  G4double CalculateXSt(G4bool oxs, G4bool xst, G4double p, G4int Z, G4int N, G4int pPDG);
111
112  // BODY
113  // Static Parameters --------------------------------------------------------------------
114  static G4int    nPartCWorld; // The#of particles for hadronization (limit of A of fragm.)
115  //--------------------------------- End of static parameters ---------------------------
116  // Working parameters
117  G4VQCrossSection* theCS;
118  G4LorentzVector EnMomConservation;                  // Residual of Energy/Momentum Cons.
119  G4int nOfNeutrons;                                  // #of neutrons in the target nucleus
120
121  // Modifires for the reaction
122  G4double Time;                                      // Time shift of the capture reaction
123  G4double EnergyDeposition;                          // Energy deposited in the reaction
124  static std::vector <G4int> ElementZ;                // Z of the element(i) in theLastCalc
125  static std::vector <G4double> ElProbInMat;          // SumProbabilityElements in Material
126  static std::vector <std::vector<G4int>*> ElIsoN;    // N of isotope(j) of Element(i)
127  static std::vector <std::vector<G4double>*> IsoProbInEl;// SumProbabIsotopes in Element i
128};
129#endif
Note: See TracBrowser for help on using the repository browser.