source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/processes/include/G4QElastic.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: 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: G4QElastic.hh,v 1.5 2010/02/16 07:53:05 mkossov Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-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// Short description: A process for CHIPS hadron-nucleus elastic scattering.
36// -------------------------------------------------------------------------------
37
38#ifndef G4QElastic_hh
39#define G4QElastic_hh
40
41// GEANT4 Headers
42#include "globals.hh"
43#include "G4ios.hh"
44#include "Randomize.hh"
45#include "G4VDiscreteProcess.hh"
46#include "G4Track.hh"
47#include "G4Step.hh"
48#include "G4ParticleTypes.hh"
49#include "G4VParticleChange.hh"
50#include "G4ParticleDefinition.hh"
51#include "G4DynamicParticle.hh"
52#include "G4ThreeVector.hh"
53#include "G4LorentzVector.hh"
54
55// CHIPS Headers
56#include "G4QPionPlusElasticCrossSection.hh"
57#include "G4QPionMinusElasticCrossSection.hh"
58#include "G4QKaonPlusElasticCrossSection.hh"
59#include "G4QKaonMinusElasticCrossSection.hh"
60#include "G4QHyperonElasticCrossSection.hh"
61#include "G4QHyperonPlusElasticCrossSection.hh"
62#include "G4QProtonElasticCrossSection.hh"
63#include "G4QNeutronElasticCrossSection.hh"
64#include "G4QAntiBaryonElasticCrossSection.hh"
65#include "G4QIsotope.hh"
66#include "G4QCHIPSWorld.hh"
67#include "G4QHadron.hh"
68#include <vector>
69
70class G4QElastic : public G4VDiscreteProcess
71{
72public:
73
74  // Constructor
75  G4QElastic(const G4String& processName ="CHIPSElasticScattering");
76
77  // Destructor
78  ~G4QElastic();
79
80  G4bool IsApplicable(const G4ParticleDefinition& particle);
81
82  G4double GetMeanFreePath(const G4Track& aTrack, G4double previousStepSize,
83                           G4ForceCondition* condition);
84  // It returns the MeanFreePath of the process for the current track :
85  // (energy, material)
86  // The previousStepSize and G4ForceCondition* are not used.
87  // This function overloads a virtual function of the base class.       
88  // It is invoked by the ProcessManager of the Particle.
89 
90
91  G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep); 
92  // It computes the final state of the process (at end of step),
93  // returned as a ParticleChange object.       
94  // This function overloads a virtual function of the base class.
95  // It is invoked by the ProcessManager of the Particle.
96
97
98  G4LorentzVector GetEnegryMomentumConservation();
99
100  G4int GetNumberOfNeutronsInTarget();
101
102private:
103
104  // Hide assignment operator as private
105  G4QElastic& operator=(const G4QElastic &right);
106
107  // Copy constructor
108  G4QElastic(const G4QElastic&);
109
110  // Calculate XS/t: oxs=true - only CS; xst=true - calculate XS, xst=false(oxs=f/t) - t/tm
111  G4double CalculateXSt(G4bool oxs, G4bool xst, 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  G4VQCrossSection* theCS;
119  G4LorentzVector EnMomConservation;                  // Residual of Energy/Momentum Cons.
120  G4int nOfNeutrons;                                  // #of neutrons in the target nucleus
121
122  // Modifires for the reaction
123  G4double Time;                                      // Time shift of the capture reaction
124  G4double EnergyDeposition;                          // Energy deposited in the reaction
125  static std::vector <G4int> ElementZ;                // Z of the element(i) in theLastCalc
126  static std::vector <G4double> ElProbInMat;          // SumProbabilityElements in Material
127  static std::vector <std::vector<G4int>*> ElIsoN;    // N of isotope(j) of Element(i)
128  static std::vector <std::vector<G4double>*> IsoProbInEl;// SumProbabIsotopes in Element i
129};
130#endif
Note: See TracBrowser for help on using the repository browser.