ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4FragmentingString.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4FragmentingString.cc
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 
44 //---------------------------------------------------------------------------------
45 
47 {
50  Ptleft=old.Ptleft;
51  Ptright=old.Ptright;
52  Pplus=old.Pplus;
53  Pminus=old.Pminus;
56  decaying=old.decaying;
57  Pstring=old.Pstring;
58  Pleft =old.Pleft;
59  Pright =old.Pright;
60 }
61 
63 {
64  if (this != &old)
65  {
68  Ptleft=old.Ptleft;
69  Ptright=old.Ptright;
70  Pplus=old.Pplus;
71  Pminus=old.Pminus;
74  decaying=old.decaying;
75  Pstring=old.Pstring;
76  Pleft =old.Pleft;
77  Pright =old.Pright;
78  }
79  return *this;
80 }
81 
82 //---------------------------------------------------------------------------------
83 
85 {
88  Ptleft=excited.GetLeftParton()->Get4Momentum().vect();
89  Ptleft.setZ(0.);
90  Ptright=excited.GetRightParton()->Get4Momentum().vect();
91  Ptright.setZ(0.);
92  G4LorentzVector P=excited.Get4Momentum();
93  Pplus =P.e() + P.pz();
94  Pminus=P.e() - P.pz();
97 
98  if (excited.GetDirection() > 0) {decaying=Left; }
99  else {decaying=Right;}
100 
101  Pleft = excited.GetLeftParton()->Get4Momentum();
102  Pright = excited.GetRightParton()->Get4Momentum();
103  Pstring= Pleft + Pright;
104 }
105 
106 //---------------------------------------------------------------------------------
107 
109  G4ParticleDefinition * newdecay,
110  const G4LorentzVector *momentum)
111 {
112  decaying=None;
113  G4LorentzVector Momentum = G4LorentzVector(momentum->vect(),momentum->e());
114  // Momentum of produced hadron
115  //G4cout<<"Had Mom "<<Momentum<<G4endl;
116  //G4cout<<"Str Mom "<<old.Pstring<<G4endl;
117  Pstring = old.Pstring - Momentum;
118  //G4cout<<"New Str Mom "<<Pstring<<" "<<Pstring.mag()<<G4endl;
119 
120  G4double StringMass = Pstring.mag();
121 
122  G4LorentzRotation toLAB(Pstring.boostVector());
123 
124  Pleft = toLAB*G4LorentzVector(0.,0., StringMass/2.,StringMass/2.);
125  Pright = toLAB*G4LorentzVector(0.,0.,-StringMass/2.,StringMass/2.);
126 
127  Ptleft =Pleft.vect(); Ptleft.setZ(0.);
128  Ptright=Pright.vect(); Ptright.setZ(0.);
129 
130  //G4cout<<"Pleft "<<Pleft<<G4endl;
131  //G4cout<<"Pright "<<Pright<<G4endl;
132  //G4cout<<"Pstring "<<Pstring<<G4endl;
133  if ( old.decaying == Left )
134  {
136  // Ptright = old.Ptright;
137  // Pright = old.Pright;
138 
139  LeftParton = newdecay;
140  // Ptleft = old.Ptleft - momentum->vect();
141  // Ptleft.setZ(0.);
142  // Pleft = old.Pleft - Momentum;
143  // Pstring=Pleft + Pright;
144 
147  decaying=Left;
148  } else if ( old.decaying == Right )
149  {
150  RightParton = newdecay;
151  // Ptright = old.Ptright - momentum->vect();
152  // Ptright.setZ(0.);
153  // Pright = old.Pright + Momentum;
154 
155  LeftParton = old.LeftParton;
156  // Ptleft = old.Ptleft;
157  // Pleft = old.Pleft;
158  // Pstring=Pleft + Pright;
159 
162  decaying=Right;
163  } else
164  {
165  throw G4HadronicException(__FILE__, __LINE__,
166  "G4FragmentingString::G4FragmentingString: no decay Direction defined");
167  }
168  Pplus = Pstring.plus(); //old.Pplus - (momentum->e() + momentum->pz());
169  Pminus = Pstring.minus(); //old.Pminus - (momentum->e() - momentum->pz());
170 }
171 
172 
173 //---------------------------------------------------------------------------------
174 
176  G4ParticleDefinition * newdecay)
177 {
178  decaying=None;
179 
180  Ptleft.setX(0.); Ptleft.setY(0.); Ptleft.setZ(0.);
181  Ptright.setX(0.); Ptright.setY(0.); Ptright.setZ(0.);
182  Pplus=0.; Pminus=0.;
184 
185  Pstring=G4LorentzVector(0.,0.,0.,0.);
186  Pleft =G4LorentzVector(0.,0.,0.,0.);
187  Pright =G4LorentzVector(0.,0.,0.,0.);
188 
189  if ( old.decaying == Left )
190  {
191  RightParton= old.RightParton;
192  LeftParton = newdecay;
193  decaying=Left;
194  } else if ( old.decaying == Right )
195  {
196  RightParton = newdecay;
197  LeftParton = old.LeftParton;
198  decaying=Right;
199  } else
200  {
201  throw G4HadronicException(__FILE__, __LINE__,
202  "G4FragmentingString::G4FragmentingString: no decay Direction defined");
203  }
204 }
205 
206 
207 //---------------------------------------------------------------------------------
208 
210 {}
211 
212 
213 //---------------------------------------------------------------------------------
214 
216 {
219  decaying=Right;
220 }
221 
222 //---------------------------------------------------------------------------------
223 
225 {
228  decaying=Left;
229 }
230 
231 //---------------------------------------------------------------------------------
232 
234 {
235  if (decaying == Left ) return +1;
236  else if (decaying == Right) return -1;
237  else throw G4HadronicException(__FILE__, __LINE__,
238  "G4FragmentingString::GetDecayDirection: decay side UNdefined!");
239  return 0;
240 }
241 
242 //---------------------------------------------------------------------------------
243 
245 {
246  return LeftParton->GetParticleSubType()== "di_quark"
247  && RightParton->GetParticleSubType()== "di_quark";
248 }
249 
250 //---------------------------------------------------------------------------------
251 
253 {
254  return theDecayParton->GetParticleSubType()== "quark";
255 }
256 
258 {
259  return theStableParton->GetParticleSubType()== "quark";
260 }
261 
262 //---------------------------------------------------------------------------------
263 
265 {
266  if (decaying == Left ) return Ptright;
267  else if (decaying == Right ) return Ptleft;
268  else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
269  return G4ThreeVector();
270 }
271 
273 {
274  if (decaying == Left ) return Ptleft;
275  else if (decaying == Right ) return Ptright;
276  else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
277  return G4ThreeVector();
278 }
279 
280 //---------------------------------------------------------------------------------
281 
283 {
284  return Pplus;
285 }
286 
288 {
289  return Pminus;
290 }
291 
293 {
294  if (decaying == Left ) return Pplus;
295  else if (decaying == Right ) return Pminus;
296  else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
297 }
298 
299 //---------------------------------------------------------------------------------
300 
302 {
304  momentum.setPz(0.5*(Pplus-Pminus));
305  momentum.setE(0.5*(Pplus+Pminus));
306  return momentum;
307 }
308 
310 {
311  // return Pplus*Pminus - (Ptleft+Ptright).mag2();
312  return Pstring.mag2();
313 }
314 
316 {
317  // return std::sqrt(this->Mass2());
318  return Pstring.mag();
319 }
320 
322 {
323  return Pplus*Pminus;
324 }
325 
327 {
328  return Pstring;
329 }
330 
332 {
333  return Pleft;
334 }
335 
337 {
338  return Pright;
339 }
340