source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/body/include/G4QString.hh @ 1340

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

update ti head

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