Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VLongitudinalStringDecay.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: G4VLongitudinalStringDecay.cc 73023 2013-08-15 09:10:04Z gcosmo $
28 //
29 // -----------------------------------------------------------------------------
30 // GEANT 4 class implementation file
31 //
32 // History: first implementation, Maxim Komogorov, 1-Jul-1998
33 // redesign Gunter Folger, August/September 2001
34 // -----------------------------------------------------------------------------
36 #include "G4PhysicalConstants.hh"
37 #include "G4SystemOfUnits.hh"
38 #include "G4ios.hh"
39 #include "Randomize.hh"
40 #include "G4FragmentingString.hh"
41 
42 #include "G4ParticleDefinition.hh"
43 #include "G4ParticleTypes.hh"
44 #include "G4ParticleChange.hh"
45 #include "G4VShortLivedParticle.hh"
47 #include "G4ParticleTable.hh"
49 #include "G4VDecayChannel.hh"
50 #include "G4DecayTable.hh"
51 
52 #include "G4DiQuarks.hh"
53 #include "G4Quarks.hh"
54 #include "G4Gluons.hh"
55 
56 //------------------------debug switches
57 //#define DEBUG_LightFragmentationTest 1
58 
59 
60 //********************************************************************************
61 // Constructors
62 
64 {
65  MassCut = 0.35*GeV;
66  ClusterMass = 0.15*GeV;
67 
68  SmoothParam = 0.9;
69  StringLoopInterrupt = 1000;
71 
72 // Changable Parameters below.
73  SigmaQT = 0.5 * GeV; // 0.5 0.1
74 
75  StrangeSuppress = 0.44; // 27 % strange quarks produced, ie. u:d:s=1:1:0.27
76  DiquarkSuppress = 0.07;
77  DiquarkBreakProb = 0.1;
78 
79  //... pspin_meson is probability to create vector meson
80  pspin_meson = 0.5;
81 
82  //... pspin_barion is probability to create 3/2 barion
83  pspin_barion = 0.5;
84 
85  //... vectorMesonMix[] is quark mixing parameters for vector mesons (Variable spin = 3)
86  vectorMesonMix.resize(6);
87  vectorMesonMix[0] = 0.5;
88  vectorMesonMix[1] = 0.0;
89  vectorMesonMix[2] = 0.5;
90  vectorMesonMix[3] = 0.0;
91  vectorMesonMix[4] = 1.0;
92  vectorMesonMix[5] = 1.0;
93 
94  //... scalarMesonMix[] is quark mixing parameters for scalar mesons (Variable spin=1)
95  scalarMesonMix.resize(6);
96  scalarMesonMix[0] = 0.5;
97  scalarMesonMix[1] = 0.25;
98  scalarMesonMix[2] = 0.5;
99  scalarMesonMix[3] = 0.25;
100  scalarMesonMix[4] = 1.0;
101  scalarMesonMix[5] = 0.5;
102 
103 // Parameters may be changed until the first fragmentation starts
104  PastInitPhase=false;
107  Kappa = 1.0 * GeV/fermi;
108 
109 
110 }
111 
112 
114  {
115  delete hadronizer;
116  }
117 
118 //=============================================================================
119 
120 // Operators
121 
122 //const & G4VLongitudinalStringDecay::operator=(const G4VLongitudinalStringDecay &)
123 // {
124 // }
125 
126 //-----------------------------------------------------------------------------
127 
128 int G4VLongitudinalStringDecay::operator==(const G4VLongitudinalStringDecay &) const
129  {
130  throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::operator== forbidden");
131  return false;
132  }
133 
134 //-------------------------------------------------------------------------------------
135 
136 int G4VLongitudinalStringDecay::operator!=(const G4VLongitudinalStringDecay &) const
137  {
138  throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::operator!= forbidden");
139  return true;
140  }
141 
142 //***********************************************************************************
143 
144 // For changing Mass Cut used for selection of very small mass strings
146 
147 //-----------------------------------------------------------------------------
148 
149 // For handling a string with very low mass
150 
152  G4ExcitedString * const string)
153 {
154  // Check string decay threshold
155 
156  G4KineticTrackVector * result=0; // return 0 when string exceeds the mass cut
157 
159 
160  G4FragmentingString aString(*string);
161 
162  if ( sqr(FragmentationMass(&aString,0,&hadrons)+MassCut) < aString.Mass2()) {
163  return 0;
164  }
165 
166 // The string mass is very low ---------------------------
167 
168  result=new G4KineticTrackVector;
169 
170  if ( hadrons.second ==0 )
171  {
172 // Substitute string by light hadron, Note that Energy is not conserved here!
173 
174 /*
175 #ifdef DEBUG_LightFragmentationTest
176  G4cout << "VlongSF Warning replacing string by single hadron " <<G4endl;
177  G4cout << hadrons.first->GetParticleName()
178  << "string .. " << string->Get4Momentum() << " "
179  << string->Get4Momentum().m() << G4endl;
180 #endif
181 */
182  G4ThreeVector Mom3 = string->Get4Momentum().vect();
183  G4LorentzVector Mom(Mom3,
184  std::sqrt(Mom3.mag2() +
185  sqr(hadrons.first->GetPDGMass())));
186  result->push_back(new G4KineticTrack(hadrons.first, 0,
187  string->GetPosition(),
188  Mom));
189  } else
190  {
191 //... string was qq--qqbar type: Build two stable hadrons,
192 
193 #ifdef DEBUG_LightFragmentationTest
194  G4cout << "VlongSF Warning replacing qq-qqbar string by TWO hadrons "
195  << hadrons.first->GetParticleName() << " / "
196  << hadrons.second->GetParticleName()
197  << "string .. " << string->Get4Momentum() << " "
198  << string->Get4Momentum().m() << G4endl;
199 #endif
200 
201  G4LorentzVector Mom1, Mom2;
202  Sample4Momentum(&Mom1, hadrons.first->GetPDGMass(),
203  &Mom2,hadrons.second->GetPDGMass(),
204  string->Get4Momentum().mag());
205 
206  result->push_back(new G4KineticTrack(hadrons.first, 0,
207  string->GetPosition(),
208  Mom1));
209  result->push_back(new G4KineticTrack(hadrons.second, 0,
210  string->GetPosition(),
211  Mom2));
212 
213  G4ThreeVector Velocity = string->Get4Momentum().boostVector();
214  result->Boost(Velocity);
215  }
216 
217  return result;
218 
219 }
220 
221 //----------------------------------------------------------------------------------------
222 
224  const G4FragmentingString * const string,
225  Pcreate build, pDefPair * pdefs )
226 {
227 
228  G4double mass;
229  static G4ThreadLocal G4bool NeedInit(true);
230  static G4ThreadLocal std::vector<double> *nomix_G4MT_TLS_ = 0 ; if (!nomix_G4MT_TLS_) nomix_G4MT_TLS_ = new std::vector<double> ; std::vector<double> &nomix = *nomix_G4MT_TLS_;
231  static G4ThreadLocal G4HadronBuilder * minMassHadronizer;
232  if ( NeedInit )
233  {
234  NeedInit = false;
235  nomix.resize(6);
236  for ( G4int i=0; i<6 ; i++ ) nomix[i]=0;
237 
238 // minMassHadronizer=new G4HadronBuilder(pspin_meson,pspin_barion,nomix,nomix);
239  minMassHadronizer=hadronizer;
240  }
241 
242  if ( build==0 ) build=&G4HadronBuilder::BuildLowSpin;
243 
244  G4ParticleDefinition *Hadron1, *Hadron2=0;
245 
246  if (!string->FourQuarkString() )
247  {
248  // spin 0 meson or spin 1/2 barion will be built
249 
250 //G4cout<<"String Left Right "<<string->GetLeftParton()<<" "<<string->GetRightParton()<<G4endl;
251  Hadron1 = (minMassHadronizer->*build)(string->GetLeftParton(),
252  string->GetRightParton());
253 //G4cout<<"Hadron1 "<<Hadron1->GetParticleName()<<G4endl;
254  mass= (Hadron1)->GetPDGMass();
255  } else
256  {
257  //... string is qq--qqbar: Build two stable hadrons,
258  //... with extra uubar or ddbar quark pair
259  G4int iflc = (G4UniformRand() < 0.5)? 1 : 2;
260  if (string->GetLeftParton()->GetPDGEncoding() < 0) iflc = -iflc;
261 
262  //... theSpin = 4; spin 3/2 baryons will be built
263  Hadron1 = (minMassHadronizer->*build)(string->GetLeftParton(),
264  FindParticle(iflc) );
265  Hadron2 = (minMassHadronizer->*build)(string->GetRightParton(),
266  FindParticle(-iflc) );
267  mass = (Hadron1)->GetPDGMass() + (Hadron2)->GetPDGMass();
268  }
269 
270  if ( pdefs != 0 )
271  { // need to return hadrons as well....
272  pdefs->first = Hadron1;
273  pdefs->second = Hadron2;
274  }
275 
276  return mass;
277 }
278 
279 //----------------------------------------------------------------------------
280 
282  {
284  if (ptr == NULL)
285  {
286  G4cout << "Particle with encoding "<<Encoding<<" does not exist!!!"<<G4endl;
287  throw G4HadronicException(__FILE__, __LINE__, "Check your particle table");
288  }
289  return ptr;
290  }
291 
292 //-----------------------------------------------------------------------------
293 // virtual void Sample4Momentum(G4LorentzVector* Mom, G4double Mass,
294 // G4LorentzVector* AntiMom, G4double AntiMass,
295 // G4double InitialMass)=0;
296 //-----------------------------------------------------------------------------
297 
298 //*********************************************************************************
299 // For decision on continue or stop string fragmentation
300 // virtual G4bool StopFragmenting(const G4FragmentingString * const string)=0;
301 // virtual G4bool IsFragmentable(const G4FragmentingString * const string)=0;
302 
303 // If a string can not fragment, make last break into 2 hadrons
304 // virtual G4bool SplitLast(G4FragmentingString * string,
305 // G4KineticTrackVector * LeftVector,
306 // G4KineticTrackVector * RightVector)=0;
307 //-----------------------------------------------------------------------------
308 //
309 // If a string fragments, do the following
310 //
311 // For transver of a string to its CMS frame
312 //-----------------------------------------------------------------------------
313 
315 {
316  G4Parton *Left=new G4Parton(*in.GetLeftParton());
317  G4Parton *Right=new G4Parton(*in.GetRightParton());
318  return new G4ExcitedString(Left,Right,in.GetDirection());
319 }
320 
321 //-----------------------------------------------------------------------------
322 
324  G4FragmentingString *string,
325  G4FragmentingString *&newString)
326 {
327 //G4cout<<"Start SplitUP"<<G4endl;
328  //... random choice of string end to use for creating the hadron (decay)
329  G4int SideOfDecay = (G4UniformRand() < 0.5)? 1: -1;
330  if (SideOfDecay < 0)
331  {
332  string->SetLeftPartonStable();
333  } else
334  {
335  string->SetRightPartonStable();
336  }
337 
338  G4ParticleDefinition *newStringEnd;
339  G4ParticleDefinition * HadronDefinition;
340  if (string->DecayIsQuark())
341  {
342  HadronDefinition= QuarkSplitup(string->GetDecayParton(), newStringEnd);
343  } else {
344  HadronDefinition= DiQuarkSplitup(string->GetDecayParton(), newStringEnd);
345  }
346 
347 //G4cout<<"New had "<<HadronDefinition->GetParticleName()<<G4endl;
348 // create new String from old, ie. keep Left and Right order, but replace decay
349 
350  newString=new G4FragmentingString(*string,newStringEnd); // To store possible
351  // quark containt of new string
352 //G4cout<<"SplitEandP "<<G4endl;
353  G4LorentzVector* HadronMomentum=SplitEandP(HadronDefinition, string, newString);
354 
355  delete newString; newString=0;
356 
357  G4KineticTrack * Hadron =0;
358  if ( HadronMomentum != 0 ) {
359 
361  Hadron = new G4KineticTrack(HadronDefinition, 0,Pos, *HadronMomentum);
362 
363  newString=new G4FragmentingString(*string,newStringEnd,
364  HadronMomentum);
365 
366  delete HadronMomentum;
367  }
368 //G4cout<<"End SplitUP"<<G4endl;
369  return Hadron;
370 }
371 
372 //--------------------------------------------------------------------------------------
373 
376  decay, G4ParticleDefinition *&created)
377 {
378  G4int IsParticle=(decay->GetPDGEncoding()>0) ? -1 : +1; // if we have a quark,
379  // we need antiquark
380  // (or diquark)
381  pDefPair QuarkPair = CreatePartonPair(IsParticle);
382  created = QuarkPair.second;
383  return hadronizer->Build(QuarkPair.first, decay);
384 
385 }
386 
387 //-----------------------------------------------------------------------------
388 
391  G4ParticleDefinition *&created)
392 {
393  //... can Diquark break or not?
395  //... Diquark break
396 
397  G4int stableQuarkEncoding = decay->GetPDGEncoding()/1000;
398  G4int decayQuarkEncoding = (decay->GetPDGEncoding()/100)%10;
399  if (G4UniformRand() < 0.5)
400  {
401  G4int Swap = stableQuarkEncoding;
402  stableQuarkEncoding = decayQuarkEncoding;
403  decayQuarkEncoding = Swap;
404  }
405 
406  G4int IsParticle=(decayQuarkEncoding>0) ? -1 : +1;
407  // if we have a quark, we need antiquark)
408  pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
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  return had;
419 // return hadronizer->Build(QuarkPair.first, decayQuark);
420 
421  } else {
422  //... Diquark does not break
423 
424  G4int IsParticle=(decay->GetPDGEncoding()>0) ? +1 : -1;
425  // if we have a diquark, we need quark)
426  pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
427  created = QuarkPair.second;
428 
429  G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay);
430  return had;
431 // return G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay);
432  }
433 }
434 
435 //-----------------------------------------------------------------------------
436 
438  {
439  return (1 + (int)(G4UniformRand()/StrangeSuppress));
440  }
441 
442 //-----------------------------------------------------------------------------
443 
445 {
446 // NeedParticle = +1 for Particle, -1 for Antiparticle
447 
448  if ( AllowDiquarks && G4UniformRand() < DiquarkSuppress )
449  {
450  // Create a Diquark - AntiDiquark pair , first in pair is anti to IsParticle
451  G4int q1 = SampleQuarkFlavor();
452  G4int q2 = SampleQuarkFlavor();
453  G4int spin = (q1 != q2 && G4UniformRand() <= 0.5)? 1 : 3;
454  // convention: quark with higher PDG number is first
455  G4int PDGcode = (std::max(q1,q2) * 1000 + std::min(q1,q2) * 100 + spin) * NeedParticle;
456  return pDefPair (FindParticle(-PDGcode),FindParticle(PDGcode));
457 
458 
459  } else {
460  // Create a Quark - AntiQuark pair, first in pair IsParticle
461  G4int PDGcode=SampleQuarkFlavor()*NeedParticle;
462  return pDefPair (FindParticle(PDGcode),FindParticle(-PDGcode));
463  }
464 
465 }
466 
467 //-----------------------------------------------------------------------------
469  {
470  G4double Pt;
471  if ( ptMax < 0 ) {
472  // sample full gaussian
473  Pt = -std::log(G4UniformRand());
474  } else {
475  // sample in limited range
476  Pt = -std::log(G4RandFlat::shoot(std::exp(-sqr(ptMax)/sqr(SigmaQT)), 1.));
477  }
478  Pt = SigmaQT * std::sqrt(Pt);
479  G4double phi = 2.*pi*G4UniformRand();
480  return G4ThreeVector(Pt * std::cos(phi),Pt * std::sin(phi),0);
481  }
482 
483 //******************************************************************************
484 
486  {
487 
488 // `yo-yo` formation time
489 // const G4double kappa = 1.0 * GeV/fermi/4.;
491  for(size_t c1 = 0; c1 < Hadrons->size(); c1++)
492  {
493  G4double SumPz = 0;
494  G4double SumE = 0;
495  for(size_t c2 = 0; c2 < c1; c2++)
496  {
497  SumPz += Hadrons->operator[](c2)->Get4Momentum().pz();
498  SumE += Hadrons->operator[](c2)->Get4Momentum().e();
499  }
500  G4double HadronE = Hadrons->operator[](c1)->Get4Momentum().e();
501  G4double HadronPz = Hadrons->operator[](c1)->Get4Momentum().pz();
502  Hadrons->operator[](c1)->SetFormationTime(
503 (theInitialStringMass - 2.*SumPz + HadronE - HadronPz)/(2.*kappa)/c_light);
504 
505  G4ThreeVector aPosition(0, 0,
506 (theInitialStringMass - 2.*SumE - HadronE + HadronPz)/(2.*kappa));
507  Hadrons->operator[](c1)->SetPosition(aPosition);
508 
509  }
510  }
511 
512 //-----------------------------------------------------------------------------
513 
515 {
516  if ( PastInitPhase ) {
517  throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetSigmaTransverseMomentum after FragmentString() not allowed");
518  } else {
519  SigmaQT = aValue;
520  }
521 }
522 
523 //----------------------------------------------------------------------------------------------------------
524 
526 {
527  if ( PastInitPhase ) {
528  throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetStrangenessSuppression after FragmentString() not allowed");
529  } else {
530  StrangeSuppress = aValue;
531  }
532 }
533 
534 //----------------------------------------------------------------------------------------------------------
535 
537 {
538  if ( PastInitPhase ) {
539  throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetDiquarkSuppression after FragmentString() not allowed");
540  } else {
541  DiquarkSuppress = aValue;
542  }
543 }
544 
545 //----------------------------------------------------------------------------------------
546 
548 {
549  if ( PastInitPhase ) {
550  throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetDiquarkBreakProbability after FragmentString() not allowed");
551  } else {
552  DiquarkBreakProb = aValue;
553  }
554 }
555 
556 //----------------------------------------------------------------------------------------------------------
557 
559 {
560  if ( PastInitPhase ) {
561  throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonProbability after FragmentString() not allowed");
562  } else {
563  pspin_meson = aValue;
564  delete hadronizer;
567  }
568 }
569 
570 //----------------------------------------------------------------------------------------------------------
571 
573 {
574  if ( PastInitPhase ) {
575  throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetSpinThreeHalfBarionProbability after FragmentString() not allowed");
576  } else {
577  pspin_barion = aValue;
578  delete hadronizer;
581  }
582 }
583 
584 //----------------------------------------------------------------------------------------------------------
585 
586 void G4VLongitudinalStringDecay::SetScalarMesonMixings(std::vector<G4double> aVector)
587 {
588  if ( PastInitPhase ) {
589  throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetScalarMesonMixings after FragmentString() not allowed");
590  } else {
591  if ( aVector.size() < 6 )
592  throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetScalarMesonMixings( argument Vector too small");
593  scalarMesonMix[0] = aVector[0];
594  scalarMesonMix[1] = aVector[1];
595  scalarMesonMix[2] = aVector[2];
596  scalarMesonMix[3] = aVector[3];
597  scalarMesonMix[4] = aVector[4];
598  scalarMesonMix[5] = aVector[5];
599  delete hadronizer;
602  }
603 }
604 
605 //----------------------------------------------------------------------------------------------------------
606 
607 void G4VLongitudinalStringDecay::SetVectorMesonMixings(std::vector<G4double> aVector)
608 {
609  if ( PastInitPhase ) {
610  throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonMixings after FragmentString() not allowed");
611  } else {
612  if ( aVector.size() < 6 )
613  throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonMixings( argument Vector too small");
614  vectorMesonMix[0] = aVector[0];
615  vectorMesonMix[1] = aVector[1];
616  vectorMesonMix[2] = aVector[2];
617  vectorMesonMix[3] = aVector[3];
618  vectorMesonMix[4] = aVector[4];
619  vectorMesonMix[5] = aVector[5];
620  delete hadronizer;
623 
624  }
625 }
626 
627 //-------------------------------------------------------------------------------------------
629 {
630  Kappa = aValue * GeV/fermi;
631 }
632 //**************************************************************************************
633 
G4ParticleDefinition * GetRightParton(void) const
ThreeVector shoot(const G4int Ap, const G4int Af)
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
G4ParticleDefinition * QuarkSplitup(G4ParticleDefinition *decay, G4ParticleDefinition *&created)
CLHEP::Hep3Vector G4ThreeVector
G4ParticleDefinition * BuildLowSpin(G4ParticleDefinition *black, G4ParticleDefinition *white)
G4ExcitedString * CPExcited(const G4ExcitedString &string)
ush Pos
Definition: deflate.h:89
void SetStrangenessSuppression(G4double aValue)
G4ParticleDefinition * Build(G4ParticleDefinition *black, G4ParticleDefinition *white)
G4Parton * GetLeftParton(void) const
G4ThreeVector SampleQuarkPt(G4double ptMax=-1.)
std::vector< G4double > vectorMesonMix
#define G4ThreadLocal
Definition: tls.hh:52
G4double FragmentationMass(const G4FragmentingString *const string, Pcreate build=0, pDefPair *pdefs=0)
virtual void SetMassCut(G4double aValue)
G4ParticleDefinition * DiQuarkSplitup(G4ParticleDefinition *decay, G4ParticleDefinition *&created)
int G4int
Definition: G4Types.hh:78
G4ParticleDefinition * GetDecayParton() const
pDefPair CreatePartonPair(G4int NeedParticle, G4bool AllowDiquarks=true)
void SetStringTensionParameter(G4double aValue)
const G4ThreeVector & GetPosition() const
G4ParticleDefinition * GetLeftParton(void) const
G4double Mass2() const
void SetDiquarkBreakProbability(G4double aValue)
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
ParticleList decay(Cluster *const c)
Carries out a cluster decay.
G4KineticTrackVector * LightFragmentationTest(const G4ExcitedString *const theString)
std::pair< G4ParticleDefinition *, G4ParticleDefinition * > pDefPair
double mag() const
bool G4bool
Definition: G4Types.hh:79
G4ParticleDefinition * FindParticle(G4int Encoding)
virtual void Sample4Momentum(G4LorentzVector *Mom, G4double Mass, G4LorentzVector *AntiMom, G4double AntiMass, G4double InitialMass)=0
void Boost(G4ThreeVector &Velocity)
G4bool FourQuarkString(void) const
std::vector< G4double > scalarMesonMix
void SetVectorMesonProbability(G4double aValue)
G4Parton * GetRightParton(void) const
static G4ParticleTable * GetParticleTable()
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void SetVectorMesonMixings(std::vector< G4double > aVector)
void SetSpinThreeHalfBarionProbability(G4double aValue)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double mag2() const
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
void SetScalarMesonMixings(std::vector< G4double > aVector)
void CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector *)
virtual G4LorentzVector * SplitEandP(G4ParticleDefinition *pHadron, G4FragmentingString *string, G4FragmentingString *newString)=0
float c_light
Definition: hepunit.py:257
tuple c1
Definition: plottest35.py:14
G4KineticTrack * Splitup(G4FragmentingString *string, G4FragmentingString *&newString)
void SetDiquarkSuppression(G4double aValue)