ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4QGSMFragmentation.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4QGSMFragmentation.cc
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 // GEANT 4 class implementation file
30 //
31 // History: first implementation, Maxim Komogorov, 10-Jul-1998
32 // -----------------------------------------------------------------------------
33 #include "G4QGSMFragmentation.hh"
34 #include "G4PhysicalConstants.hh"
35 #include "G4SystemOfUnits.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
45 
46 // Class G4QGSMFragmentation
47 //****************************************************************************************
48 
50 {
51  MassCut = 0.35*GeV;
52 
53  SetStrangenessSuppression((1.0 - 0.16)/2.);
54 
55  // For the time being, set to 0.0 the probabilities for c-cbar and b-bbar creation.
56  SetProbCCbar(0.0); //(0.0033); // According to O.I. Piskunova Yad. Fiz. 56 (1993) 1094
57  SetProbBBbar(0.0); //(5.0e-5); // According to O.I. Piskunova Yad. Fiz. 56 (1993) 1094
58 
61 
62  SetMinMasses();
63 
64  arho = 0.5; // alpha_rho0
65  aphi = 0.0; // alpha_fi
66  aJPs =-2.2; // alpha_J/Psi
67  aUps =-8.0; // alpha_Y ??? O. Piskunova Yad. Phys. 56 (1993) 1094.
68 
69  aksi =-1.0;
70  alft = 0.5; // 2 * alpha'_R *<Pt^2>
71 
72  an = -0.5 ;
73  ala = -0.75; // an - arho/2 + aphi/2
74  alaC = an - arho/2.0 + aJPs/2.0;
75  alaB = an - arho/2.0 + aUps/2.0;
76  aXi = 0.0; // ??
77  aXiC = 0.0; // ??
78  aXiB = 0.0; // ??
79  aXiCC = 0.0; // ??
80  aXiCB = 0.0; // ??
81  aXiBB = 0.0; // ??
82 
83  SetFFq2q();
84  SetFFq2qq();
85  SetFFqq2q();
86  SetFFqq2qq();
87  // d u s c b
88  G4int Index[5][5] = { { 0, 1, 2, 3, 4 }, // d
89  { 1, 5, 6, 7, 8 }, // u
90  { 2, 6, 9, 10, 11 }, // s
91  { 3, 7, 10, 12, 13 }, // c
92  { 4, 8, 11, 13, 14 } }; // b
93  for (G4int i = 0; i < 5; i++ ) {
94  for ( G4int j = 0; j < 5; j++ ) {
95  IndexDiQ[i][j] = Index[i][j];
96  }
97  };
98 }
99 
101 {}
102 
103 //----------------------------------------------------------------------------------------------------------
104 
106 {
107 
108  G4FragmentingString aString(theString);
109  SetMinimalStringMass(&aString);
110 
111  #ifdef debug_QGSMfragmentation
112  G4cout<<G4endl<<"QGSM StringFragm: String Mass "
113  <<theString.Get4Momentum().mag()<<" Pz "
114  <<theString.Get4Momentum().pz()
115  <<"------------------------------------"<<G4endl;
116  G4cout<<"String ends Direct "<<theString.GetLeftParton()->GetPDGcode()<<" "
117  <<theString.GetRightParton()->GetPDGcode()<<" "
118  <<theString.GetDirection()<< G4endl;
119  G4cout<<"Left mom "<<theString.GetLeftParton()->Get4Momentum()<<G4endl;
120  G4cout<<"Right mom "<<theString.GetRightParton()->Get4Momentum()<<G4endl;
121  G4cout<<"Check for Fragmentation "<<G4endl;
122  #endif
123 
124  // Can no longer modify Parameters for Fragmentation.
125  PastInitPhase=true;
126 
127  // Check if string has enough mass to fragment...
128  G4KineticTrackVector * LeftVector=NULL;
129 
130  if ( !IsItFragmentable(&aString) ) {
131  LeftVector=ProduceOneHadron(&theString);
132 
133  #ifdef debug_QGSMfragmentation
134  if ( LeftVector != 0 ) G4cout<<"Non fragmentable - the string is converted to one hadron "<<G4endl;
135  #endif
136 
137  if ( LeftVector == nullptr ) LeftVector = new G4KineticTrackVector;
138  return LeftVector;
139  }
140 
141  #ifdef debug_QGSMfragmentation
142  G4cout<<"The string will be fragmented. "<<G4endl;
143  #endif
144 
145  LeftVector = new G4KineticTrackVector;
146  G4KineticTrackVector * RightVector=new G4KineticTrackVector;
147 
148  G4ExcitedString *theStringInCMS=CopyExcited(theString);
149  G4LorentzRotation toCms=theStringInCMS->TransformToAlignedCms();
150 
151  G4bool success=false, inner_sucess=true;
152  G4int attempt=0;
153  while ( !success && attempt++ < StringLoopInterrupt ) /* Loop checking, 07.08.2015, A.Ribon */
154  {
155  #ifdef debug_QGSMfragmentation
156  G4cout<<"Loop_toFrag "<<theStringInCMS->GetLeftParton()->GetPDGcode()<<" "
157  <<theStringInCMS->GetRightParton()->GetPDGcode()<<" "
158  <<theStringInCMS->GetDirection()<< G4endl;
159  #endif
160 
161  G4FragmentingString *currentString=new G4FragmentingString(*theStringInCMS);
162 
163  std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
164  LeftVector->clear();
165  std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
166  RightVector->clear();
167 
168  inner_sucess=true; // set false on failure..
169  const G4int maxNumberOfLoops = 1000;
170  G4int loopCounter = -1;
171  while (! StopFragmenting(currentString) && ++loopCounter < maxNumberOfLoops ) /* Loop checking, 07.08.2015, A.Ribon */
172  { // Split current string into hadron + new string
173 
174  #ifdef debug_QGSMfragmentation
175  G4cout<<"The string can fragment. "<<G4endl;;
176  #endif
177  G4FragmentingString *newString=0; // used as output from SplitUp...
178  G4KineticTrack * Hadron=Splitup(currentString,newString);
179 
180  if ( Hadron != 0 )
181  {
182  #ifdef debug_QGSMfragmentation
183  G4cout<<"Hadron prod at fragm. "<<Hadron->GetDefinition()->GetParticleName()<<G4endl;
184  #endif
185  // To close the production of hadrons at fragmentation stage
186  if ( currentString->GetDecayDirection() > 0 )
187  LeftVector->push_back(Hadron);
188  else
189  RightVector->push_back(Hadron);
190 
191  delete currentString;
192  currentString=newString;
193 
194  } else {
195 
196  #ifdef debug_QGSMfragmentation
197  G4cout<<"abandon ... start from the beginning ---------------"<<G4endl;
198  #endif
199 
200  // Abandon ... start from the beginning
201  if (newString) delete newString;
202  inner_sucess=false;
203  break;
204  }
205  }
206  if ( loopCounter >= maxNumberOfLoops ) {
207  inner_sucess=false;
208  }
209 
210  // Split current string into 2 final Hadrons
211  #ifdef debug_QGSMfragmentation
212  G4cout<<"Split remaining string into 2 final hadrons."<<G4endl;
213  #endif
214  // To the close production of hadrons at last string decay
215  if ( inner_sucess &&
216  SplitLast(currentString,LeftVector, RightVector) )
217  {
218  success=true;
219  }
220  delete currentString;
221  }
222 
223  delete theStringInCMS;
224 
225  if ( ! success )
226  {
227  std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
228  LeftVector->clear();
229  std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
230  delete RightVector;
231  return LeftVector;
232  }
233 
234  // Join Left- and RightVector into LeftVector in correct order.
235  while(!RightVector->empty()) /* Loop checking, 07.08.2015, A.Ribon */
236  {
237  LeftVector->push_back(RightVector->back());
238  RightVector->erase(RightVector->end()-1);
239  }
240  delete RightVector;
241 
242  CalculateHadronTimePosition(theString.Get4Momentum().mag(), LeftVector);
243 
244  G4LorentzRotation toObserverFrame(toCms.inverse());
245 
246  for (size_t C1 = 0; C1 < LeftVector->size(); C1++)
247  {
248  G4KineticTrack* Hadron = LeftVector->operator[](C1);
249  G4LorentzVector Momentum = Hadron->Get4Momentum();
250  Momentum = toObserverFrame*Momentum;
251  Hadron->Set4Momentum(Momentum);
252  G4LorentzVector Coordinate(Hadron->GetPosition(), Hadron->GetFormationTime());
253  Momentum = toObserverFrame*Coordinate;
254  Hadron->SetFormationTime(Momentum.e());
255  G4ThreeVector aPosition(Momentum.vect());
256  Hadron->SetPosition(theString.GetPosition()+aPosition);
257  }
258  return LeftVector;
259 }
260 
261 //----------------------------------------------------------------------------------------------------------
262 
264 {
265  return sqr( PossibleHadronMass(string) + MassCut ) < string->Mass2();
266 }
267 
268 //----------------------------------------------------------------------------------------------------------
269 
271 {
272  SetMinimalStringMass(string);
273  if ( MinimalStringMass < 0.0 ) return true;
274 
275  if (string->IsAFourQuarkString())
276  {
277  return G4UniformRand() < G4Exp(-0.005*(string->Mass() - MinimalStringMass));
278  } else {
279  G4bool Result = G4UniformRand() <
280  G4Exp(-0.66e-6*(string->Mass()*string->Mass() - MinimalStringMass*MinimalStringMass));
281  // G4bool Result = string->Mass() < MinimalStringMass + 150.*MeV*G4UniformRand(); // a'la LUND
282 
283  #ifdef debug_QGSMfragmentation
284  G4cout<<"StopFragmenting MinimalStringMass string->Mass() "<<MinimalStringMass<<" "<<string->Mass()<<G4endl;
285  G4cout<<"StopFragmenting - Yes/No "<<Result<<G4endl;
286  #endif
287  return Result;
288  }
289 }
290 
291 //-----------------------------------------------------------------------------
292 
294  G4FragmentingString *&newString )
295 {
296  #ifdef debug_QGSMfragmentation
297  G4cout<<G4endl;
298  G4cout<<"Start SplitUP (G4VLongitudinalStringDecay) ========================="<<G4endl;
299  G4cout<<"String partons: " <<string->GetLeftParton()->GetPDGEncoding()<<" "
300  <<string->GetRightParton()->GetPDGEncoding()<<" "
301  <<"Direction " <<string->GetDecayDirection()<<G4endl;
302  #endif
303 
304  //... random choice of string end to use for creating the hadron (decay)
305  G4int SideOfDecay = (G4UniformRand() < 0.5)? 1: -1;
306  if (SideOfDecay < 0)
307  {
308  string->SetLeftPartonStable();
309  } else
310  {
311  string->SetRightPartonStable();
312  }
313 
314  G4ParticleDefinition *newStringEnd;
315  G4ParticleDefinition * HadronDefinition;
316  if (string->DecayIsQuark())
317  {
318  G4double ProbDqADq = GetDiquarkSuppress();
319 
320  G4int NumberOfpossibleBaryons = 2;
321 
322  if (string->GetLeftParton()->GetParticleSubType() != "quark") NumberOfpossibleBaryons++;
323  if (string->GetRightParton()->GetParticleSubType() != "quark") NumberOfpossibleBaryons++;
324 
325  G4double ActualProb = ProbDqADq ;
326  ActualProb *= (1.0-G4Exp(2.0*(1.0 - string->Mass()/(NumberOfpossibleBaryons*1400.0))));
327 
328  SetDiquarkSuppression(ActualProb);
329 
330  HadronDefinition= QuarkSplitup(string->GetDecayParton(), newStringEnd);
331 
332  SetDiquarkSuppression(ProbDqADq);
333  } else {
334  HadronDefinition= DiQuarkSplitup(string->GetDecayParton(), newStringEnd);
335  }
336 
337  if ( HadronDefinition == NULL ) return NULL;
338 
339  #ifdef debug_QGSMfragmentation
340  G4cout<<"The parton "<<string->GetDecayParton()->GetPDGEncoding()<<" "
341  <<" produces hadron "<<HadronDefinition->GetParticleName()
342  <<" and is transformed to "<<newStringEnd->GetPDGEncoding()<<G4endl;
343  G4cout<<"The side of the string decay Left/Right (1/-1) "<<SideOfDecay<<G4endl;
344  #endif
345  // create new String from old, ie. keep Left and Right order, but replace decay
346 
347  newString=new G4FragmentingString(*string,newStringEnd); // To store possible
348  // quark containt of new string
349 
350  #ifdef debug_QGSMfragmentation
351  G4cout<<"An attempt to determine its energy (SplitEandP)"<<G4endl;
352  #endif
353  G4LorentzVector* HadronMomentum=SplitEandP(HadronDefinition, string, newString);
354 
355  delete newString; newString=0;
356 
357  G4KineticTrack * Hadron =0;
358  if ( HadronMomentum != 0 ) {
359 
360  #ifdef debug_QGSMfragmentation
361  G4cout<<"The attempt was successful"<<G4endl;
362  #endif
364  Hadron = new G4KineticTrack(HadronDefinition, 0,Pos, *HadronMomentum);
365 
366  newString=new G4FragmentingString(*string,newStringEnd,HadronMomentum);
367 
368  delete HadronMomentum;
369  }
370  else
371  {
372  #ifdef debug_QGSMfragmentation
373  G4cout<<"The attempt was not successful !!!"<<G4endl;
374  #endif
375  }
376 
377  #ifdef debug_VStringDecay
378  G4cout<<"End SplitUP (G4VLongitudinalStringDecay) ====================="<<G4endl;
379  #endif
380 
381  return Hadron;
382 }
383 
384 //-----------------------------------------------------------------------------
385 
387  G4ParticleDefinition *&created )
388 {
389  //... can Diquark break or not?
390  if (G4UniformRand() < DiquarkBreakProb ) //... Diquark break
391  {
392  G4int stableQuarkEncoding = decay->GetPDGEncoding()/1000;
393  G4int decayQuarkEncoding = (decay->GetPDGEncoding()/100)%10;
394 
395  if (G4UniformRand() < 0.5)
396  {
397  G4int Swap = stableQuarkEncoding;
398  stableQuarkEncoding = decayQuarkEncoding;
399  decayQuarkEncoding = Swap;
400  }
401 
402  G4int IsParticle=(decayQuarkEncoding>0) ? -1 : +1; // if we have a quark, we need antiquark
403 
404  G4double StrSup=GetStrangeSuppress();
405  SetStrangenessSuppression((1.0 - 0.07)/2.);
406  pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
408 
409  //... Build new Diquark
410  G4int QuarkEncoding=QuarkPair.second->GetPDGEncoding();
411  G4int i10 = std::max(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding));
412  G4int i20 = std::min(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding));
413  G4int spin = (i10 != i20 && G4UniformRand() <= 0.5)? 1 : 3;
414  G4int NewDecayEncoding = -1*IsParticle*(i10 * 1000 + i20 * 100 + spin);
415  created = FindParticle(NewDecayEncoding);
416  G4ParticleDefinition * decayQuark=FindParticle(decayQuarkEncoding);
417  G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decayQuark);
418 
419  DecayQuark = decayQuarkEncoding;
420  NewQuark = QuarkPair.first->GetPDGEncoding();
421 
422  return had;
423 
424  } else { //... Diquark does not break
425 
426  G4int IsParticle=(decay->GetPDGEncoding()>0) ? +1 : -1; // if we have a diquark, we need quark)
427  G4double StrSup=GetStrangeSuppress(); // for changing s-sbar production
428  SetStrangenessSuppression((1.0 - 0.07)/2.);
429  pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
431 
432  created = QuarkPair.second;
433 
434  DecayQuark = decay->GetPDGEncoding();
435  NewQuark = created->GetPDGEncoding();
436 
437  G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay);
438  return had;
439  }
440 }
441 
442 //-----------------------------------------------------------------------------------------
443 
445  G4FragmentingString * string,
447 {
448  G4double HadronMass = pHadron->GetPDGMass();
449 
450  SetMinimalStringMass(NewString);
451 
452  if ( MinimalStringMass < 0.0 ) return nullptr;
453 
454  #ifdef debug_QGSMfragmentation
455  G4cout<<"G4QGSMFragmentation::SplitEandP "<<pHadron->GetParticleName()<<G4endl;
456  G4cout<<"String 4 mom, String M "<<string->Get4Momentum()<<" "<<string->Mass()<<G4endl;
457  G4cout<<"HadM MinimalStringMassLeft StringM hM+sM "<<HadronMass<<" "<<MinimalStringMass<<" "
458  <<string->Mass()<<" "<<HadronMass+MinimalStringMass<<G4endl;
459  #endif
460 
461  if (HadronMass + MinimalStringMass > string->Mass())
462  {
463  #ifdef debug_QGSMfragmentation
464  G4cout<<"Mass of the string is not sufficient to produce the hadron!"<<G4endl;
465  #endif
466  return 0;
467  } // have to start all over!
468 
469  // calculate and assign hadron transverse momentum component HadronPx andHadronPy
470  G4double StringMT2 = string->MassT2();
471  G4double StringMT = std::sqrt(StringMT2);
472 
473  G4LorentzVector String4Momentum = string->Get4Momentum();
474  String4Momentum.setPz(0.);
475  G4ThreeVector StringPt = String4Momentum.vect();
476 
477  G4ThreeVector HadronPt , RemSysPt;
478  G4double HadronMassT2, ResidualMassT2;
479 
480  //... sample Pt of the hadron
481  G4int attempt=0;
482  do
483  {
484  attempt++; if (attempt > StringLoopInterrupt) return 0;
485 
486  HadronPt =SampleQuarkPt() + string->DecayPt();
487  HadronPt.setZ(0);
488  RemSysPt = StringPt - HadronPt;
489 
490  HadronMassT2 = sqr(HadronMass) + HadronPt.mag2();
491  ResidualMassT2=sqr(MinimalStringMass) + RemSysPt.mag2();
492 
493  } while (std::sqrt(HadronMassT2) + std::sqrt(ResidualMassT2) > StringMT); /* Loop checking, 07.08.2015, A.Ribon */
494 
495  //... sample z to define hadron longitudinal momentum and energy
496  //... but first check the available phase space
497 
498  G4double Pz2 = (sqr(StringMT2 - HadronMassT2 - ResidualMassT2) -
499  4*HadronMassT2 * ResidualMassT2)/4./StringMT2;
500 
501  if ( Pz2 < 0 ) {return 0;} // have to start all over!
502 
503  //... then compute allowed z region z_min <= z <= z_max
504 
505  G4double Pz = std::sqrt(Pz2);
506  G4double zMin = (std::sqrt(HadronMassT2+Pz2) - Pz)/std::sqrt(StringMT2);
507  G4double zMax = (std::sqrt(HadronMassT2+Pz2) + Pz)/std::sqrt(StringMT2);
508 
509  if (zMin >= zMax) return 0; // have to start all over!
510 
511  G4double z = GetLightConeZ(zMin, zMax,
512  string->GetDecayParton()->GetPDGEncoding(), pHadron,
513  HadronPt.x(), HadronPt.y());
514 
515  //... now compute hadron longitudinal momentum and energy
516  // longitudinal hadron momentum component HadronPz
517 
518  HadronPt.setZ( 0.5* string->GetDecayDirection() *
519  (z * string->LightConeDecay() -
520  HadronMassT2/(z * string->LightConeDecay())) );
521  G4double HadronE = 0.5* (z * string->LightConeDecay() +
522  HadronMassT2/(z * string->LightConeDecay()) );
523 
524  G4LorentzVector * a4Momentum= new G4LorentzVector(HadronPt,HadronE);
525 
526  #ifdef debug_QGSMfragmentation
527  G4cout<<"string->GetDecayDirection() string->LightConeDecay() "
528  <<string->GetDecayDirection()<<" "<<string->LightConeDecay()<<G4endl;
529  G4cout<<"HadronPt,HadronE "<<HadronPt<<" "<<HadronE<<G4endl;
530  G4cout<<"Out of QGSM SplitEandP "<<G4endl;
531  #endif
532 
533  return a4Momentum;
534 }
535 
536 //----------------------------------------------------------------------------------------------------------
537 
539  G4ParticleDefinition* /* pHadron */, G4double , G4double )
540 {
541  #ifdef debug_QGSMfragmentation
542  G4cout<<"GetLightConeZ zmin zmax Parton pHadron "<<zmin<<" "<<zmax<<" "<< /* PartonEncoding */
543  <<" "<</* pHadron->GetParticleName() */ <<G4endl;
544  #endif
545 
546  G4double z(0.);
547  G4int DiQold(0), DiQnew(0);
548  G4double d1(-1.0), d2(0.);
549  G4double invD1(0.),invD2(0.), r1(0.),r2(0.),r12(0.);
550 
551  G4int absDecayQuarkCode = std::abs( DecayQuark );
552  G4int absNewQuarkCode = std::abs( NewQuark );
553 
554  G4int q1(0), q2(0);
555  // q1 = absDecayQuarkCode/1000; q2 = (absDecayQuarkCode % 1000)/100;
556 
557  G4int qA(0), qB(0);
558  // qA = absNewQuarkCode/1000; qB = (absNewQuarkCode % 1000)/100;
559 
560  if ( (absDecayQuarkCode < 6) && (absNewQuarkCode < 6) ) {
561  d1 = FFq2q[absDecayQuarkCode-1][absNewQuarkCode-1][0]; d2 = FFq2q[absDecayQuarkCode-1][absNewQuarkCode-1][1];
562  }
563 
564  if ( (absDecayQuarkCode < 6) && (absNewQuarkCode > 6) ) {
565  qA = absNewQuarkCode/1000; qB = (absNewQuarkCode % 1000)/100; DiQnew = IndexDiQ[qA-1][qB-1];
566  d1 = FFq2qq[absDecayQuarkCode-1][DiQnew][0]; d2 = FFq2q[absDecayQuarkCode-1][DiQnew][1];
567  }
568 
569  if ( (absDecayQuarkCode > 6) && (absNewQuarkCode < 6) ) {
570  q1 = absDecayQuarkCode/1000; q2 = (absDecayQuarkCode % 1000)/100; DiQold = IndexDiQ[q1-1][q2-1];
571  d1 = FFqq2q[DiQold][absNewQuarkCode-1][0]; d2 = FFqq2q[DiQold][absNewQuarkCode-1][1];
572  }
573 
574  if ( d1 < 0. ) {
575  q1 = absDecayQuarkCode/1000; q2 = (absDecayQuarkCode % 1000)/100; DiQold = IndexDiQ[q1-1][q2-1];
576  d1 = FFqq2qq[DiQold][absNewQuarkCode-1][0]; d2 = FFqq2qq[DiQold][absNewQuarkCode-1][1];
577  }
578 
579  d1+=1.0; d2+=1.0;
580 
581  invD1=1./d1; invD2=1./d2;
582 
583  const G4int maxNumberOfLoops = 10000;
584  G4int loopCounter = 0;
585  do // Jong's algorithm
586  {
589  r12=r1+r2;
590  z=r1/r12;
591  } while ( ( (r12 > 1.0) || !((zmin <= z)&&(z <= zmax))) &&
592  ++loopCounter < maxNumberOfLoops ); /* Loop checking, 07.08.2015, A.Ribon */
593 
594  if ( loopCounter >= maxNumberOfLoops ) {
595  z = 0.5*(zmin + zmax); // Just a value between zmin and zmax, no physics considerations at all!
596  }
597 
598  return z;
599 }
600 
601 //-----------------------------------------------------------------------------------------
602 
604  G4KineticTrackVector * LeftVector,
605  G4KineticTrackVector * RightVector)
606 {
607  //... perform last cluster decay
608 
609  G4ThreeVector ClusterVel =string->Get4Momentum().boostVector();
610  G4double ResidualMass =string->Mass();
611 
612  #ifdef debug_QGSMfragmentation
613  G4cout<<"Split last-----------------------------------------"<<G4endl;
614  G4cout<<"StrMass "<<ResidualMass<<" q's "
615  <<string->GetLeftParton()->GetParticleName()<<" "
616  <<string->GetRightParton()->GetParticleName()<<G4endl;
617  #endif
618 
619  G4ParticleDefinition *LeftHadron = nullptr;
620  G4ParticleDefinition *RightHadron = nullptr;
621  const G4int maxNumberOfLoops = 1000;
622  G4int loopCounter = 0;
623  G4bool isOK = false;
624 
625  G4double LeftHadronMass(0.); G4double RightHadronMass(0.);
626  do
627  {
628  LeftHadronMass = -MaxMass; RightHadronMass = -MaxMass;
629 
630  G4ParticleDefinition * quark = nullptr;
631  string->SetLeftPartonStable(); // to query quark contents..
632 
633  if (string->DecayIsQuark() && string->StableIsQuark() )
634  {
635  //... there are quarks on cluster ends
636 
637  G4int IsParticle=(string->GetLeftParton()->GetPDGEncoding()>0) ? -1 : +1;
638  // if we have a quark, we need antiquark or diquark
639 
640  pDefPair QuarkPair = CreatePartonPair(IsParticle);
641  quark = QuarkPair.second;
642 
643  LeftHadron= hadronizer->BuildLowSpin(QuarkPair.first, string->GetLeftParton());
644  } else {
645  //... there is a Diquark on cluster ends
646  G4int IsParticle;
647  if ( string->StableIsQuark() ) {
648  IsParticle=(string->GetLeftParton()->GetPDGEncoding()>0) ? -1 : +1;
649  } else {
650  IsParticle=(string->GetLeftParton()->GetPDGEncoding()>0) ? +1 : -1;
651  }
652 
653  //G4double ProbSaS = 1.0 - 2.0 * GetStrangeSuppress();
654  //G4double ActualProb = ProbSaS * 1.4;
655  //SetStrangenessSuppression((1.0-ActualProb)/2.0);
656 
657  pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
658  //SetStrangenessSuppression((1.0-ProbSaS)/2.0);
659  quark = QuarkPair.second;
660  LeftHadron=hadronizer->BuildLowSpin(QuarkPair.first, string->GetLeftParton());
661  }
662 
663  if ( LeftHadron != nullptr ) {
664  RightHadron = hadronizer->BuildLowSpin(string->GetRightParton(), quark);
665 
666  if ( RightHadron != nullptr ) {
667 
668  LeftHadronMass = LeftHadron->GetPDGMass();
669  RightHadronMass = RightHadron->GetPDGMass();
670  isOK = (ResidualMass > LeftHadronMass + RightHadronMass);
671  }
672  }
673  ++loopCounter;
674  if ( loopCounter >= maxNumberOfLoops ) {
675  return false;
676  }
677  //... repeat procedure, if mass of cluster is too low to produce hadrons
678  //... ClusterMassCut = 0.15*GeV model parameter
679  }
680  while (isOK == false);
681  /* Loop checking, 07.08.2015, A.Ribon */
682 
683  //... compute hadron momenta and energies
684  G4LorentzVector LeftMom, RightMom;
686  Sample4Momentum(&LeftMom , LeftHadron->GetPDGMass() ,
687  &RightMom, RightHadron->GetPDGMass(), ResidualMass);
688  LeftMom.boost(ClusterVel);
689  RightMom.boost(ClusterVel);
690 
691  #ifdef debug_QGSMfragmentation
692  G4cout<<LeftHadron->GetParticleName()<<" "<<RightHadron->GetParticleName()<<G4endl;
693  G4cout<<"Left Hadrom P M "<<LeftMom<<" "<<LeftMom.mag()<<G4endl;
694  G4cout<<"Right Hadrom P M "<<RightMom<<" "<<RightMom.mag()<<G4endl;
695  #endif
696 
697  LeftVector->push_back(new G4KineticTrack(LeftHadron, 0, Pos, LeftMom));
698  RightVector->push_back(new G4KineticTrack(RightHadron, 0, Pos, RightMom));
699 
700  return true;
701 }
702 
703 //----------------------------------------------------------------------------------------------------------
704 
706  G4LorentzVector* AntiMom, G4double AntiMass, G4double InitialMass)
707 {
708  G4double r_val = sqr(InitialMass*InitialMass - Mass*Mass - AntiMass*AntiMass) - sqr(2.*Mass*AntiMass);
709  G4double Pabs = (r_val > 0.)? std::sqrt(r_val)/(2.*InitialMass) : 0;
710 
711  //... sample unit vector
712  G4double pz = 1. - 2.*G4UniformRand();
713  G4double st = std::sqrt(1. - pz * pz)*Pabs;
714  G4double phi = 2.*pi*G4UniformRand();
715  G4double px = st*std::cos(phi);
716  G4double py = st*std::sin(phi);
717  pz *= Pabs;
718 
719  Mom->setPx(px); Mom->setPy(py); Mom->setPz(pz);
720  Mom->setE(std::sqrt(Pabs*Pabs + Mass*Mass));
721 
722  AntiMom->setPx(-px); AntiMom->setPy(-py); AntiMom->setPz(-pz);
723  AntiMom->setE (std::sqrt(Pabs*Pabs + AntiMass*AntiMass));
724 }
725 
726 //----------------------------------------------------------------------------------------------------------
727 
728 void G4QGSMFragmentation::SetFFq2q() // q-> q' + Meson (q anti q')
729 {
730  for (G4int i=0; i < 5; i++) {
731  FFq2q[i][0][0] = 2.0 ; FFq2q[i][0][1] = -arho + alft; // q->d + (q dbar) Pi0, Eta, Eta', Rho0, omega
732  FFq2q[i][1][0] = 2.0 ; FFq2q[i][1][1] = -arho + alft; // q->u + (q ubar) Pi-, Rho-
733  FFq2q[i][2][0] = 1.0 ; FFq2q[i][2][1] = -aphi + alft; // q->s + (q sbar) K0, K*0
734  FFq2q[i][3][0] = 1.0 ; FFq2q[i][3][1] = -aJPs + alft; // q->c + (q+cbar) D-, D*-
735  FFq2q[i][4][0] = 1.0 ; FFq2q[i][4][1] = -aUps + alft; // q->b + (q bbar) EtaB, Upsilon
736  }
737 }
738 
739 //----------------------------------------------------------------------------------------------------------
740 
741 void G4QGSMFragmentation::SetFFq2qq() // q-> anti (q1'q2') + Baryon (q + q1 + q2)
742 {
743  for (G4int i=0; i < 5; i++) {
744  FFq2qq[i][ 0][0] = 0.0 ; FFq2qq[i][ 0][1] = arho - 2.0*an + alft; // q->dd bar + (q dd)
745  FFq2qq[i][ 1][0] = 0.0 ; FFq2qq[i][ 1][1] = arho - 2.0*an + alft; // q->ud bar + (q ud)
746  FFq2qq[i][ 2][0] = 0.0 ; FFq2qq[i][ 2][1] = arho - 2.0*ala + alft; // q->sd bar + (q sd)
747  FFq2qq[i][ 3][0] = 0.0 ; FFq2qq[i][ 3][1] = arho - 2.0*alaC + alft; // q->cd bar + (q cd)
748  FFq2qq[i][ 4][0] = 0.0 ; FFq2qq[i][ 4][1] = arho - 2.0*alaB + alft; // q->bd bar + (q bd)
749  FFq2qq[i][ 5][0] = 0.0 ; FFq2qq[i][ 5][1] = arho - 2.0*an + alft; // q->uu bar + (q uu)
750  FFq2qq[i][ 6][0] = 0.0 ; FFq2qq[i][ 6][1] = arho - 2.0*ala + alft; // q->su bar + (q su)
751  FFq2qq[i][ 7][0] = 0.0 ; FFq2qq[i][ 7][1] = arho - 2.0*alaC + alft; // q->cu bar + (q cu)
752  FFq2qq[i][ 8][0] = 0.0 ; FFq2qq[i][ 8][1] = arho - 2.0*alaB + alft; // q->bu bar + (q bu)
753  FFq2qq[i][ 9][0] = 0.0 ; FFq2qq[i][ 9][1] = arho - 2.0*aXi + alft; // q->ss bar + (q ss)
754  FFq2qq[i][10][0] = 0.0 ; FFq2qq[i][10][1] = arho - 2.0*aXiC + alft; // q->cs bar + (q cs)
755  FFq2qq[i][11][0] = 0.0 ; FFq2qq[i][11][1] = arho - 2.0*aXiB + alft; // q->bs bar + (q bc)
756  FFq2qq[i][12][0] = 0.0 ; FFq2qq[i][12][1] = arho - 2.0*aXiCC + alft; // q->cc bar + (q cc)
757  FFq2qq[i][13][0] = 0.0 ; FFq2qq[i][13][1] = arho - 2.0*aXiCB + alft; // q->bc bar + (q bc)
758  FFq2qq[i][14][0] = 0.0 ; FFq2qq[i][14][1] = arho - 2.0*aXiBB + alft; // q->bb bar + (q bb)
759  }
760 }
761 
762 //----------------------------------------------------------------------------------------------------------
763 
764 void G4QGSMFragmentation::SetFFqq2q() // q1q2-> anti(q') + Baryon (q1 + q2 + q')
765 {
766  for (G4int i=0; i < 15; i++) {
767  FFqq2q[i][0][0] = 2.0*(arho - an); FFqq2q[i][0][1] = -arho + alft;
768  FFqq2q[i][1][0] = 2.0*(arho - an); FFqq2q[i][1][1] = -arho + alft;
769  FFqq2q[i][2][0] = 2.0*(arho - an); FFqq2q[i][2][1] = -aphi + alft;
770  FFqq2q[i][3][0] = 2.0*(arho - an); FFqq2q[i][3][1] = -aJPs + alft;
771  FFqq2q[i][4][0] = 2.0*(arho - an); FFqq2q[i][4][1] = -aUps + alft;
772  }
773 }
774 
775 //----------------------------------------------------------------------------------------------------------
776 
777 void G4QGSMFragmentation::SetFFqq2qq() // q1(q2)-> q'(q2) + Meson(q1 anti q')
778 {
779  for (G4int i=0; i < 15; i++) {
780  FFqq2qq[i][0][0] = 0. ; FFqq2qq[i][0][1] = 2.0*arho - 2.0*an -arho + alft; // dd -> dd + Pi0 (d d bar)
781  FFqq2qq[i][1][0] = 0. ; FFqq2qq[i][1][1] = 2.0*arho - 2.0*an -arho + alft; // dd -> ud + Pi- (d u bar)
782  FFqq2qq[i][2][0] = 0. ; FFqq2qq[i][2][1] = 2.0*arho - 2.0*an -aphi + alft; // dd -> sd + K0 (d s bar)
783  FFqq2qq[i][3][0] = 0. ; FFqq2qq[i][3][1] = 2.0*arho - 2.0*an -aJPs + alft; // dd -> cd + D- (d c bar)
784  FFqq2qq[i][4][0] = 0. ; FFqq2qq[i][4][1] = 2.0*arho - 2.0*an -aUps + alft; // dd -> bd + B0 (d b bar)
785  }
786 }
787