Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QNucleus.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//
27// $Id$
28//
29// ---------------- G4QNucleus ----------------
30// by Mikhail Kossov, Sept 1999.
31// class header for the nuclei and nuclear environment of the CHIPS Model
32// -----------------------------------------------------------------------
33// Short description: a class describing properties of nuclei, which
34// are necessary for the CHIPS Model.
35// -----------------------------------------------------------------------
36
37#ifndef G4QNucleus_h
38#define G4QNucleus_h 1
39
40#include <utility>
41#include <vector>
43
44#include "globals.hh"
45#include "G4RandomDirection.hh"
46#include "G4QCandidateVector.hh"
47#include "G4QHadronVector.hh"
48#include "G4LorentzRotation.hh"
49#include "G4QChipolino.hh"
50
51class G4QNucleus : public G4QHadron
52{
53public:
54 G4QNucleus(); // Default Constructor
55 G4QNucleus(G4int nucPDG); // At Rest PDG-Constructor
56 G4QNucleus(G4LorentzVector p, G4int nucPDG); // Full PDG-Constructor
57 G4QNucleus(G4QContent nucQC); // At Rest QuarkCont-Constructor
58 G4QNucleus(G4QContent nucQC, G4LorentzVector p); // Full QuarkCont-Constructor
59 G4QNucleus(G4int z, G4int n, G4int s=0); // At Rest ZNS-Constructor
60 G4QNucleus(G4int z, G4int n, G4int s, G4LorentzVector p);// Full ZNS-Constructor
61 G4QNucleus(G4QNucleus* right, G4bool cop3D = false); // Copy Constructor by pointer
62 G4QNucleus(const G4QNucleus &right, G4bool cop3D=false); // Copy Constructor by value
63 ~G4QNucleus(); // Public Destructor
64 // Overloaded Operators
65 const G4QNucleus& operator=(const G4QNucleus& right);
66 G4bool operator==(const G4QNucleus &right) const {return this==&right;}
67 G4bool operator!=(const G4QNucleus &right) const {return this!=&right;}
68 // Specific Selectors
69 G4int GetPDG() const {return 90000000+1000*(1000*S+Z)+N;}// PDG Code of Nucleus
70 G4int GetZ() const {return Z;} // Get a#of protons
71 G4int GetN() const {return N;} // Get a#of neutrons
72 G4int GetS() const {return S;} // Get a#of lambdas
73 G4int GetA() const {return Z+N+S;} // Get A of the nucleus
74 G4int GetDZ() const {return dZ;} // Get a#of protons in dense region
75 G4int GetDN() const {return dN;} // Get a#of neutrons in dense region
76 G4int GetDS() const {return dS;} // Get a#of lambdas in dense region
77 G4int GetDA() const {return dZ+dN+dS;} // Get A of the dense part of nucleus
78 G4int GetMaxClust() const {return maxClust;} // Get Max BarNum of Clusters
79 G4double GetProbability(G4int bn=0) const {return probVect[bn];} // clust(BarN)probabil
80 G4double GetMZNS() const {return GetQPDG().GetNuclMass(Z,N,S);} // not H or Q
81 G4double GetTbIntegral(); // Calculate the integral of T(b)
82 G4double GetGSMass() const {return GetQPDG().GetMass();}//Nucleus GSMass (not Hadron)
83 G4QContent GetQCZNS() const // Get ZNS quark content of Nucleus
84 {
85 if(S>=0) return G4QContent(Z+N+N+S,Z+Z+N+S,S,0,0,0);
86 else return G4QContent(Z+N+N+S,Z+Z+N+S,0,0,0,-S);
87 }
88 G4int GetNDefMesonC() const{return nDefMesonC;}; // max#of predefed mesonCandidates
89 G4int GetNDefBaryonC()const{return nDefBaryonC;};// max#of predefed baryonCandidates
90 G4double GetDensity(const G4ThreeVector&aPos) {return rho0*GetRelativeDensity(aPos);}
91 G4double GetRho0() {return rho0;} // One nucleon prob-density
92 G4double GetRelativeDensity(const G4ThreeVector& aPosition); // Densyty/rho0
93 G4double GetRelWSDensity(const G4double& r) // Wood-Saxon rho/rho0(r)
94 {return 1./(1.+std::exp((r-radius)/WoodSaxonSurf));}
95 G4double GetRelOMDensity(const G4double& r2){return std::exp(-r2/radius);} // OscModelRelDens
96 G4double GetRadius(const G4double maxRelativeDenisty=0.5); // Radius of %ofDensity
97 G4double GetOuterRadius(); // Get radius of the most far nucleon
98 G4double GetDeriv(const G4ThreeVector& point); // Derivitive of density
99 G4double GetFermiMomentum(G4double density); // Returns modul of FermyMomentum(dens)
101 {
102 //G4cout<<"G4QNucleus::GetNextNucleon: cN="<<currentNucleon<<", A="<<GetA()<<G4endl;
103 return (currentNucleon>=0&&currentNucleon<GetA()) ? theNucleons[currentNucleon++] : 0;
104 }
105 void SubtractNucleon(G4QHadron* pNucleon); // Subtract the nucleon from the 3D Nucleus
106 void DeleteNucleons(); // Deletes all residual nucleons
108 {
109 G4LorentzVector sum(0.,0.,0.,0.);
110 for(unsigned i=0; i<theNucleons.size(); i++) sum += theNucleons[i]->Get4Momentum();
111 sum.setE(std::sqrt(sqr(GetGSMass())+sum.v().mag2())); // Energy is corrected !
112 return sum;
113 }
114 std::vector<G4double> const* GetBThickness() {return &Tb;} // T(b) function, step .1 fm
115
116 // Specific Modifiers
117 G4bool EvaporateBaryon(G4QHadron* h1,G4QHadron* h2); // Evaporate Baryon from Nucleus
118 void EvaporateNucleus(G4QHadron* hA, G4QHadronVector* oHV);// Evaporate Nucleus
119 //void DecayBaryon(G4QHadron* dB, G4QHadronVector* oHV); // gamma+N or Delt->N+Pi @@later
120 void DecayDibaryon(G4QHadron* dB, G4QHadronVector* oHV); // deuteron is kept
121 void DecayAntiDibaryon(G4QHadron* dB, G4QHadronVector* oHV);// antiDeuteron is kept
122 void DecayIsonucleus(G4QHadron* dB, G4QHadronVector* oHV); // nP+(Pi+) or nN+(Pi-)
123 void DecayMultyBaryon(G4QHadron* dB, G4QHadronVector* oHV);// A*p, A*n or A*L
124 void DecayAntiStrange(G4QHadron* dB, G4QHadronVector* oHV);// nuclei with K+/K0
125 void DecayAlphaBar(G4QHadron* dB, G4QHadronVector* oHV); // alpha+p or alpha+n
126 void DecayAlphaDiN(G4QHadron* dB, G4QHadronVector* oHV); // alpha+p+p
127 void DecayAlphaAlpha(G4QHadron* dB, G4QHadronVector* oHV); // alpha+alpha
128 G4int SplitBaryon(); // Is it possible to split baryon/alpha
129 G4int HadrToNucPDG(G4int hPDG); // Converts hadronic PDGCode to nuclear
130 G4int NucToHadrPDG(G4int nPDG); // Converts nuclear PDGCode to hadronic
131 G4bool Split2Baryons(); // Is it possible to split two baryons?
132 void ActivateBThickness(); // Calculate T(b) for nucleus (db=.1fm)
133 G4double GetBThickness(G4double b); // Calculates T(b)
134 G4double GetThickness(G4double b); // Calculates T(b)/rho(0)
135 void InitByPDG(G4int newPDG); // Init existing nucleus by new PDG
136 void InitByQC(G4QContent newQC) // Init existing nucleus by new QCont
137 {G4int PDG=G4QPDGCode(newQC).GetPDGCode(); InitByPDG(PDG);}
138 void IncProbability(G4int bn); // Add one cluster to probability
139 void Increase(G4int PDG, G4LorentzVector LV = G4LorentzVector(0.,0.,0.,0.));
140 void Increase(G4QContent QC, G4LorentzVector LV = G4LorentzVector(0.,0.,0.,0.));
141 void Reduce(G4int PDG); // Reduce Nucleus by PDG fragment
143 void SetMaxClust(G4int maxC){maxClust=maxC;}// Set Max BarNum of Clusters
144 void InitCandidateVector(G4QCandidateVector& theQCandidates,
145 G4int nM=45, G4int nB=72, G4int nC=117);
146 void PrepareCandidates(G4QCandidateVector& theQCandidates, G4bool piF=false, G4bool
147 gaF=false, G4LorentzVector LV=G4LorentzVector(0.,0.,0.,0.));
148 G4int UpdateClusters(G4bool din); // Return a#of clusters & calc.probab's
149 G4QNucleus operator+=(const G4QNucleus& rhs); // Add a cluster to the nucleus
150 G4QNucleus operator-=(const G4QNucleus& rhs); // Subtract a cluster from a nucleus
151 G4QNucleus operator*=(const G4int& rhs); // Multiplication of the Nucleus
152 G4bool StartLoop(); // returns size of theNucleons (cN=0)
153 G4bool ReduceSum(G4ThreeVector* vectors, G4ThreeVector sum);// Reduce zero-sum of vectors
154 void SimpleSumReduction(G4ThreeVector* vectors, G4ThreeVector sum); // Reduce zero-V-sum
155 void DoLorentzBoost(const G4LorentzVector& theBoost) // Boost nucleons by 4-vector
156 {
157 theMomentum.boost(theBoost);
158 for(unsigned i=0; i<theNucleons.size(); i++) theNucleons[i]->Boost(theBoost);
159 }
160 void DoLorentzRotation(const G4LorentzRotation& theLoRot) // Lorentz Rotate nucleons
161 {
162 theMomentum=theLoRot*theMomentum;
163 for(unsigned i=0; i<theNucleons.size(); i++) theNucleons[i]->LorentzRotate(theLoRot);
164 }
165 void DoLorentzBoost(const G4ThreeVector& theBeta)// Boost nucleons by v/c
166 {
167 theMomentum.boost(theBeta);
168 for(unsigned i=0; i<theNucleons.size(); i++) theNucleons[i]->Boost(theBeta);
169 }
171 void DoLorentzContraction(const G4ThreeVector& theBeta); // Lorentz Contraction by v/c
172 void DoTranslation(const G4ThreeVector& theShift); // Used only in G4QFragmentation
173
174 // Static functions
175 static void SetParameters(G4double fN=.1,G4double fD=.05, G4double cP=4., G4double mR=1.,
176 G4double nD=.8*CLHEP::fermi);
177
178 // Specific General Functions
179 G4int RandomizeBinom(G4double p,G4int N); // Randomize according to Binomial Law
180 G4double CoulombBarGen(const G4double& rZ, const G4double& rA, const G4double& cZ,
181 const G4double& cA); // CoulombBarrier in MeV (General)
182 G4double CoulombBarrier(const G4double& cZ=1, const G4double& cA=1, G4double dZ=0.,
183 G4double dA=0.); // CoulombBarrier in MeV
184 G4double FissionCoulombBarrier(const G4double& cZ, const G4double& cA, G4double dZ=0.,
185 G4double dA=0.); // Fission CoulombBarrier in MeV
186 G4double BindingEnergy(const G4double& cZ=0, const G4double& cA=0, G4double dZ=0.,
187 G4double dA=0.);
188 G4double CoulBarPenProb(const G4double& CB, const G4double& E, const G4int& C,
189 const G4int& B);
190 std::pair<G4double, G4double> ChooseImpactXandY(G4double maxImpact); // Randomize bbar
191 void ChooseNucleons(); // Initializes 3D Nucleons
192 void ChoosePositions(); // Initializes positions of 3D nucleons
193 void ChooseFermiMomenta(); // Initializes FermyMoms of 3D nucleons
194 void InitDensity(); // Initializes density distribution
195 void Init3D(); // automatically starts the LOOP
196private:
197 // Specific Encapsulated Functions
198 void SetZNSQC(G4int z, G4int n, G4int s); // Set QC, using Z,N,S
199 G4QNucleus GetThis() const {return G4QNucleus(Z,N,S);} // @@ Check for memory leak
200
201// Body
202private:
203 // Static Parameters
204 static const G4int nDefMesonC =45;
205 static const G4int nDefBaryonC=72;
206 //
207 static G4double freeNuc; // probability of the quasi-free baryon on surface
208 static G4double freeDib; // probability of the quasi-free dibaryon on surface
209 static G4double clustProb; // clusterization probability in dense region
210 static G4double mediRatio; // relative vacuum hadronization probability
211 static G4double nucleonDistance;// Distance between nucleons (0.8 fm)
212 static G4double WoodSaxonSurf; // Surface parameter of Wood-Saxon density (0.545 fm)
213 // The basic
214 G4int Z; // Z of the Nucleus
215 G4int N; // N of the Nucleus
216 G4int S; // S of the Nucleus
217 // The secondaries
218 G4int dZ; // Z of the dense region of the nucleus
219 G4int dN; // N of the dense region of the nucleus
220 G4int dS; // S of the dense region of the nucleus
221 G4int maxClust; // Baryon Number of the last calculated cluster
222 G4double probVect[256]; // Cluster probability ("a#of issues" can be real) Vector
223 // 3D
224 G4QHadronVector theNucleons; // Vector of nucleons of which the Nucleus consists of
225 G4int currentNucleon; // Current nucleon for the NextNucleon (? M.K.)
226 G4double rho0; // Normalazation density
227 G4double radius; // Nuclear radius
228 std::vector<G4double> Tb; // T(b) function with step .1 fm (@@ make .1 a parameter)
229 G4bool TbActive; // Flag that the T(b) is activated
230 G4bool RhoActive; // Flag that the Density is activated
231};
232
233std::ostream& operator<<(std::ostream& lhs, G4QNucleus& rhs);
234std::ostream& operator<<(std::ostream& lhs, const G4QNucleus& rhs);
235
236#endif
CLHEP::HepLorentzVector G4LorentzVector
std::vector< G4QCandidate * > G4QCandidateVector
std::vector< G4QHadron * > G4QHadronVector
std::ostream & operator<<(std::ostream &lhs, G4QNucleus &rhs)
Definition: G4QNucleus.cc:357
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
double mag2() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
Hep3Vector v() const
G4LorentzVector Get4Momentum() const
Definition: G4QHadron.hh:79
void Boost(const G4LorentzVector &theBoost)
Definition: G4QHadron.cc:1293
void LorentzRotate(const G4LorentzRotation &rotation)
Definition: G4QHadron.hh:112
G4LorentzVector theMomentum
Definition: G4QHadron.hh:143
void Set4Momentum(const G4LorentzVector &aMom)
Definition: G4QHadron.hh:187
G4QPDGCode GetQPDG() const
Definition: G4QHadron.hh:172
G4double BindingEnergy(const G4double &cZ=0, const G4double &cA=0, G4double dZ=0., G4double dA=0.)
Definition: G4QNucleus.cc:3418
G4int GetNDefMesonC() const
Definition: G4QNucleus.hh:88
void SimpleSumReduction(G4ThreeVector *vectors, G4ThreeVector sum)
Definition: G4QNucleus.cc:3853
void InitByPDG(G4int newPDG)
Definition: G4QNucleus.cc:371
G4int GetN() const
Definition: G4QNucleus.hh:71
G4int GetZ() const
Definition: G4QNucleus.hh:70
G4double GetRelativeDensity(const G4ThreeVector &aPosition)
Definition: G4QNucleus.cc:3757
void DecayAlphaBar(G4QHadron *dB, G4QHadronVector *oHV)
Definition: G4QNucleus.cc:7414
G4double GetProbability(G4int bn=0) const
Definition: G4QNucleus.hh:79
std::pair< G4double, G4double > ChooseImpactXandY(G4double maxImpact)
Definition: G4QNucleus.cc:3551
G4int GetPDG() const
Definition: G4QNucleus.hh:69
G4double GetRelOMDensity(const G4double &r2)
Definition: G4QNucleus.hh:95
G4bool operator!=(const G4QNucleus &right) const
Definition: G4QNucleus.hh:67
G4double GetRelWSDensity(const G4double &r)
Definition: G4QNucleus.hh:93
G4double GetOuterRadius()
Definition: G4QNucleus.cc:3947
G4double CoulombBarGen(const G4double &rZ, const G4double &rA, const G4double &cZ, const G4double &cA)
Definition: G4QNucleus.cc:3358
G4int RandomizeBinom(G4double p, G4int N)
Definition: G4QNucleus.cc:3015
G4bool StartLoop()
Definition: G4QNucleus.cc:3982
G4QNucleus operator-=(const G4QNucleus &rhs)
Definition: G4QNucleus.cc:4081
void SubtractNucleon(G4QHadron *pNucleon)
Definition: G4QNucleus.cc:612
G4int HadrToNucPDG(G4int hPDG)
Definition: G4QNucleus.cc:4117
void InitCandidateVector(G4QCandidateVector &theQCandidates, G4int nM=45, G4int nB=72, G4int nC=117)
Definition: G4QNucleus.cc:3037
void DecayAntiDibaryon(G4QHadron *dB, G4QHadronVector *oHV)
Definition: G4QNucleus.cc:6591
void DoLorentzBoost(const G4LorentzVector &theBoost)
Definition: G4QNucleus.hh:155
void Increase(G4int PDG, G4LorentzVector LV=G4LorentzVector(0., 0., 0., 0.))
Definition: G4QNucleus.cc:729
void SetMaxClust(G4int maxC)
Definition: G4QNucleus.hh:143
G4int GetA() const
Definition: G4QNucleus.hh:73
G4double GetDeriv(const G4ThreeVector &point)
Definition: G4QNucleus.cc:3735
void DoLorentzRotation(const G4LorentzRotation &theLoRot)
Definition: G4QNucleus.hh:160
void DecayIsonucleus(G4QHadron *dB, G4QHadronVector *oHV)
Definition: G4QNucleus.cc:5925
G4double GetMZNS() const
Definition: G4QNucleus.hh:80
G4int GetS() const
Definition: G4QNucleus.hh:72
void DecayAntiStrange(G4QHadron *dB, G4QHadronVector *oHV)
Definition: G4QNucleus.cc:6896
void DecayAlphaAlpha(G4QHadron *dB, G4QHadronVector *oHV)
Definition: G4QNucleus.cc:7665
void DecayMultyBaryon(G4QHadron *dB, G4QHadronVector *oHV)
Definition: G4QNucleus.cc:7190
void Reduce(G4int PDG)
Definition: G4QNucleus.cc:705
G4double GetRadius(const G4double maxRelativeDenisty=0.5)
Definition: G4QNucleus.cc:3746
G4QContent GetQCZNS() const
Definition: G4QNucleus.hh:83
const G4QNucleus & operator=(const G4QNucleus &right)
Definition: G4QNucleus.cc:307
G4int GetDA() const
Definition: G4QNucleus.hh:77
void CalculateMass()
Definition: G4QNucleus.hh:142
G4bool Split2Baryons()
Definition: G4QNucleus.cc:863
void IncProbability(G4int bn)
G4int GetMaxClust() const
Definition: G4QNucleus.hh:78
void ChooseFermiMomenta()
Definition: G4QNucleus.cc:3773
G4double CoulBarPenProb(const G4double &CB, const G4double &E, const G4int &C, const G4int &B)
Definition: G4QNucleus.cc:3441
G4double GetTbIntegral()
Definition: G4QNucleus.cc:4018
G4bool operator==(const G4QNucleus &right) const
Definition: G4QNucleus.hh:66
G4int GetDN() const
Definition: G4QNucleus.hh:75
G4int GetDS() const
Definition: G4QNucleus.hh:76
void DeleteNucleons()
Definition: G4QNucleus.cc:697
G4int GetNDefBaryonC() const
Definition: G4QNucleus.hh:89
void ChooseNucleons()
Definition: G4QNucleus.cc:3570
void DecayDibaryon(G4QHadron *dB, G4QHadronVector *oHV)
Definition: G4QNucleus.cc:6286
static void SetParameters(G4double fN=.1, G4double fD=.05, G4double cP=4., G4double mR=1., G4double nD=.8 *CLHEP::fermi)
Definition: G4QNucleus.cc:347
void PrepareCandidates(G4QCandidateVector &theQCandidates, G4bool piF=false, G4bool gaF=false, G4LorentzVector LV=G4LorentzVector(0., 0., 0., 0.))
Definition: G4QNucleus.cc:3115
void ActivateBThickness()
Definition: G4QNucleus.cc:3991
G4int NucToHadrPDG(G4int nPDG)
Definition: G4QNucleus.cc:4148
G4bool ReduceSum(G4ThreeVector *vectors, G4ThreeVector sum)
Definition: G4QNucleus.cc:3861
void ChoosePositions()
Definition: G4QNucleus.cc:3599
void DecayAlphaDiN(G4QHadron *dB, G4QHadronVector *oHV)
Definition: G4QNucleus.cc:7331
G4double GetThickness(G4double b)
Definition: G4QNucleus.cc:4048
void DoLorentzBoost(const G4ThreeVector &theBeta)
Definition: G4QNucleus.hh:165
G4LorentzVector GetNucleons4Momentum()
Definition: G4QNucleus.hh:107
void Init3D()
Definition: G4QNucleus.cc:3910
G4QHadron * GetNextNucleon()
Definition: G4QNucleus.hh:100
G4double GetRho0()
Definition: G4QNucleus.hh:91
void InitDensity()
Definition: G4QNucleus.cc:3692
std::vector< G4double > const * GetBThickness()
Definition: G4QNucleus.hh:114
G4bool EvaporateBaryon(G4QHadron *h1, G4QHadron *h2)
Definition: G4QNucleus.cc:951
G4double GetFermiMomentum(G4double density)
Definition: G4QNucleus.cc:3765
G4int UpdateClusters(G4bool din)
Definition: G4QNucleus.cc:417
void InitByQC(G4QContent newQC)
Definition: G4QNucleus.hh:136
G4double FissionCoulombBarrier(const G4double &cZ, const G4double &cA, G4double dZ=0., G4double dA=0.)
Definition: G4QNucleus.cc:3403
G4QNucleus operator+=(const G4QNucleus &rhs)
Definition: G4QNucleus.cc:4063
G4int SplitBaryon()
Definition: G4QNucleus.cc:774
G4int GetDZ() const
Definition: G4QNucleus.hh:74
G4QNucleus operator*=(const G4int &rhs)
Definition: G4QNucleus.cc:4099
G4double GetDensity(const G4ThreeVector &aPos)
Definition: G4QNucleus.hh:90
void DoLorentzContraction(const G4LorentzVector &B)
Definition: G4QNucleus.hh:170
void DoTranslation(const G4ThreeVector &theShift)
Definition: G4QNucleus.cc:3974
G4double CoulombBarrier(const G4double &cZ=1, const G4double &cA=1, G4double dZ=0., G4double dA=0.)
Definition: G4QNucleus.cc:3386
void EvaporateNucleus(G4QHadron *hA, G4QHadronVector *oHV)
Definition: G4QNucleus.cc:4171
G4double GetGSMass() const
Definition: G4QNucleus.hh:82
G4int GetPDGCode() const
Definition: G4QPDGCode.hh:326
G4double GetMass()
Definition: G4QPDGCode.cc:693
G4double GetNuclMass(G4int Z, G4int N, G4int S)
Definition: G4QPDGCode.cc:766
T sqr(const T &x)
Definition: templates.hh:145
#define const
Definition: zconf.h:118