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

Last change on this file since 819 was 819, checked in by garnier, 16 years ago

import all except CVS

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