Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
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 // $Id: G4QuarkExchange.cc 99348 2016-09-19 08:39:04Z vuzhinsk $
28 // ------------------------------------------------------------
29 // GEANT 4 class implemetation file
30 //
31 // ---------------- G4QuarkExchange --------------
32 // by V. Uzhinsky, October 2016.
33 // QuarkExchange is used by strings models.
34 // Take a projectile and a target.
35 //Simulate Q exchange with excitation of projectile or target.
36 // ------------------------------------------------------------
37 
38 #include "G4QuarkExchange.hh"
39 #include "globals.hh"
40 #include "G4PhysicalConstants.hh"
41 #include "G4SystemOfUnits.hh"
42 #include "Randomize.hh"
43 #include "G4LorentzRotation.hh"
44 #include "G4ThreeVector.hh"
45 #include "G4ParticleDefinition.hh"
46 #include "G4VSplitableHadron.hh"
47 #include "G4ExcitedString.hh"
48 
49 #include "G4Log.hh"
50 #include "G4Pow.hh"
51 
52 //#define debugQuarkExchange
53 
55 
57 
60 {
61  #ifdef debugQuarkExchange
62  G4cout<<G4endl<<"G4QuarkExchange::ExciteParticipants"<<G4endl;
63  #endif
64 
65  G4LorentzVector Pprojectile = projectile->Get4Momentum();
66  G4double Mprojectile = projectile->GetDefinition()->GetPDGMass();
67  G4double Mprojectile2 = sqr(Mprojectile);
68 
69  G4LorentzVector Ptarget = target->Get4Momentum();
70  G4double Mtarget = target->GetDefinition()->GetPDGMass();
71  G4double Mtarget2 = sqr(Mtarget);
72 
73  #ifdef debugQuarkExchange
74  G4cout<<"Proj Targ "<<projectile->GetDefinition()->GetPDGEncoding()<<" "<<target->GetDefinition()->GetPDGEncoding()<<G4endl;
75  G4cout<<"Proj. 4-Mom "<<Pprojectile<<" "<<Pprojectile.mag()<<G4endl
76  <<"Targ. 4-Mom "<<Ptarget <<" "<<Ptarget.mag() <<G4endl;
77  #endif
78 
79  G4LorentzVector Psum=Pprojectile+Ptarget;
80  G4double SqrtS=Psum.mag();
81  G4double S =Psum.mag2();
82 
83  #ifdef debugQuarkExchange
84  G4cout<<"SS Mpr Mtr SqrtS-Mprojectile-Mtarget "<<SqrtS<<" "<<Mprojectile<<" "<<Mtarget
85  <<" "<<SqrtS-Mprojectile-Mtarget<<G4endl;
86  #endif
87  if(SqrtS-Mprojectile-Mtarget <= 250.0*MeV) {
88  #ifdef debugQuarkExchange
89  G4cerr<<"Energy is too small for quark exchange!"<<G4endl;
90  G4cerr<<"Projectile: "<<projectile->GetDefinition()->GetPDGEncoding()<<" "
91  <<Pprojectile<<" "<<Pprojectile.mag()<<G4endl;
92  G4cerr<<"Target: "<<target->GetDefinition()->GetPDGEncoding()<<" "
93  <<Ptarget<<" "<<Ptarget.mag()<<G4endl;
94  G4cerr<<"sqrt(S) = "<<SqrtS<<" Mp + Mt = "<<Pprojectile.mag()+Ptarget.mag()<<G4endl;
95  #endif
96  return true;
97  }
98 
99  G4LorentzRotation toCms(-1*Psum.boostVector());
100 
101  G4LorentzVector Ptmp=toCms*Pprojectile;
102 
103  if ( Ptmp.pz() <= 0. )
104  {
105  // "String" moving backwards in CMS, abort collision !!
106  // G4cout << " abort Collision!! " << G4endl;
107  return false;
108  }
109 
110  toCms.rotateZ(-1*Ptmp.phi());
111  toCms.rotateY(-1*Ptmp.theta());
112 
113  G4LorentzRotation toLab(toCms.inverse());
114 
115  Pprojectile.transform(toCms);
116  Ptarget.transform(toCms);
117 
118  #ifdef debugQuarkExchange
119  G4cout << "Pprojectile in CMS " << Pprojectile << G4endl;
120  G4cout << "Ptarget in CMS " << Ptarget << G4endl;
121  #endif
122  G4double maxPtSquare=sqr(Ptarget.pz());
123 
124 // G4double ProjectileMinDiffrMass(0.), TargetMinDiffrMass(0.);
125  G4double ProjectileMinDiffrMass = Pprojectile.mag()/GeV;
126  G4double TargetMinDiffrMass = Ptarget.mag()/GeV;
127 
128  G4double AveragePt2(0.);
129 
130  G4int PDGcode=projectile->GetDefinition()->GetPDGEncoding();
131  G4int absPDGcode=std::abs(PDGcode);
132 
133  G4bool ProjectileDiffraction = true;
134 
135  if( absPDGcode > 1000 ) { ProjectileDiffraction = G4UniformRand() <= 0.5; }
136  if( (absPDGcode == 211) || (absPDGcode == 111) ) { ProjectileDiffraction = G4UniformRand() <= 0.66; }
137  if( (absPDGcode == 321) || (absPDGcode == 311) ||
138  ( PDGcode == 130) || ( PDGcode == 310) ) { ProjectileDiffraction = G4UniformRand() <= 0.5; } // Uzhi ???
139 
140  if ( ProjectileDiffraction ) {
141  if( absPDGcode > 1000 ) //------Projectile is baryon --------
142  {
143  ProjectileMinDiffrMass = 1.16; // GeV
144  AveragePt2 = 0.3; // GeV^2
145  }
146  else if( absPDGcode == 211 || absPDGcode == 111) //------Projectile is Pion -----------
147  {
148  ProjectileMinDiffrMass = 1.0; // GeV
149  AveragePt2 = 0.3; // GeV^2
150  }
151  else if( absPDGcode == 321 || absPDGcode == 130 || absPDGcode == 310) //Projectile is Kaon
152  {
153  ProjectileMinDiffrMass = 1.1; // GeV
154  AveragePt2 = 0.3; // GeV^2
155  }
156  else if( absPDGcode == 22) //------Projectile is Gamma -----------
157  {
158  ProjectileMinDiffrMass = 0.25; // GeV
159  AveragePt2 = 0.36; // GeV^2
160  }
161  else //------Projectile is undefined, Nucleon assumed
162  {
163  ProjectileMinDiffrMass = 1.1; // GeV
164  AveragePt2 = 0.3; // GeV^2
165  };
166 
167  ProjectileMinDiffrMass = ProjectileMinDiffrMass * GeV;
168  Mprojectile2=sqr(ProjectileMinDiffrMass);
169 
170  if(G4UniformRand() <= 0.5) TargetMinDiffrMass += 0.22;
171  TargetMinDiffrMass *= GeV;
172  Mtarget2 = sqr( TargetMinDiffrMass) ;
173  }
174  else
175  {
176  if(G4UniformRand() <= 0.5) ProjectileMinDiffrMass += 0.22;
177  ProjectileMinDiffrMass *=GeV;
178  Mprojectile2=sqr(ProjectileMinDiffrMass);
179 
180  TargetMinDiffrMass = 1.16*GeV; // For target nucleon
181  Mtarget2 = sqr( TargetMinDiffrMass) ;
182  AveragePt2 = 0.3; // GeV^2
183  } // end of if ( ProjectileDiffraction )
184 
185  AveragePt2 = AveragePt2 * GeV*GeV; // Uzhi 6 Oct. 2016
186 
187 if( SqrtS - (ProjectileMinDiffrMass+TargetMinDiffrMass) < 220* MeV ) return false; // Uzhi Dec. 2017
188 //-----------------------
189  G4double Pt2, PZcms, PZcms2;
190  G4double ProjMassT2, ProjMassT;
191  G4double TargMassT2, TargMassT;
192  G4double PMinusMin, PMinusMax, sqrtPMinusMin, sqrtPMinusMax;
193  //G4double PPlusMin , PPlusMax;
194  G4double TPlusMin, TPlusMax, sqrtTPlusMin, sqrtTPlusMax;
195  G4double PMinusNew, PPlusNew, TPlusNew(0.), TMinusNew;
196 
197  G4LorentzVector Qmomentum;
198  G4double Qminus, Qplus;
199 
200  G4double x(0.), y(0.);
201  G4int whilecount=0;
202  do {
203  whilecount++;
204 
205  if (whilecount > 1000 )
206  {
207  Qmomentum=G4LorentzVector(0.,0.,0.,0.);
208  return false; // Ignore this interaction
209  }
210 
211  // Generate pt
212  Qmomentum=G4LorentzVector(GaussianPt(AveragePt2,maxPtSquare),0);
213 
214  Pt2 = G4ThreeVector( Qmomentum.vect() ).mag2();
215  ProjMassT2 = Mprojectile2 + Pt2;
216  ProjMassT = std::sqrt( ProjMassT2 );
217  TargMassT2 = Mtarget2 + Pt2;
218  TargMassT = std::sqrt( TargMassT2 );
219 
220  #ifdef debugQuarkExchange
221  G4cout<<"whilecount Pt2 ProjMassT TargMassT SqrtS S ProjectileDiffraction"<<G4endl;
222  G4cout<<whilecount<<" "<<Pt2<<" "<<ProjMassT<<" "<<TargMassT<<" "<<SqrtS<<" "<<S<<" "<<ProjectileDiffraction<<G4endl;
223  #endif
224 
225  if ( SqrtS < ProjMassT + TargMassT + 220.0*MeV ) continue;
226 
227  PZcms2 = ( S*S + ProjMassT2*ProjMassT2 + TargMassT2*TargMassT2
228  - 2.0*S*ProjMassT2 - 2.0*S*TargMassT2 - 2.0*ProjMassT2*TargMassT2 ) / 4.0 / S;
229 
230  if ( PZcms2 < 0 ) continue;
231 
232  PZcms = std::sqrt( PZcms2 );
233 
234  if ( ProjectileDiffraction )
235  {// The projectile will fragment, the target will saved.
236  PMinusMin = std::sqrt( ProjMassT2 + PZcms2 ) - PZcms;
237  PMinusMax = SqrtS - TargMassT;
238  sqrtPMinusMin = std::sqrt(PMinusMin); sqrtPMinusMax = std::sqrt(PMinusMax);
239 //===========
240 
241  if( absPDGcode > 1000 )
242  {
243 /* // Uzhi Dec. 2017
244  while(true)
245  {
246  x=PMinusMax-(PMinusMax-PMinusMin)*G4UniformRand();
247  y=G4UniformRand();
248 
249  if( y < G4Pow::GetInstance()->powA(1.0-x/PMinusMax,3.) ) break;
250 //if( y < G4Pow::GetInstance()->powA(x/PMinusMax,3.) ) break; // Uzhi
251  }
252 */ // Uzhi Dec. 2017
253  PMinusNew = PMinusMax * (1.0 - (1.0 - PMinusMin/PMinusMax)
254  * G4Pow::GetInstance()->powA(G4UniformRand(),0.3333) ); // Uzhi Dec. 2017
255 
256  } else if( (absPDGcode == 211) || (absPDGcode == 111) )
257  {
258  while(true)
259  {
260  x=sqrtPMinusMax-(sqrtPMinusMax-sqrtPMinusMin)*G4UniformRand();
261  y=G4UniformRand();
262  if( y < 1.0-0.7 * x/sqrtPMinusMax ) break; // 0.7 for Pi Uzhi
263  }
264  PMinusNew = sqr(x);
265  } else if( (absPDGcode == 321) || (absPDGcode == 311) ||
266  ( PDGcode == 130) || ( PDGcode == 310) )
267  { // For K-mesons it must be found !!! Uzhi
268  while(true)
269  {
270  x=sqrtPMinusMax-(sqrtPMinusMax-sqrtPMinusMin)*G4UniformRand();
271  y=G4UniformRand();
272  if( y < 1.0-0.7 * x/sqrtPMinusMax ) break;
273  }
274  PMinusNew = sqr(x);
275  } else
276  {
277 /* // Uzhi Dec. 2017
278  while(true)
279  {
280  x=PMinusMax-(PMinusMax-PMinusMin)*G4UniformRand();
281  y=G4UniformRand();
282  if( y < G4Pow::GetInstance()->powA(1.0-x/PMinusMax,3.) ) break;
283  }
284 */ // Uzhi Dec. 2017
285  PMinusNew = PMinusMax * (1.0 - (1.0 - PMinusMin/PMinusMax)
286  * G4Pow::GetInstance()->powA(G4UniformRand(),0.3333) ); // Uzhi Dec. 2017
287 
288  };
289 
290  TMinusNew = SqrtS - PMinusNew;
291 
292  Qminus = Ptarget.minus() - TMinusNew;
293  TPlusNew = TargMassT2 / TMinusNew;
294  Qplus = Ptarget.plus() - TPlusNew;
295 
296  } else
297 //-------------------------------------------------------------------------------
298  {// The target will fragment, the projectile will saved.
299  TPlusMin = std::sqrt( TargMassT2 + PZcms2 ) - PZcms;
300  TPlusMax = SqrtS - ProjMassT;
301  sqrtTPlusMin = std::sqrt(TPlusMin); sqrtTPlusMax = std::sqrt(TPlusMax);
302 
303  if( absPDGcode > 1000 )
304  {
305 /* // Uzhi Dec. 2017
306  while(true)
307  {
308  x=TPlusMax-(TPlusMax-TPlusMin)*G4UniformRand();
309  y=G4UniformRand();
310  if( y < G4Pow::GetInstance()->powA(1.0-x/TPlusMax,3.) ) break;
311 //if( y < G4Pow::GetInstance()->powA(x/TPlusMax,3.) ) break; // Uzhi
312  }
313 */ // Uzhi Dec. 2017
314  TPlusNew = TPlusMax * (1.0 - (1.0 - TPlusMin/TPlusMax)
315  * G4Pow::GetInstance()->powA(G4UniformRand(),0.3333) ); // Uzhi Dec. 2017
316 
317 
318  } else if( (absPDGcode == 211) || (absPDGcode == 111) )
319  {
320  while(true)
321  {
322  x=sqrtTPlusMax-(sqrtTPlusMax-sqrtTPlusMin)*G4UniformRand();
323  y=G4UniformRand();
324  if( y < 1.0-0.7 * x/sqrtTPlusMax ) break; // 0.7 for Pi Uzhi
325  }
326  TPlusNew = sqr(x);
327  } else if( (absPDGcode == 321) || (absPDGcode == 311) ||
328  ( PDGcode == 130) || ( PDGcode == 310) )
329  { // For K-mesons it must be found !!! Uzhi
330  while(true)
331  {
332  x=sqrtTPlusMax-(sqrtTPlusMax-sqrtTPlusMin)*G4UniformRand();
333  y=G4UniformRand();
334  if( y < 1.0-0.7 * x/sqrtTPlusMax ) break;
335  }
336  } else
337  {
338 /* // Uzhi Dec. 2017
339  while(true)
340  {
341  x=TPlusMax-(TPlusMax-TPlusMin)*G4UniformRand();
342  y=G4UniformRand();
343  if( y < G4Pow::GetInstance()->powA(1.0-x/TPlusMax,3.) ) break;
344  }
345  TPlusNew = sqr(x);
346 */ // Uzhi Dec. 2017
347  TPlusNew = TPlusMax * (1.0 - (1.0 - TPlusMin/TPlusMax)
348  * G4Pow::GetInstance()->powA(G4UniformRand(),0.3333) ); // Uzhi Dec. 2017
349  };
350 
351  PPlusNew = SqrtS - TPlusNew;
352 
353  Qplus = PPlusNew - Pprojectile.plus();
354  PMinusNew = ProjMassT2 / PPlusNew;
355  Qminus = PMinusNew - Pprojectile.minus();
356  }
357 
358 
359  Qmomentum.setPz( (Qplus - Qminus)/2 );
360  Qmomentum.setE( (Qplus + Qminus)/2 );
361 
362  #ifdef debugQuarkExchange
363  G4cout<<"ProjectileDiffraction (Pprojectile + Qmomentum).mag2() Mprojectile2"<<G4endl;
364  G4cout<<ProjectileDiffraction<<" "<<( Pprojectile + Qmomentum ).mag2()<<" "<< Mprojectile2<<G4endl;
365  G4cout<<"TargetDiffraction (Ptarget - Qmomentum).mag2() Mtarget2"<<G4endl;
366  G4cout<<!ProjectileDiffraction<<" "<<( Ptarget - Qmomentum ).mag2()<<" "<< Mtarget2<<G4endl;
367  #endif
368 
369  } while ( ( ProjectileDiffraction&&( Pprojectile + Qmomentum ).mag2() < Mprojectile2 ) ||
370  (!ProjectileDiffraction&&( Ptarget - Qmomentum ).mag2() < Mtarget2 ) );
371  // Repeat the sampling because there was not any excitation
372 
373 
374  Pprojectile += Qmomentum;
375 
376  Ptarget -= Qmomentum;
377 
378  // Transform back and update SplitableHadron Participant.
379  Pprojectile.transform(toLab);
380  Ptarget.transform(toLab);
381 
382  #ifdef debugQuarkExchange
383  G4cout << "Pprojectile in Lab. " << Pprojectile << G4endl;
384  G4cout << "Ptarget in Lab. " << Ptarget << G4endl;
385  G4cout << "G4QuarkExchange: Projectile mass " << Pprojectile.mag() << G4endl;
386  G4cout << "G4QuarkExchange: Target mass " << Ptarget.mag() << G4endl;
387  #endif
388 
389  target->Set4Momentum(Ptarget);
390  projectile->Set4Momentum(Pprojectile);
391 
392 //=================================== Quark exchange ================================
393  projectile->SplitUp();
394  target->SplitUp();
395 
396  G4Parton* PrQuark = projectile->GetNextParton();
397  G4Parton* TrQuark = target->GetNextParton();
398 
399 //G4cout<<"Pr quark "<<PrQuark->GetDefinition()->GetParticleName()<<G4endl;
400 //G4cout<<"Tr quark "<<TrQuark->GetDefinition()->GetParticleName()<<G4endl;
401 
402 
403  G4ParticleDefinition * Tmp = PrQuark->GetDefinition();
404  PrQuark->SetDefinition(TrQuark->GetDefinition());
405  TrQuark->SetDefinition(Tmp);
406 
407 //============================================
408  return true;
409 }
410 
411 
412 // --------- private methods ----------------------
413 
415 { // @@ this method is used in FTFModel as well. Should go somewhere common!
416 
417  G4double pt2;
418 
419  const G4int maxNumberOfLoops = 1000;
420  G4int loopCounter = 0;
421  do {
422  pt2=-widthSquare * G4Log( G4UniformRand() );
423  } while ( ( pt2 > maxPtSquare) && ++loopCounter < maxNumberOfLoops ); /* Loop checking, 07.08.2015, A.Ribon */
424  if ( loopCounter >= maxNumberOfLoops ) {
425  pt2 = 0.99*maxPtSquare; // Just an acceptable value, without any physics consideration.
426  }
427 
428  pt2=std::sqrt(pt2);
429 
430  G4double phi=G4UniformRand() * twopi;
431 
432  return G4ThreeVector (pt2*std::cos(phi), pt2*std::sin(phi), 0.);
433 }
434 
435 
436 
437 
438 
Float_t x
Definition: compare.C:6
double mag2() const
const G4ParticleDefinition * GetDefinition() const
HepLorentzVector & transform(const HepRotation &)
CLHEP::Hep3Vector G4ThreeVector
const XML_Char * target
Definition: expat.h:268
G4ParticleDefinition * GetDefinition()
Definition: G4Parton.hh:161
static constexpr double MeV
Definition: G4SIunits.hh:214
#define G4endl
Definition: G4ios.hh:61
Float_t y
Definition: compare.C:6
void SetDefinition(G4ParticleDefinition *aDefinition)
Definition: G4Parton.hh:166
G4ThreeVector GaussianPt(G4double widthSquare, G4double maxPtSquare) const
virtual G4Parton * GetNextParton()=0
G4double GetPDGMass() const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
double S(double temp)
G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner) const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
static G4Pow * GetInstance()
Definition: G4Pow.cc:57
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:242
#define G4UniformRand()
Definition: Randomize.hh:53
static constexpr double twopi
Definition: G4SIunits.hh:76
HepLorentzVector & rotateY(double)
virtual void SplitUp()=0
G4GLOB_DLL std::ostream G4cerr
HepLorentzVector & rotateZ(double)
double mag() const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
Hep3Vector boostVector() const
const G4LorentzVector & Get4Momentum() const
Hep3Vector vect() const
T sqr(const T &x)
Definition: templates.hh:145
CLHEP::HepLorentzVector G4LorentzVector
static constexpr double GeV
Definition: G4SIunits.hh:217
HepLorentzRotation & transform(const HepBoost &b)
void Set4Momentum(const G4LorentzVector &a4Momentum)