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

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

maj sur la beta de geant 4.9.3

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