Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/parton_string/hadronization/src/G4FragmentingString.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 
 46 G4FragmentingString::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   theStableParton=old.theStableParton;
 55   theDecayParton=old.theDecayParton;
 56   decaying=old.decaying;
 57         Pstring=old.Pstring;
 58         Pleft  =old.Pleft;
 59         Pright =old.Pright;
 60 }
 61 
 62 G4FragmentingString & G4FragmentingString::operator =(const G4FragmentingString &old)
 63 {
 64    if (this != &old)
 65    {
 66      LeftParton=old.LeftParton;
 67      RightParton=old.RightParton;
 68      Ptleft=old.Ptleft;
 69      Ptright=old.Ptright;
 70      Pplus=old.Pplus;
 71      Pminus=old.Pminus;
 72      theStableParton=old.theStableParton;
 73      theDecayParton=old.theDecayParton;
 74      decaying=old.decaying;
 75      Pstring=old.Pstring;
 76      Pleft  =old.Pleft;
 77      Pright =old.Pright;
 78    }
 79    return *this;
 80 }
 81 
 82 //---------------------------------------------------------------------------------
 83 
 84 G4FragmentingString::G4FragmentingString(const G4ExcitedString &excited)
 85 {
 86   LeftParton=excited.GetLeftParton()->GetDefinition();
 87   RightParton=excited.GetRightParton()->GetDefinition();
 88   Ptleft=excited.GetLeftParton()->Get4Momentum().vect();
 89   Ptleft.setZ(0.);
 90   Ptright=excited.GetRightParton()->Get4Momentum().vect();
 91   Ptright.setZ(0.);
 92   theStableParton=0;
 93   theDecayParton=0;
 94 
 95         if (excited.GetDirection() > 0) {decaying=Left; }
 96         else                            {decaying=Right;}
 97 
 98         Pleft  = excited.GetLeftParton()->Get4Momentum();
 99         Pright = excited.GetRightParton()->Get4Momentum();
100         Pstring= Pleft + Pright;
101 
102         Pplus  = Pstring.plus();
103         Pminus = Pstring.minus(); 
104 }
105 
106 //---------------------------------------------------------------------------------
107 
108 G4FragmentingString::G4FragmentingString(const G4FragmentingString &old,
109            G4ParticleDefinition * newdecay,
110            const G4LorentzVector *momentum)
111 {
112   decaying=None;
113   // Momentum of produced hadron
114         G4LorentzVector Momentum = G4LorentzVector(momentum->vect(),momentum->e());
115              
116         if ( old.decaying == Left )
117         {
118                 RightParton= old.RightParton;
119                 Ptright    = old.Ptright;
120                 Pright     = old.Pright;
121 
122                 LeftParton = newdecay;
123                 Ptleft     = old.Ptleft - momentum->vect();
124                 Ptleft.setZ(0.);
125                 Pleft = old.Pleft - Momentum;
126 
127                 Pstring = Pleft + Pright;
128                 Pplus   = Pstring.plus();
129                 Pminus  = Pstring.minus();
130 
131                 theDecayParton=GetLeftParton();
132                 theStableParton=GetRightParton();
133                 decaying = Left;
134         } else if ( old.decaying == Right )
135         {
136           RightParton = newdecay;
137                 Ptright     = old.Ptright - momentum->vect();
138                 Ptright.setZ(0.);
139                 Pright = old.Pright - Momentum;
140 
141                 LeftParton  = old.LeftParton;
142                 Ptleft      = old.Ptleft;
143                 Pleft  = old.Pleft;
144 
145                 Pstring = Pleft + Pright;
146                 Pplus   = Pstring.plus();
147                 Pminus  = Pstring.minus();
148 
149                 theDecayParton=GetRightParton();
150                 theStableParton=GetLeftParton();
151                 decaying = Right;
152         } else
153   {
154     throw G4HadronicException(__FILE__, __LINE__, 
155                   "G4FragmentingString::G4FragmentingString: no decay Direction defined");
156   }
157 }
158 
159 
160 //---------------------------------------------------------------------------------
161 
162 G4FragmentingString::G4FragmentingString(const G4FragmentingString &old,  
163            G4ParticleDefinition * newdecay) 
164 {                                                                         
165   decaying=None;                                                    
166 
167         Ptleft.setX(0.);  Ptleft.setY(0.);  Ptleft.setZ(0.);
168         Ptright.setX(0.); Ptright.setY(0.); Ptright.setZ(0.);
169         Pplus=0.; Pminus=0.;                               
170         theStableParton=0; theDecayParton=0;              
171 
172         Pstring=G4LorentzVector(0.,0.,0.,0.);
173         Pleft  =G4LorentzVector(0.,0.,0.,0.);
174         Pright =G4LorentzVector(0.,0.,0.,0.);
175 
176   if ( old.decaying == Left )                                       
177   {                                                                 
178     RightParton= old.RightParton;                             
179     LeftParton = newdecay;                                    
180                 decaying=Left;
181   } else if ( old.decaying == Right )                               
182   {                                                                 
183     RightParton = newdecay;                                   
184     LeftParton  = old.LeftParton;                             
185                 decaying=Right;
186   } else                                                            
187   {
188     throw G4HadronicException(__FILE__, __LINE__, 
189                   "G4FragmentingString::G4FragmentingString: no decay Direction defined");
190   }
191 }
192 
193 
194 //---------------------------------------------------------------------------------
195 
196 G4FragmentingString::~G4FragmentingString()
197 {}
198 
199 //---------------------------------------------------------------------------------
200 
201 void G4FragmentingString::SetLeftPartonStable()
202 {
203      theStableParton=GetLeftParton();
204      theDecayParton=GetRightParton();
205      decaying=Right;
206 }
207 
208 //---------------------------------------------------------------------------------
209 
210 void G4FragmentingString::SetRightPartonStable()
211 {
212      theStableParton=GetRightParton();
213      theDecayParton=GetLeftParton();
214      decaying=Left;
215 }
216 
217 //---------------------------------------------------------------------------------
218 
219 G4int G4FragmentingString::GetDecayDirection() const
220 {
221   if      (decaying == Left ) return +1;
222   else if (decaying == Right) return -1;
223   else throw G4HadronicException(__FILE__, __LINE__, 
224                "G4FragmentingString::GetDecayDirection: decay side UNdefined!");
225   return 0;
226 }
227  
228 //---------------------------------------------------------------------------------
229 
230 G4bool G4FragmentingString::IsAFourQuarkString() const
231 {
232   return   LeftParton->GetParticleSubType()== "di_quark" 
233        && RightParton->GetParticleSubType()== "di_quark";
234 }
235 
236 //---------------------------------------------------------------------------------
237 
238 G4bool G4FragmentingString::DecayIsQuark()
239 {
240   return theDecayParton->GetParticleSubType()== "quark";
241 }
242 
243 G4bool G4FragmentingString::StableIsQuark()
244 {
245   return theStableParton->GetParticleSubType()== "quark";
246 }
247 
248 //---------------------------------------------------------------------------------
249 
250 G4ThreeVector G4FragmentingString::StablePt()
251 {
252   if (decaying == Left ) return Ptright;
253   else if (decaying == Right ) return Ptleft;
254   else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
255   return G4ThreeVector();
256 }
257 
258 G4ThreeVector G4FragmentingString::DecayPt()
259 {
260   if (decaying == Left ) return Ptleft;
261   else if (decaying == Right ) return Ptright;
262   else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
263   return G4ThreeVector();
264 }
265 
266 //---------------------------------------------------------------------------------
267 
268 G4double G4FragmentingString::LightConePlus()
269 {
270   return Pplus;
271 }
272 
273 G4double G4FragmentingString::LightConeMinus()
274 {
275   return Pminus;
276 }
277 
278 G4double G4FragmentingString::LightConeDecay()
279 {
280   if (decaying == Left ) return Pplus;
281   else if (decaying == Right ) return Pminus;
282   else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!");
283 }
284 
285 //---------------------------------------------------------------------------------
286 
287 G4LorentzVector G4FragmentingString::Get4Momentum() const
288 {
289         return Pstring;
290 }
291 
292 G4double G4FragmentingString::Mass2() const
293 {
294   return Pstring.mag2();
295 }
296 
297 G4double G4FragmentingString::Mass() const
298 {
299   return Pstring.mag();
300 }
301 
302 G4double G4FragmentingString::MassT2() const
303 {
304   return Pplus*Pminus;
305 }
306 
307 G4LorentzVector G4FragmentingString::GetPstring()
308 { 
309   return Pstring;
310 }
311 
312 G4LorentzVector G4FragmentingString::GetPleft()
313 { 
314   return Pleft;
315 }
316 
317 G4LorentzVector G4FragmentingString::GetPright()
318 { 
319   return Pright;
320 }
321 
322 G4LorentzRotation G4FragmentingString::TransformToAlignedCms()
323 {
324      G4LorentzVector momentum = Pstring;
325      G4LorentzRotation toAlignedCms(-1*momentum.boostVector());
326      momentum = toAlignedCms * Pleft;
327 
328      toAlignedCms.rotateZ(-1*momentum.phi());
329      toAlignedCms.rotateY(-1*momentum.theta());
330 
331      Pleft   *= toAlignedCms;
332      Pright  *= toAlignedCms;
333      Pstring *= toAlignedCms;
334 
335      Ptleft  = G4ThreeVector(Pleft.vect());
336      Ptleft.setZ(0.);
337      Ptright = G4ThreeVector(Pright.vect());
338      Pplus  = Pstring.plus();
339      Pminus = Pstring.minus();
340 
341      return toAlignedCms;
342 }
343