source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/body/include/G4QString.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: 10.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//
28// $Id: G4QString.hh,v 1.3 2009/02/23 09:49:24 mkossov Exp $
29// GEANT4 tag $Name: geant4-09-03-beta-cand-01 $
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 "G4QHadronBuilder.hh"
63#include "G4QPartonPair.hh"
64#include "G4QPartonVector.hh"
65#include <algorithm>
66
67class G4QString 
68{
69 public:
70
71  enum {PROJECTILE  = 1, TARGET  = -1}; // The same as in quark-pair (@@ ? M.K.)
72
73  G4QString(); // formal creation of the string with future excitation
74  G4QString(G4QParton* Color,G4QParton* Gluon, G4QParton* AntiColor, G4int Dir=PROJECTILE);
75  G4QString(G4QParton* Col,G4QParton* AntiCol, G4int Dir=PROJECTILE);
76  G4QString(const G4QString &right);
77
78  ~G4QString();
79
80  // Selectors
81  G4int operator==(const G4QString &right) const {return this == &right;}
82  G4int operator!=(const G4QString &right) const {return this != &right;}
83  const G4ThreeVector& GetPosition() const       {return thePosition;}
84  const G4QPartonVector* GetPartonList() const   {return &thePartons;}
85  G4QParton* GetGluon() const                    {return thePartons[1];}
86  G4QParton* GetGluon(G4int GluonPos) const      {return thePartons[1 + GluonPos];}
87  G4QParton* GetLeftParton() const               {return *thePartons.begin();}
88  G4QParton* GetRightParton() const              {return *(thePartons.end()-1);}
89  G4bool     IsItKinkyString() const             {return (thePartons.size() > 2);}
90  G4int      GetDirection() const                {return theDirection;}
91  G4bool     IsExcited() const                   {return true;} // Can be absolete? @@ M.K.
92  G4LorentzVector Get4Momentum() const;
93  //G4QHadron* GetAsQHadron(); //@@ Create a QHadron, corresponding to not excited string
94  G4QParton* GetColorParton() const;
95  G4QParton* GetAntiColorParton() const;
96  G4QParton* GetStableParton() const{ return theStableParton;} // stable at the moment
97  G4QParton* GetDecayParton() const{return theDecayParton;}// now involved in fragmentation
98  G4int GetDecayDirection() const;
99  G4bool FourQuarkString() const; // Checks that the string is qq-(qq)bar string
100  G4bool DecayIsQuark() const {return theDecayParton->GetParticleSubType()=="quark";}
101  G4bool StableIsQuark() const {return theStableParton->GetParticleSubType()=="quark";}
102  G4ThreeVector StablePt(); // Get Pt of the stable quark
103  G4ThreeVector DecayPt();  // Get Pt of the decaying quark
104  G4double LightConePlus(){return Pplus;}
105  G4double LightConeMinus() {return Pminus;}
106  G4double LightConeDecay();
107  G4LorentzVector GetFragmentation4Mom() const; //!! Instead of Get$Momentum in FragmString
108  G4double Mass2() const {return Pplus*Pminus-(Ptleft+Ptright).mag2();}
109  G4double Mass() const{return std::sqrt(Mass2());}
110  G4bool StopFragmentation() {return FragmentationMass(&G4QHadronBuilder::BuildHighSpin)
111                                                 + MassCut > GetFragmentation4Mom().mag();}
112  G4bool IsFragmentable() {return FragmentationMass() + MassCut < Mass();}
113  G4ThreeVector SampleQuarkPt();
114  G4bool SplitLast(G4QString* string, G4QHadronVector* LeftHV, G4QHadronVector* RightHV);
115
116  void Sample4Momentum(G4LorentzVector* Mom, G4double M, G4LorentzVector* aMom,
117                       G4double aM, G4double InitialMass);
118  // Modifiers
119  void SetPosition(const G4ThreeVector& aPosition) {thePosition= aPosition;}
120  void SetDirection(G4int dir)                     {theDirection=dir;}
121  void LorentzRotate(const G4LorentzRotation& rotation);
122  void InsertParton(G4QParton* aParton, const G4QParton* addafter = NULL);
123  void Boost(G4ThreeVector& Velocity);
124  G4LorentzRotation TransformToCenterOfMass();
125  G4LorentzRotation TransformToAlignedCms();
126  void DiffString(G4QHadron* aHadron, G4bool isProjectile);
127  void ExciteString(G4QParton* Col,G4QParton* AntiCol, G4int Dir);
128  G4QHadronVector* FragmentString(G4bool QL); // Fragment String using QGSM=true/LUND=false
129  G4QString* CPExcited(); // Creates a string, using only CMS end-partons of the string
130  G4QHadronVector* LightFragmentationTest();
131  G4double FragmentationMass(G4QHcreate build=0, G4QHadronPair* pdefs=0);
132  void CalculateHadronTimePosition(G4double theInitialStringMass, G4QHadronVector*);
133  void SetLeftPartonStable();
134  void SetRightPartonStable();
135  void UpdateString(G4QParton* decay, const G4LorentzVector* mom);
136  G4QString(G4QParton* newdecay, const G4LorentzVector* momentum);
137  G4QHadron* Splitup(G4bool QL);  // Split Hadron & update the string, QGSM:true,Lund:false
138  G4LorentzVector* SplitEandP(G4QHadron* pHadron, G4bool QL); // QGSM:QL=true,Lund:QL=false
139  G4QParton* CreateParton(G4int PDGcode)                   {return new G4QParton(PDGcode);}
140  G4QPartonPair CreatePartonPair(G4int NeedParticle, G4bool AllowDiquarks=true);
141  G4QHadron* QuarkSplitup(G4QParton* decay, G4QParton* &created);
142  G4QHadron* DiQuarkSplitup(G4QParton* decay, G4QParton* &created);
143  G4int SampleQuarkFlavor() {return (1+G4int(G4UniformRand()/StrangeSuppress));} // ? M.K.
144
145  // Static functions
146  static void SetParameters(G4double mCut, G4double clustM, G4double sigQT, G4double DQSup,
147   G4double DQBU, G4double smPar, G4double SSup, G4double SigPt, G4int SLmax, G4int CLmax);
148
149 private:
150  // Private functions
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 ClusterMass;       // minimum cluster mass for the fragmentation
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  static G4int StringLoopInterrupt;  // String fragmentation LOOP limit
168  static G4int ClusterLoopInterrupt; // Cluster fragmentation LOOP limit
169
170  // Body
171  G4int         theDirection;        // must be 1 or -1 (PROJECTILE or TARGET)
172  G4ThreeVector thePosition;         // Defined by the first quark position
173  G4QPartonVector thePartons;        // would like initial capacity for 3 Partons (? M.K.)
174  G4QHadronBuilder* hadronizer;      // Hadronizer of hodrons out of partons
175  G4ThreeVector Ptleft,Ptright;      // Pt (px,py) for partons (pz ignored!)
176  G4double Pplus, Pminus;            // p-, p+ of string, Plus is assigned to Left!
177  G4QParton* theStableParton;        // Parton on the stable side of the string
178  G4QParton* theDecayParton;         // Parton on the decaying part of the string
179  enum DecaySide {None, Left, Right};// @@ To have two         @@ Leav   :  1=Left
180  DecaySide decaying;                // @@   it's too much     @@  only  :  0=Unknown
181  G4int     SideOfDecay;             // @@     of a good thing @@   one! : -1=Right
182};
183
184#endif
Note: See TracBrowser for help on using the repository browser.