Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4FragmentingString.cc
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 
29 // ------------------------------------------------------------
30 // GEANT 4 class implementation file
31 //
32 // ---------------- G4FragmentingString ----------------
33 // by Gunter Folger, September 2001.
34 // class for an excited string used in Fragmention
35 // ------------------------------------------------------------
36 
37 
38 // G4FragmentingString
39 #include "G4FragmentingString.hh"
40 #include "G4ExcitedString.hh"
41 
42 //---------------------------------------------------------------------------------
43 
45 {
46  LeftParton=old.LeftParton;
47  RightParton=old.RightParton;
48  Ptleft=old.Ptleft;
49  Ptright=old.Ptright;
50  Pplus=old.Pplus;
51  Pminus=old.Pminus;
52  theStableParton=old.theStableParton;
53  theDecayParton=old.theDecayParton;
54  decaying=old.decaying;
55 }
56 
58 {
59  if (this != &old)
60  {
61  LeftParton=old.LeftParton;
62  RightParton=old.RightParton;
63  Ptleft=old.Ptleft;
64  Ptright=old.Ptright;
65  Pplus=old.Pplus;
66  Pminus=old.Pminus;
67  theStableParton=old.theStableParton;
68  theDecayParton=old.theDecayParton;
69  decaying=old.decaying;
70  }
71  return *this;
72 }
73 
74 //---------------------------------------------------------------------------------
75 
77 {
78  LeftParton=excited.GetLeftParton()->GetDefinition();
79  RightParton=excited.GetRightParton()->GetDefinition();
80  Ptleft=excited.GetLeftParton()->Get4Momentum().vect();
81  Ptleft.setZ(0.);
82  Ptright=excited.GetRightParton()->Get4Momentum().vect();
83  Ptright.setZ(0.);
84  G4LorentzVector P=excited.Get4Momentum();
85  Pplus =P.e() + P.pz();
86  Pminus=P.e() - P.pz();
87  theStableParton=0;
88  theDecayParton=0;
89  //decaying=None;
90  if(excited.GetDirection() > 0) {decaying=Left; }
91  else {decaying=Right;}
92 }
93 
94 //---------------------------------------------------------------------------------
95 
97  G4ParticleDefinition * newdecay,
98  const G4LorentzVector *momentum)
99 {
100  decaying=None;
101  if ( old.decaying == Left )
102  {
103  RightParton= old.RightParton;
104  Ptright = old.Ptright;
105  LeftParton = newdecay;
106  Ptleft = old.Ptleft - momentum->vect();
107  Ptleft.setZ(0.);
108  theDecayParton=GetLeftParton();
109  theStableParton=GetRightParton();
110  decaying=Left;
111  } else if ( old.decaying == Right )
112  {
113  RightParton = newdecay;
114  Ptright = old.Ptright - momentum->vect();
115  Ptright.setZ(0.);
116  LeftParton = old.LeftParton;
117  Ptleft = old.Ptleft;
118  theDecayParton=GetRightParton();
119  theStableParton=GetLeftParton();
120  decaying=Right;
121  } else
122  {
123  throw G4HadronicException(__FILE__, __LINE__,
124  "G4FragmentingString::G4FragmentingString: no decay Direction defined");
125  }
126  Pplus = old.Pplus - (momentum->e() + momentum->pz());
127  Pminus = old.Pminus - (momentum->e() - momentum->pz());
128 
129  //G4double Eold=0.5 * (old.Pplus + old.Pminus);
130  //G4double Enew=0.5 * (Pplus + Pminus);
131 }
132 
133 //---------------------------------------------------------------------------------
134 
136  G4ParticleDefinition * newdecay)
137 {
138  decaying=None;
139 
140  Ptleft.setX(0.); Ptleft.setY(0.); Ptleft.setZ(0.);
141  Ptright.setX(0.); Ptright.setY(0.); Ptright.setZ(0.);
142  Pplus=0.; Pminus=0.;
143  theStableParton=0; theDecayParton=0;
144 
145  if ( old.decaying == Left )
146  {
147  RightParton= old.RightParton;
148  LeftParton = newdecay;
149  decaying=Left;
150  } else if ( old.decaying == Right )
151  {
152  RightParton = newdecay;
153  LeftParton = old.LeftParton;
154  decaying=Right;
155  } else
156  {
157  throw G4HadronicException(__FILE__, __LINE__,
158  "G4FragmentingString::G4FragmentingString: no decay Direction defined");
159  }
160 }
161 
162 //---------------------------------------------------------------------------------
163 
165 {}
166 
167 //---------------------------------------------------------------------------------
168 
170 {
171  theStableParton=GetLeftParton();
172  theDecayParton=GetRightParton();
173  decaying=Right;
174 }
175 
176 //---------------------------------------------------------------------------------
177 
179 {
180  theStableParton=GetRightParton();
181  theDecayParton=GetLeftParton();
182  decaying=Left;
183 }
184 
185 //---------------------------------------------------------------------------------
186 
188 {
189  if (decaying == Left ) return +1;
190  else if (decaying == Right) return -1;
191  else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::GetDecayDirection: decay side UNdefined!");
192  return 0;
193 }
194 
195 //---------------------------------------------------------------------------------
196 
198 {
199  return LeftParton->GetParticleSubType()== "di_quark"
200  && RightParton->GetParticleSubType()== "di_quark";
201 }
202 
203 //---------------------------------------------------------------------------------
204 
206 {
207  return theDecayParton->GetParticleSubType()== "quark";
208 }
209 
211 {
212  return theStableParton->GetParticleSubType()== "quark";
213 }
214 
215 //---------------------------------------------------------------------------------
216 
218 {
219  if (decaying == Left ) return Ptright;
220  else if (decaying == Right ) return Ptleft;
221  else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
222  return G4ThreeVector();
223 }
224 
226 {
227  if (decaying == Left ) return Ptleft;
228  else if (decaying == Right ) return Ptright;
229  else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
230  return G4ThreeVector();
231 }
232 
233 //---------------------------------------------------------------------------------
234 
236 {
237  return Pplus;
238 }
239 
241 {
242  return Pminus;
243 }
244 
246 {
247  if (decaying == Left ) return Pplus;
248  else if (decaying == Right ) return Pminus;
249  else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
250 }
251 
252 //---------------------------------------------------------------------------------
253 
255 {
256  G4LorentzVector momentum(Ptleft+Ptright,0);
257  momentum.setPz(0.5*(Pplus-Pminus));
258  momentum.setE(0.5*(Pplus+Pminus));
259  return momentum;
260 }
261 
263 {
264  return Pplus*Pminus - (Ptleft+Ptright).mag2();
265 }
266 
268 {
269  return std::sqrt(this->Mass2());
270 }
271 
273 {
274  return Pplus*Pminus;
275 }
276 
G4FragmentingString(const G4FragmentingString &right)
G4ParticleDefinition * GetRightParton(void) const
CLHEP::Hep3Vector G4ThreeVector
const G4LorentzVector & Get4Momentum() const
Definition: G4Parton.hh:143
G4Parton * GetLeftParton(void) const
G4double Mass() const
G4ParticleDefinition * GetDefinition()
Definition: G4Parton.hh:161
const G4String & GetParticleSubType() const
int G4int
Definition: G4Types.hh:78
void setY(double)
void setZ(double)
void setX(double)
static double P[]
G4ParticleDefinition * GetLeftParton(void) const
G4double Mass2() const
Hep3Vector vect() const
G4LorentzVector Get4Momentum() const
bool G4bool
Definition: G4Types.hh:79
G4int GetDecayDirection() const
G4bool FourQuarkString(void) const
G4LorentzVector Get4Momentum() const
G4Parton * GetRightParton(void) const
G4FragmentingString & operator=(const G4FragmentingString &)
G4double MassT2() const
double pz() const
G4int GetDirection(void) const
double G4double
Definition: G4Types.hh:76