source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/processes/include/G4QLowEnergy.hh @ 1315

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

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 6.1 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: G4QLowEnergy.hh,v 1.2 2010/01/14 11:24:36 mkossov Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-cand-01 $
28//
29//      ---------------- G4QLowEnergy header ----------------
30//                 by Mikhail Kossov, Aug 2007.
31//  Header of G4QLowEnergy class (A+A) of the CHIPS Simulation Branch in GEANT4
32// -------------------------------------------------------------------------------
33// This is a unique CHIPS class for the Ion-Ion Low Energy Inelastic Interactions
34// -------------------------------------------------------------------------------
35// At present (Aug-07) it is not tested.
36// The normalization is based on the temporary G4QIonIonCrossSection class
37// -------------------------------------------------------------------------------
38// Short description: This is a fast low energy algorithm for the
39// inelastic interactions of nucleons and nuclei (ions) with nuclei.
40// This is a fase-space algorithm, but not quark level. Provides
41// nuclear fragments upto alpha only. Never was tumed (but can be).
42// ---------------------------------------------------------------
43
44#ifndef G4QLowEnergy_hh
45#define G4QLowEnergy_hh
46
47// GEANT4 Headers
48#include "globals.hh"
49#include "G4ios.hh"
50#include "Randomize.hh"
51#include "G4VDiscreteProcess.hh"
52#include "G4Track.hh"
53#include "G4Step.hh"
54#include "G4ParticleTypes.hh"
55#include "G4ParticleTable.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 "G4QNucleus.hh"
64#include "G4QIonIonCrossSection.hh"
65#include "G4QProtonNuclearCrossSection.hh"
66#include "G4QProtonElasticCrossSection.hh"
67#include "G4QNeutronElasticCrossSection.hh"
68#include "G4QIsotope.hh"
69#include "G4QCHIPSWorld.hh"
70#include "G4QHadronVector.hh"
71#include <vector>
72
73class G4QLowEnergy : public G4VDiscreteProcess
74{
75public:
76
77  // Constructor
78  G4QLowEnergy(const G4String& processName ="CHIPS_LowEnergyIonIonInelastic");
79
80  // Destructor
81  ~G4QLowEnergy();
82
83  G4bool IsApplicable(const G4ParticleDefinition& particle);
84
85  G4double GetMeanFreePath(const G4Track& aTrack, G4double previousStepSize,
86                           G4ForceCondition* condition);
87  // It returns the MeanFreePath of the process for the current track :
88  // (energy, material)
89  // The previousStepSize and G4ForceCondition* are not used.
90  // This function overloads a virtual function of the base class.       
91  // It is invoked by the ProcessManager of the Particle.
92 
93
94  G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep); 
95  // It computes the final state of the process (at end of step),
96  // returned as a ParticleChange object.       
97  // This function overloads a virtual function of the base class.
98  // It is invoked by the ProcessManager of the Particle.
99
100  G4LorentzVector GetEnegryMomentumConservation();
101
102  G4int GetNumberOfNeutronsInTarget();
103
104  void SwitchOnEvaporation() {evaporate=true;}   // Evaporation instead of gammas (default)
105
106  void SwitchOffEvaporation() {evaporate=false;} // Gammas instead of evaporation
107
108private:
109
110  // Hide assignment operator as private
111  G4QLowEnergy& operator=(const G4QLowEnergy &right);
112
113  // Copy constructor
114  G4QLowEnergy(const G4QLowEnergy&);
115
116  // Calculate Cross-Section of the Diffraction Reaction (p is in GeV @@ units)
117  G4double CalculateXS(G4double p, G4int Z, G4int N, G4int pPDG);
118
119  // BODY
120  // Static Parameters --------------------------------------------------------------------
121  static G4int    nPartCWorld; // The#of particles for hadronization (limit of A of fragm.)
122  //--------------------------------- End of static parameters ---------------------------
123  // Working parameters
124  G4bool evaporate;
125  G4VQCrossSection* theCS;
126  G4LorentzVector EnMomConservation;                  // Residual of Energy/Momentum Cons.
127  G4int nOfNeutrons;                                  // #of neutrons in the target nucleus
128
129  // Modifires for the reaction
130  G4double Time;                                      // Time shift of the capture reaction
131  G4double EnergyDeposition;                          // Energy deposited in the reaction
132  static std::vector <G4int> ElementZ;                // Z of the element(i) in theLastCalc
133  static std::vector <G4double> ElProbInMat;          // SumProbabilityElements in Material
134  static std::vector <std::vector<G4int>*> ElIsoN;    // N of isotope(j) of Element(i)
135  static std::vector <std::vector<G4double>*> IsoProbInEl;// SumProbabIsotopes in Element i
136};
137#endif
Note: See TracBrowser for help on using the repository browser.