Geant4  10.02.p02
G4QGSMFragmentation.cc
Go to the documentation of this file.
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 // $Id: G4QGSMFragmentation.cc 92025 2015-08-13 14:17:05Z gcosmo $
28 //
29 // -----------------------------------------------------------------------------
30 // GEANT 4 class implementation file
31 //
32 // History: first implementation, Maxim Komogorov, 10-Jul-1998
33 // -----------------------------------------------------------------------------
34 #include "G4QGSMFragmentation.hh"
35 #include "G4PhysicalConstants.hh"
36 #include "Randomize.hh"
37 #include "G4ios.hh"
38 #include "G4FragmentingString.hh"
39 #include "G4DiQuarks.hh"
40 #include "G4Quarks.hh"
41 
42 #include "G4Pow.hh"
43 
44 //#define debug_QGSMfragmentation // Uzhi Oct. 2014
45 
46 // Class G4QGSMFragmentation
47 //****************************************************************************************
48 
50 arho(0.5), aphi(0.), an(-0.5), ala(-0.75), aksi(-1.), alft(0.5)
51  {
52  SetStrangenessSuppression(0.41); // 0.47 0.447 Uzhi 27.09.2014 0.43 last 0.425
53  SetDiquarkSuppression(0.25); // 0.087 std 0.07 Uzhi 0.25 Last
54  SetDiquarkBreakProbability(0.4); // 0.05 std 0.1 Uzhi 27.09.2014
55  }
56 
58  {
59  }
60 
61 //----------------------------------------------------------------------------------------------------------
62 
64 {
65 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
66  G4cout<<G4endl<<"QGSM StringFragm: String Mass "
67  <<theString.Get4Momentum().mag()<<" Pz "
68  <<theString.Get4Momentum().pz()
69  <<"------------------------------------"<<G4endl;
70  G4cout<<"String ends Direct "<<theString.GetLeftParton()->GetPDGcode()<<" "
71  <<theString.GetRightParton()->GetPDGcode()<<" "
72  <<theString.GetDirection()<< G4endl;
73  G4cout<<"Left mom "<<theString.GetLeftParton()->Get4Momentum()<<G4endl;
74  G4cout<<"Right mom "<<theString.GetRightParton()->Get4Momentum()<<G4endl;
75  G4cout<<"Check for Fragmentation "<<G4endl;
76 #endif
77 
78 // Can no longer modify Parameters for Fragmentation.
79  PastInitPhase=true;
80 
81 // check if string has enough mass to fragment...
82 
83  G4KineticTrackVector * LeftVector=LightFragmentationTest(&theString);
84 
85 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
86  if ( LeftVector != 0 ) G4cout<<"Non fragmentable - the string is converted to one hadron "<<G4endl;
87 #endif
88 
89  if ( LeftVector != 0 ) return LeftVector;
90 
91 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
92  G4cout<<"The string will be fragmented. "<<G4endl;
93 #endif
94 
95  LeftVector = new G4KineticTrackVector;
97 
98 // this should work but its only a semi deep copy. %GF G4ExcitedString theStringInCMS(theString);
99  G4ExcitedString *theStringInCMS=CPExcited(theString);
100  G4LorentzRotation toCms=theStringInCMS->TransformToAlignedCms();
101 
102  G4bool success=false, inner_sucess=true;
103  G4int attempt=0;
104  while ( !success && attempt++ < StringLoopInterrupt ) /* Loop checking, 07.08.2015, A.Ribon */
105  {
106 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
107  G4cout<<"Loop_toFrag "<<theStringInCMS->GetLeftParton()->GetPDGcode()<<" "
108  <<theStringInCMS->GetRightParton()->GetPDGcode()<<" "
109  <<theStringInCMS->GetDirection()<< G4endl;
110 #endif
111 
112  G4FragmentingString *currentString=new G4FragmentingString(*theStringInCMS);
113 
114  std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
115  LeftVector->clear();
116  std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
117  RightVector->clear();
118 
119  inner_sucess=true; // set false on failure..
120  const G4int maxNumberOfLoops = 1000;
121  G4int loopCounter = -1;
122  while (! StopFragmenting(currentString) && ++loopCounter < maxNumberOfLoops ) /* Loop checking, 07.08.2015, A.Ribon */
123  { // Split current string into hadron + new string
124 
125 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
126  G4cout<<"The string can fragment. "<<G4endl;;
127 #endif
128  G4FragmentingString *newString=0; // used as output from SplitUp...
129  G4KineticTrack * Hadron=Splitup(currentString,newString);
130 
131  if ( Hadron != 0 ) // && IsFragmentable(newString)) // Closed by Uzhi, Oct. 2014
132  {
133 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
134  G4cout<<"Hadron prod at fragm. "<<Hadron->GetDefinition()->GetParticleName()<<G4endl;
135 #endif
136  if ( currentString->GetDecayDirection() > 0 )
137  LeftVector->push_back(Hadron);
138  else
139  RightVector->push_back(Hadron);
140 
141  delete currentString;
142  currentString=newString;
143 
144  } else {
145 
146 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
147  G4cout<<"abandon ... start from the beginning ---------------"<<G4endl;
148 #endif
149 
150  // abandon ... start from the beginning
151  if (newString) delete newString;
152  inner_sucess=false;
153  break;
154  }
155  }
156  if ( loopCounter >= maxNumberOfLoops ) {
157  inner_sucess=false;
158  }
159 
160  // Split current string into 2 final Hadrons
161 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
162  G4cout<<"Split remaining string into 2 final hadrons."<<G4endl;
163 #endif
164 
165  if ( inner_sucess && //true) // Uzhi -- No Last Splitting
166  SplitLast(currentString,LeftVector, RightVector) )
167  {
168  success=true;
169  }
170  delete currentString;
171  }
172 
173  delete theStringInCMS;
174 
175  if ( ! success )
176  {
177  std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
178  LeftVector->clear();
179  std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
180  delete RightVector;
181  return LeftVector;
182  }
183 
184  // Join Left- and RightVector into LeftVector in correct order.
185  while(!RightVector->empty()) /* Loop checking, 07.08.2015, A.Ribon */
186  {
187  LeftVector->push_back(RightVector->back());
188  RightVector->erase(RightVector->end()-1);
189  }
190  delete RightVector;
191 
192  CalculateHadronTimePosition(theString.Get4Momentum().mag(), LeftVector);
193 
194  G4LorentzRotation toObserverFrame(toCms.inverse());
195 
196  for(size_t C1 = 0; C1 < LeftVector->size(); C1++)
197  {
198  G4KineticTrack* Hadron = LeftVector->operator[](C1);
199  G4LorentzVector Momentum = Hadron->Get4Momentum();
200  Momentum = toObserverFrame*Momentum;
201  Hadron->Set4Momentum(Momentum);
202  G4LorentzVector Coordinate(Hadron->GetPosition(), Hadron->GetFormationTime());
203  Momentum = toObserverFrame*Coordinate;
204  Hadron->SetFormationTime(Momentum.e());
205  G4ThreeVector aPosition(Momentum.vect());
206  Hadron->SetPosition(theString.GetPosition()+aPosition);
207  }
208  return LeftVector;
209 
210 
211 
212 }
213 
214 //----------------------------------------------------------------------------------------------------------
215 
218 {
219 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
220  G4cout<<"GetLightConeZ zmin zmax Parton pHadron "<<zmin<<" "<<zmax<<" "<<PartonEncoding<<" "<<pHadron->GetParticleName()<<G4endl;
221 #endif
222  G4double z;
223  G4double d1, d2, yf;
224  G4double invD1(0.),invD2(0.), r1(0.),r2(0.),r12(0.);
225 
226  G4int absCode = std::abs( PartonEncoding );
227  G4int absHadronCode=std::abs(pHadron->GetPDGEncoding());
228 
229  G4int q1, q2, q3;
230  q1 = absHadronCode/1000; q2 = (absHadronCode % 1000)/100; q3 = (absHadronCode % 100)/10;
231 
232  G4bool StrangeHadron = (q1 == 3) || (q2 == 3) || (q3 == 3);
233 
234  if (absCode < 10)
235  { // A quark fragmentation ----------------------------
236  if(absCode == 1 || absCode == 2)
237  {
238  if(absHadronCode < 1000)
239  { // Meson produced
240  if( !StrangeHadron ) {d1=2.0; d2 = -arho + alft;}
241  else {d1=1.0; d2 = -aphi + alft;}
242  } else
243  { // Baryon produced
244  if( !StrangeHadron ) {d1=0.0; d2 = arho - 2.0*an + alft;}
245  else {d1=0.0; d2 = 2.0*arho - 2.0*an - aphi + alft;}
246  }
247  }
248  else if(absCode == 3)
249  {
250  if(absHadronCode < 1000){d1=1.0 - aphi; d2 = -arho + alft;} // Meson produced s->K + u/d
251  else {d1=1.0 - aphi; d2 = arho - 2.0*an + alft;} // Baryon produced
252 
253  } else throw G4HadronicException(__FILE__, __LINE__, "Unknown PDGencoding in G4QGSMFragmentation::G4LightConeZ");
254 
255 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
256  G4cout<<"d1 d2 "<<d1<<" "<<d2<<G4endl;
257 #endif
258 
259  d1+=1.0; d2+=1.0;
260 
261  invD1=1./d1; invD2=1./d2;
262 
263  const G4int maxNumberOfLoops = 10000;
264  G4int loopCounter = 0;
265  do
266  {
267  r1=G4Pow::GetInstance()->powA(G4UniformRand(),invD1);
268  r2=G4Pow::GetInstance()->powA(G4UniformRand(),invD2);
269  r12=r1+r2;
270  z=r1/r12;
271  } while( ( (r12 > 1.0) || !((zmin <= z)&&(z <= zmax))) && ++loopCounter < maxNumberOfLoops ); /* Loop checking, 07.08.2015, A.Ribon */
272  if ( loopCounter >= maxNumberOfLoops ) {
273  z = 0.5*(zmin + zmax); // Just a value between zmin and zmax, no physics considerations at all!
274  }
275 
276  return z;
277  }
278  else
279  { // A di-quark fragmentation -------------------------
280  if(absCode == 1103 || absCode == 2101 ||
281  absCode == 2203 || absCode == 2103)
282  {
283  if(absHadronCode < 1000) // Meson production
284  {
285  if( !StrangeHadron ) {d1=1.0; d2= arho - 2.0*an + alft;}
286  else {d1=1.0; d2 = 2.*arho - 2.0*an - aphi + alft;}
287  } else // Baryon production
288  {
289  if( !StrangeHadron ) {d1=2.0*(arho - an); d2= -arho + alft;}
290  else {d1=2.0*(arho - an); d2 =-aphi + alft;}
291  }
292 
293 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
294  G4cout<<"d1 d2 "<<d1<<" "<<d2<<G4endl;
295 #endif
296 
297  d1+=1.0; d2+=1.0;
298  invD1=1./d1; invD2=1./d2;
299 
300  const G4int maxNumberOfLoops = 10000;
301  G4int loopCounter = 0;
302  do
303  {
304  r1=G4Pow::GetInstance()->powA(G4UniformRand(),invD1);
305  r2=G4Pow::GetInstance()->powA(G4UniformRand(),invD2);
306  r12=r1+r2;
307  z=r1/r12;
308  } while( ( (r12 > 1.0) || !((zmin <= z)&&(z <= zmax))) && ++loopCounter < maxNumberOfLoops ); /* Loop checking, 07.08.2015, A.Ribon */
309  if ( loopCounter >= maxNumberOfLoops ) {
310  z = 0.5*(zmin + zmax); // Just a value between zmin and zmax, no physics considerations at all!
311  }
312 
313  return z;
314  }
315  else if(absCode == 3101 || absCode == 3103 || // For strange d-quarks
316  absCode == 3201 || absCode == 3203)
317  {
318  d2 = (alft - (2.*ala - arho));
319 
320  }
321  else
322  {
323  d2 = (alft - (2.*aksi - arho));
324  }
325 
326  const G4int maxNumberOfLoops = 1000;
327  G4int loopCounter = 0;
328  do
329  {
330  z = zmin + G4UniformRand() * (zmax - zmin);
331  d1 = (1. - z);
332  yf = G4Pow::GetInstance()->powA(d1, d2);
333  }
334  while( (G4UniformRand() > yf) && ++loopCounter < maxNumberOfLoops ); /* Loop checking, 07.08.2015, A.Ribon */
335  }
336 
337  return z;
338 }
339 //-----------------------------------------------------------------------------------------
340 
342  G4FragmentingString * string, // Uzhi
343  G4FragmentingString * NewString) // Uzhi Oct. 2014
344 {
345  G4double HadronMass = pHadron->GetPDGMass();
346 
347 // G4double MinimalStringMass= FragmentationMass(NewString,&G4HadronBuilder::BuildHighSpin);
348  G4double MinimalStringMass=
349  FragmentationMass(NewString,&G4HadronBuilder::Build); // Uzhi 03.06.2015
350 // FragmentationMass(NewString,&G4HadronBuilder::BuildLowSpin); // Uzhi 03.06.2015
351 // Uzhi 03.06.2015 It would be well to sample randomly HighSpin
352 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
353  G4cout<<"G4QGSMFragmentation::SplitEandP "<<pHadron->GetParticleName()<<G4endl;
354  G4cout<<"String 4 mom, String M "<<string->Get4Momentum()<<" "<<string->Mass()<<G4endl;
355  G4cout<<"HadM MinimalStringMassLeft StringM hM+sM "<<HadronMass<<" "<<MinimalStringMass<<" "
356  <<string->Mass()<<" "<<HadronMass+MinimalStringMass<<G4endl;
357 #endif
358 
359  if(HadronMass + MinimalStringMass > string->Mass())
360  {
361 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
362  G4cout<<"Mass of the string is not sufficient to produce the hadron!"<<G4endl;
363 #endif
364  return 0;
365  }// have to start all over!
366 
367  // calculate and assign hadron transverse momentum component HadronPx andHadronPy
368  G4double StringMT2 = string->MassT2();
369  G4double StringMT = std::sqrt(StringMT2);
370 
371  G4LorentzVector String4Momentum = string->Get4Momentum();
372  String4Momentum.setPz(0.);
373  G4ThreeVector StringPt = String4Momentum.vect();
374 
375  G4ThreeVector HadronPt , RemSysPt;
376  G4double HadronMassT2, ResidualMassT2;
377 
378  //... sample Pt of the hadron
379  G4int attempt=0;
380  do
381  {
382  attempt++; if(attempt > StringLoopInterrupt) return 0;
383 
384  HadronPt =SampleQuarkPt() + string->DecayPt();
385  HadronPt.setZ(0);
386  RemSysPt = StringPt - HadronPt;
387 
388  HadronMassT2 = sqr(HadronMass) + HadronPt.mag2();
389  ResidualMassT2=sqr(MinimalStringMass) + RemSysPt.mag2();
390 
391  } while(std::sqrt(HadronMassT2) + std::sqrt(ResidualMassT2) > StringMT); /* Loop checking, 07.08.2015, A.Ribon */
392 
393  //... sample z to define hadron longitudinal momentum and energy
394  //... but first check the available phase space
395 
396  G4double Pz2 = (sqr(StringMT2 - HadronMassT2 - ResidualMassT2) -
397  4*HadronMassT2 * ResidualMassT2)/4./StringMT2;
398 
399  if(Pz2 < 0 ) {return 0;} // have to start all over!
400 
401  //... then compute allowed z region z_min <= z <= z_max
402 
403  G4double Pz = std::sqrt(Pz2);
404  G4double zMin = (std::sqrt(HadronMassT2+Pz2) - Pz)/std::sqrt(StringMT2);
405  G4double zMax = (std::sqrt(HadronMassT2+Pz2) + Pz)/std::sqrt(StringMT2);
406 
407 /* close by Uzhi, Oct. 2014
408  G4double DecayQuarkMass2 = sqr(string->GetDecayParton()->GetPDGMass());
409  G4double HadronMass2T = sqr(HadronMass) + HadronPt.mag2();
410 
411  if (DecayQuarkMass2 + HadronMass2T >= SmoothParam*(string->Mass2()) )
412  return 0; // have to start all over!
413 
414  //... then compute allowed z region z_min <= z <= z_max
415 
416 // G4double zMin = HadronMass2T/(string->Mass2());
417 // G4double zMax = 1. - DecayQuarkMass2/(string->Mass2());
418 */
419  if (zMin >= zMax) return 0; // have to start all over!
420 
421  G4double z = GetLightConeZ(zMin, zMax,
422  string->GetDecayParton()->GetPDGEncoding(), pHadron,
423  HadronPt.x(), HadronPt.y());
424 
425  //... now compute hadron longitudinal momentum and energy
426  // longitudinal hadron momentum component HadronPz
427 
428  HadronPt.setZ(0.5* string->GetDecayDirection() *
429  (z * string->LightConeDecay() -
430  HadronMassT2/(z * string->LightConeDecay())));
431  G4double HadronE = 0.5* (z * string->LightConeDecay() +
432  HadronMassT2/(z * string->LightConeDecay()));
433 
434  G4LorentzVector * a4Momentum= new G4LorentzVector(HadronPt,HadronE);
435 
436 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
437  G4cout<<"string->GetDecayDirection() string->LightConeDecay() "
438  <<string->GetDecayDirection()<<" "<<string->LightConeDecay()<<G4endl;
439  G4cout<<"HadronPt,HadronE "<<HadronPt<<" "<<HadronE<<G4endl;
440 // G4cout<<"String4Momentum "<<String4Momentum<<G4endl;
441 //G4int Uzhi; G4cin>>Uzhi;
442  G4cout<<"Out of QGSM SplitEandP "<<G4endl;
443 #endif
444 
445  return a4Momentum;
446 }
447 
448 
449 //-----------------------------------------------------------------------------------------
450 
452  G4KineticTrackVector * LeftVector,
453  G4KineticTrackVector * RightVector)
454 {
455  //... perform last cluster decay
456 
457  G4ThreeVector ClusterVel =string->Get4Momentum().boostVector();
458  G4double ResidualMass =string->Mass();
459 
460 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
461  G4cout<<"Split last-----------------------------------------"<<G4endl;
462  G4cout<<"StrMass "<<ResidualMass<<" q's "
463  <<string->GetLeftParton()->GetParticleName()<<" "
464  <<string->GetRightParton()->GetParticleName()<<G4endl;
465 #endif
466 
467  G4double ClusterMassCut = ClusterMass; // Taken from G4VLongitudinalStringDecay
468  G4int cClusterInterrupt = 0;
469  G4ParticleDefinition * LeftHadron, * RightHadron;
470  const G4int maxNumberOfLoops = 1000;
471  G4int loopCounter = 0;
472  do
473  {
474  if (cClusterInterrupt++ >= ClusterLoopInterrupt)
475  {
476  return false;
477  }
478 
479  G4ParticleDefinition * quark = NULL;
480  string->SetLeftPartonStable(); // to query quark contents..
481 
482  if (string->DecayIsQuark() && string->StableIsQuark() )
483  {
484  //... there are quarks on cluster ends
485  LeftHadron= QuarkSplitup(string->GetLeftParton(), quark);
486  } else {
487  //... there is a Diquark on cluster ends
488  G4int IsParticle;
489  if ( string->StableIsQuark() ) {
490  IsParticle=(string->GetLeftParton()->GetPDGEncoding()>0) ? -1 : +1;
491  } else {
492  IsParticle=(string->GetLeftParton()->GetPDGEncoding()>0) ? +1 : -1;
493  }
494 
495  pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
496  quark = QuarkPair.second;
497  LeftHadron=hadronizer->Build(QuarkPair.first, string->GetLeftParton());
498  }
499 
500  RightHadron = hadronizer->Build(string->GetRightParton(), quark);
501 
502  //... repeat procedure, if mass of cluster is too low to produce hadrons
503  //... ClusterMassCut = 0.15*GeV model parameter
504  if ( quark->GetParticleSubType()== "quark" ) {ClusterMassCut = 0.;}
505  else {ClusterMassCut = ClusterMass;}
506  }
507  while ( (ResidualMass <= LeftHadron->GetPDGMass() + RightHadron->GetPDGMass() + ClusterMassCut)
508  && ++loopCounter < maxNumberOfLoops ); /* Loop checking, 07.08.2015, A.Ribon */
509 
510  if ( loopCounter >= maxNumberOfLoops ) {
511  return false;
512  }
513 
514  //... compute hadron momenta and energies
515  G4LorentzVector LeftMom, RightMom;
516  G4ThreeVector Pos;
517  Sample4Momentum(&LeftMom , LeftHadron->GetPDGMass() ,
518  &RightMom, RightHadron->GetPDGMass(), ResidualMass);
519  LeftMom.boost(ClusterVel);
520  RightMom.boost(ClusterVel);
521 
522 #ifdef debug_QGSMfragmentation // Uzhi Oct. 2014
523  G4cout<<LeftHadron->GetParticleName()<<" "<<RightHadron->GetParticleName()<<G4endl;
524  G4cout<<"Left Hadrom P M "<<LeftMom<<" "<<LeftMom.mag()<<G4endl;
525  G4cout<<"Right Hadrom P M "<<RightMom<<" "<<RightMom.mag()<<G4endl;
526 #endif
527 
528  LeftVector->push_back(new G4KineticTrack(LeftHadron, 0, Pos, LeftMom));
529  RightVector->push_back(new G4KineticTrack(RightHadron, 0, Pos, RightMom));
530 
531  return true;
532 
533 }
534 
535 //----------------------------------------------------------------------------------------------------------
536 
538 {
539  return sqr(FragmentationMass(string)+MassCut) <
540  string->Mass2();
541 }
542 
543 //----------------------------------------------------------------------------------------------------------
544 
546 {
547  return
549  string->Get4Momentum().mag2();
550 }
551 
552 //----------------------------------------------------------------------------------------------------------
553 
554 //----------------------------------------------------------------------------------------------------------
555 
557  G4LorentzVector* AntiMom, G4double AntiMass, G4double InitialMass)
558  {
559  G4double r_val = sqr(InitialMass*InitialMass - Mass*Mass - AntiMass*AntiMass) - sqr(2.*Mass*AntiMass);
560  G4double Pabs = (r_val > 0.)? std::sqrt(r_val)/(2.*InitialMass) : 0;
561 
562  //... sample unit vector
563  G4double pz = 1. - 2.*G4UniformRand();
564  G4double st = std::sqrt(1. - pz * pz)*Pabs;
565  G4double phi = 2.*pi*G4UniformRand();
566  G4double px = st*std::cos(phi);
567  G4double py = st*std::sin(phi);
568  pz *= Pabs;
569 
570  Mom->setPx(px); Mom->setPy(py); Mom->setPz(pz);
571  Mom->setE(std::sqrt(Pabs*Pabs + Mass*Mass));
572 
573  AntiMom->setPx(-px); AntiMom->setPy(-py); AntiMom->setPz(-pz);
574  AntiMom->setE (std::sqrt(Pabs*Pabs + AntiMass*AntiMass));
575  }
576 
577 
578 //*********************************************************************************************
579 // Uzhi June 2014 Insert from G4ExcitedStringDecay.cc
580 //-----------------------------------------------------------------------------
581 
584  G4ParticleDefinition *&created)
585 {
586  //... can Diquark break or not?
588  //... Diquark break
589 
590  G4int stableQuarkEncoding = decay->GetPDGEncoding()/1000;
591  G4int decayQuarkEncoding = (decay->GetPDGEncoding()/100)%10;
592 
593  if (G4UniformRand() < 0.5)
594  {
595  G4int Swap = stableQuarkEncoding;
596  stableQuarkEncoding = decayQuarkEncoding;
597  decayQuarkEncoding = Swap;
598  }
599 
600  G4int IsParticle=(decayQuarkEncoding>0) ? -1 : +1;
601  // if we have a quark, we need antiquark)
602 
603  G4double StrSup=GetStrangeSuppress(); // for changing s-sbar production, Uzhi Oct. 2014
604  StrangeSuppress=0.41; // was 0.47
605  pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
606  StrangeSuppress=StrSup;
607 
608  //... Build new Diquark
609  G4int QuarkEncoding=QuarkPair.second->GetPDGEncoding();
610  G4int i10 = std::max(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding));
611  G4int i20 = std::min(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding));
612  G4int spin = (i10 != i20 && G4UniformRand() <= 0.5)? 1 : 3;
613  G4int NewDecayEncoding = -1*IsParticle*(i10 * 1000 + i20 * 100 + spin);
614  created = FindParticle(NewDecayEncoding);
615  G4ParticleDefinition * decayQuark=FindParticle(decayQuarkEncoding);
616  G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decayQuark);
617 
618  return had;
619 // return hadronizer->Build(QuarkPair.first, decayQuark);
620 
621  } else {
622  //... Diquark does not break
623 
624  G4int IsParticle=(decay->GetPDGEncoding()>0) ? +1 : -1;
625  // if we have a diquark, we need quark)
626  G4double StrSup=GetStrangeSuppress(); // for changing s-sbar production, Uzhi Oct. 2014
627  StrangeSuppress=0.41; //0.41; 0.47
628  pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
629  StrangeSuppress=StrSup;
630 
631  created = QuarkPair.second;
632 
633  G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay);
634  return had;
635 // return G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay);
636  }
637 }
638 // Uzhi June 2014 End of the inserting
639 
G4ParticleDefinition * BuildHighSpin(G4ParticleDefinition *black, G4ParticleDefinition *white)
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4ParticleDefinition * GetRightParton(void) const
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
static const G4double d1
G4ParticleDefinition * QuarkSplitup(G4ParticleDefinition *decay, G4ParticleDefinition *&created)
G4int GetPDGcode() const
Definition: G4Parton.hh:124
const double C1
virtual G4bool StopFragmenting(const G4FragmentingString *const string)
CLHEP::Hep3Vector G4ThreeVector
virtual G4bool SplitLast(G4FragmentingString *string, G4KineticTrackVector *LeftVector, G4KineticTrackVector *RightVector)
CLHEP::HepLorentzRotation G4LorentzRotation
const G4ThreeVector & GetPosition() const
G4ExcitedString * CPExcited(const G4ExcitedString &string)
G4double z
Definition: TRTMaterials.hh:39
void SetFormationTime(G4double aFormationTime)
void SetStrangenessSuppression(G4double aValue)
const G4LorentzVector & Get4Momentum() const
Definition: G4Parton.hh:140
G4ParticleDefinition * Build(G4ParticleDefinition *black, G4ParticleDefinition *white)
virtual G4bool IsFragmentable(const G4FragmentingString *const string)
G4Parton * GetLeftParton(void) const
G4double Mass() const
G4ThreeVector SampleQuarkPt(G4double ptMax=-1.)
virtual G4KineticTrackVector * FragmentString(const G4ExcitedString &theString)
const G4String & GetParticleSubType() const
G4double FragmentationMass(const G4FragmentingString *const string, Pcreate build=0, pDefPair *pdefs=0)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4ParticleDefinition * GetDecayParton() const
pDefPair CreatePartonPair(G4int NeedParticle, G4bool AllowDiquarks=true)
#define NewString(str)
const G4ThreeVector & GetPosition() const
G4ParticleDefinition * GetLeftParton(void) const
void SetDiquarkBreakProbability(G4double aValue)
virtual G4LorentzVector * SplitEandP(G4ParticleDefinition *pHadron, G4FragmentingString *string, G4FragmentingString *newString)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
ParticleList decay(Cluster *const c)
Carries out a cluster decay.
G4KineticTrackVector * LightFragmentationTest(const G4ExcitedString *const theString)
G4LorentzVector Get4Momentum() const
G4double GetFormationTime() const
std::pair< G4ParticleDefinition *, G4ParticleDefinition * > pDefPair
bool G4bool
Definition: G4Types.hh:79
G4LorentzRotation TransformToAlignedCms()
G4ParticleDefinition * FindParticle(G4int Encoding)
void SetPosition(const G4ThreeVector aPosition)
G4int GetDecayDirection() const
void Set4Momentum(const G4LorentzVector &a4Momentum)
virtual G4ParticleDefinition * DiQuarkSplitup(G4ParticleDefinition *decay, G4ParticleDefinition *&created)
G4Parton * GetRightParton(void) const
G4double GetPDGMass() const
static const double pi
Definition: G4SIunits.hh:74
virtual void Sample4Momentum(G4LorentzVector *Mom, G4double Mass, G4LorentzVector *AntiMom, G4double AntiMass, G4double InitialMass)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4int GetDirection(void) const
#define G4endl
Definition: G4ios.hh:61
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
virtual G4double GetLightConeZ(G4double zmin, G4double zmax, G4int PartonEncoding, G4ParticleDefinition *pHadron, G4double Px, G4double Py)
const G4LorentzVector & Get4Momentum() const
static const G4double d2
void CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector *)
const G4ParticleDefinition * GetDefinition() const
G4KineticTrack * Splitup(G4FragmentingString *string, G4FragmentingString *&newString)
void SetDiquarkSuppression(G4double aValue)
CLHEP::HepLorentzVector G4LorentzVector