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

Last change on this file since 1358 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

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