Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/cross_sections/src/G4ComponentGGNuclNuclXsc.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 // 24.11.08 V. Grichine - first implementation
 27 //
 28 // 04.09.18 V. Ivantchenko Major revision of interfaces and implementation
 29 // 27.05.19 V. Ivantchenko Removed obsolete methods and members 
 30 
 31 #include "G4ComponentGGNuclNuclXsc.hh"
 32 
 33 #include "G4PhysicalConstants.hh"
 34 #include "G4SystemOfUnits.hh"
 35 #include "G4NucleiProperties.hh"
 36 #include "G4ParticleDefinition.hh"
 37 #include "G4HadronNucleonXsc.hh"
 38 #include "G4ComponentGGHadronNucleusXsc.hh" 
 39 #include "G4NuclearRadii.hh"
 40 #include "G4Pow.hh"
 41 
 42 namespace
 43 {
 44   const G4double inve = 1./CLHEP::eplus;
 45 }
 46 
 47 G4ComponentGGNuclNuclXsc::G4ComponentGGNuclNuclXsc() 
 48  : G4VComponentCrossSection("Glauber-Gribov Nucl-nucl")
 49 {
 50   theProton   = G4Proton::Proton();
 51   theNeutron  = G4Neutron::Neutron();
 52   theLambda  = G4Lambda::Lambda();
 53   fHNXsc = new G4HadronNucleonXsc();
 54   fHadrNucl = new G4ComponentGGHadronNucleusXsc(); 
 55 }
 56 
 57 G4ComponentGGNuclNuclXsc::~G4ComponentGGNuclNuclXsc()
 58 {
 59   delete fHNXsc;
 60 }
 61 
 62 //////////////////////////////////////////////////////////////////////
 63 
 64 G4double G4ComponentGGNuclNuclXsc::GetTotalElementCrossSection(
 65          const G4ParticleDefinition* aParticle, G4double kinEnergy, 
 66    G4int Z, G4double A)
 67 {
 68   ComputeCrossSections(aParticle, kinEnergy, Z, G4lrint(A));
 69   return fTotalXsc;
 70 }
 71 
 72 ////////////////////////////////////////////////////////////////////
 73 
 74 G4double G4ComponentGGNuclNuclXsc::GetTotalIsotopeCrossSection(
 75          const G4ParticleDefinition* aParticle, G4double kinEnergy,
 76    G4int Z, G4int A)
 77 {
 78   ComputeCrossSections(aParticle, kinEnergy, Z, A);
 79   return fTotalXsc;
 80 }
 81 
 82 /////////////////////////////////////////////////////////////////////
 83 
 84 G4double G4ComponentGGNuclNuclXsc::GetInelasticElementCrossSection(
 85          const G4ParticleDefinition* aParticle, G4double kinEnergy, 
 86    G4int Z, G4double A)
 87 {
 88   ComputeCrossSections(aParticle, kinEnergy, Z, G4lrint(A));
 89   return fInelasticXsc;
 90 }
 91 
 92 ////////////////////////////////////////////////////////////////////
 93 
 94 G4double G4ComponentGGNuclNuclXsc::GetInelasticIsotopeCrossSection(
 95          const G4ParticleDefinition* aParticle, G4double kinEnergy, 
 96    G4int Z, G4int A)
 97 {
 98   ComputeCrossSections(aParticle, kinEnergy, Z, A);
 99   return fInelasticXsc;
100 }
101 
102 //////////////////////////////////////////////////////////////////
103 
104 G4double G4ComponentGGNuclNuclXsc::GetElasticElementCrossSection(
105          const G4ParticleDefinition* aParticle, G4double kinEnergy, 
106    G4int Z, G4double A)
107 {
108   ComputeCrossSections(aParticle, kinEnergy, Z, G4lrint(A));
109   return fElasticXsc;
110 }
111 
112 ///////////////////////////////////////////////////////////////////
113 
114 G4double G4ComponentGGNuclNuclXsc::GetElasticIsotopeCrossSection(
115          const G4ParticleDefinition* aParticle, G4double kinEnergy, 
116    G4int Z, G4int A)
117 {
118   ComputeCrossSections(aParticle, kinEnergy, Z, A);
119   return fElasticXsc;
120 }
121 
122 ////////////////////////////////////////////////////////////////
123  
124 G4double G4ComponentGGNuclNuclXsc::ComputeQuasiElasticRatio(
125          const G4ParticleDefinition* aParticle, G4double kinEnergy, 
126    G4int Z, G4int A)
127 {
128   ComputeCrossSections(aParticle, kinEnergy, Z, A);
129   return (fInelasticXsc > fProductionXsc) 
130     ? (fInelasticXsc - fProductionXsc)/fInelasticXsc : 0.0;
131 }
132 
133 //////////////////////////////////////////////////////////////////////
134  
135 void G4ComponentGGNuclNuclXsc::BuildPhysicsTable(const G4ParticleDefinition&)
136 {}
137 
138 //////////////////////////////////////////////////////////////////////
139 
140 void G4ComponentGGNuclNuclXsc::DumpPhysicsTable(const G4ParticleDefinition&)
141 {
142   G4cout << "G4ComponentGGNuclNuclXsc: uses Glauber-Gribov formula" << G4endl;
143 }
144 
145 //////////////////////////////////////////////////////////////////////
146 
147 void G4ComponentGGNuclNuclXsc::Description(std::ostream& outFile) const
148 {
149   outFile << "G4ComponentGGNuclNuclXsc calculates total, inelastic and\n"
150           << "elastic cross sections for nucleus-nucleus collisions using\n"
151           << "the Glauber model with Gribov corrections.  It is valid for\n"
152           << "all incident energies above 100 keV./n"
153     << "For the hydrogen target G4HadronNucleonXsc class is used.\n";
154 }
155 
156 ///////////////////////////////////////////////////////////////////////////////
157 //
158 // Calculates total and inelastic Xsc, derives elastic as total - inelastic 
159 // accordong to Glauber model with Gribov correction calculated in the dipole 
160 // approximation on light cone. Gaussian density of point-like nucleons helps 
161 // to calculate rest integrals of the model. [1] B.Z. Kopeliovich,
162 // nucl-th/0306044 + simplification above
163 
164 void G4ComponentGGNuclNuclXsc::ComputeCrossSections(
165      const G4ParticleDefinition* aParticle, G4double kinEnergy, 
166      G4int Z, G4int A)
167 {
168   // check cache
169   if(aParticle == fParticle && fZ == Z && fA == A && kinEnergy == fEnergy)
170     { return; }
171   fParticle = aParticle;
172   fZ = Z;
173   fA = A;
174   fEnergy = kinEnergy;
175   G4Pow* pG4Pow = G4Pow::GetInstance();
176   
177   G4int pZ = G4lrint(aParticle->GetPDGCharge()*inve);
178   G4int pA = aParticle->GetBaryonNumber();
179   G4int pL = aParticle->GetNumberOfLambdasInHypernucleus();
180   G4bool pHN = aParticle->IsHypernucleus();
181   G4double cHN(0.88);
182 
183   // hydrogen
184   if(1 == Z && 1 == A) {
185     G4double e = kinEnergy*CLHEP::proton_mass_c2/aParticle->GetPDGMass();
186     fHadrNucl->ComputeCrossSections( theProton, e, pZ, pA, pL );
187     fTotalXsc = fHadrNucl->GetTotalGlauberGribovXsc(); 
188     fElasticXsc = fHadrNucl->GetElasticGlauberGribovXsc(); 
189     fInelasticXsc = fHadrNucl->GetInelasticGlauberGribovXsc(); 
190     fProductionXsc = fHadrNucl->GetProductionGlauberGribovXsc(); 
191     fDiffractionXsc = fHadrNucl->GetDiffractionGlauberGribovXsc(); 
192     return;
193   }
194   static const G4double cofInelastic = 2.4;
195   static const G4double cofTotal = 2.0;
196 
197   G4double pTkin = kinEnergy/(G4double)pA;
198 
199   G4int pN = pA - pZ;
200   G4int tN = A - Z;
201 
202   G4double tR = G4NuclearRadii::Radius(Z, A);  
203   G4double pR = G4NuclearRadii::Radius(pZ, pA);
204   
205   if(pHN)
206     pR *= std::sqrt( pG4Pow->Z23( pA - pL ) + cHN*pG4Pow->Z23( pL ) )/pG4Pow->Z13(pA);
207   
208   G4double cB = ComputeCoulombBarier(aParticle, kinEnergy, Z, A, pR, tR);
209 
210   if ( cB > 0. ) 
211   {
212     G4double sigma = (pZ*Z+pN*tN)*fHNXsc->HadronNucleonXscNS(theProton, theProton, pTkin);
213     if(pHN) sigma += pL*A*fHNXsc->HadronNucleonXsc(theLambda, theProton, pTkin);
214     G4double ppInXsc = fHNXsc->GetInelasticHadronNucleonXsc();
215 
216     sigma += (pZ*tN+pN*Z)*fHNXsc->HadronNucleonXscNS(theNeutron, theProton, pTkin);
217     G4double npInXsc = fHNXsc->GetInelasticHadronNucleonXsc();
218 
219     G4double nucleusSquare = cofTotal*CLHEP::pi*( pR*pR + tR*tR ); // basically 2piRR
220 
221     G4double ratio= sigma/nucleusSquare;
222     fTotalXsc     = nucleusSquare*G4Log( 1. + ratio )*cB;
223     fInelasticXsc = nucleusSquare*G4Log( 1. + cofInelastic*ratio )*cB/cofInelastic;
224     fElasticXsc   = std::max(fTotalXsc - fInelasticXsc, 0.0);
225 
226     G4double difratio = ratio/(1.+ratio);
227     fDiffractionXsc = 0.5*nucleusSquare*( difratio - G4Log( 1. + difratio ) );
228 
229     G4double xratio= ((pZ*Z+pN*tN)*ppInXsc + (pZ*tN+pN*Z)*npInXsc)/nucleusSquare;
230     fProductionXsc = nucleusSquare*G4Log( 1. + cofInelastic*xratio)*cB/cofInelastic;
231     fProductionXsc = std::min(fProductionXsc, fInelasticXsc);
232   }
233   else
234   {
235     fInelasticXsc = fTotalXsc = fElasticXsc = fProductionXsc = fDiffractionXsc = 0.;
236   }
237 }
238 
239 ///////////////////////////////////////////////////////////////////////////////
240 
241 G4double G4ComponentGGNuclNuclXsc::ComputeCoulombBarier(
242                      const G4ParticleDefinition* aParticle,
243          G4double pTkin, G4int Z, G4int A,
244          G4double pR, G4double tR)
245 {
246   G4int pZ = aParticle->GetPDGCharge()*inve;
247   G4double pM = aParticle->GetPDGMass();
248   G4double tM = G4NucleiProperties::GetNuclearMass(A, Z); 
249   G4double pElab = pTkin + pM;
250   G4double totEcm = std::sqrt(pM*pM + tM*tM + 2.*pElab*tM);
251   G4double totTcm = totEcm - pM -tM;
252 
253   // 0.5 defines shape of Cross section correction
254   // at cB = totTcm it become zero
255   static const G4double qfact = 0.5*CLHEP::elm_coupling; 
256   G4double bC = qfact*pZ*Z/(pR + tR);
257 
258   G4double ratio = (totTcm <= bC) ? 0. : 1. - bC/totTcm;
259 
260 #ifdef G4VERBOSE
261   if (GetVerboseLevel() > 1) { 
262     G4cout << "G4ComponentGGNuclNuclXsc::ComputeCoulombBarier(..)=" <<ratio
263      << "; pTkin(GeV)=" << pTkin/CLHEP::MeV
264      << " totTcm= " << totTcm/CLHEP::MeV<< "; bC=" << bC/CLHEP::MeV
265      << G4endl;
266   }
267 #endif
268   return ratio;
269 }
270 
271 //////////////////////////////////////////////////////////////////////////
272 //
273 // Return single-diffraction/inelastic cross-section ratio
274 
275 G4double G4ComponentGGNuclNuclXsc::GetRatioSD(
276          const G4DynamicParticle* aParticle, G4double tA, G4double tZ)
277 {
278   ComputeCrossSections(aParticle->GetDefinition(), 
279                        aParticle->GetKineticEnergy(), 
280            G4lrint(tZ), G4lrint(tA));
281 
282   return (fInelasticXsc > 0.0) ? fDiffractionXsc/fInelasticXsc : 0.0;
283 }
284 
285 //////////////////////////////////////////////////////////////////////////
286 //
287 // Return quasi-elastic/inelastic cross-section ratio
288 
289 G4double G4ComponentGGNuclNuclXsc::GetRatioQE(
290          const G4DynamicParticle* aParticle, G4double tA, G4double tZ)
291 {
292   ComputeCrossSections(aParticle->GetDefinition(), 
293                        aParticle->GetKineticEnergy(), 
294            G4lrint(tZ), G4lrint(tA));
295 
296   return (fInelasticXsc > 0.0) ? 1.0 - fProductionXsc/fInelasticXsc : 0.0;
297 }
298 
299 ///////////////////////////////////////////////////////////////////////////////
300