Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LundStringFragmentation Class Reference

#include <G4LundStringFragmentation.hh>

Inheritance diagram for G4LundStringFragmentation:
Collaboration diagram for G4LundStringFragmentation:

Public Member Functions

 G4LundStringFragmentation ()
 
virtual ~G4LundStringFragmentation ()
 
virtual G4KineticTrackVectorFragmentString (const G4ExcitedString &theString)
 
- Public Member Functions inherited from G4VLongitudinalStringDecay
 G4VLongitudinalStringDecay ()
 
virtual ~G4VLongitudinalStringDecay ()
 
G4int SampleQuarkFlavor (void)
 
G4ThreeVector SampleQuarkPt (G4double ptMax=-1.)
 
G4KineticTrackVectorDecayResonans (G4KineticTrackVector *aHadrons)
 
void SetSigmaTransverseMomentum (G4double aQT)
 
void SetStrangenessSuppression (G4double aValue)
 
void SetDiquarkSuppression (G4double aValue)
 
void SetDiquarkBreakProbability (G4double aValue)
 
void SetVectorMesonProbability (G4double aValue)
 
void SetSpinThreeHalfBarionProbability (G4double aValue)
 
void SetScalarMesonMixings (std::vector< G4double > aVector)
 
void SetVectorMesonMixings (std::vector< G4double > aVector)
 
void SetStringTensionParameter (G4double aValue)
 

Additional Inherited Members

- Protected Types inherited from G4VLongitudinalStringDecay
typedef std::pair
< G4ParticleDefinition
*, G4ParticleDefinition * > 
pDefPair
 
typedef G4ParticleDefinition
*(G4HadronBuilder::* 
Pcreate )(G4ParticleDefinition *, G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VLongitudinalStringDecay
virtual void SetMassCut (G4double aValue)
 
G4KineticTrackVectorLightFragmentationTest (const G4ExcitedString *const theString)
 
G4double FragmentationMass (const G4FragmentingString *const string, Pcreate build=0, pDefPair *pdefs=0)
 
G4ParticleDefinitionFindParticle (G4int Encoding)
 
G4ExcitedStringCPExcited (const G4ExcitedString &string)
 
G4KineticTrackSplitup (G4FragmentingString *string, G4FragmentingString *&newString)
 
G4ParticleDefinitionQuarkSplitup (G4ParticleDefinition *decay, G4ParticleDefinition *&created)
 
pDefPair CreatePartonPair (G4int NeedParticle, G4bool AllowDiquarks=true)
 
void CalculateHadronTimePosition (G4double theInitialStringMass, G4KineticTrackVector *)
 
void ConstructParticle ()
 
G4ParticleDefinitionCreateHadron (G4int id1, G4int id2, G4bool theGivenSpin, G4int theSpin)
 
G4double GetDiquarkSuppress ()
 
G4double GetDiquarkBreakProb ()
 
G4double GetStrangeSuppress ()
 
G4double GetClusterMass ()
 
G4int GetClusterLoopInterrupt ()
 
G4double GetStringTensionParameter ()
 
- Protected Attributes inherited from G4VLongitudinalStringDecay
G4double MassCut
 
G4double ClusterMass
 
G4double SigmaQT
 
G4double DiquarkSuppress
 
G4double DiquarkBreakProb
 
G4double SmoothParam
 
G4double StrangeSuppress
 
G4int StringLoopInterrupt
 
G4int ClusterLoopInterrupt
 
G4HadronBuilderhadronizer
 
G4double pspin_meson
 
G4double pspin_barion
 
std::vector< G4doublevectorMesonMix
 
std::vector< G4doublescalarMesonMix
 
G4bool PastInitPhase
 
G4double Kappa
 

Detailed Description

Definition at line 42 of file G4LundStringFragmentation.hh.

Constructor & Destructor Documentation

G4LundStringFragmentation::G4LundStringFragmentation ( )

Definition at line 51 of file G4LundStringFragmentation.cc.

52 {
53  // ------ For estimation of a minimal string mass ---------------
54  Mass_of_light_quark =140.*MeV;
55  Mass_of_heavy_quark =500.*MeV;
56  Mass_of_string_junction=720.*MeV;
57  // ------ An estimated minimal string mass ----------------------
58  MinimalStringMass = 0.;
59  MinimalStringMass2 = 0.;
60  // ------ Minimal invariant mass used at a string fragmentation -
61  WminLUND = 0.45*GeV; //0.23*GeV; // Uzhi 0.7 -> 0.23 3.8.10 //0.8 1.5
62  // ------ Smooth parameter used at a string fragmentation for ---
63  // ------ smearing sharp mass cut-off ---------------------------
64  SmoothParam = 0.2;
65 
68  SetStrangenessSuppression(0.46); //(0.447);
70 
71  // For treating of small string decays
72  for(G4int i=0; i<3; i++)
73  { for(G4int j=0; j<3; j++)
74  { for(G4int k=0; k<6; k++)
75  {
76  Meson[i][j][k]=0; MesonWeight[i][j][k]=0.;
77  }
78  }
79  }
80  //--------------------------
81  Meson[0][0][0]=111; // dbar-d Pi0
82  MesonWeight[0][0][0]=(1.-pspin_meson)*(1.-scalarMesonMix[0]);
83 
84  Meson[0][0][1]=221; // dbar-d Eta
85  MesonWeight[0][0][1]=(1.-pspin_meson)*(scalarMesonMix[0]-scalarMesonMix[1]);
86 
87  Meson[0][0][2]=331; // dbar-d EtaPrime
88  MesonWeight[0][0][2]=(1.-pspin_meson)*(scalarMesonMix[1]);
89 
90  Meson[0][0][3]=113; // dbar-d Rho0
91  MesonWeight[0][0][3]=pspin_meson*(1.-vectorMesonMix[0]);
92 
93  Meson[0][0][4]=223; // dbar-d Omega
94  MesonWeight[0][0][4]=pspin_meson*(vectorMesonMix[0]);
95  //--------------------------
96 
97  Meson[0][1][0]=211; // dbar-u Pi+
98  MesonWeight[0][1][0]=(1.-pspin_meson);
99 
100  Meson[0][1][1]=213; // dbar-u Rho+
101  MesonWeight[0][1][1]=pspin_meson;
102  //--------------------------
103 
104  Meson[0][2][0]=311; // dbar-s K0bar
105  MesonWeight[0][2][0]=(1.-pspin_meson);
106 
107  Meson[0][2][1]=313; // dbar-s K*0bar
108  MesonWeight[0][2][1]=pspin_meson;
109  //--------------------------
110 
111  Meson[1][0][0]=211; // ubar-d Pi-
112  MesonWeight[1][0][0]=(1.-pspin_meson);
113 
114  Meson[1][0][1]=213; // ubar-d Rho-
115  MesonWeight[1][0][1]=pspin_meson;
116  //--------------------------
117 
118  Meson[1][1][0]=111; // ubar-u Pi0
119  MesonWeight[1][1][0]=(1.-pspin_meson)*(1.-scalarMesonMix[0]);
120 
121  Meson[1][1][1]=221; // ubar-u Eta
122  MesonWeight[1][1][1]=(1.-pspin_meson)*(scalarMesonMix[0]-scalarMesonMix[1]);
123 
124  Meson[1][1][2]=331; // ubar-u EtaPrime
125  MesonWeight[1][1][2]=(1.-pspin_meson)*(scalarMesonMix[1]);
126 
127  Meson[1][1][3]=113; // ubar-u Rho0
128  MesonWeight[1][1][3]=pspin_meson*(1.-vectorMesonMix[0]);
129 
130  Meson[1][1][4]=223; // ubar-u Omega
131  //MesonWeight[1][1][4]=pspin_meson*(scalarMesonMix[0]);
132  MesonWeight[1][1][4]=pspin_meson*(vectorMesonMix[0]); // Uzhi 2015 scalar -> vector
133  //--------------------------
134 
135  Meson[1][2][0]=321; // ubar-s K-
136  MesonWeight[1][2][0]=(1.-pspin_meson);
137 
138  Meson[1][2][1]=323; // ubar-s K*-bar -
139  MesonWeight[1][2][1]=pspin_meson;
140  //--------------------------
141 
142  Meson[2][0][0]=311; // sbar-d K0
143  MesonWeight[2][0][0]=(1.-pspin_meson);
144 
145  Meson[2][0][1]=313; // sbar-d K*0
146  MesonWeight[2][0][1]=pspin_meson;
147  //--------------------------
148 
149  Meson[2][1][0]=321; // sbar-u K+
150  MesonWeight[2][1][0]=(1.-pspin_meson);
151 
152  Meson[2][1][1]=323; // sbar-u K*+
153  MesonWeight[2][1][1]=pspin_meson;
154  //--------------------------
155 
156  Meson[2][2][0]=221; // sbar-s Eta
157  MesonWeight[2][2][0]=(1.-pspin_meson)*(1.-scalarMesonMix[5]);
158 
159  Meson[2][2][1]=331; // sbar-s EtaPrime
160  MesonWeight[2][2][1]=(1.-pspin_meson)*(1.-scalarMesonMix[5]);
161 
162  Meson[2][2][3]=333; // sbar-s EtaPrime
163  MesonWeight[2][2][3]=pspin_meson*(vectorMesonMix[5]);
164  //--------------------------
165 
166  for(G4int i=0; i<3; i++)
167  { for(G4int j=0; j<3; j++)
168  { for(G4int k=0; k<3; k++)
169  { for(G4int l=0; l<4; l++)
170  { Baryon[i][j][k][l]=0; BaryonWeight[i][j][k][l]=0.;}
171  }
172  }
173  }
174 
175  G4double pspin_barion_in=pspin_barion;
176  //pspin_barion=0.75;
177  //---------------------------------------
178 
179  Baryon[0][0][0][0]=1114; // Delta-
180  BaryonWeight[0][0][0][0]=1.;
181  //---------------------------------------
182 
183  Baryon[0][0][1][0]=2112; // neutron
184  BaryonWeight[0][0][1][0]=1.-pspin_barion;
185 
186  Baryon[0][0][1][1]=2114; // Delta0
187  BaryonWeight[0][0][1][1]=pspin_barion;
188  //---------------------------------------
189 
190  Baryon[0][0][2][0]=3112; // Sigma-
191  BaryonWeight[0][0][2][0]=1.-pspin_barion;
192 
193  Baryon[0][0][2][1]=3114; // Sigma*-
194  BaryonWeight[0][0][2][1]=pspin_barion;
195  //---------------------------------------
196 
197  Baryon[0][1][0][0]=2112; // neutron
198  BaryonWeight[0][1][0][0]=1.-pspin_barion;
199 
200  Baryon[0][1][0][1]=2114; // Delta0
201  BaryonWeight[0][1][0][1]=pspin_barion;
202  //---------------------------------------
203 
204  Baryon[0][1][1][0]=2212; // proton
205  BaryonWeight[0][1][1][0]=1.-pspin_barion;
206 
207  Baryon[0][1][1][1]=2214; // Delta+
208  BaryonWeight[0][1][1][1]=pspin_barion;
209  //---------------------------------------
210 
211  Baryon[0][1][2][0]=3122; // Lambda
212  BaryonWeight[0][1][2][0]=(1.-pspin_barion)*0.5;
213 
214  Baryon[0][1][2][1]=3212; // Sigma0
215  BaryonWeight[0][1][2][1]=(1.-pspin_barion)*0.5;
216 
217  Baryon[0][1][2][2]=3214; // Sigma*0
218  BaryonWeight[0][1][2][2]=pspin_barion;
219  //---------------------------------------
220 
221  Baryon[0][2][0][0]=3112; // Sigma-
222  BaryonWeight[0][2][0][0]=1.-pspin_barion;
223 
224  Baryon[0][2][0][1]=3114; // Sigma*-
225  BaryonWeight[0][2][0][1]=pspin_barion;
226  //---------------------------------------
227 
228  Baryon[0][2][1][0]=3122; // Lambda
229  BaryonWeight[0][2][1][0]=(1.-pspin_barion)*0.5;
230 
231  Baryon[0][2][1][1]=3212; // Sigma0
232  BaryonWeight[0][2][1][1]=(1.-pspin_barion)*0.5;
233 
234  Baryon[0][2][1][2]=3214; // Sigma*0
235  BaryonWeight[0][2][1][2]=pspin_barion;
236  //---------------------------------------
237 
238  Baryon[0][2][2][0]=3312; // Theta-
239  BaryonWeight[0][2][2][0]=1.-pspin_barion;
240 
241  Baryon[0][2][2][1]=3314; // Theta*-
242  BaryonWeight[0][2][2][1]=pspin_barion;
243  //---------------------------------------
244 
245  Baryon[1][0][0][0]=2112; // neutron
246  BaryonWeight[1][0][0][0]=1.-pspin_barion;
247 
248  Baryon[1][0][0][1]=2114; // Delta0
249  BaryonWeight[1][0][0][1]=pspin_barion;
250  //---------------------------------------
251 
252  Baryon[1][0][1][0]=2212; // proton
253  BaryonWeight[1][0][1][0]=1.-pspin_barion;
254 
255  Baryon[1][0][1][1]=2214; // Delta+
256  BaryonWeight[1][0][1][1]=pspin_barion;
257  //---------------------------------------
258 
259  Baryon[1][0][2][0]=3122; // Lambda
260  BaryonWeight[1][0][2][0]=(1.-pspin_barion)*0.5;
261 
262  Baryon[1][0][2][1]=3212; // Sigma0
263  BaryonWeight[1][0][2][1]=(1.-pspin_barion)*0.5;
264 
265  Baryon[1][0][2][2]=3214; // Sigma*0
266  BaryonWeight[1][0][2][2]=pspin_barion;
267  //---------------------------------------
268 
269  Baryon[1][1][0][0]=2212; // proton
270  BaryonWeight[1][1][0][0]=1.-pspin_barion;
271 
272  Baryon[1][1][0][1]=2214; // Delta+
273  BaryonWeight[1][1][0][1]=pspin_barion;
274  //---------------------------------------
275 
276  Baryon[1][1][1][0]=2224; // Delta++
277  BaryonWeight[1][1][1][0]=1.;
278  //---------------------------------------
279 
280  Baryon[1][1][2][0]=3222; // Sigma+
281  BaryonWeight[1][1][2][0]=1.-pspin_barion;
282 
283  Baryon[1][1][2][1]=3224; // Sigma*+
284  BaryonWeight[1][1][2][1]=pspin_barion;
285  //---------------------------------------
286 
287  Baryon[1][2][0][0]=3122; // Lambda
288  BaryonWeight[1][2][0][0]=(1.-pspin_barion)*0.5;
289 
290  Baryon[1][2][0][1]=3212; // Sigma0
291  BaryonWeight[1][2][0][1]=(1.-pspin_barion)*0.5;
292 
293  Baryon[1][2][0][2]=3214; // Sigma*0
294  BaryonWeight[1][2][0][2]=pspin_barion;
295  //---------------------------------------
296 
297  Baryon[1][2][1][0]=3222; // Sigma+
298  BaryonWeight[1][2][1][0]=1.-pspin_barion;
299 
300  Baryon[1][2][1][1]=3224; // Sigma*+
301  BaryonWeight[1][2][1][1]=pspin_barion;
302  //---------------------------------------
303 
304  Baryon[1][2][2][0]=3322; // Theta0
305  BaryonWeight[1][2][2][0]=1.-pspin_barion;
306 
307  Baryon[1][2][2][1]=3324; // Theta*0
308  BaryonWeight[1][2][2][1]=pspin_barion;
309  //---------------------------------------
310 
311  Baryon[2][0][0][0]=3112; // Sigma-
312  BaryonWeight[2][0][0][0]=1.-pspin_barion;
313 
314  Baryon[2][0][0][1]=3114; // Sigma*-
315  BaryonWeight[2][0][0][1]=pspin_barion;
316  //---------------------------------------
317 
318  Baryon[2][0][1][0]=3122; // Lambda
319  BaryonWeight[2][0][1][0]=(1.-pspin_barion)*0.5;
320 
321  Baryon[2][0][1][1]=3212; // Sigma0
322  BaryonWeight[2][0][1][1]=(1.-pspin_barion)*0.5;
323 
324  Baryon[2][0][1][2]=3214; // Sigma*0
325  BaryonWeight[2][0][1][2]=pspin_barion;
326  //---------------------------------------
327 
328  Baryon[2][0][2][0]=3312; // Sigma-
329  BaryonWeight[2][0][2][0]=1.-pspin_barion;
330 
331  Baryon[2][0][2][1]=3314; // Sigma*-
332  BaryonWeight[2][0][2][1]=pspin_barion;
333  //---------------------------------------
334 
335  Baryon[2][1][0][0]=3122; // Lambda
336  BaryonWeight[2][1][0][0]=(1.-pspin_barion)*0.5;
337 
338  Baryon[2][1][0][1]=3212; // Sigma0
339  BaryonWeight[2][1][0][1]=(1.-pspin_barion)*0.5;
340 
341  Baryon[2][1][0][2]=3214; // Sigma*0
342  BaryonWeight[2][1][0][2]=pspin_barion;
343  //---------------------------------------
344 
345  Baryon[2][1][1][0]=3222; // Sigma+
346  BaryonWeight[2][1][1][0]=1.-pspin_barion;
347 
348  Baryon[2][1][1][1]=3224; // Sigma*+
349  BaryonWeight[2][1][1][1]=pspin_barion;
350  //---------------------------------------
351 
352  Baryon[2][1][2][0]=3322; // Theta0
353  BaryonWeight[2][1][2][0]=1.-pspin_barion;
354 
355  Baryon[2][1][2][1]=3324; // Theta*0
356  BaryonWeight[2][1][2][1]=pspin_barion;
357  //---------------------------------------
358 
359  Baryon[2][2][0][0]=3312; // Theta-
360  BaryonWeight[2][2][0][0]=1.-pspin_barion;
361 
362  Baryon[2][2][0][1]=3314; // Theta*-
363  BaryonWeight[2][2][0][1]=pspin_barion;
364  //---------------------------------------
365 
366  Baryon[2][2][1][0]=3322; // Theta0
367  BaryonWeight[2][2][1][0]=1.-pspin_barion;
368 
369  Baryon[2][2][1][1]=3324; // Theta*0
370  BaryonWeight[2][2][1][1]=pspin_barion;
371  //---------------------------------------
372 
373  Baryon[2][2][2][0]=3334; // Omega
374  BaryonWeight[2][2][2][0]=1.;
375  //---------------------------------------
376 
377  pspin_barion=pspin_barion_in;
378  /*
379  for(G4int i=0; i<3; i++)
380  { for(G4int j=0; j<3; j++)
381  { for(G4int k=0; k<3; k++)
382  { for(G4int l=0; l<4; l++)
383  { G4cout<<i<<" "<<j<<" "<<k<<" "<<l<<" "<<Baryon[i][j][k][l]<<G4endl;}
384  }
385  }
386  }
387  G4int Uzhi;
388  G4cin>>Uzhi;
389  */
390 
392  Prob_QQbar[0]=StrangeSuppress; // Probability of ddbar production
393  Prob_QQbar[1]=StrangeSuppress; // Probability of uubar production
394  Prob_QQbar[2]=1.0-2.*StrangeSuppress; // Probability of ssbar production
395  SetStrangenessSuppression(0.46); //(0.447);
396 
397  //A.R. 25-Jul-2012 : Coverity fix.
398  for ( G4int i=0 ; i<35 ; i++ ) {
399  FS_LeftHadron[i] = 0;
400  FS_RightHadron[i] = 0;
401  FS_Weight[i] = 0.0;
402  }
403  NumberOf_FS = 0;
404 
405 }
void SetStrangenessSuppression(G4double aValue)
std::vector< G4double > vectorMesonMix
int G4int
Definition: G4Types.hh:78
void SetStringTensionParameter(G4double aValue)
void SetDiquarkBreakProbability(G4double aValue)
std::vector< G4double > scalarMesonMix
static constexpr double GeV
Definition: G4SIunits.hh:217
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
void SetDiquarkSuppression(G4double aValue)

Here is the call graph for this function:

G4LundStringFragmentation::~G4LundStringFragmentation ( )
virtual

Definition at line 408 of file G4LundStringFragmentation.cc.

409 {}

Member Function Documentation

G4KineticTrackVector * G4LundStringFragmentation::FragmentString ( const G4ExcitedString theString)
virtual

Implements G4VLongitudinalStringDecay.

Definition at line 512 of file G4LundStringFragmentation.cc.

513 {
514  // Can no longer modify Parameters for Fragmentation.
515  PastInitPhase=true;
516 
517  SetMassCut(160.*MeV); // For LightFragmentationTest it is required that no one pi-meson can be produced
518 
519  G4FragmentingString aString(theString);
520  SetMinimalStringMass(&aString);
521 
522  #ifdef debug_LUNDfragmentation
523  G4cout<<G4endl<<"LUND StringFragm: String Mass " <<theString.Get4Momentum().mag()<<" Pz "
524  <<theString.Get4Momentum().pz()
525  <<"------------------------------------"<<G4endl;
526  G4cout<<"String ends Direct "<<theString.GetLeftParton()->GetPDGcode()<<" "
527  <<theString.GetRightParton()->GetPDGcode()<<" "
528  <<theString.GetDirection()<< G4endl;
529  G4cout<<"Left mom "<<theString.GetLeftParton()->Get4Momentum()<<G4endl;
530  G4cout<<"Right mom "<<theString.GetRightParton()->Get4Momentum()<<G4endl;
531  G4cout<<"Check for Fragmentation "<<G4endl;
532  #endif
533 
534  G4KineticTrackVector * LeftVector(0);
535 
536  //if(!IsFragmentable(&aString)) // produce 1 hadron True ===============
537  if(!aString.FourQuarkString() && !IsFragmentable(&aString))
538  {
539  #ifdef debug_LUNDfragmentation
540  G4cout<<"Non fragmentable - the string is converted to one hadron "<<G4endl;
541  #endif
542 
543  SetMassCut(10000.*MeV);
544  LeftVector=LightFragmentationTest(&theString);
545  SetMassCut(160.*MeV);
546 
547  LeftVector->operator[](0)->SetFormationTime(theString.GetTimeOfCreation());
548  LeftVector->operator[](0)->SetPosition(theString.GetPosition());
549 
550  if(LeftVector->size() > 1)
551  {
552  // 2 hadrons created from qq-qqbar are stored
553  LeftVector->operator[](1)->SetFormationTime(theString.GetTimeOfCreation());
554  LeftVector->operator[](1)->SetPosition(theString.GetPosition());
555  }
556  return LeftVector;
557  }
558 
559  #ifdef debug_LUNDfragmentation
560  G4cout<<"The string will be fragmented. "<<G4endl;
561  #endif
562 
563  // The string can fragment. At least two particles can be produced.
564  LeftVector =new G4KineticTrackVector;
565  G4KineticTrackVector * RightVector=new G4KineticTrackVector;
566 
567  G4ExcitedString *theStringInCMS=CPExcited(theString);
568 
569  #ifdef debug_LUNDfragmentation
570  G4cout<<"CMS Left mom "<<theStringInCMS->GetLeftParton()->Get4Momentum()<<G4endl;
571  G4cout<<"CMS Right mom "<<theStringInCMS->GetRightParton()->Get4Momentum()<<G4endl;
572  #endif
573 
574  G4LorentzRotation toCms=theStringInCMS->TransformToAlignedCms();
575 
576  #ifdef debug_LUNDfragmentation
577  G4cout<<"aligCMS Left mom "<<theStringInCMS->GetLeftParton()->Get4Momentum()<<G4endl;
578  G4cout<<"aligCMS Right mom "<<theStringInCMS->GetRightParton()->Get4Momentum()<<G4endl;
579  G4cout<<G4endl<<"LUND StringFragm: String Mass " <<theStringInCMS->Get4Momentum().mag()<<" Pz Lab "
580  <<theStringInCMS->Get4Momentum().pz()
581  <<"------------------------------------"<<G4endl;
582  G4cout<<"String ends and Direction "<<theStringInCMS->GetLeftParton()->GetPDGcode()<<" "
583  <<theStringInCMS->GetRightParton()->GetPDGcode()<<" "
584  <<theStringInCMS->GetDirection()<< G4endl;
585  #endif
586 
587  G4bool success = Loop_toFragmentString(theStringInCMS, LeftVector, RightVector);
588 
589  delete theStringInCMS;
590 
591  if ( ! success )
592  {
593  std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
594  LeftVector->clear();
595  std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
596  delete RightVector;
597  return LeftVector;
598  }
599 
600  // Join Left- and RightVector into LeftVector in correct order.
601  while(!RightVector->empty()) /* Loop checking, 07.08.2015, A.Ribon */
602  {
603  LeftVector->push_back(RightVector->back());
604  RightVector->erase(RightVector->end()-1);
605  }
606  delete RightVector;
607 
608  CalculateHadronTimePosition(theString.Get4Momentum().mag(), LeftVector);
609 
610  G4LorentzRotation toObserverFrame(toCms.inverse());
611 
612  G4double TimeOftheStringCreation=theString.GetTimeOfCreation();
613  G4ThreeVector PositionOftheStringCreation(theString.GetPosition());
614 
615  for(size_t C1 = 0; C1 < LeftVector->size(); C1++)
616  {
617  G4KineticTrack* Hadron = LeftVector->operator[](C1);
618  G4LorentzVector Momentum = Hadron->Get4Momentum();
619  //G4cout<<"Hadron "<<Hadron->GetDefinition()->GetParticleName()<<" "<<Momentum<<G4endl;
620  Momentum = toObserverFrame*Momentum;
621  Hadron->Set4Momentum(Momentum);
622 
623  G4LorentzVector Coordinate(Hadron->GetPosition(), Hadron->GetFormationTime());
624  Momentum = toObserverFrame*Coordinate;
625  Hadron->SetFormationTime(TimeOftheStringCreation + Momentum.e() - fermi/c_light);
626  G4ThreeVector aPosition(Momentum.vect());
627  Hadron->SetPosition(PositionOftheStringCreation+aPosition);
628  }
629 
630  return LeftVector;
631 }
G4int GetPDGcode() const
Definition: G4Parton.hh:127
const double C1
const G4ThreeVector & GetPosition() const
G4ExcitedString * CPExcited(const G4ExcitedString &string)
void SetFormationTime(G4double aFormationTime)
const G4LorentzVector & Get4Momentum() const
Definition: G4Parton.hh:143
G4Parton * GetLeftParton(void) const
virtual void SetMassCut(G4double aValue)
const G4ThreeVector & GetPosition() const
G4GLOB_DLL std::ostream G4cout
G4KineticTrackVector * LightFragmentationTest(const G4ExcitedString *const theString)
G4LorentzVector Get4Momentum() const
G4double GetFormationTime() const
double mag() const
bool G4bool
Definition: G4Types.hh:79
G4LorentzRotation TransformToAlignedCms()
void SetPosition(const G4ThreeVector aPosition)
void Set4Momentum(const G4LorentzVector &a4Momentum)
G4Parton * GetRightParton(void) const
G4double GetTimeOfCreation() const
double pz() const
G4int GetDirection(void) const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
HepLorentzRotation inverse() const
double G4double
Definition: G4Types.hh:76
static constexpr double fermi
Definition: G4SIunits.hh:103
const G4LorentzVector & Get4Momentum() const
void CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector *)
float c_light
Definition: hepunit.py:257

Here is the call graph for this function:


The documentation for this class was generated from the following files: