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