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

Last change on this file since 836 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 9.3 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.2 2006/12/12 11:02:22 mkossov Exp $
29// GEANT4 tag $Name:  $
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
46#include "G4ios.hh"
47#include "globals.hh"
48#include "G4ThreeVector.hh"
49#include "G4LorentzVector.hh"
50#include "G4LorentzRotation.hh"
51#include "G4QHadronVector.hh"
52#include "G4QHadronBuilder.hh"
53#include "G4QPartonPair.hh"
54#include "G4QPartonVector.hh"
55#include <algorithm>
56
57class G4QString 
58{
59 public:
60
61                enum {PROJECTILE  = 1, TARGET  = -1}; // The same as in quark-pair (@@ ? M.K.)
62
63  G4QString(); // formal creation of the string with future excitation
64  G4QString(G4QParton* Color,G4QParton* Gluon, G4QParton* AntiColor, G4int Dir=PROJECTILE);
65  G4QString(G4QParton* Col,G4QParton* AntiCol, G4int Dir=PROJECTILE);
66  G4QString(const G4QString &right);
67
68  ~G4QString();
69
70  // Selectors
71  G4int operator==(const G4QString &right) const {return this == &right;}
72  G4int operator!=(const G4QString &right) const {return this != &right;}
73  const G4ThreeVector& GetPosition() const       {return thePosition;}
74  const G4QPartonVector* GetPartonList() const   {return &thePartons;}
75  G4QParton* GetGluon() const                    {return thePartons[1];}
76  G4QParton* GetGluon(G4int GluonPos) const      {return thePartons[1 + GluonPos];}
77  G4QParton* GetLeftParton() const               {return *thePartons.begin();}
78  G4QParton* GetRightParton() const              {return *(thePartons.end()-1);}
79  G4bool     IsItKinkyString() const             {return (thePartons.size() > 2);}
80  G4int      GetDirection() const                {return theDirection;}
81  G4bool     IsExcited() const                   {return true;} // Can be absolete? @@ M.K.
82  G4LorentzVector Get4Momentum() const;
83  //G4QHadron* GetAsQHadron(); //@@ Create a QHadron, corresponding to not excited string
84  G4QParton* GetColorParton() const;
85  G4QParton* GetAntiColorParton() const;
86  G4QParton* GetStableParton() const{ return theStableParton;} // stable at the moment
87  G4QParton* GetDecayParton() const{return theDecayParton;}// now involved in fragmentation
88  G4int GetDecayDirection() const;
89  G4bool FourQuarkString() const; // Checks that the string is qq-(qq)bar string
90  G4bool DecayIsQuark() const {return theDecayParton->GetParticleSubType()=="quark";}
91  G4bool StableIsQuark() const {return theStableParton->GetParticleSubType()=="quark";}
92  G4ThreeVector StablePt(); // Get Pt of the stable quark
93  G4ThreeVector DecayPt();  // Get Pt of the decaying quark
94  G4double LightConePlus(){return Pplus;}
95  G4double LightConeMinus() {return Pminus;}
96  G4double LightConeDecay();
97  G4LorentzVector GetFragmentation4Mom() const; //!! Instead of Get$Momentum in FragmString
98  G4double Mass2() const {return Pplus*Pminus-(Ptleft+Ptright).mag2();}
99  G4double Mass() const{return std::sqrt(Mass2());}
100  G4bool StopFragmentation() {return FragmentationMass(&G4QHadronBuilder::BuildHighSpin)
101                                                 + MassCut > GetFragmentation4Mom().mag();}
102  G4bool IsFragmentable() {return FragmentationMass() + MassCut < Mass();}
103  G4ThreeVector SampleQuarkPt();
104  G4bool SplitLast(G4QString* string, G4QHadronVector* LeftHV, G4QHadronVector* RightHV);
105
106  void Sample4Momentum(G4LorentzVector* Mom, G4double M, G4LorentzVector* aMom,
107                       G4double aM, G4double InitialMass);
108  // Modifiers
109  void SetPosition(const G4ThreeVector& aPosition) {thePosition= aPosition;}
110  void SetDirection(G4int dir)                     {theDirection=dir;}
111  void LorentzRotate(const G4LorentzRotation& rotation);
112  void InsertParton(G4QParton* aParton, const G4QParton* addafter = NULL);
113  void Boost(G4ThreeVector& Velocity);
114  G4LorentzRotation TransformToCenterOfMass();
115  G4LorentzRotation TransformToAlignedCms();
116  void DiffString(G4QHadron* aHadron, G4bool isProjectile);
117  void ExciteString(G4QParton* Col,G4QParton* AntiCol, G4int Dir);
118  G4QHadronVector* FragmentString(G4bool QL); // Fragment String using QGSM=true/LUND=false
119  G4QString* CPExcited(); // Creates a string, using only CMS end-partons of the string
120  G4QHadronVector* LightFragmentationTest();
121  G4double FragmentationMass(G4QHcreate build=0, G4QHadronPair* pdefs=0);
122  void CalculateHadronTimePosition(G4double theInitialStringMass, G4QHadronVector*);
123  void SetLeftPartonStable();
124  void SetRightPartonStable();
125  void UpdateString(G4QParton* decay, const G4LorentzVector* mom);
126  G4QString(G4QParton* newdecay, const G4LorentzVector* momentum);
127  G4QHadron* Splitup(G4bool QL);  // Split Hadron & update the string, QGSM:true,Lund:false
128  G4LorentzVector* SplitEandP(G4QHadron* pHadron, G4bool QL); // QGSM:QL=true,Lund:QL=false
129  G4QParton* CreateParton(G4int PDGcode)                   {return new G4QParton(PDGcode);}
130  G4QPartonPair CreatePartonPair(G4int NeedParticle, G4bool AllowDiquarks=true);
131  G4QHadron* QuarkSplitup(G4QParton* decay, G4QParton* &created);
132  G4QHadron* DiQuarkSplitup(G4QParton* decay, G4QParton* &created);
133  G4int SampleQuarkFlavor() {return (1+G4int(G4UniformRand()/StrangeSuppress));} // ? M.K.
134
135  // Static functions
136  static void SetParameters(G4double mCut, G4double clustM, G4double sigQT, G4double DQSup,
137   G4double DQBU, G4double smPar, G4double SSup, G4double SigPt, G4int SLmax, G4int CLmax);
138
139 private:
140  // Private functions
141  G4ThreeVector GaussianPt(G4double widthSquare, G4double maxPtSquare) const;
142  G4double GetLundLightConeZ(G4double zmin, G4double zmax, G4int PartonEncoding,
143                             G4QHadron* pHadron, G4double Px, G4double Py);
144  G4double GetQGSMLightConeZ(G4double zmin, G4double zmax, G4int PartonEncoding,
145                             G4QHadron* pHadron, G4double Px, G4double Py);
146
147  // Static parameters
148  // Parameters of Longitudinal String Fragmentation
149  static G4double MassCut;           // minimum mass cut for the string
150  static G4double ClusterMass;       // minimum cluster mass for the fragmentation
151  static G4double SigmaQT;           // quark transverse momentum distribution parameter
152  static G4double DiquarkSuppress;   // is Diquark suppression parameter 
153  static G4double DiquarkBreakProb;  // is Diquark breaking probability
154  static G4double SmoothParam;       // QGS model parameter
155  static G4double StrangeSuppress;   // Strangeness suppression parameter
156         static G4double widthOfPtSquare;         // width^2 of pt for string excitation
157  static G4int StringLoopInterrupt;  // String fragmentation LOOP limit
158  static G4int ClusterLoopInterrupt; // Cluster fragmentation LOOP limit
159
160  // Body
161  G4int         theDirection;  // must be 1 or -1 (PROJECTILE or TARGET)
162  G4ThreeVector thePosition;   // Defined by the first quark position
163  G4QPartonVector thePartons;  // would like initial capacity for 3 Partons (? M.K.)
164  G4QHadronBuilder* hadronizer;// Hadronizer of hodrons out of partons
165  G4ThreeVector Ptleft,Ptright;// Pt (px,py) for partons (pz ignored!)
166  G4double Pplus, Pminus;      // p-, p+ of string, Plus is assigned to Left!
167  G4QParton* theStableParton;  // Parton on the stable side of the string
168  G4QParton* theDecayParton;   // Parton on the decaying part of the string
169  enum DecaySide {None, Left, Right}; // @@ To have two         @@ Leav   :  1=Left
170  DecaySide decaying;                 // @@   it's too much     @@  only  :  0=Unknown
171  G4int     SideOfDecay;              // @@     of a good thing @@   one! : -1=Right
172};
173
174#endif
Note: See TracBrowser for help on using the repository browser.