ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4QuarkExchange.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4QuarkExchange.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 // GEANT 4 class implemetation file
29 //
30 // ---------------- G4QuarkExchange --------------
31 // by V. Uzhinsky, October 2016.
32 // QuarkExchange is used by strings models.
33 // Take a projectile and a target.
34 //Simulate Q exchange with excitation of projectile or target.
35 // ------------------------------------------------------------
36 
37 #include "G4QuarkExchange.hh"
38 #include "globals.hh"
39 #include "G4PhysicalConstants.hh"
40 #include "G4SystemOfUnits.hh"
41 #include "Randomize.hh"
42 #include "G4LorentzRotation.hh"
43 #include "G4ThreeVector.hh"
44 #include "G4ParticleDefinition.hh"
45 #include "G4VSplitableHadron.hh"
46 #include "G4ExcitedString.hh"
47 
48 #include "G4Log.hh"
49 #include "G4Pow.hh"
50 
51 //#define debugQuarkExchange
52 
54 
56 
59 {
60  #ifdef debugQuarkExchange
61  G4cout<<G4endl<<"G4QuarkExchange::ExciteParticipants"<<G4endl;
62  #endif
63 
64  G4LorentzVector Pprojectile = projectile->Get4Momentum();
65  G4double Mprojectile = projectile->GetDefinition()->GetPDGMass();
66  G4double Mprojectile2 = sqr(Mprojectile);
67 
68  G4LorentzVector Ptarget = target->Get4Momentum();
69  G4double Mtarget = target->GetDefinition()->GetPDGMass();
70  G4double Mtarget2 = sqr(Mtarget);
71 
72  #ifdef debugQuarkExchange
73  G4cout<<"Proj Targ "<<projectile->GetDefinition()->GetPDGEncoding()<<" "<<target->GetDefinition()->GetPDGEncoding()<<G4endl;
74  G4cout<<"Proj. 4-Mom "<<Pprojectile<<" "<<Pprojectile.mag()<<G4endl
75  <<"Targ. 4-Mom "<<Ptarget <<" "<<Ptarget.mag() <<G4endl;
76  #endif
77 
78  G4LorentzVector Psum=Pprojectile+Ptarget;
79  G4double SqrtS=Psum.mag();
80  G4double S =Psum.mag2();
81 
82  #ifdef debugQuarkExchange
83  G4cout<<"SS Mpr Mtr SqrtS-Mprojectile-Mtarget "<<SqrtS<<" "<<Mprojectile<<" "<<Mtarget
84  <<" "<<SqrtS-Mprojectile-Mtarget<<G4endl;
85  #endif
86  if (SqrtS-Mprojectile-Mtarget <= 250.0*MeV) {
87  #ifdef debugQuarkExchange
88  G4cerr<<"Energy is too small for quark exchange!"<<G4endl;
89  G4cerr<<"Projectile: "<<projectile->GetDefinition()->GetPDGEncoding()<<" "
90  <<Pprojectile<<" "<<Pprojectile.mag()<<G4endl;
91  G4cerr<<"Target: "<<target->GetDefinition()->GetPDGEncoding()<<" "
92  <<Ptarget<<" "<<Ptarget.mag()<<G4endl;
93  G4cerr<<"sqrt(S) = "<<SqrtS<<" Mp + Mt = "<<Pprojectile.mag()+Ptarget.mag()<<G4endl;
94  #endif
95  return true;
96  }
97 
98  G4LorentzRotation toCms(-1*Psum.boostVector());
99 
100  G4LorentzVector Ptmp=toCms*Pprojectile;
101 
102  if ( Ptmp.pz() <= 0. )
103  {
104  // "String" moving backwards in CMS, abort collision !!
105  // G4cout << " abort Collision!! " << G4endl;
106  return false;
107  }
108 
109  toCms.rotateZ(-1*Ptmp.phi());
110  toCms.rotateY(-1*Ptmp.theta());
111 
112  G4LorentzRotation toLab(toCms.inverse());
113 
114  Pprojectile.transform(toCms);
115  Ptarget.transform(toCms);
116 
117  #ifdef debugQuarkExchange
118  G4cout << "Pprojectile in CMS " << Pprojectile << G4endl;
119  G4cout << "Ptarget in CMS " << Ptarget << G4endl;
120  #endif
121  G4double maxPtSquare=sqr(Ptarget.pz());
122 
123  G4double ProjectileMinDiffrMass = Pprojectile.mag()/GeV;
124  G4double TargetMinDiffrMass = Ptarget.mag()/GeV;
125 
126  G4double AveragePt2(0.);
127 
128  G4int PDGcode=projectile->GetDefinition()->GetPDGEncoding();
129  G4int absPDGcode=std::abs(PDGcode);
130 
131  G4bool ProjectileDiffraction = true;
132 
133  if ( absPDGcode > 1000 ) { ProjectileDiffraction = G4UniformRand() <= 0.5; }
134  if ( (absPDGcode == 211) || (absPDGcode == 111) ) { ProjectileDiffraction = G4UniformRand() <= 0.66; }
135  if ( (absPDGcode == 321) || (absPDGcode == 311) ||
136  ( PDGcode == 130) || ( PDGcode == 310) ) { ProjectileDiffraction = G4UniformRand() <= 0.5; }
137 
138  if ( ProjectileDiffraction ) {
139  if ( absPDGcode > 1000 ) //------Projectile is baryon --------
140  {
141  ProjectileMinDiffrMass = 1.16; // GeV
142  AveragePt2 = 0.3; // GeV^2
143  }
144  else if( absPDGcode == 211 || absPDGcode == 111) //------Projectile is Pion -----------
145  {
146  ProjectileMinDiffrMass = 1.0; // GeV
147  AveragePt2 = 0.3; // GeV^2
148  }
149  else if( absPDGcode == 321 || absPDGcode == 130 || absPDGcode == 310) //Projectile is Kaon
150  {
151  ProjectileMinDiffrMass = 1.1; // GeV
152  AveragePt2 = 0.3; // GeV^2
153  }
154  else if( absPDGcode == 22) //------Projectile is Gamma -----------
155  {
156  ProjectileMinDiffrMass = 0.25; // GeV
157  AveragePt2 = 0.36; // GeV^2
158  }
159  else //------Projectile is undefined, Nucleon assumed
160  {
161  ProjectileMinDiffrMass = 1.1; // GeV
162  AveragePt2 = 0.3; // GeV^2
163  };
164 
165  ProjectileMinDiffrMass = ProjectileMinDiffrMass * GeV;
166  Mprojectile2=sqr(ProjectileMinDiffrMass);
167 
168  if (G4UniformRand() <= 0.5) TargetMinDiffrMass += 0.22;
169  TargetMinDiffrMass *= GeV;
170  Mtarget2 = sqr( TargetMinDiffrMass) ;
171  }
172  else
173  {
174  if (G4UniformRand() <= 0.5) ProjectileMinDiffrMass += 0.22;
175  ProjectileMinDiffrMass *=GeV;
176  Mprojectile2=sqr(ProjectileMinDiffrMass);
177 
178  TargetMinDiffrMass = 1.16*GeV; // For target nucleon
179  Mtarget2 = sqr( TargetMinDiffrMass) ;
180  AveragePt2 = 0.3; // GeV^2
181  } // end of if ( ProjectileDiffraction )
182 
183  AveragePt2 = AveragePt2 * GeV*GeV;
184 
185  if ( SqrtS - (ProjectileMinDiffrMass+TargetMinDiffrMass) < 220* MeV ) return false;
186 
187  //-----------------------
188  G4double Pt2, PZcms, PZcms2;
189  G4double ProjMassT2, ProjMassT;
190  G4double TargMassT2, TargMassT;
191  G4double PMinusMin, PMinusMax, sqrtPMinusMin, sqrtPMinusMax;
192  G4double TPlusMin, TPlusMax, sqrtTPlusMin, sqrtTPlusMax;
193  G4double PMinusNew, PPlusNew, TPlusNew(0.), TMinusNew;
194 
195  G4LorentzVector Qmomentum;
196  G4double Qminus, Qplus;
197 
198  G4double x(0.), y(0.);
199  G4int whilecount=0;
200  do {
201  whilecount++;
202 
203  if (whilecount > 1000 )
204  {
205  Qmomentum=G4LorentzVector(0.,0.,0.,0.);
206  return false; // Ignore this interaction
207  }
208 
209  // Generate pt
210  Qmomentum=G4LorentzVector(GaussianPt(AveragePt2,maxPtSquare),0);
211 
212  Pt2 = G4ThreeVector( Qmomentum.vect() ).mag2();
213  ProjMassT2 = Mprojectile2 + Pt2;
214  ProjMassT = std::sqrt( ProjMassT2 );
215  TargMassT2 = Mtarget2 + Pt2;
216  TargMassT = std::sqrt( TargMassT2 );
217 
218  #ifdef debugQuarkExchange
219  G4cout<<"whilecount Pt2 ProjMassT TargMassT SqrtS S ProjectileDiffraction"<<G4endl;
220  G4cout<<whilecount<<" "<<Pt2<<" "<<ProjMassT<<" "<<TargMassT<<" "<<SqrtS<<" "<<S<<" "<<ProjectileDiffraction<<G4endl;
221  #endif
222 
223  if ( SqrtS < ProjMassT + TargMassT + 220.0*MeV ) continue;
224 
225  PZcms2 = ( S*S + ProjMassT2*ProjMassT2 + TargMassT2*TargMassT2
226  - 2.0*S*ProjMassT2 - 2.0*S*TargMassT2 - 2.0*ProjMassT2*TargMassT2 ) / 4.0 / S;
227 
228  if ( PZcms2 < 0 ) continue;
229 
230  PZcms = std::sqrt( PZcms2 );
231 
232  if ( ProjectileDiffraction )
233  { // The projectile will fragment, the target will saved.
234  PMinusMin = std::sqrt( ProjMassT2 + PZcms2 ) - PZcms;
235  PMinusMax = SqrtS - TargMassT;
236  sqrtPMinusMin = std::sqrt(PMinusMin); sqrtPMinusMax = std::sqrt(PMinusMax);
237 
238  if ( absPDGcode > 1000 )
239  {
240  PMinusNew = PMinusMax * (1.0 - (1.0 - PMinusMin/PMinusMax)
241  * G4Pow::GetInstance()->powA(G4UniformRand(),0.3333) );
242  }
243  else if ( (absPDGcode == 211) || (absPDGcode == 111) )
244  {
245  while (true)
246  {
247  x=sqrtPMinusMax-(sqrtPMinusMax-sqrtPMinusMin)*G4UniformRand();
248  y=G4UniformRand();
249  if ( y < 1.0-0.7 * x/sqrtPMinusMax ) break;
250  }
251  PMinusNew = sqr(x);
252  }
253  else if ( (absPDGcode == 321) || (absPDGcode == 311) ||
254  ( PDGcode == 130) || ( PDGcode == 310) )
255  { // For K-mesons it must be found !!! Uzhi
256  while (true)
257  {
258  x=sqrtPMinusMax-(sqrtPMinusMax-sqrtPMinusMin)*G4UniformRand();
259  y=G4UniformRand();
260  if ( y < 1.0-0.7 * x/sqrtPMinusMax ) break;
261  }
262  PMinusNew = sqr(x);
263  }
264  else
265  {
266  PMinusNew = PMinusMax * (1.0 - (1.0 - PMinusMin/PMinusMax)
267  * G4Pow::GetInstance()->powA(G4UniformRand(),0.3333) );
268  };
269 
270  TMinusNew = SqrtS - PMinusNew;
271 
272  Qminus = Ptarget.minus() - TMinusNew;
273  TPlusNew = TargMassT2 / TMinusNew;
274  Qplus = Ptarget.plus() - TPlusNew;
275 
276  }
277  else
278  { // The target will fragment, the projectile will saved.
279  TPlusMin = std::sqrt( TargMassT2 + PZcms2 ) - PZcms;
280  TPlusMax = SqrtS - ProjMassT;
281  sqrtTPlusMin = std::sqrt(TPlusMin); sqrtTPlusMax = std::sqrt(TPlusMax);
282 
283  if ( absPDGcode > 1000 )
284  {
285  TPlusNew = TPlusMax * (1.0 - (1.0 - TPlusMin/TPlusMax)
286  * G4Pow::GetInstance()->powA(G4UniformRand(),0.3333) );
287  }
288  else if ( (absPDGcode == 211) || (absPDGcode == 111) )
289  {
290  while (true)
291  {
292  x=sqrtTPlusMax-(sqrtTPlusMax-sqrtTPlusMin)*G4UniformRand();
293  y=G4UniformRand();
294  if ( y < 1.0-0.7 * x/sqrtTPlusMax ) break;
295  }
296  TPlusNew = sqr(x);
297  }
298  else if ( (absPDGcode == 321) || (absPDGcode == 311) ||
299  ( PDGcode == 130) || ( PDGcode == 310) )
300  { // For K-mesons it must be found !!! Uzhi
301  while (true)
302  {
303  x=sqrtTPlusMax-(sqrtTPlusMax-sqrtTPlusMin)*G4UniformRand();
304  y=G4UniformRand();
305  if ( y < 1.0-0.7 * x/sqrtTPlusMax ) break;
306  }
307  }
308  else
309  {
310  TPlusNew = TPlusMax * (1.0 - (1.0 - TPlusMin/TPlusMax)
311  * G4Pow::GetInstance()->powA(G4UniformRand(),0.3333) );
312  };
313 
314  PPlusNew = SqrtS - TPlusNew;
315 
316  Qplus = PPlusNew - Pprojectile.plus();
317  PMinusNew = ProjMassT2 / PPlusNew;
318  Qminus = PMinusNew - Pprojectile.minus();
319  }
320 
321  Qmomentum.setPz( (Qplus - Qminus)/2 );
322  Qmomentum.setE( (Qplus + Qminus)/2 );
323 
324  #ifdef debugQuarkExchange
325  G4cout<<"ProjectileDiffraction (Pprojectile + Qmomentum).mag2() Mprojectile2"<<G4endl;
326  G4cout<<ProjectileDiffraction<<" "<<( Pprojectile + Qmomentum ).mag2()<<" "<< Mprojectile2<<G4endl;
327  G4cout<<"TargetDiffraction (Ptarget - Qmomentum).mag2() Mtarget2"<<G4endl;
328  G4cout<<!ProjectileDiffraction<<" "<<( Ptarget - Qmomentum ).mag2()<<" "<< Mtarget2<<G4endl;
329  #endif
330 
331  } while ( ( ProjectileDiffraction&&( Pprojectile + Qmomentum ).mag2() < Mprojectile2 ) ||
332  (!ProjectileDiffraction&&( Ptarget - Qmomentum ).mag2() < Mtarget2 ) );
333  // Repeat the sampling because there was not any excitation
334 
335  Pprojectile += Qmomentum;
336 
337  Ptarget -= Qmomentum;
338 
339  // Transform back and update SplitableHadron Participant.
340  Pprojectile.transform(toLab);
341  Ptarget.transform(toLab);
342 
343  #ifdef debugQuarkExchange
344  G4cout << "Pprojectile in Lab. " << Pprojectile << G4endl;
345  G4cout << "Ptarget in Lab. " << Ptarget << G4endl;
346  G4cout << "G4QuarkExchange: Projectile mass " << Pprojectile.mag() << G4endl;
347  G4cout << "G4QuarkExchange: Target mass " << Ptarget.mag() << G4endl;
348  #endif
349 
350  target->Set4Momentum(Ptarget);
351  projectile->Set4Momentum(Pprojectile);
352 
353  //=================================== Quark exchange ================================
354  projectile->SplitUp();
355  target->SplitUp();
356 
357  G4Parton* PrQuark = projectile->GetNextParton();
358  G4Parton* TrQuark = target->GetNextParton();
359 
360  G4ParticleDefinition * Tmp = PrQuark->GetDefinition();
361  PrQuark->SetDefinition(TrQuark->GetDefinition());
362  TrQuark->SetDefinition(Tmp);
363 
364  return true;
365 }
366 
367 
368 // --------- private methods ----------------------
369 
371 { // @@ this method is used in FTFModel as well. Should go somewhere common!
372 
373  G4double pt2;
374 
375  const G4int maxNumberOfLoops = 1000;
376  G4int loopCounter = 0;
377  do {
378  pt2=-widthSquare * G4Log( G4UniformRand() );
379  } while ( ( pt2 > maxPtSquare) && ++loopCounter < maxNumberOfLoops ); /* Loop checking, 07.08.2015, A.Ribon */
380  if ( loopCounter >= maxNumberOfLoops ) {
381  pt2 = 0.99*maxPtSquare; // Just an acceptable value, without any physics consideration.
382  }
383 
384  pt2=std::sqrt(pt2);
385 
387 
388  return G4ThreeVector (pt2*std::cos(phi), pt2*std::sin(phi), 0.);
389 }
390