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