source: trunk/source/processes/hadronic/models/parton_string/hadronization/src/G4FragmentingString.cc@ 880

Last change on this file since 880 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 6.9 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
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
46G4FragmentingString::G4FragmentingString(const G4FragmentingString &old)
47{
48 LeftParton=old.LeftParton;
49 RightParton=old.RightParton;
50 Ptleft=old.Ptleft;
51 Ptright=old.Ptright;
52 Pplus=old.Pplus;
53 Pminus=old.Pminus;
54 decaying=old.decaying;
55}
56
57//---------------------------------------------------------------------------------
58
59G4FragmentingString::G4FragmentingString(const G4ExcitedString &excited)
60{
61 LeftParton=excited.GetLeftParton()->GetDefinition();
62 RightParton=excited.GetRightParton()->GetDefinition();
63 Ptleft=excited.GetLeftParton()->Get4Momentum().vect();
64 Ptleft.setZ(0.);
65 Ptright=excited.GetRightParton()->Get4Momentum().vect();
66 Ptright.setZ(0.);
67 G4LorentzVector P=excited.Get4Momentum();
68 Pplus =P.e() + P.pz();
69 Pminus=P.e() - P.pz();
70 decaying=None;
71}
72
73//---------------------------------------------------------------------------------
74
75G4FragmentingString::G4FragmentingString(const G4FragmentingString &old,
76 G4ParticleDefinition * newdecay,
77 const G4LorentzVector *momentum)
78{
79 decaying=None;
80 if ( old.decaying == Left )
81 {
82 RightParton= old.RightParton;
83 Ptright = old.Ptright;
84 LeftParton = newdecay;
85 Ptleft = old.Ptleft - momentum->vect();
86 Ptleft.setZ(0.);
87 } else if ( old.decaying == Right )
88 {
89 RightParton = newdecay;
90 Ptright = old.Ptright - momentum->vect();
91 Ptright.setZ(0.);
92 LeftParton = old.LeftParton;
93 Ptleft = old.Ptleft;
94 } else
95 {
96 throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::G4FragmentingString: no decay Direction defined");
97 }
98 Pplus = old.Pplus - (momentum->e() + momentum->pz());
99 Pminus = old.Pminus - (momentum->e() - momentum->pz());
100
101 //G4double Eold=0.5 * (old.Pplus + old.Pminus);
102 //G4double Enew=0.5 * (Pplus + Pminus);
103}
104
105
106//---------------------------------------------------------------------------------
107
108G4FragmentingString::~G4FragmentingString()
109{}
110
111
112//---------------------------------------------------------------------------------
113
114void G4FragmentingString::SetLeftPartonStable()
115{
116 theStableParton=GetLeftParton();
117 theDecayParton=GetRightParton();
118 decaying=Right;
119}
120
121//---------------------------------------------------------------------------------
122
123void G4FragmentingString::SetRightPartonStable()
124{
125 theStableParton=GetRightParton();
126 theDecayParton=GetLeftParton();
127 decaying=Left;
128}
129
130//---------------------------------------------------------------------------------
131
132G4int G4FragmentingString::GetDecayDirection() const
133{
134 if (decaying == Left ) return +1;
135 else if (decaying == Right) return -1;
136 else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::GetDecayDirection: decay side UNdefined!");
137 return 0;
138}
139
140//---------------------------------------------------------------------------------
141
142G4bool G4FragmentingString::FourQuarkString() const
143{
144 return LeftParton->GetParticleSubType()== "di_quark"
145 && RightParton->GetParticleSubType()== "di_quark";
146}
147
148//---------------------------------------------------------------------------------
149
150G4bool G4FragmentingString::DecayIsQuark()
151{
152 return theDecayParton->GetParticleSubType()== "quark";
153}
154
155G4bool G4FragmentingString::StableIsQuark()
156{
157 return theStableParton->GetParticleSubType()== "quark";
158}
159
160//---------------------------------------------------------------------------------
161
162G4ThreeVector G4FragmentingString::StablePt()
163{
164 if (decaying == Left ) return Ptright;
165 else if (decaying == Right ) return Ptleft;
166 else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
167 return G4ThreeVector();
168}
169
170G4ThreeVector G4FragmentingString::DecayPt()
171{
172 if (decaying == Left ) return Ptleft;
173 else if (decaying == Right ) return Ptright;
174 else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
175 return G4ThreeVector();
176}
177
178//---------------------------------------------------------------------------------
179
180G4double G4FragmentingString::LightConePlus()
181{
182 return Pplus;
183}
184
185G4double G4FragmentingString::LightConeMinus()
186{
187 return Pminus;
188}
189
190G4double G4FragmentingString::LightConeDecay()
191{
192 if (decaying == Left ) return Pplus;
193 else if (decaying == Right ) return Pminus;
194 else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
195 return 0;
196}
197
198//---------------------------------------------------------------------------------
199
200G4LorentzVector G4FragmentingString::Get4Momentum() const
201{
202 G4LorentzVector momentum(Ptleft+Ptright,0);
203 momentum.setPz(0.5*(Pplus-Pminus));
204 momentum.setE(0.5*(Pplus+Pminus));
205 return momentum;
206}
207
208G4double G4FragmentingString::Mass2() const
209{
210 return Pplus*Pminus - (Ptleft+Ptright).mag2();
211}
212
213G4double G4FragmentingString::Mass() const
214{
215 return std::sqrt(this->Mass2());
216}
Note: See TracBrowser for help on using the repository browser.