source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/interface/include/G4QIonIonElastic.hh @ 967

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

fichier ajoutes

File size: 5.6 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: G4QIonIonElastic.hh,v 1.2 2008/10/02 21:10:07 dennis Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29//      ---------------- G4QIonIonElastic header ----------------
30//                 by Mikhail Kossov, December 2006.
31//  Header of G4QIonIonElastic class (A+A) of the CHIPS Simulation Branch in GEANT4
32// -------------------------------------------------------------------------------
33// This is a unique CHIPS class for the Ion-Ion Elastic Scattering Prosesses
34// -------------------------------------------------------------------------------
35// @@ This class is on the testing level @@
36// -------------------------------------------------------------------------------
37// ****************************************************************************************
38// ********** This CLASS is temporary moved from the photolepton_hadron directory *********
39// ****************************************************************************************
40
41#ifndef G4QIonIonElastic_hh
42#define G4QIonIonElastic_hh
43
44// GEANT4 Headers
45#include "globals.hh"
46#include "G4ios.hh"
47#include "Randomize.hh"
48#include "G4VDiscreteProcess.hh"
49#include "G4Track.hh"
50#include "G4Step.hh"
51#include "G4ParticleTypes.hh"
52#include "G4VParticleChange.hh"
53#include "G4ParticleDefinition.hh"
54#include "G4DynamicParticle.hh"
55#include "G4NucleiPropertiesTable.hh"
56#include "G4ThreeVector.hh"
57#include "G4LorentzVector.hh"
58#include "G4HadronicProcessType.hh"
59
60// CHIPS Headers
61#include "G4QIonIonCrossSection.hh"
62#include "G4QElasticCrossSection.hh"
63#include "G4QIsotope.hh"
64#include "G4QPDGToG4Particle.hh"
65#include "G4QCHIPSWorld.hh"
66#include "G4QHadron.hh"
67#include <vector>
68
69class G4QIonIonElastic : public G4VDiscreteProcess
70{
71public:
72
73  // Constructor
74  G4QIonIonElastic(const G4String& processName ="CHIPS_IonIonElasticScattering");
75
76  // Destructor
77  ~G4QIonIonElastic();
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(){return EnMomConservation;}
98
99  G4int GetNumberOfNeutronsInTarget(){return nOfNeutrons;}
100
101private:
102
103  // Hide assignment operator as private
104  G4QIonIonElastic& operator=(const G4QIonIonElastic &right);
105
106  // Copy constructor
107  G4QIonIonElastic(const G4QIonIonElastic&);
108
109                // BODY
110  // Static Parameters --------------------------------------------------------------------
111  static G4int    nPartCWorld; // The#of particles for hadronization (limit of A of fragm.)
112  //--------------------------------- End of static parameters ---------------------------
113  // Working parameters
114  G4VQCrossSection* theCS;
115  G4LorentzVector EnMomConservation;                  // Residual of Energy/Momentum Cons.
116  G4int nOfNeutrons;                                  // #of neutrons in the target nucleus
117
118  // Modifires for the reaction
119  G4double Time;                                      // Time shift of the capture reaction
120  G4double EnergyDeposition;                          // Energy deposited in the reaction
121  static std::vector <G4int> ElementZ;                // Z of the element(i) in theLastCalc
122  static std::vector <G4double> ElProbInMat;          // SumProbabilityElements in Material
123  static std::vector <std::vector<G4int>*> ElIsoN;    // N of isotope(j) of Element(i)
124  static std::vector <std::vector<G4double>*> IsoProbInEl;// SumProbabIsotopes in Element i
125};
126#endif
127
Note: See TracBrowser for help on using the repository browser.