source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/body/include/G4QHadron.hh @ 1315

Last change on this file since 1315 was 1228, checked in by garnier, 14 years ago

update geant4.9.3 tag

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