File indexing completed on 2025-07-05 08:55:10
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029 #ifndef G4FragmentingString_h
0030 #define G4FragmentingString_h 1
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040 #include "G4ios.hh"
0041 #include "globals.hh"
0042 #include "G4ThreeVector.hh"
0043 #include "G4LorentzVector.hh"
0044 #include "G4LorentzRotation.hh"
0045 #include "G4ParticleDefinition.hh"
0046
0047 class G4ExcitedString;
0048
0049 class G4FragmentingString
0050 {
0051 public:
0052 G4FragmentingString(const G4FragmentingString &right);
0053 G4FragmentingString(const G4ExcitedString &excited);
0054 G4FragmentingString(const G4FragmentingString &old,
0055 G4ParticleDefinition * newdecay,
0056 const G4LorentzVector *momentum);
0057 G4FragmentingString(const G4FragmentingString &old,
0058 G4ParticleDefinition * newdecay);
0059
0060 ~G4FragmentingString();
0061
0062 G4FragmentingString& operator=(const G4FragmentingString &);
0063 G4bool operator==(const G4FragmentingString &right) const;
0064
0065 G4bool operator!=(const G4FragmentingString &right) const;
0066
0067 G4LorentzVector Get4Momentum() const;
0068
0069 G4ThreeVector StablePt();
0070 G4ThreeVector DecayPt();
0071
0072 G4double LightConePlus();
0073 G4double LightConeMinus();
0074 G4double LightConeDecay();
0075
0076 G4double Mass() const;
0077 G4double Mass2() const;
0078 G4double MassT2() const;
0079
0080 G4ParticleDefinition* GetLeftParton(void) const;
0081 G4ParticleDefinition* GetRightParton(void) const;
0082
0083 G4ParticleDefinition* GetStableParton() const;
0084 G4ParticleDefinition* GetDecayParton() const;
0085
0086 void SetLeftPartonStable();
0087 void SetRightPartonStable();
0088
0089 G4int GetDecayDirection() const;
0090
0091 G4bool DecayIsQuark();
0092 G4bool StableIsQuark();
0093 G4bool IsAFourQuarkString(void) const;
0094
0095 G4LorentzVector GetPstring();
0096 G4LorentzVector GetPleft();
0097 void SetPleft(G4LorentzVector a4momentum);
0098 G4LorentzVector GetPright();
0099 void SetPright(G4LorentzVector a4momentum);
0100 void LorentzRotate(const G4LorentzRotation & rotation);
0101 G4LorentzRotation TransformToCenterOfMass();
0102 G4LorentzRotation TransformToAlignedCms();
0103 void Boost(G4ThreeVector& Velocity);
0104
0105 private:
0106 G4ParticleDefinition *LeftParton, *RightParton;
0107 G4ThreeVector Ptleft,Ptright;
0108 G4double Pplus, Pminus;
0109
0110 G4ParticleDefinition * theStableParton, * theDecayParton;
0111
0112 G4LorentzVector Pstring, Pleft, Pright;
0113 enum DecaySide { None, Left, Right };
0114 DecaySide decaying;
0115 };
0116
0117 inline
0118 G4bool G4FragmentingString::operator==(const G4FragmentingString &right) const
0119 {
0120 return this == &right;
0121 }
0122
0123 inline
0124 G4bool G4FragmentingString::operator!=(const G4FragmentingString &right) const
0125 {
0126 return this != &right;
0127 }
0128
0129
0130 inline
0131 G4ParticleDefinition * G4FragmentingString::GetStableParton() const
0132 {
0133 return theStableParton;
0134 }
0135
0136 inline
0137 G4ParticleDefinition * G4FragmentingString::GetDecayParton() const
0138 {
0139 return theDecayParton;
0140 }
0141
0142 inline
0143 G4ParticleDefinition* G4FragmentingString::GetLeftParton(void) const
0144 {
0145 return LeftParton;
0146 }
0147
0148 inline
0149 G4ParticleDefinition* G4FragmentingString::GetRightParton(void) const
0150 {
0151 return RightParton;
0152 }
0153
0154
0155 inline
0156 void G4FragmentingString::LorentzRotate(const G4LorentzRotation & rotation)
0157 {
0158 SetPleft(rotation*Pleft);
0159 SetPright(rotation*Pright);
0160 Pstring = Pleft+Pright;
0161 Ptleft =Pleft.vect(); Ptleft.setZ(0.);
0162 Ptright=Pright.vect(); Ptright.setZ(0.);
0163 Pplus =Pstring.plus();
0164 Pminus=Pstring.minus();
0165 }
0166
0167 inline
0168 G4LorentzRotation G4FragmentingString::TransformToCenterOfMass()
0169 {
0170 G4LorentzVector momentum=Pstring;
0171 G4LorentzRotation toCMS(-1*momentum.boostVector());
0172
0173 Pleft *= toCMS;
0174 Pright *= toCMS;
0175 Pstring *= toCMS;
0176 Ptleft =Pleft.vect(); Ptleft.setZ(0.);
0177 Ptright=Pright.vect(); Ptright.setZ(0.);
0178 Pplus =Pstring.plus();
0179 Pminus=Pstring.minus();
0180 return toCMS;
0181 }
0182
0183 inline
0184 void G4FragmentingString::SetPleft(G4LorentzVector a4momentum)
0185 {
0186 Pleft = a4momentum;
0187 Ptleft = Pleft.vect(); Ptleft.setZ(0.);
0188 Pstring = Pleft + Pright;
0189 Pplus = Pstring.plus();
0190 Pminus = Pstring.minus();
0191 }
0192
0193 inline
0194 void G4FragmentingString::SetPright(G4LorentzVector a4momentum)
0195 {
0196 Pright = a4momentum;
0197 Ptright = Pright.vect(); Ptright.setZ(0.);
0198 Pstring = Pleft + Pright;
0199 Pplus = Pstring.plus();
0200 Pminus = Pstring.minus();
0201 }
0202
0203 #endif
0204