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

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

update processes

File size: 5.7 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.4 2008/10/02 21:10:07 dennis Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
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 Interaction Prosesses
34// -------------------------------------------------------------------------------
35// At present (Aug-07) it is not tested.
36// The normalization is based on the temporary G4QIonIonCrossSection class
37// -------------------------------------------------------------------------------
38
39#ifndef G4QLowEnergy_hh
40#define G4QLowEnergy_hh
41
42// GEANT4 Headers
43#include "globals.hh"
44#include "G4ios.hh"
45#include "Randomize.hh"
46#include "G4VDiscreteProcess.hh"
47#include "G4Track.hh"
48#include "G4Step.hh"
49#include "G4ParticleTypes.hh"
50#include "G4ParticleTable.hh"
51#include "G4VParticleChange.hh"
52#include "G4ParticleDefinition.hh"
53#include "G4DynamicParticle.hh"
54#include "G4NucleiPropertiesTable.hh"
55#include "G4ThreeVector.hh"
56#include "G4LorentzVector.hh"
57#include "G4HadronicProcessType.hh"
58
59// CHIPS Headers
60#include "G4QNucleus.hh"
61#include "G4QIonIonCrossSection.hh"
62#include "G4QIsotope.hh"
63#include "G4QPDGToG4Particle.hh"
64#include "G4QCHIPSWorld.hh"
65#include "G4QHadronVector.hh"
66#include <vector>
67
68class G4QLowEnergy : public G4VDiscreteProcess
69{
70public:
71
72  // Constructor
73  G4QLowEnergy(const G4String& processName ="CHIPS_LowEnergyIonIonInelastic");
74
75  // Destructor
76  ~G4QLowEnergy();
77
78  G4bool IsApplicable(const G4ParticleDefinition& particle);
79
80  G4double GetMeanFreePath(const G4Track& aTrack, G4double previousStepSize,
81                           G4ForceCondition* condition);
82  // It returns the MeanFreePath of the process for the current track :
83  // (energy, material)
84  // The previousStepSize and G4ForceCondition* are not used.
85  // This function overloads a virtual function of the base class.                   
86  // It is invoked by the ProcessManager of the Particle.
87 
88
89  G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep); 
90  // It computes the final state of the process (at end of step),
91  // returned as a ParticleChange object.                           
92  // This function overloads a virtual function of the base class.
93  // It is invoked by the ProcessManager of the Particle.
94
95  G4LorentzVector GetEnegryMomentumConservation();
96
97  G4int GetNumberOfNeutronsInTarget();
98
99  void SwitchOnEvaporation() {evaporate=true;}   // Evaporation instead of gammas (default)
100
101  void SwitchOffEvaporation() {evaporate=false;} // Gammas instead of evaporation
102
103private:
104
105  // Hide assignment operator as private
106  G4QLowEnergy& operator=(const G4QLowEnergy &right);
107
108  // Copy constructor
109  G4QLowEnergy(const G4QLowEnergy&);
110
111  // Calculate Cross-Section of the Diffraction Reaction (p is in GeV @@ units)
112  G4double CalculateXS(G4double p, G4int Z, G4int N, G4int pPDG);
113
114                // BODY
115  // Static Parameters --------------------------------------------------------------------
116  static G4int    nPartCWorld; // The#of particles for hadronization (limit of A of fragm.)
117  //--------------------------------- End of static parameters ---------------------------
118  // Working parameters
119  G4bool evaporate;
120  G4VQCrossSection* theCS;
121  G4LorentzVector EnMomConservation;                  // Residual of Energy/Momentum Cons.
122  G4int nOfNeutrons;                                  // #of neutrons in the target nucleus
123
124  // Modifires for the reaction
125  G4double Time;                                      // Time shift of the capture reaction
126  G4double EnergyDeposition;                          // Energy deposited in the reaction
127  static std::vector <G4int> ElementZ;                // Z of the element(i) in theLastCalc
128  static std::vector <G4double> ElProbInMat;          // SumProbabilityElements in Material
129  static std::vector <std::vector<G4int>*> ElIsoN;    // N of isotope(j) of Element(i)
130  static std::vector <std::vector<G4double>*> IsoProbInEl;// SumProbabIsotopes in Element i
131};
132#endif
Note: See TracBrowser for help on using the repository browser.