Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QString.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//
28// $Id$
29
30#ifndef G4QString_h
31#define G4QString_h 1
32
33// ------------------------------------------------------------
34// GEANT 4 class header file
35//
36// ---------------- G4QString ----------------
37// by Mikhail Kosov, October 2006.
38// class for an excited string used by Parton String Models
39// For comparison mirror member functions are taken from G4 classes:
40// G4FragmentingString
41// G4ExcitedStringDecay
42//
43// ------------------------------------------------------------
44// Short description: If partons from the G4QPartonPair are close in
45// rapidity, they create Quasmons, but if they are far in the rapidity
46// space, they can not interact directly. Say the bottom parton (quark)
47// has rapidity 0, and the top parton (antiquark) has rapidity 8, then
48// the top quark splits in two by radiating gluon, and each part has
49// rapidity 4, then the gluon splits in quark-antiquark pair (rapidity
50// 2 each), and then the quark gadiates anothe gluon and reachs rapidity
51// 1. Now it can interact with the bottom antiquark, creating a Quasmon
52// or a hadron. The intermediate partons is the string ladder.
53// ---------------------------------------------------------------------
54
55#include "G4ios.hh"
56#include "globals.hh"
57#include "G4ThreeVector.hh"
58#include "G4LorentzVector.hh"
59#include "G4LorentzRotation.hh"
60#include "G4QHadronVector.hh"
61#include "G4QPartonPair.hh"
62#include "G4QPartonVector.hh"
63#include <algorithm>
64
66{
67 public:
68
69 enum {PROJECTILE = 1, TARGET = -1}; // The same as in quark-pair (@@ ? M.K.)
70
71 G4QString(); // formal creation of the string with future excitation
72 G4QString(G4QParton* Color,G4QParton* Gluon, G4QParton* AntiColor, G4int Dir=PROJECTILE);
73 G4QString(G4QParton* Col, G4QParton* AntiCol, G4int Dir=PROJECTILE);
74 G4QString(G4QPartonPair* ColAntiCol);
75 G4QString(const G4QString &right);
76 G4QString(const G4QString &old, G4QParton* newdecay, const G4LorentzVector *momentum);
77 G4QString(const G4QString &old, G4QParton* newdecay);
78 G4QString(G4QParton* newdecay, const G4LorentzVector* momentum);
79
80 ~G4QString();
81
82 // Selectors
83 G4int operator==(const G4QString &right) const {return this == &right;}
84 G4int operator!=(const G4QString &right) const {return this != &right;}
85 const G4ThreeVector& GetPosition() const {return thePosition;}
86 const G4QPartonVector* GetPartonList() const {return &thePartons;}
87 G4QParton* GetLeftParton() const {return *thePartons.begin();}
88 G4QParton* GetRightParton() const {return *(thePartons.end()-1);}
89 G4int GetDirection() const {return theDirection;}
92 G4int GetCharge() const {return GetQC().GetCharge();}
96 G4bool DecayIsQuark() const {return theDecayParton->GetType()==1;}
97 G4bool StableIsQuark() const {return theStableParton->GetType()==1;}
98 G4ThreeVector DecayPt(); // Get Pt of the decaying quark @@ Called once
99 G4double Mass2() const { return Pplus*Pminus-(Ptleft+Ptright).mag2();}
100 G4double Mass() const // @@ Very dangerous! USE ONLY FORE THE LIGHT CONE ALGORITHM !!
101 {
102 G4double mass2=Mass2();
103 if(mass2>0) return std::sqrt(Mass2());
104 else return 0.; // @@ Make Warning
105 }
107 {
108 G4LorentzVector mom(Ptleft+Ptright, 0.5*(Pplus+Pminus));
109 mom.setPz(0.5*(Pplus-Pminus));
110 return FragmentationMass(1) + MassCut > mom.mag();
111 }
112 G4bool IsFragmentable() {return FragmentationMass() + MassCut < Mass();} // @@ Mass() ?
113 G4ThreeVector SampleQuarkPt(); // @@ Called once
114 G4QHadron* CreateHadron(G4QParton* black, G4QParton* white);
117
118 // Modifiers
119 void SetPosition(const G4ThreeVector& aPosition){thePosition= aPosition;}
120 void SetDirection(G4int dir) {if(dir==1 || dir==-1) theDirection=dir;}
121 void SetLeftParton(G4QParton* LP) {thePartons[0]=LP;} // !! Not deleting the substituty
122 void SetRightParton(G4QParton* RP) {thePartons.pop_back(); thePartons.push_back(RP);}
123 void KillString() {theDirection=0;} // @@ Can be absolete
124 void LorentzRotate(const G4LorentzRotation& rotation);
125 //void InsertParton(G4QParton* aParton, const G4QParton* addafter = NULL);
126 void Boost(G4ThreeVector& Velocity);
128 //void DiffString(G4QHadron* aHadron, G4bool isProjectile); @@ Mast be used!!
129 void ExciteString(G4QParton* Col,G4QParton* AntiCol, G4int Dir);
130 G4QHadronVector* FragmentString(G4bool QL); // Fragment String using QGSM=true/LUND=false
132 G4double FragmentationMass(G4int HighSpin = 0, G4QHadronPair* pdefs = 0);
133 void SetLeftPartonStable();
136 G4LorentzVector* SplitEandP(G4QHadron* pHadron, G4bool QL); // QGSM:QL=true,Lund:QL=false
137 G4QPartonPair CreatePartonPair(G4int NeedParticle, G4bool AllowDiquarks=true);
138 G4QHadron* QuarkSplitup(G4QParton* decay, G4QParton* &created);
139 G4QHadron* DiQuarkSplitup(G4QParton* decay, G4QParton* &created);
140 G4int SampleQuarkFlavor() {return (1+G4int(G4UniformRand()/StrangeSuppress));} // ? M.K.
141
142 // Static functions
143 static void SetParameters(G4double mCut, G4double sigQT, G4double DQSup, G4double DQBU,
144 G4double smPar, G4double SSup, G4double SigPt);
145
146 private:
147 enum Spin {SpinZero=1, SpinHalf=2, SpinOne=3, SpinThreeHalf=4};
148 // Private functions
149 G4QHadron* CreateMeson(G4QParton* black, G4QParton* white, Spin spin);
150 G4QHadron* CreateBaryon(G4QParton* black,G4QParton* white, Spin spin);
151 G4ThreeVector GaussianPt(G4double widthSquare, G4double maxPtSquare) const;
152 G4double GetLundLightConeZ(G4double zmin, G4double zmax, G4int PartonEncoding,
153 G4QHadron* pHadron, G4double Px, G4double Py);
154 G4double GetQGSMLightConeZ(G4double zmin, G4double zmax, G4int PartonEncoding,
155 G4QHadron* pHadron, G4double Px, G4double Py);
156
157 // Static parameters
158 // Parameters of Longitudinal String Fragmentation
159 static G4double MassCut; // minimum mass cut for the string
160 static G4double SigmaQT; // quark transverse momentum distribution parameter
161 static G4double DiquarkSuppress; // is Diquark suppression parameter
162 static G4double DiquarkBreakProb; // is Diquark breaking probability
163 static G4double SmoothParam; // QGS model parameter
164 static G4double StrangeSuppress; // Strangeness suppression parameter
165 static G4double widthOfPtSquare; // width^2 of pt for string excitation
166
167 // Body
168 G4int theDirection; // must be 1 (PROJECTILE) or -1 (TARGET), 0 - DEAD
169 G4ThreeVector thePosition; // Defined by the first quark position
170 G4QPartonVector thePartons; // Partons on the ends of the string @@ Use PartonPair
171 G4ThreeVector Ptleft,Ptright; // Pt (px,py) for partons (pz ignored!)
172 G4double Pplus, Pminus; // p-, p+ of string, Plus is assigned to Left!
173 G4QParton* theStableParton; // Parton on the stable side of the string
174 G4QParton* theDecayParton; // Parton on the decaying part of the string
175 enum DecaySide {None, Left, Right};// @@ To have two @@ Leav : 1=Left
176 DecaySide decaying; // @@ it's too much @@ only : 0=Unknown
177 G4int SideOfDecay; // @@ of a good thing @@ one! : -1=Right
178};
179
180#endif
std::vector< G4QHadron * > G4QHadronVector
std::pair< G4QHadron *, G4QHadron * > G4QHadronPair
Definition: G4QHadron.hh:165
std::vector< G4QParton * > G4QPartonVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4UniformRand()
Definition: Randomize.hh:53
G4int GetCharge() const
Definition: G4QContent.cc:1159
G4int GetBaryonNumber() const
Definition: G4QContent.cc:1182
G4int GetStrangeness() const
Definition: G4QContent.hh:184
G4QContent GetQC() const
Definition: G4QParton.hh:82
const G4int & GetType() const
Definition: G4QParton.hh:88
G4QHadron * CreateHadron(G4QParton *black, G4QParton *white)
Definition: G4QString.cc:1223
G4int operator!=(const G4QString &right) const
Definition: G4QString.hh:84
G4bool StopFragmentation()
Definition: G4QString.hh:106
G4QHadron * CreateLowSpinHadron(G4QParton *black, G4QParton *white)
Definition: G4QString.cc:1256
G4QHadron * Splitup(G4bool QL)
Definition: G4QString.cc:877
G4ThreeVector DecayPt()
Definition: G4QString.cc:864
G4bool IsFragmentable()
Definition: G4QString.hh:112
G4QHadronVector * FragmentString(G4bool QL)
Definition: G4QString.cc:348
void LorentzRotate(const G4LorentzRotation &rotation)
Definition: G4QString.cc:131
G4QParton * GetLeftParton() const
Definition: G4QString.hh:87
G4int GetCharge() const
Definition: G4QString.hh:92
void SetLeftPartonStable()
Definition: G4QString.cc:826
G4int GetDirection() const
Definition: G4QString.hh:89
void SetLeftParton(G4QParton *LP)
Definition: G4QString.hh:121
G4LorentzVector * SplitEandP(G4QHadron *pHadron, G4bool QL)
Definition: G4QString.cc:959
G4QPartonPair CreatePartonPair(G4int NeedParticle, G4bool AllowDiquarks=true)
Definition: G4QString.cc:1072
G4int GetStrangeness() const
Definition: G4QString.hh:94
G4QString(G4QParton *newdecay, const G4LorentzVector *momentum)
G4LorentzRotation TransformToAlignedCms()
void SetRightPartonStable()
Definition: G4QString.cc:833
void SetDirection(G4int dir)
Definition: G4QString.hh:120
void SetRightParton(G4QParton *RP)
Definition: G4QString.hh:122
void SetPosition(const G4ThreeVector &aPosition)
Definition: G4QString.hh:119
static void SetParameters(G4double mCut, G4double sigQT, G4double DQSup, G4double DQBU, G4double smPar, G4double SSup, G4double SigPt)
Definition: G4QString.cc:181
G4bool StableIsQuark() const
Definition: G4QString.hh:97
G4LorentzVector Get4Momentum() const
Definition: G4QString.cc:124
G4int SampleQuarkFlavor()
Definition: G4QString.hh:140
G4QHadron * QuarkSplitup(G4QParton *decay, G4QParton *&created)
Definition: G4QString.cc:1007
G4QHadron * CreateHighSpinHadron(G4QParton *black, G4QParton *white)
Definition: G4QString.cc:1281
G4ThreeVector SampleQuarkPt()
Definition: G4QString.cc:1000
G4int GetBaryonNumber() const
Definition: G4QString.hh:93
G4QParton * GetRightParton() const
Definition: G4QString.hh:88
G4double Mass() const
Definition: G4QString.hh:100
G4double FragmentationMass(G4int HighSpin=0, G4QHadronPair *pdefs=0)
Definition: G4QString.cc:769
G4int GetDecayDirection() const
Definition: G4QString.cc:840
G4QContent GetQC() const
Definition: G4QString.hh:91
void KillString()
Definition: G4QString.hh:123
G4int operator==(const G4QString &right) const
Definition: G4QString.hh:83
const G4ThreeVector & GetPosition() const
Definition: G4QString.hh:85
G4bool DecayIsQuark() const
Definition: G4QString.hh:96
G4QString(const G4QString &old, G4QParton *newdecay, const G4LorentzVector *momentum)
G4double Mass2() const
Definition: G4QString.hh:99
G4QHadron * DiQuarkSplitup(G4QParton *decay, G4QParton *&created)
Definition: G4QString.cc:1022
G4QString(const G4QString &old, G4QParton *newdecay)
const G4QPartonVector * GetPartonList() const
Definition: G4QString.hh:86
@ PROJECTILE
Definition: G4QString.hh:69
void Boost(G4ThreeVector &Velocity)
Definition: G4QString.cc:152
G4QHadronVector * LightFragmentationTest()
Definition: G4QString.cc:702
void ExciteString(G4QParton *Col, G4QParton *AntiCol, G4int Dir)
Definition: G4QString.cc:254
#define const
Definition: zconf.h:118