Geant4_10
G4QGSMSplitableHadron.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 #include "G4QGSMSplitableHadron.hh"
27 #include "G4PhysicalConstants.hh"
28 #include "G4SystemOfUnits.hh"
29 #include "G4ParticleTable.hh"
30 #include "G4PionPlus.hh"
31 #include "G4PionMinus.hh"
32 #include "G4Gamma.hh"
33 #include "G4PionZero.hh"
34 #include "G4KaonPlus.hh"
35 #include "G4KaonMinus.hh"
36 
37 // based on prototype by Maxim Komogorov
38 // Splitting into methods, and centralizing of model parameters HPW Feb 1999
39 // restructuring HPW Feb 1999
40 // fixing bug in the sampling of 'x', HPW Feb 1999
41 // fixing bug in sampling pz, HPW Feb 1999.
42 // Code now also good for p-nucleus scattering (before only p-p), HPW Feb 1999.
43 // Using Parton more directly, HPW Feb 1999.
44 // Shortening the algorithm for sampling x, HPW Feb 1999.
45 // sampling of x replaced by formula, taking X_min into account in the correlated sampling. HPW, Feb 1999.
46 // logic much clearer now. HPW Feb 1999
47 // Removed the ordering problem. No Direction needed in selection of valence quark types. HPW Mar'99.
48 // Fixing p-t distributions for scattering of nuclei.
49 // Separating out parameters.
50 
51 void G4QGSMSplitableHadron::InitParameters()
52 {
53  // changing rapidity distribution for all
54  alpha = -0.5; // Note that this number is still assumed in the algorithm
55  // needs to be generalized.
56  // changing rapidity distribution for projectile like
57  beta = 2.5;// Note that this number is still assumed in the algorithm
58  // needs to be generalized.
59  theMinPz = 0.5*G4PionMinus::PionMinus()->GetPDGMass();
60  // theMinPz = 0.1*G4PionMinus::PionMinus()->GetPDGMass();
61  // theMinPz = G4PionMinus::PionMinus()->GetPDGMass();
62  // as low as possible, otherwise, we have unphysical boundary conditions in the sampling.
63  StrangeSuppress = 0.48;
64  sigmaPt = 0.*GeV; // widens eta slightly, if increased to 1.7,
65  // but Maxim's algorithm breaks energy conservation
66  // to be revised.
67  widthOfPtSquare = 0.01*GeV*GeV;
68  Direction = FALSE;
69  minTransverseMass = 1*keV;
70 }
71 
73 {
74  InitParameters();
75 }
76 
78 :G4VSplitableHadron(aPrimary)
79 {
80  InitParameters();
81  Direction = aDirection;
82 }
83 
84 
86 : G4VSplitableHadron(aPrimary)
87 {
88  InitParameters();
89 }
90 
92 : G4VSplitableHadron(aNucleon)
93 {
94  InitParameters();
95 }
96 
98 : G4VSplitableHadron(aNucleon)
99 {
100  InitParameters();
101  Direction = aDirection;
102 }
103 
105 
106 
107 
108 //**************************************************************************************************************************
109 
111 {
112  if (IsSplit()) return;
113  Splitting();
114  if (Color.size()!=0) return;
115  if (GetSoftCollisionCount() == 0)
116  {
117  DiffractiveSplitUp();
118  }
119  else
120  {
121  SoftSplitUp();
122  }
123 }
124 
125 void G4QGSMSplitableHadron::DiffractiveSplitUp()
126 {
127  // take the particle definitions and get the partons HPW
128  G4Parton * Left = NULL;
129  G4Parton * Right = NULL;
130  GetValenceQuarkFlavors(GetDefinition(), Left, Right);
131  Left->SetPosition(GetPosition());
132  Right->SetPosition(GetPosition());
133 
134  G4LorentzVector HadronMom = Get4Momentum();
135  //std::cout << "DSU 1 - "<<HadronMom<<std::endl;
136 
137  // momenta of string ends
138  G4double pt2 = HadronMom.perp2();
139  G4double transverseMass2 = HadronMom.plus()*HadronMom.minus();
140  G4double maxAvailMomentum2 = sqr(std::sqrt(transverseMass2) - std::sqrt(pt2));
141  G4ThreeVector pt(minTransverseMass, minTransverseMass, 0);
142  if(maxAvailMomentum2/widthOfPtSquare>0.01) pt = GaussianPt(widthOfPtSquare, maxAvailMomentum2);
143  //std::cout << "DSU 1.1 - "<< maxAvailMomentum2<< pt <<std::endl;
144 
145  G4LorentzVector LeftMom(pt, 0.);
146  G4LorentzVector RightMom;
147  RightMom.setPx(HadronMom.px() - pt.x());
148  RightMom.setPy(HadronMom.py() - pt.y());
149  //std::cout << "DSU 2 - "<<RightMom<<" "<< LeftMom <<std::endl;
150 
151  G4double Local1 = HadronMom.minus() + (RightMom.perp2() - LeftMom.perp2())/HadronMom.plus();
152  G4double Local2 = std::sqrt(std::max(0., sqr(Local1) - 4.*RightMom.perp2()*HadronMom.minus()/HadronMom.plus()));
153  //std::cout << "DSU 3 - "<< Local1 <<" "<< Local2 <<std::endl;
154  if (Direction) Local2 = -Local2;
155  G4double RightMinus = 0.5*(Local1 + Local2);
156  G4double LeftMinus = HadronMom.minus() - RightMinus;
157  //std::cout << "DSU 4 - "<< RightMinus <<" "<< LeftMinus << " "<<HadronMom.minus() <<std::endl;
158 
159  G4double LeftPlus = LeftMom.perp2()/LeftMinus;
160  G4double RightPlus = HadronMom.plus() - LeftPlus;
161  //std::cout << "DSU 5 - "<< RightPlus <<" "<< LeftPlus <<std::endl;
162  LeftMom.setPz(0.5*(LeftPlus - LeftMinus));
163  LeftMom.setE (0.5*(LeftPlus + LeftMinus));
164  RightMom.setPz(0.5*(RightPlus - RightMinus));
165  RightMom.setE (0.5*(RightPlus + RightMinus));
166  //std::cout << "DSU 6 - "<< LeftMom <<" "<< RightMom <<std::endl;
167  Left->Set4Momentum(LeftMom);
168  Right->Set4Momentum(RightMom);
169  Color.push_back(Left);
170  AntiColor.push_back(Right);
171 }
172 
173 
174 void G4QGSMSplitableHadron::SoftSplitUp()
175 {
176  //... sample transversal momenta for sea and valence quarks
177  G4double phi, pts;
178  G4double SumPy = 0.;
179  G4double SumPx = 0.;
181  G4int nSeaPair = GetSoftCollisionCount()-1;
182 
183  // here the condition,to ensure viability of splitting, also in cases
184  // where difractive excitation occured together with soft scattering.
185  // G4double LightConeMomentum = (Direction)? Get4Momentum().plus() : Get4Momentum().minus();
186  // G4double Xmin = theMinPz/LightConeMomentum;
187  G4double Xmin = theMinPz/( Get4Momentum().e() - GetDefinition()->GetPDGMass() );
188  while(Xmin>=1-(2*nSeaPair+1)*Xmin) Xmin*=0.95;
189 
190  G4int aSeaPair;
191  for (aSeaPair = 0; aSeaPair < nSeaPair; aSeaPair++)
192  {
193  // choose quark flavour, d:u:s = 1:1:(1/StrangeSuppress-2)
194 
195  G4int aPDGCode = 1 + (G4int)(G4UniformRand()/StrangeSuppress);
196 
197  // BuildSeaQuark() determines quark spin, isospin and colour
198  // via parton-constructor G4Parton(aPDGCode)
199 
200  G4Parton * aParton = BuildSeaQuark(false, aPDGCode, nSeaPair);
201 
202  // G4cerr << "G4QGSMSplitableHadron::SoftSplitUp()" << G4endl;
203 
204  // G4cerr << "Parton 1: "
205  // << " PDGcode: " << aPDGCode
206  // << " - Name: " << aParton->GetDefinition()->GetParticleName()
207  // << " - Type: " << aParton->GetDefinition()->GetParticleType()
208  // << " - Spin-3: " << aParton->GetSpinZ()
209  // << " - Colour: " << aParton->GetColour() << G4endl;
210 
211  // save colour a spin-3 for anti-quark
212 
213  G4int firstPartonColour = aParton->GetColour();
214  G4double firstPartonSpinZ = aParton->GetSpinZ();
215 
216  SumPx += aParton->Get4Momentum().px();
217  SumPy += aParton->Get4Momentum().py();
218  Color.push_back(aParton);
219 
220  // create anti-quark
221 
222  aParton = BuildSeaQuark(true, aPDGCode, nSeaPair);
223  aParton->SetSpinZ(-firstPartonSpinZ);
224  aParton->SetColour(-firstPartonColour);
225 
226  // G4cerr << "Parton 2: "
227  // << " PDGcode: " << -aPDGCode
228  // << " - Name: " << aParton->GetDefinition()->GetParticleName()
229  // << " - Type: " << aParton->GetDefinition()->GetParticleType()
230  // << " - Spin-3: " << aParton->GetSpinZ()
231  // << " - Colour: " << aParton->GetColour() << G4endl;
232  // G4cerr << "------------" << G4endl;
233 
234  SumPx += aParton->Get4Momentum().px();
235  SumPy += aParton->Get4Momentum().py();
236  AntiColor.push_back(aParton);
237  }
238  // Valence quark
239  G4Parton* pColorParton = NULL;
240  G4Parton* pAntiColorParton = NULL;
241  GetValenceQuarkFlavors(GetDefinition(), pColorParton, pAntiColorParton);
242  G4int ColorEncoding = pColorParton->GetPDGcode();
243 
244  pts = sigmaPt*std::sqrt(-std::log(G4UniformRand()));
245  phi = 2.*pi*G4UniformRand();
246  G4double Px = pts*std::cos(phi);
247  G4double Py = pts*std::sin(phi);
248  SumPx += Px;
249  SumPy += Py;
250 
251  if (ColorEncoding < 0) // use particle definition
252  {
253  G4LorentzVector ColorMom(-SumPx, -SumPy, 0, 0);
254  pColorParton->Set4Momentum(ColorMom);
255  G4LorentzVector AntiColorMom(Px, Py, 0, 0);
256  pAntiColorParton->Set4Momentum(AntiColorMom);
257  }
258  else
259  {
260  G4LorentzVector ColorMom(Px, Py, 0, 0);
261  pColorParton->Set4Momentum(ColorMom);
262  G4LorentzVector AntiColorMom(-SumPx, -SumPy, 0, 0);
263  pAntiColorParton->Set4Momentum(AntiColorMom);
264  }
265  Color.push_back(pColorParton);
266  AntiColor.push_back(pAntiColorParton);
267 
268  // Sample X
269  G4int nAttempt = 0;
270  G4double SumX = 0;
271  G4double aBeta = beta;
272  G4double ColorX, AntiColorX;
273  if (GetDefinition() == G4PionMinus::PionMinusDefinition()) aBeta = 1.;
274  if (GetDefinition() == G4Gamma::GammaDefinition()) aBeta = 1.;
275  if (GetDefinition() == G4PionPlus::PionPlusDefinition()) aBeta = 1.;
276  if (GetDefinition() == G4PionZero::PionZeroDefinition()) aBeta = 1.;
277  if (GetDefinition() == G4KaonPlus::KaonPlusDefinition()) aBeta = 0.;
278  if (GetDefinition() == G4KaonMinus::KaonMinusDefinition()) aBeta = 0.;
279  do
280  {
281  SumX = 0;
282  nAttempt++;
283  G4int NumberOfUnsampledSeaQuarks = 2*nSeaPair;
284  ColorX = SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta);
285  Color.back()->SetX(SumX = ColorX);// this is the valenz quark.
286  for(G4int aPair = 0; aPair < nSeaPair; aPair++)
287  {
288  NumberOfUnsampledSeaQuarks--;
289  ColorX = SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta);
290  Color[aPair]->SetX(ColorX);
291  SumX += ColorX;
292  NumberOfUnsampledSeaQuarks--;
293  AntiColorX = SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta);
294  AntiColor[aPair]->SetX(AntiColorX); // the 'sea' partons
295  SumX += AntiColorX;
296  if (1. - SumX <= Xmin) break;
297  }
298  }
299  while (1. - SumX <= Xmin);
300 
301  (*(AntiColor.end()-1))->SetX(1. - SumX); // the di-quark takes the rest, then go to momentum
302  G4double lightCone = ((!Direction) ? Get4Momentum().minus() : Get4Momentum().plus());
303  G4double lightCone2 = ((!Direction) ? Get4Momentum().plus() : Get4Momentum().minus());
304  for(aSeaPair = 0; aSeaPair < nSeaPair+1; aSeaPair++)
305  {
306  G4Parton* aParton = Color[aSeaPair];
307  aParton->DefineMomentumInZ(lightCone, lightCone2, Direction);
308 
309  aParton = AntiColor[aSeaPair];
310  aParton->DefineMomentumInZ(lightCone, lightCone2, Direction);
311  }
312  return;
313 }
314 
315 void G4QGSMSplitableHadron::GetValenceQuarkFlavors(const G4ParticleDefinition * aPart, G4Parton *& Parton1, G4Parton *& Parton2)
316 {
317  // Note! convention aEnd = q or (qq)bar and bEnd = qbar or qq.
318  G4int aEnd;
319  G4int bEnd;
320  G4int HadronEncoding = aPart->GetPDGEncoding();
321  if (aPart->GetBaryonNumber() == 0)
322  {
323  theMesonSplitter.SplitMeson(HadronEncoding, &aEnd, &bEnd);
324  }
325  else
326  {
327  theBaryonSplitter.SplitBarion(HadronEncoding, &aEnd, &bEnd);
328  }
329 
330  Parton1 = new G4Parton(aEnd);
331  Parton1->SetPosition(GetPosition());
332 
333  // G4cerr << "G4QGSMSplitableHadron::GetValenceQuarkFlavors()" << G4endl;
334  // G4cerr << "Parton 1: "
335  // << " PDGcode: " << aEnd
336  // << " - Name: " << Parton1->GetDefinition()->GetParticleName()
337  // << " - Type: " << Parton1->GetDefinition()->GetParticleType()
338  // << " - Spin-3: " << Parton1->GetSpinZ()
339  // << " - Colour: " << Parton1->GetColour() << G4endl;
340 
341  Parton2 = new G4Parton(bEnd);
342  Parton2->SetPosition(GetPosition());
343 
344  // G4cerr << "Parton 2: "
345  // << " PDGcode: " << bEnd
346  // << " - Name: " << Parton2->GetDefinition()->GetParticleName()
347  // << " - Type: " << Parton2->GetDefinition()->GetParticleType()
348  // << " - Spin-3: " << Parton2->GetSpinZ()
349  // << " - Colour: " << Parton2->GetColour() << G4endl;
350  // G4cerr << "... now checking for color and spin conservation - yielding: " << G4endl;
351 
352  // colour of parton 1 choosen at random by G4Parton(aEnd)
353  // colour of parton 2 is the opposite:
354 
355  Parton2->SetColour(-(Parton1->GetColour()));
356 
357  // isospin-3 of both partons is handled by G4Parton(PDGCode)
358 
359  // spin-3 of parton 1 and 2 choosen at random by G4Parton(aEnd)
360  // spin-3 of parton 2 may be constrained by spin of original particle:
361 
362  if ( std::abs(Parton1->GetSpinZ() + Parton2->GetSpinZ()) > aPart->GetPDGSpin())
363  {
364  Parton2->SetSpinZ(-(Parton2->GetSpinZ()));
365  }
366 
367  // G4cerr << "Parton 2: "
368  // << " PDGcode: " << bEnd
369  // << " - Name: " << Parton2->GetDefinition()->GetParticleName()
370  // << " - Type: " << Parton2->GetDefinition()->GetParticleType()
371  // << " - Spin-3: " << Parton2->GetSpinZ()
372  // << " - Colour: " << Parton2->GetColour() << G4endl;
373  // G4cerr << "------------" << G4endl;
374 
375 }
376 
377 
378 G4ThreeVector G4QGSMSplitableHadron::GaussianPt(G4double widthSquare, G4double maxPtSquare)
379 {
380  G4double R;
381  while((R = -widthSquare*std::log(G4UniformRand())) > maxPtSquare) {;}
382  R = std::sqrt(R);
383  G4double phi = twopi*G4UniformRand();
384  return G4ThreeVector (R*std::cos(phi), R*std::sin(phi), 0.);
385 }
386 
387 G4Parton * G4QGSMSplitableHadron::
388 BuildSeaQuark(G4bool isAntiQuark, G4int aPDGCode, G4int /* nSeaPair*/)
389 {
390  if (isAntiQuark) aPDGCode*=-1;
391  G4Parton* result = new G4Parton(aPDGCode);
392  result->SetPosition(GetPosition());
393  G4ThreeVector aPtVector = GaussianPt(sigmaPt, DBL_MAX);
394  G4LorentzVector a4Momentum(aPtVector, 0);
395  result->Set4Momentum(a4Momentum);
396  return result;
397 }
398 
399 G4double G4QGSMSplitableHadron::
400 SampleX(G4double anXmin, G4int nSea, G4int totalSea, G4double aBeta)
401 {
403  G4double x1, x2;
404  G4double ymax = 0;
405  for(G4int ii=1; ii<100; ii++)
406  {
407  G4double y = std::pow(1./G4double(ii), alpha);
408  y *= std::pow( std::pow(1-anXmin-totalSea*anXmin, alpha+1) - std::pow(anXmin, alpha+1), nSea);
409  y *= std::pow(1-anXmin-totalSea*anXmin, aBeta+1) - std::pow(anXmin, aBeta+1);
410  if(y>ymax) ymax = y;
411  }
412  G4double y;
413  G4double xMax=1-(totalSea+1)*anXmin;
414  if(anXmin > xMax)
415  {
416  G4cout << "anXmin = "<<anXmin<<" nSea = "<<nSea<<" totalSea = "<< totalSea<<G4endl;
417  throw G4HadronicException(__FILE__, __LINE__, "G4QGSMSplitableHadron - Fatal: Cannot sample parton densities under these constraints.");
418  }
419  do
420  {
421  x1 = G4RandFlat::shoot(anXmin, xMax);
422  y = std::pow(x1, alpha);
423  y *= std::pow( std::pow(1-x1-totalSea*anXmin, alpha+1) - std::pow(anXmin, alpha+1), nSea);
424  y *= std::pow(1-x1-totalSea*anXmin, aBeta+1) - std::pow(anXmin, aBeta+1);
425  x2 = ymax*G4UniformRand();
426  }
427  while(x2>y);
428  result = x1;
429  return result;
430 }
ThreeVector shoot(const G4int Ap, const G4int Af)
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
G4int GetPDGcode() const
Definition: G4Parton.hh:124
Double_t x2[nxs]
Definition: Style.C:19
CLHEP::Hep3Vector G4ThreeVector
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
ush Pos
Definition: deflate.h:89
const G4LorentzVector & Get4Momentum() const
Definition: G4Parton.hh:140
G4double G4NeutronHPJENDLHEData::G4double result
void Set4Momentum(const G4LorentzVector &aMomentum)
Definition: G4Parton.hh:145
int G4int
Definition: G4Types.hh:78
static G4PionZero * PionZeroDefinition()
Definition: G4PionZero.cc:99
G4ParticleDefinition * GetDefinition() const
Double_t y
Definition: plot.C:279
void DefineMomentumInZ(G4double aLightConeMomentum, G4bool aDirection)
Definition: G4Parton.cc:142
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
double py() const
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
bool G4bool
Definition: G4Types.hh:79
#define FALSE
Definition: globals.hh:52
Double_t x1[nxs]
Definition: Style.C:18
double px() const
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
void SetColour(G4int aColour)
Definition: G4Parton.hh:88
TMarker * pt
Definition: egs.C:22
double minus() const
const G4LorentzVector & Get4Momentum() const
void SetSpinZ(G4double aSpinZ)
Definition: G4Parton.hh:94
G4double GetSpinZ()
Definition: G4Parton.hh:95
G4double GetPDGMass() const
G4int GetColour()
Definition: G4Parton.hh:89
Color
Definition: test07.cc:36
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
G4bool SplitMeson(G4int PDGcode, G4int *aEnd, G4int *bEnd)
double perp2() const
const G4ThreeVector & GetPosition() const
G4double GetPDGSpin() const
#define G4endl
Definition: G4ios.hh:61
T sqr(const T &x)
Definition: templates.hh:145
double plus() const
double G4double
Definition: G4Types.hh:76
int Xmin
#define DBL_MAX
Definition: templates.hh:83
void SetPosition(const G4ThreeVector &aPosition)
Definition: G4Parton.hh:134
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81
G4bool SplitBarion(G4int PDGCode, G4int *q_or_qqbar, G4int *qbar_or_qq)