[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 | // |
---|
| 27 | // $Id: G4QHadron.hh,v 1.36 2008/01/09 09:37:24 gcosmo Exp $ |
---|
| 28 | // GEANT4 tag $Name: $ |
---|
| 29 | // |
---|
| 30 | // ---------------- G4QHadron ---------------- |
---|
| 31 | // by Mikhail Kossov, Sept 1999. |
---|
| 32 | // class header for Hadrons generated by the CHIPS Model |
---|
| 33 | // ------------------------------------------------------ |
---|
| 34 | |
---|
| 35 | #ifndef G4QHadron_h |
---|
| 36 | #define G4QHadron_h 1 |
---|
| 37 | |
---|
| 38 | #include "globals.hh" |
---|
| 39 | #include "G4ThreeVector.hh" |
---|
| 40 | #include "G4LorentzVector.hh" |
---|
| 41 | #include "Randomize.hh" |
---|
| 42 | #include "G4QParticle.hh" |
---|
| 43 | #include "G4QPartonVector.hh" |
---|
| 44 | #include <list> |
---|
| 45 | |
---|
| 46 | class G4QHadron |
---|
| 47 | { |
---|
| 48 | public: |
---|
| 49 | // Constructors |
---|
| 50 | G4QHadron(); // Default Constructor |
---|
| 51 | G4QHadron(G4LorentzVector p); // Kinematical Constructor |
---|
| 52 | G4QHadron(G4int PDGcode, G4LorentzVector p=G4LorentzVector(0.,0.,0.,0.));//CHIPS-W Hadron |
---|
| 53 | G4QHadron(G4QPDGCode QPDG, G4LorentzVector p=G4LorentzVector(0.,0.,0.,0.));//CHIPS-W Had. |
---|
| 54 | G4QHadron(G4QContent QC, G4LorentzVector p=G4LorentzVector(0.,0.,0.,0.));//QC C-W Hadron |
---|
| 55 | G4QHadron(G4int PDG, G4double m, G4QContent QC); // Constructor for Chipolino or Quasmon |
---|
| 56 | G4QHadron(G4QPDGCode QPDG, G4double m, G4QContent QC);// Constr. for Chipolino or Quasmon |
---|
| 57 | G4QHadron(G4int PDG, G4LorentzVector p, G4QContent QC);// Constr for Chipolino or Quasmon |
---|
| 58 | G4QHadron(G4QPDGCode QPDG, G4LorentzVector p, G4QContent QC);// Con. for Chipo or Quasmon |
---|
| 59 | G4QHadron(G4QParticle* pPart, G4double maxM); // Constructor for Res with RANDOM mass |
---|
| 60 | G4QHadron(const G4QHadron& right); // Copy constructor by object |
---|
| 61 | G4QHadron(const G4QHadron* right); // Copy constructor by pointer |
---|
| 62 | G4QHadron(const G4QHadron* right, G4int ColC, G4ThreeVector Pos, G4LorentzVector Mom); |
---|
| 63 | virtual ~G4QHadron(); // Destructor |
---|
| 64 | // Operators |
---|
| 65 | const G4QHadron& operator=(const G4QHadron& right); |
---|
| 66 | G4bool operator==(const G4QHadron& right) const; |
---|
| 67 | G4bool operator!=(const G4QHadron& right) const; |
---|
| 68 | // Selectors |
---|
| 69 | G4int GetPDGCode() const; // Get PDG code of the Hadron |
---|
| 70 | G4int GetQCode() const; // Get Q code of the Hadron |
---|
| 71 | G4QPDGCode GetQPDG() const; // Get QPDG of the Hadron |
---|
| 72 | G4double GetSpin() const{return .5*(GetPDGCode()%10-1);} |
---|
| 73 | G4LorentzVector Get4Momentum() const; // Get 4-mom of the Hadron |
---|
| 74 | G4QContent GetQC() const; // Get private quark content |
---|
| 75 | G4double GetMass() const; // Get a mass of the Hadron |
---|
| 76 | G4double GetMass2() const; // Get an m^2 value for the Hadron |
---|
| 77 | G4double GetWidth() const; // Get Width of Hadron |
---|
| 78 | G4int GetNFragments() const; // Get a#of Fragments of this Hadron |
---|
| 79 | G4int GetCharge() const; // Get Charge of the Hadron |
---|
| 80 | G4int GetStrangeness() const; // Get Strangeness of the Hadron |
---|
| 81 | G4int GetBaryonNumber() const; // Get Baryon Number of the Hadron |
---|
| 82 | const G4ThreeVector& GetPosition() const; // Get hadron coordinates |
---|
| 83 | G4int GetSoftCollisionCount(); // Get QGS Counter of collisions |
---|
| 84 | G4bool IsSplit() {return isSplit;} // Check that hadron has been split |
---|
| 85 | G4double GetBindingEnergy() {return bindE;}// Returns binding E in NucMatter |
---|
| 86 | G4double GetFormationTime() {return formTime;} // Returns formation time |
---|
| 87 | |
---|
| 88 | // Modifiers |
---|
| 89 | void SetQPDG(const G4QPDGCode& QPDG); // Set QPDG of the Hadron |
---|
| 90 | void Set4Momentum(const G4LorentzVector& aMom); // Set 4-mom of the Hadron |
---|
| 91 | void SetQC(const G4QContent& newQC); // Set new private quark content |
---|
| 92 | void SetNFragments(const G4int& nf); // Set a#of Fragments of this Hadron |
---|
| 93 | void NegPDGCode(); // Change a sign of the PDG code |
---|
| 94 | void MakeAntiHadron(); // Make AntiHadron of this Hadron |
---|
| 95 | void SetPosition(const G4ThreeVector& aPosition); // Set coordinates of hadron position |
---|
| 96 | void IncrementCollisionCount(G4int aCount); // Increnment counter of collisions |
---|
| 97 | void SetCollisionCount(G4int aCount); // Set counter of QGSM collisions |
---|
| 98 | void Splitting() {isSplit = true;} // Put Up a flag that splitting is done |
---|
| 99 | void SplitUp(); // Make QGSM String Splitting of Hadron |
---|
| 100 | G4QParton* GetNextParton(); // Next Parton in a string |
---|
| 101 | G4QParton* GetNextAntiParton(); // Next Anti-Parton in a string |
---|
| 102 | void SetBindingEnergy(G4double aBindE){bindE=aBindE;}// Set Binding E in Nuclear Matter |
---|
| 103 | void Boost(const G4LorentzVector& theBoost); // Boosts hadron's 4-Momentum using 4M |
---|
| 104 | void Boost(const G4ThreeVector& B){theMomentum.boost(B);} // Boosts 4-Momentum using v/c |
---|
| 105 | void SetFormationTime(G4double fT){formTime=fT;} // Defines formationTime for the Hadron |
---|
| 106 | |
---|
| 107 | // General |
---|
| 108 | G4double RandomizeMass(G4QParticle* pPart, G4double maxM); // Randomize a mass value |
---|
| 109 | G4bool TestRealNeutral(); |
---|
| 110 | G4bool DecayIn2(G4LorentzVector& f4Mom, G4LorentzVector& s4Mom); |
---|
| 111 | G4bool CorMDecayIn2(G4double corM, G4LorentzVector& fr4Mom);// This(newMass corM)+fr4Mom |
---|
| 112 | G4bool CorEDecayIn2(G4double corE, G4LorentzVector& fr4Mom);// This(E+=cE,P)+f(fE-=cE,fP) |
---|
| 113 | G4bool RelDecayIn2(G4LorentzVector& f4Mom, G4LorentzVector& s4Mom, G4LorentzVector& dir, |
---|
| 114 | G4double maxCost = 1., G4double minCost = -1.); |
---|
| 115 | G4bool CopDecayIn2(G4LorentzVector& f4Mom, G4LorentzVector& s4Mom, G4LorentzVector& dir, |
---|
| 116 | G4double cop); |
---|
| 117 | G4bool DecayIn3(G4LorentzVector& f4Mom, G4LorentzVector& s4Mom, G4LorentzVector& t4Mom); |
---|
| 118 | G4bool RelDecayIn3(G4LorentzVector& fh4M, G4LorentzVector& sh4M, G4LorentzVector& th4Mom, |
---|
| 119 | G4LorentzVector& dir, G4double maxCost = 1., G4double minCost = -1.); |
---|
| 120 | G4bool CopDecayIn3(G4LorentzVector& fh4M, G4LorentzVector& sh4M, G4LorentzVector& th4Mom, |
---|
| 121 | G4LorentzVector& dir, G4double cosp); |
---|
| 122 | void Init3D(); // Initializes 3D nucleus with (Pos,4M)nucleons |
---|
| 123 | |
---|
| 124 | private: |
---|
| 125 | // Private methods |
---|
| 126 | void DefineQC(G4int PDGCode); |
---|
| 127 | G4QParton* BuildSeaQuark(G4bool isAntiQuark, G4int aPDGCode); |
---|
| 128 | G4double SampleX(G4double anXmin, G4int nSea, G4int theTotalSea, G4double aBeta); |
---|
| 129 | void GetValenceQuarkFlavors(G4QParton* &Part1,G4QParton* &Part2); |
---|
| 130 | G4ThreeVector GaussianPt(G4double widthSquare, G4double maxPtSquare); |
---|
| 131 | G4bool SplitMeson(G4int PDGcode, G4int* aEnd, G4int* bEnd); |
---|
| 132 | G4bool SplitBaryon(G4int PDGcode, G4int* aEnd, G4int* bEnd); |
---|
| 133 | |
---|
| 134 | private: |
---|
| 135 | // Static Parameters of QGSM Splitting |
---|
| 136 | static G4double alpha; // changing rapidity distribution for all |
---|
| 137 | static G4double beta; // changing rapidity distribution for projectile region |
---|
| 138 | static G4double theMinPz; // Can be from 14 to 140 MeV |
---|
| 139 | static G4double StrangeSuppress; // ? M.K. |
---|
| 140 | static G4double sigmaPt; // Can be 0 |
---|
| 141 | static G4double widthOfPtSquare; // ? M.K. |
---|
| 142 | static G4double minTransverseMass;// ? M.K. |
---|
| 143 | // Body |
---|
| 144 | G4QPDGCode theQPDG; // Instance of QPDG for the Hadron |
---|
| 145 | G4LorentzVector theMomentum; // The 4-mom of Hadron |
---|
| 146 | G4QContent valQ; // QC (@@ for Quasmon and Chipolino?) |
---|
| 147 | G4int nFragm; // 0 - stable, N - decayed in N part's |
---|
| 148 | // Body of Splitable Hadron and Nuclear Nucleon |
---|
| 149 | G4ThreeVector thePosition; // Coordinates of Hadron position |
---|
| 150 | G4int theCollisionCount; // ? |
---|
| 151 | G4bool isSplit; // Flag, that splitting was done |
---|
| 152 | G4bool Direction; // FALSE=target, TRUE=projectile |
---|
| 153 | std::list<G4QParton*> Color; // container for quarks & anti-diquarks |
---|
| 154 | std::list<G4QParton*> AntiColor; // container for anti-quarks & diquarks |
---|
| 155 | G4double bindE; // Binding energy in nuclear matter |
---|
| 156 | G4double formTime; // Formation time for the hadron |
---|
| 157 | }; |
---|
| 158 | |
---|
| 159 | typedef std::pair<G4QHadron*, G4QHadron*> G4QHadronPair; |
---|
| 160 | |
---|
| 161 | inline G4bool G4QHadron::operator==(const G4QHadron &rhs) const {return this==&rhs;} |
---|
| 162 | inline G4bool G4QHadron::operator!=(const G4QHadron &rhs) const {return this!=&rhs;} |
---|
| 163 | |
---|
| 164 | inline G4int G4QHadron::GetPDGCode() const {return theQPDG.GetPDGCode();} |
---|
| 165 | inline G4int G4QHadron::GetQCode() const {return theQPDG.GetQCode();} |
---|
| 166 | inline G4QPDGCode G4QHadron::GetQPDG() const {return theQPDG;} |
---|
| 167 | inline G4QContent G4QHadron::GetQC() const {return valQ;} |
---|
| 168 | inline G4LorentzVector G4QHadron::Get4Momentum() const {return theMomentum;} |
---|
| 169 | inline G4int G4QHadron::GetNFragments() const {return nFragm;} |
---|
| 170 | //@@ This is an example how to make other inline selectors for the 4-Momentum of the Hadron |
---|
| 171 | inline G4double G4QHadron::GetMass() const {return theMomentum.m();} |
---|
| 172 | inline G4double G4QHadron::GetMass2() const {return theMomentum.m2();} |
---|
| 173 | //@@ This is an example how to make other inline selectors for the Hadron |
---|
| 174 | inline G4int G4QHadron::GetCharge() const {return valQ.GetCharge();} |
---|
| 175 | inline G4int G4QHadron::GetStrangeness() const {return valQ.GetStrangeness();} |
---|
| 176 | inline G4int G4QHadron::GetBaryonNumber() const {return valQ.GetBaryonNumber();} |
---|
| 177 | inline const G4ThreeVector& G4QHadron::GetPosition() const {return thePosition;} |
---|
| 178 | inline G4int G4QHadron::GetSoftCollisionCount() {return theCollisionCount;} |
---|
| 179 | |
---|
| 180 | inline void G4QHadron::MakeAntiHadron() {if(TestRealNeutral()) NegPDGCode();} |
---|
| 181 | inline void G4QHadron::SetQC(const G4QContent& newQC) {valQ=newQC;} |
---|
| 182 | inline void G4QHadron::Set4Momentum(const G4LorentzVector& aMom) {theMomentum=aMom;} |
---|
| 183 | inline void G4QHadron::SetNFragments(const G4int& nf) {nFragm=nf;} |
---|
| 184 | inline void G4QHadron::SetPosition(const G4ThreeVector& position) {thePosition=position;} |
---|
| 185 | inline void G4QHadron::IncrementCollisionCount(G4int aCount) {theCollisionCount+=aCount;} |
---|
| 186 | inline void G4QHadron::SetCollisionCount(G4int aCount) {theCollisionCount =aCount;} |
---|
| 187 | |
---|
| 188 | inline void G4QHadron::NegPDGCode() {theQPDG.NegPDGCode(); valQ.Anti();} |
---|
| 189 | inline G4bool G4QHadron::TestRealNeutral() { return theQPDG.TestRealNeutral();} |
---|
| 190 | |
---|
| 191 | inline G4QParton* G4QHadron::GetNextParton() |
---|
| 192 | { |
---|
| 193 | if(Color.size()==0) return 0; |
---|
| 194 | G4QParton* result = Color.back(); |
---|
| 195 | Color.pop_back(); |
---|
| 196 | return result; |
---|
| 197 | } |
---|
| 198 | |
---|
| 199 | inline G4QParton* G4QHadron::GetNextAntiParton() |
---|
| 200 | { |
---|
| 201 | if(AntiColor.size() == 0) return 0; |
---|
| 202 | G4QParton* result = AntiColor.front(); |
---|
| 203 | AntiColor.pop_front(); |
---|
| 204 | return result; |
---|
| 205 | } |
---|
| 206 | #endif |
---|
| 207 | |
---|