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

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

update processes

File size: 10.9 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: G4QCollision.hh,v 1.11 2008/10/02 21:10:07 dennis Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29//      ---------------- G4QCollision header ----------------
30//                 by Mikhail Kossov, December 2003.
31//  Header of G4QCollision class (mu-,pi-,K-) of the CHIPS Simulation Branch in GEANT4
32// -------------------------------------------------------------------------------
33// This is a unique CHIPS class for the Hadron-Nuclear Inelastic Interaction Prosesses.
34// -------------------------------------------------------------------------------
35// At present (Dec.04) only pi+/-, K+/- proton, neutron, antiproton and antineutron
36// collisions with protons are implemented, which are fundamental for the in matter
37// simulation of hadronic reactions. The interactions of the same particles with
38// nuclei are planned only. The collisions of nuclei with nuclei are possible...
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), which can be simulated only by the SU(6) QUIPS (QUark Invariant
44// Phase Space) model which is an expantion of the CHIPS.-December 2003.M.Kossov.-
45// -------------------------------------------------------------------------------
46// Algorithms: the interactions in CHIPS are described by the quark exchange (QE) process.
47// The first step is the low energy quark exchange. If as a result of the QE one or
48// 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 fiill 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// in the another Quasmon, while the possibility exist to create more residual Quasmons
57// instead of one - one per each target-quark+projectile-antiquark(diquark) pair. This
58// possibility is tuned by the Drell-Yan pair production process. If the target (or
59// pojectile) are nuclei, 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 crated in nuclear matter, they are fragmented by
62// the G4QEnvironment class with the subsequent Quark Exchange nuclear fragmentation.
63// This is the planned scenario.- December 2004.Mikhail Kossov.-
64// --------------------------------------------------------------------------------
65// ****************************************************************************************
66// ********* This HEADER is temporary moved from the photolepton_hadron directory *********
67// ******* DO NOT MAKE ANY CHANGE! With time it'll move back to photolepton...(M.K.) ******
68// ****************************************************************************************
69
70#ifndef G4QCollision_hh
71#define G4QCollision_hh
72
73// GEANT4 Headers
74#include "globals.hh"
75#include "G4ios.hh"
76#include "Randomize.hh"
77#include "G4VDiscreteProcess.hh"
78#include "G4Track.hh"
79#include "G4Step.hh"
80#include "G4ParticleTypes.hh"
81#include "G4VParticleChange.hh"
82#include "G4ParticleDefinition.hh"
83#include "G4DynamicParticle.hh"
84#include "G4NucleiPropertiesTable.hh"
85#include "G4ThreeVector.hh"
86#include "G4LorentzVector.hh"
87#include "G4HadronicProcessType.hh"
88
89// CHIPS Headers
90#include "G4QEnvironment.hh"
91#include "G4VQCrossSection.hh"
92#include "G4QIsotope.hh"
93#include "G4QProtonNuclearCrossSection.hh"
94#include "G4QPhotonNuclearCrossSection.hh"
95#include "G4QElectronNuclearCrossSection.hh"
96#include "G4QMuonNuclearCrossSection.hh"
97#include "G4QTauNuclearCrossSection.hh"
98#include "G4QNuMuNuclearCrossSection.hh"
99#include "G4QANuMuNuclearCrossSection.hh"
100#include "G4QNuENuclearCrossSection.hh"
101#include "G4QANuENuclearCrossSection.hh"
102#include "G4QNuNuNuclearCrossSection.hh"
103#include "G4QANuANuNuclearCrossSection.hh"
104//#include "G4QuasmonString.hh"
105#include "G4QuasiFreeRatios.hh"
106#include "G4QPDGToG4Particle.hh"
107#include <vector>
108
109class G4QCollision : public G4VDiscreteProcess
110{
111public:
112
113  // Constructor
114  G4QCollision(const G4String& processName ="CHIPSNuclearCollision");
115
116  // Destructor
117  ~G4QCollision();
118
119  G4bool IsApplicable(const G4ParticleDefinition& particle);
120
121  G4double GetMeanFreePath(const G4Track& aTrack, G4double previousStepSize,
122                           G4ForceCondition* condition);
123  // It returns the MeanFreePath of the process for the current track :
124  // (energy, material)
125  // The previousStepSize and G4ForceCondition* are not used.
126  // This function overloads a virtual function of the base class.                   
127  // It is invoked by the ProcessManager of the Particle.
128 
129
130  G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep); 
131  // It computes the final state of the process (at end of step),
132  // returned as a ParticleChange object.                           
133  // This function overloads a virtual function of the base class.
134  // It is invoked by the ProcessManager of the Particle.
135
136  // Fake void functions
137  void SetPhysicsTableBining(G4double, G4double, G4int) {;}
138  void BuildPhysicsTable(const G4ParticleDefinition&) {;}
139  void PrintInfoDefinition() {;}
140
141  // Internal Energy-Momentum Residual
142  G4LorentzVector GetEnegryMomentumConservation();
143
144  // Number of neutrons in the target nucleus (primary)
145  G4int GetNumberOfNeutronsInTarget();
146
147  // Static functions ---------------------------------------------------------------------
148  static void SetManual();
149  static void SetStandard();
150  static void SetParameters(G4double temper=180., G4double ssin2g=.1, G4double etaetap=.3,
151                            G4double fN=0., G4double fD=0., G4double cP=1., G4double mR=1.,
152                            G4int npCHIPSWorld=234, G4double solAn=.5, G4bool efFlag=false,
153                            G4double piTh=141.4,G4double mpi2=20000.,G4double dinum=1880.);
154  static void SetPhotNucBias(G4double phnB=1.);
155  static void SetWeakNucBias(G4double ccnB=1.);
156  //--- End of static member functions ----------------------------------------------------
157
158  G4double GetPhotNucBias(){return photNucBias;}
159  G4double GetWeakNucBias(){return weakNucBias;}
160
161private:
162
163  // Hide assignment operator as private
164  G4QCollision& operator=(const G4QCollision &right);
165
166  // Copy constructor
167  G4QCollision(const G4QCollision&);
168
169  // Random direction in two dimentions pair(first=sin(phi), second=cos(phi))
170  std::pair<G4double,G4double> Random2DDirection();
171
172                // BODY
173  // Static Parameters --------------------------------------------------------------------
174  static G4bool   manualFlag;  // If false then standard parameters are used
175  static G4int    nPartCWorld; // The#of particles for hadronization (limit of A of fragm.)
176  // -> Parameters of the G4Quasmon class:
177  static G4double Temperature; // Quasmon Temperature
178  static G4double SSin2Gluons; // Percent of ssbar sea in a constituen gluon
179  static G4double EtaEtaprime; // Part of eta-prime in all etas
180  // -> Parameters of the G4QNucleus class:
181  static G4double freeNuc;     // probability of the quasi-free baryon on surface
182  static G4double freeDib;     // probability of the quasi-free dibaryon on surface
183  static G4double clustProb;   // clusterization probability in dense region
184  static G4double mediRatio;   // relative vacuum hadronization probability
185  // -> Parameters of the G4QEnvironment class:
186  static G4bool   EnergyFlux;  // Flag for Energy Flux use instead of Multy Quasmon
187  static G4double SolidAngle;  // Part of Solid Angle to capture secondaries(@@A-dep)
188  static G4double PiPrThresh;  // Pion Production Threshold for gammas
189  static G4double M2ShiftVir;  // Shift for M2=-Q2=m_pi^2 of the virtual gamma
190  static G4double DiNuclMass;  // Double Nucleon Mass for virtual normalization
191  // -> Biasing parameters:
192  static G4double photNucBias; // Biasing parameter for photo-($e,mu,tau)Nuclear reactions
193  static G4double weakNucBias; // Biasing parameter for Charged Currents (nu,mu) reactions
194  //--------------------------------- End of static parameters ---------------------------
195  // Working parameters
196  G4VQCrossSection* theCS;
197  G4LorentzVector EnMomConservation;                  // Residual of Energy/Momentum Cons.
198  G4int nOfNeutrons;                                  // #of neutrons in the target nucleus
199
200  // Modifires for the reaction
201  G4double Time;                                      // Time shift of the capture reaction
202  G4double EnergyDeposition;                          // Energy deposited in the reaction
203  static std::vector <G4int> ElementZ;                // Z of the element(i) in theLastCalc
204  static std::vector <G4double> ElProbInMat;          // SumProbabilityElements in Material
205  static std::vector <std::vector<G4int>*> ElIsoN;    // N of isotope(j) of Element(i)
206  static std::vector <std::vector<G4double>*> IsoProbInEl;// SumProbabIsotopes in Element i
207};
208#endif
209
Note: See TracBrowser for help on using the repository browser.