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

Last change on this file since 1196 was 962, checked in by garnier, 15 years ago

update processes

File size: 8.5 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//G4cout<<" Left "<<G4endl;
83//G4cout<<"Pt right "<<Ptright<<G4endl;
84//G4cout<<"Pt left  "<<Ptleft <<G4endl;
85                RightParton= old.RightParton;
86                Ptright    = old.Ptright;
87                LeftParton = newdecay;
88                Ptleft     = old.Ptleft - momentum->vect();
89                Ptleft.setZ(0.);
90//G4cout<<"Pt right "<<Ptright<<G4endl;
91//G4cout<<"Pt left  "<<Ptleft <<G4endl;
92        } else if ( old.decaying == Right )
93        {
94//G4cout<<" Right "<<G4endl;
95//G4cout<<"Pt right "<<Ptright<<G4endl;
96//G4cout<<"Pt left  "<<Ptleft <<G4endl;
97                RightParton = newdecay;
98                Ptright     = old.Ptright - momentum->vect();
99                Ptright.setZ(0.);
100                LeftParton  = old.LeftParton;
101                Ptleft      = old.Ptleft;
102//G4cout<<"Pt right "<<Ptright<<G4endl;
103//G4cout<<"Pt left  "<<Ptleft <<G4endl;
104        } else
105        {
106                throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::G4FragmentingString: no decay Direction defined");
107        }
108        Pplus  = old.Pplus  - (momentum->e() + momentum->pz());
109        Pminus = old.Pminus - (momentum->e() - momentum->pz());
110       
111        //G4double Eold=0.5 * (old.Pplus + old.Pminus);
112        //G4double Enew=0.5 * (Pplus + Pminus);
113}
114
115
116//---------------------------------------------------------------------------------
117
118G4FragmentingString::G4FragmentingString(const G4FragmentingString &old,  // Uzhi
119                                         G4ParticleDefinition * newdecay) // Uzhi
120{                                                                         // Uzhi
121        decaying=None;                                                    // Uzhi
122        if ( old.decaying == Left )                                       // Uzhi
123        {                                                                 // Uzhi
124                RightParton= old.RightParton;                             // Uzhi
125                LeftParton = newdecay;                                    // Uzhi
126        } else if ( old.decaying == Right )                               // Uzhi
127        {                                                                 // Uzhi
128                RightParton = newdecay;                                   // Uzhi
129                LeftParton  = old.LeftParton;                             // Uzhi
130        } else                                                            // Uzhi
131        {
132                throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::G4FragmentingString: no decay Direction defined");
133        }
134}
135
136
137//---------------------------------------------------------------------------------
138
139G4FragmentingString::~G4FragmentingString()
140{}
141
142
143//---------------------------------------------------------------------------------
144
145void G4FragmentingString::SetLeftPartonStable()
146{
147     theStableParton=GetLeftParton();
148     theDecayParton=GetRightParton();
149     decaying=Right;
150}
151
152//---------------------------------------------------------------------------------
153
154void G4FragmentingString::SetRightPartonStable()
155{
156     theStableParton=GetRightParton();
157     theDecayParton=GetLeftParton();
158     decaying=Left;
159}
160
161//---------------------------------------------------------------------------------
162
163G4int G4FragmentingString::GetDecayDirection() const
164{
165        if      (decaying == Left ) return +1;
166        else if (decaying == Right) return -1;
167        else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::GetDecayDirection: decay side UNdefined!");
168        return 0;
169}
170 
171//---------------------------------------------------------------------------------
172
173G4bool G4FragmentingString::FourQuarkString() const
174{
175        return   LeftParton->GetParticleSubType()== "di_quark" 
176             && RightParton->GetParticleSubType()== "di_quark";
177}
178
179//---------------------------------------------------------------------------------
180
181G4bool G4FragmentingString::DecayIsQuark()
182{
183        return theDecayParton->GetParticleSubType()== "quark";
184}
185
186G4bool G4FragmentingString::StableIsQuark()
187{
188        return theStableParton->GetParticleSubType()== "quark";
189}
190
191//---------------------------------------------------------------------------------
192
193G4ThreeVector G4FragmentingString::StablePt()
194{
195        if (decaying == Left ) return Ptright;
196        else if (decaying == Right ) return Ptleft;
197        else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
198        return G4ThreeVector();
199}
200
201G4ThreeVector G4FragmentingString::DecayPt()
202{
203        if (decaying == Left ) return Ptleft;
204        else if (decaying == Right ) return Ptright;
205        else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
206        return G4ThreeVector();
207}
208
209//---------------------------------------------------------------------------------
210
211G4double G4FragmentingString::LightConePlus()
212{
213        return Pplus;
214}
215
216G4double G4FragmentingString::LightConeMinus()
217{
218        return Pminus;
219}
220
221G4double G4FragmentingString::LightConeDecay()
222{
223        if (decaying == Left ) return Pplus;
224        else if (decaying == Right ) return Pminus;
225        else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
226        return 0;
227}
228
229//---------------------------------------------------------------------------------
230
231G4LorentzVector G4FragmentingString::Get4Momentum() const
232{
233        G4LorentzVector momentum(Ptleft+Ptright,0);
234        momentum.setPz(0.5*(Pplus-Pminus));
235        momentum.setE(0.5*(Pplus+Pminus));
236        return momentum;
237}
238
239G4double G4FragmentingString::Mass2() const
240{
241        return Pplus*Pminus - (Ptleft+Ptright).mag2();
242}
243
244G4double G4FragmentingString::Mass() const
245{
246        return std::sqrt(this->Mass2());
247}
248
249G4double G4FragmentingString::MassT2() const
250{
251        return Pplus*Pminus;
252}
Note: See TracBrowser for help on using the repository browser.