Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4QInelastic.hh
Go to the documentation of this file.
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$
27//
28// ---------------- G4QInelastic header ----------------
29// by Mikhail Kossov, December 2003.
30// Header of G4QInelastic class (all particles) of the CHIPS Physics Pachage in GEANT4
31// ------------------------------------------------------------------------------------
32// This is a unique CHIPS class for the Hadron-Nuclear Inelastic Interaction Prosesses.
33// ------------------------------------------------------------------------------------
34// At present (Dec.04) only pi+/-, K+/- proton, neutron, antiproton and antineutron
35// collisions with protons are implemented, which are fundamental for the in matter
36// simulation of hadronic reactions. The interactions of the same particles with
37// nuclei are implemented only for the low energy (below 1 GeV) nucleon-nuclear
38// reactions only. The collisions of nuclei with nuclei are planned for the near future.
39// The simulation is based on the G4QuasmonString class, which extends the CHIPS model
40// to the highest energyes, implementing the Quasmon string with the
41// String->Quasmons->Hadrons scenario of the quark-gluon string fragmentation
42// --> CHIPS is a SU(3) event generator, so it does not include reactions with the
43// heavy c,b,t-quarks, which can be simulated only by the SU(6) QUIPS (QUark Invariant
44// Phase Space) model which is an expantion of the CHIPS.- May 2009, M.Kossov.-
45// --------------------------------------------------------------------------------------
46// Algorithms: the vacuum interactions in CHIPS are described by the quark exchange (QE)
47// process. The first step is the low energy quark exchange. If as a result of the QE one
48// or both secondary hadrons are below the pi0 threshold (roughly) they are pushed to the
49// Ground State (GS) value(s). The excited (above the pi0 production threshold) hadronic
50// state is considered as a Quasmon, which is filled in the G4QuasmonVector of the
51// G4QuasmonString class. On the second step all G4Quasmons are decayed by the
52// G4Quasmon class and fill the G4QHadronVector output. If the exchange quark is too far
53// in the rapidity space (a parameter of the G4QuasmonString class) from any of the quarks
54// of the other hadron it creates a string with the nearest in the rapidity space quark.
55// This string is converted into a Quasmon. This forces the coalescence of the residuals
56// to create another Quasmon, while the possibility exists to create more residual
57// Quasmons instead of one - one per each target-quark+projectile-antiquark(diquark) pair.
58// This possibility is tuned by the Drell-Yan pair production process. If the target (or
59// pojectile) is a nucleus, then the Quasmons are created not only in vacuum, where they
60// can be fragmented by the G4Quasmon class, but in nuclear matter of the residual target
61// (or projectile). If the Quasmons are created in nuclear matter, they are fragmented by
62// the G4QEnvironment class with the subsequent Quark Exchange nuclear fragmentation.
63// This is the present general scenario.- May 2009, Mikhail Kossov.-
64// --------------------------------------------------------------------------------
65// ****************************************************************************************
66// This Header is a part of the CHIPS Physics Package (author: M. Kosov)
67// ****************************************************************************************
68// Short description: This is a universal class for the incoherent (inelastic)
69// nuclear interactions within the framework of the CHIPS model.
70// ---------------------------------------------------------------------------
71
72#ifndef G4QInelastic_hh
73#define G4QInelastic_hh
74
75// GEANT4 Headers
76#include "globals.hh"
77#include "G4ios.hh"
78#include "Randomize.hh"
79#include "G4QThd.hh"
80#include "G4VDiscreteProcess.hh"
81#include "G4Track.hh"
82#include "G4Step.hh"
83#include "G4ParticleTypes.hh"
84#include "G4VParticleChange.hh"
86#include "G4DynamicParticle.hh"
87#include "G4ThreeVector.hh"
88#include "G4LorentzVector.hh"
89#include "G4RandomDirection.hh"
90
91// CHIPS Headers
92#include "G4QEnvironment.hh"
93#include "G4VQCrossSection.hh"
94#include "G4QIsotope.hh"
117#include "G4QIonIonCollision.hh"
118#include "G4QFragmentation.hh"
119#include "G4QuasiFreeRatios.hh"
120#include "G4QPDGToG4Particle.hh"
121
123{
124public:
125
126 // Constructor
127 G4QInelastic(const G4String& processName ="CHIPS_Inelastic");
128
129 // Destructor
131
133
134 G4double GetMeanFreePath(const G4Track& aTrack, G4double previousStepSize,
136 // It returns the MeanFreePath of the process for the current track :
137 // (energy, material)
138 // The previousStepSize and G4ForceCondition* are not used.
139 // This function overloads a virtual function of the base class.
140 // It is invoked by the ProcessManager of the Particle.
141
142
143 G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep);
144 // It computes the final state of the process (at end of step),
145 // returned as a ParticleChange object.
146 // This function overloads a virtual function of the base class.
147 // It is invoked by the ProcessManager of the Particle.
148
149 // Fake void functions
153
154 // Internal Energy-Momentum Residual
156
157 // Number of neutrons in the target nucleus (primary)
159
160 // Static functions ---------------------------------------------------------------------
161 static void SetManual();
162 static void SetStandard();
163 static void SetParameters(G4double temper=180., G4double ssin2g=.1, G4double etaetap=.3,
164 G4double fN=0., G4double fD=0., G4double cP=1., G4double mR=1.,
165 G4int npCHIPSWorld=234, G4double solAn=.5, G4bool efFlag=false,
166 G4double piTh=141.4,G4double mpi2=20000.,G4double dinum=1880.);
167 static void SetPhotNucBias(G4double phnB=1.);
168 static void SetWeakNucBias(G4double ccnB=1.);
169 //--- End of static member functions ----------------------------------------------------
170
171 G4double GetPhotNucBias(){return photNucBias;}
172 G4double GetWeakNucBias(){return weakNucBias;}
173
174private:
175
176 // Hide assignment operator as private
177 G4QInelastic& operator=(const G4QInelastic &right);
178
179 // Copy constructor
181
182 // Random direction in two dimentions pair(first=sin(phi), second=cos(phi))
183 std::pair<G4double,G4double> Random2DDirection();
184
185 // BODY
186 // Static Parameters --------------------------------------------------------------------
187 static G4bool manualFlag; // If false then standard parameters are used
188 static G4int nPartCWorld; // The#of particles for hadronization (limit of A of fragm.)
189 // -> Parameters of the G4Quasmon class:
190 static G4double Temperature; // Quasmon Temperature
191 static G4double SSin2Gluons; // Percent of ssbar sea in a constituen gluon
192 static G4double EtaEtaprime; // Part of eta-prime in all etas
193 // -> Parameters of the G4QNucleus class:
194 static G4double freeNuc; // probability of the quasi-free baryon on surface
195 static G4double freeDib; // probability of the quasi-free dibaryon on surface
196 static G4double clustProb; // clusterization probability in dense region
197 static G4double mediRatio; // relative vacuum hadronization probability
198 // -> Parameters of the G4QEnvironment class:
199 static G4bool EnergyFlux; // Flag for Energy Flux use instead of Multy Quasmon
200 static G4double SolidAngle; // Part of Solid Angle to capture secondaries(@@A-dep)
201 static G4double PiPrThresh; // Pion Production Threshold for gammas
202 static G4double M2ShiftVir; // Shift for M2=-Q2=m_pi^2 of the virtual gamma
203 static G4double DiNuclMass; // Double Nucleon Mass for virtual normalization
204 // -> Biasing parameters:
205 static G4double photNucBias; // Biasing parameter for photo-($e,mu,tau)Nuclear reactions
206 static G4double weakNucBias; // Biasing parameter for Charged Currents (nu,mu) reactions
207 //--------------------------------- End of static parameters ---------------------------
208 // Working parameters
209 G4VQCrossSection* theCS;
210 G4LorentzVector EnMomConservation; // Residual of Energy/Momentum Cons.
211 G4int nOfNeutrons; // #of neutrons in the target nucleus
212
213 // Modifires for the reaction
214 G4double Time; // Time shift of the capture reaction
215 G4double EnergyDeposition; // Energy deposited in the reaction
216 static std::vector <G4int> ElementZ; // Z of the element(i) in theLastCalc
217 static std::vector <G4double> ElProbInMat; // SumProbabilityElements in Material
218 static std::vector <std::vector<G4int>*> ElIsoN; // N of isotope(j) of Element(i)
219 static std::vector <std::vector<G4double>*> IsoProbInEl;// SumProbabIsotopes in Element i
220};
221#endif
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4double GetWeakNucBias()
static void SetStandard()
Definition: G4QInelastic.cc:95
G4int GetNumberOfNeutronsInTarget()
void BuildPhysicsTable(const G4ParticleDefinition &)
static void SetManual()
Definition: G4QInelastic.cc:94
void PrintInfoDefinition()
static void SetPhotNucBias(G4double phnB=1.)
void SetPhysicsTableBining(G4double, G4double, G4int)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
static void SetParameters(G4double temper=180., G4double ssin2g=.1, G4double etaetap=.3, G4double fN=0., G4double fD=0., G4double cP=1., G4double mR=1., G4int npCHIPSWorld=234, G4double solAn=.5, G4bool efFlag=false, G4double piTh=141.4, G4double mpi2=20000., G4double dinum=1880.)
Definition: G4QInelastic.cc:98
G4double GetPhotNucBias()
static void SetWeakNucBias(G4double ccnB=1.)
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
G4bool IsApplicable(const G4ParticleDefinition &particle)
G4LorentzVector GetEnegryMomentumConservation()
Definition: G4Step.hh:78