[819] | 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 | // |
---|
[1055] | 27 | // $Id: G4QNucleus.hh,v 1.34 2009/02/23 09:49:24 mkossov Exp $ |
---|
| 28 | // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ |
---|
[819] | 29 | // |
---|
| 30 | // ---------------- G4QNucleus ---------------- |
---|
| 31 | // by Mikhail Kossov, Sept 1999. |
---|
[1055] | 32 | // class header for the nuclei and nuclear environment of the CHIPS Model |
---|
| 33 | // ----------------------------------------------------------------------- |
---|
| 34 | // Short description: a class describing properties of nuclei, which |
---|
| 35 | // are necessary for the CHIPS Model. |
---|
| 36 | // ----------------------------------------------------------------------- |
---|
[819] | 37 | |
---|
| 38 | #ifndef G4QNucleus_h |
---|
| 39 | #define G4QNucleus_h 1 |
---|
| 40 | |
---|
| 41 | #include "G4QCandidateVector.hh" |
---|
| 42 | #include "G4QHadronVector.hh" |
---|
| 43 | #include "G4QChipolino.hh" |
---|
| 44 | #include <utility> |
---|
| 45 | #include <vector> |
---|
| 46 | #include "globals.hh" |
---|
| 47 | |
---|
| 48 | class G4QNucleus : public G4QHadron |
---|
| 49 | { |
---|
| 50 | public: |
---|
| 51 | G4QNucleus(); // Default Constructor |
---|
| 52 | G4QNucleus(G4int nucPDG); // At Rest PDG-Constructor |
---|
| 53 | G4QNucleus(G4LorentzVector p, G4int nucPDG); // Full PDG-Constructor |
---|
| 54 | G4QNucleus(G4QContent nucQC); // At Rest QuarkCont-Constructor |
---|
| 55 | G4QNucleus(G4QContent nucQC, G4LorentzVector p); // Full QuarkCont-Constructor |
---|
| 56 | G4QNucleus(G4int z, G4int n, G4int s=0); // At Rest ZNS-Constructor |
---|
| 57 | G4QNucleus(G4int z, G4int n, G4int s, G4LorentzVector p);// Full ZNS-Constructor |
---|
| 58 | //G4QNucleus(const G4QNucleus& right); // Copy Constructor by value |
---|
| 59 | G4QNucleus(G4QNucleus* right); // Copy Constructor by pointer |
---|
| 60 | ~G4QNucleus(); // Public Destructor |
---|
| 61 | // Overloaded Operators |
---|
| 62 | const G4QNucleus& operator=(const G4QNucleus& right); |
---|
| 63 | G4bool operator==(const G4QNucleus &right) const {return this==&right;} |
---|
| 64 | G4bool operator!=(const G4QNucleus &right) const {return this!=&right;} |
---|
| 65 | // Specific Selectors |
---|
| 66 | G4int GetPDG() const {return 90000000+1000*(1000*S+Z)+N;}// PDG Code of Nucleus |
---|
| 67 | G4int GetZ() const {return Z;} // Get a#of protons |
---|
| 68 | G4int GetN() const {return N;} // Get a#of neutrons |
---|
| 69 | G4int GetS() const {return S;} // Get a#of lambdas |
---|
| 70 | G4int GetA() const {return Z+N+S;} // Get A of the nucleus |
---|
| 71 | G4int GetDZ() const {return dZ;} // Get a#of protons in dense region |
---|
| 72 | G4int GetDN() const {return dN;} // Get a#of neutrons in dense region |
---|
| 73 | G4int GetDS() const {return dS;} // Get a#of lambdas in dense region |
---|
| 74 | G4int GetDA() const {return dZ+dN+dS;} // Get A of the dense part of nucleus |
---|
| 75 | G4int GetMaxClust() const {return maxClust;} // Get Max BarNum of Clusters |
---|
| 76 | G4double GetProbability(G4int bn=0) const {return probVect[bn];} // clust(BarN)probabil |
---|
| 77 | G4double GetMZNS() const {return GetQPDG().GetNuclMass(Z,N,S);} // not H or Q |
---|
| 78 | G4double GetGSMass() const {return GetQPDG().GetMass();}//Nucleus GSMass (not Hadron) |
---|
| 79 | G4QContent GetQCZNS() const // Get ZNS quark content of Nucleus |
---|
| 80 | { |
---|
| 81 | if(S>=0) return G4QContent(Z+N+N+S,Z+Z+N+S,S,0,0,0); |
---|
| 82 | else return G4QContent(Z+N+N+S,Z+Z+N+S,0,0,0,-S); |
---|
| 83 | } |
---|
| 84 | G4int GetNDefMesonC() const{return nDefMesonC;}; // max#of predefed mesonCandidates |
---|
| 85 | G4int GetNDefBaryonC()const{return nDefBaryonC;};// max#of predefed baryonCandidates |
---|
| 86 | std::pair<G4double, G4double> RefetchImpactXandY() const {return theImpactParameter;} |
---|
| 87 | G4double GetDensity(const G4ThreeVector&aPos) {return rho0*GetRelativeDensity(aPos);} |
---|
| 88 | G4double GetRelativeDensity(const G4ThreeVector& aPosition); // Densyty/rho0 |
---|
| 89 | G4double GetRadius(const G4double maxRelativeDenisty=0.5); // Radius of %ofDensity |
---|
| 90 | G4double GetOuterRadius(); // Get radius of the most far nucleon |
---|
| 91 | G4double GetDeriv(const G4ThreeVector& point); // Derivitive of density |
---|
| 92 | G4double GetFermiMomentum(G4double density); // Returns modul of FermyMomentum(dens) |
---|
| 93 | G4ThreeVector Get3DFermiMomentum(G4double density, G4double maxMom=-1.) |
---|
| 94 | { |
---|
| 95 | if(maxMom<0) maxMom=GetFermiMomentum(density); |
---|
| 96 | return maxMom*RandomUnitSphere(); |
---|
| 97 | } |
---|
| 98 | G4QHadron* GetNextNucleon() |
---|
| 99 | {return (currentNucleon>=0&¤tNucleon<GetA()) ? theNucleons[currentNucleon++] :0;} |
---|
| 100 | //std::vector<G4double>* GetBThickness() const {return Tb;} // T(b) function, step .1 fm |
---|
| 101 | std::vector<G4double> const* GetBThickness() {return &Tb;} // T(b) function, step .1 fm |
---|
| 102 | |
---|
| 103 | // Specific Modifiers |
---|
| 104 | G4bool EvaporateBaryon(G4QHadron* h1,G4QHadron* h2); // Evaporate Baryon from Nucleus |
---|
| 105 | void EvaporateNucleus(G4QHadron* hA, G4QHadronVector* oHV);// Evaporate Nucleus |
---|
| 106 | //void DecayBaryon(G4QHadron* dB, G4QHadronVector* oHV); // gamma+N or Delt->N+Pi @@later |
---|
| 107 | void DecayDibaryon(G4QHadron* dB, G4QHadronVector* oHV); // deuteron is kept |
---|
| 108 | void DecayIsonucleus(G4QHadron* dB, G4QHadronVector* oHV); // nP+(Pi+) or nN+(Pi-) |
---|
| 109 | void DecayMultyBaryon(G4QHadron* dB, G4QHadronVector* oHV);// A*p, A*n or A*L |
---|
| 110 | void DecayAntiStrange(G4QHadron* dB, G4QHadronVector* oHV);// nuclei with K+/K0 |
---|
| 111 | void DecayAlphaBar(G4QHadron* dB, G4QHadronVector* oHV); // alpha+p or alpha+n |
---|
| 112 | void DecayAlphaDiN(G4QHadron* dB, G4QHadronVector* oHV); // alpha+p+p |
---|
| 113 | void DecayAlphaAlpha(G4QHadron* dB, G4QHadronVector* oHV); // alpha+alpha |
---|
| 114 | G4int SplitBaryon(); // Is it possible to split baryon/alpha |
---|
| 115 | G4int HadrToNucPDG(G4int hPDG); // Converts hadronic PDGCode to nuclear |
---|
| 116 | G4int NucToHadrPDG(G4int nPDG); // Converts nuclear PDGCode to hadronic |
---|
| 117 | G4bool Split2Baryons(); // Is it possible to split two baryons? |
---|
| 118 | void ActivateBThickness(); // Calculate T(b) for nucleus (db=.1fm) |
---|
| 119 | void InitByPDG(G4int newPDG); // Init existing nucleus by new PDG |
---|
| 120 | void InitByQC(G4QContent newQC) // Init existing nucleus by new QCont |
---|
| 121 | {G4int PDG=G4QPDGCode(newQC).GetPDGCode(); InitByPDG(PDG);} |
---|
| 122 | void IncProbability(G4int bn); // Add one cluster to probability |
---|
| 123 | void Increase(G4int PDG, G4LorentzVector LV = G4LorentzVector(0.,0.,0.,0.)); |
---|
| 124 | void Increase(G4QContent QC, G4LorentzVector LV = G4LorentzVector(0.,0.,0.,0.)); |
---|
| 125 | void Reduce(G4int PDG); // Reduce Nucleus by PDG fragment |
---|
| 126 | void CalculateMass() {Set4Momentum(G4LorentzVector(0.,0.,0.,GetGSMass()));} |
---|
| 127 | void SetMaxClust(G4int maxC){maxClust=maxC;}// Set Max BarNum of Clusters |
---|
| 128 | void InitCandidateVector(G4QCandidateVector& theQCandidates, |
---|
| 129 | G4int nM=45, G4int nB=72, G4int nC=117); |
---|
| 130 | void PrepareCandidates(G4QCandidateVector& theQCandidates, G4bool piF=false, G4bool |
---|
| 131 | gaF=false, G4LorentzVector LV=G4LorentzVector(0.,0.,0.,0.)); |
---|
| 132 | G4int UpdateClusters(G4bool din); // Return a#of clusters & calc.probab's |
---|
| 133 | G4QNucleus operator+=(const G4QNucleus& rhs); // Add a cluster to the nucleus |
---|
| 134 | G4QNucleus operator-=(const G4QNucleus& rhs); // Subtract a cluster from a nucleus |
---|
| 135 | G4QNucleus operator*=(const G4int& rhs); // Multiplication of the Nucleus |
---|
| 136 | G4bool StartLoop(); // returns size of theNucleons (cN=0) |
---|
| 137 | G4bool ReduceSum(G4ThreeVector* momentum, G4double*); // Reduce momentum nonconservation |
---|
| 138 | void DoLorentzBoost(const G4LorentzVector& theBoost); // Boost nucleons by 4-vector |
---|
| 139 | void DoLorentzBoost(const G4ThreeVector& theBeta);// Boost nucleons by v/c |
---|
| 140 | void DoLorentzContraction(const G4LorentzVector&B){DoLorentzContraction(B.vect()/B.e());} |
---|
| 141 | void DoLorentzContraction(const G4ThreeVector& theBeta); // Lorentz Contraction by v/c |
---|
| 142 | void DoTranslation(const G4ThreeVector& theShift); // Used only in GHAD-TFT |
---|
| 143 | |
---|
| 144 | // Static functions |
---|
| 145 | static void SetParameters(G4double fN=.1,G4double fD=.05, G4double cP=4., G4double mR=1., |
---|
| 146 | G4double nD=.8*fermi); |
---|
| 147 | |
---|
| 148 | // Specific General Functions |
---|
| 149 | G4ThreeVector RandomUnitSphere(); // Randomize position inside UnitSphere |
---|
| 150 | G4int RandomizeBinom(G4double p,G4int N); // Randomize according to Binomial Law |
---|
| 151 | G4double CoulombBarrier(const G4double& cZ=1, const G4double& cA=1, G4double dZ=0., |
---|
| 152 | G4double dA=0.); // CoulombBarrier in MeV |
---|
| 153 | G4double FissionCoulombBarrier(const G4double& cZ, const G4double& cA, G4double dZ=0., |
---|
| 154 | G4double dA=0.); // Fission CoulombBarrier in MeV |
---|
| 155 | G4double BindingEnergy(const G4double& cZ=0, const G4double& cA=0, G4double dZ=0., |
---|
| 156 | G4double dA=0.); |
---|
| 157 | G4double CoulBarPenProb(const G4double& CB, const G4double& E, const G4int& C, |
---|
| 158 | const G4int& B); |
---|
| 159 | std::pair<G4double, G4double> ChooseImpactXandY(G4double maxImpact); // Randomize bbar |
---|
| 160 | void ChooseNucleons(); // Initializes 3D Nucleons |
---|
| 161 | void ChoosePositions(); // Initializes positions of 3D nucleons |
---|
| 162 | void ChooseFermiMomenta(); // Initializes FermyMoms of 3D nucleons |
---|
| 163 | void InitDensity(); // Initializes density distribution |
---|
| 164 | void Init3D(); // automatically starts the LOOP |
---|
| 165 | private: |
---|
| 166 | // Specific Encapsulated Functions |
---|
| 167 | void SetZNSQC(G4int z, G4int n, G4int s); // Set QC, using Z,N,S |
---|
| 168 | G4QNucleus GetThis() const {return G4QNucleus(Z,N,S);} // @@ Check for memory leak |
---|
| 169 | |
---|
| 170 | // Body |
---|
| 171 | private: |
---|
| 172 | // Static Parameters |
---|
| 173 | static const G4int nDefMesonC =45; |
---|
| 174 | static const G4int nDefBaryonC=72; |
---|
| 175 | // |
---|
| 176 | static G4double freeNuc; // probability of the quasi-free baryon on surface |
---|
| 177 | static G4double freeDib; // probability of the quasi-free dibaryon on surface |
---|
| 178 | static G4double clustProb; // clusterization probability in dense region |
---|
| 179 | static G4double mediRatio; // relative vacuum hadronization probability |
---|
| 180 | static G4double nucleonDistance;// Distance between nucleons (0.8 fm) |
---|
| 181 | // The basic |
---|
| 182 | G4int Z; // Z of the Nucleus |
---|
| 183 | G4int N; // N of the Nucleus |
---|
| 184 | G4int S; // S of the Nucleus |
---|
| 185 | // The secondaries |
---|
| 186 | G4int dZ; // Z of the dense region of the nucleus |
---|
| 187 | G4int dN; // N of the dense region of the nucleus |
---|
| 188 | G4int dS; // S of the dense region of the nucleus |
---|
| 189 | G4int maxClust; // Baryon Number of the last calculated cluster |
---|
| 190 | G4double probVect[256]; // Cluster probability ("a#of issues" can be real) Vector |
---|
| 191 | // 3D |
---|
| 192 | std::pair<G4double, G4double> theImpactParameter; // 2D impact parameter vector bbar |
---|
| 193 | G4QHadronVector theNucleons; // Vector of nucleons of which Nucleus consists of |
---|
| 194 | G4int currentNucleon; // Current nucleon for the NextNucleon (? M.K.) |
---|
| 195 | G4double rho0; // Normalazation density |
---|
| 196 | G4double radius; // Nuclear radius |
---|
| 197 | //std::vector<G4double>* Tb; // T(b) function with step .1 fm (@@ make .1 a parameter) |
---|
| 198 | std::vector<G4double> Tb; // T(b) function with step .1 fm (@@ make .1 a parameter) |
---|
| 199 | }; |
---|
| 200 | |
---|
| 201 | std::ostream& operator<<(std::ostream& lhs, G4QNucleus& rhs); |
---|
| 202 | std::ostream& operator<<(std::ostream& lhs, const G4QNucleus& rhs); |
---|
| 203 | |
---|
| 204 | |
---|
| 205 | #endif |
---|