source: trunk/source/processes/hadronic/cross_sections/src/G4GlauberGribovCrossSection.cc@ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 15 years ago

update ti head

File size: 52.8 KB
Line 
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// author: V. Grichine
27//
28// 17.07.06 V. Grichine - first implementation
29// 22.01.07 V.Ivanchenko - add interface with Z and A
30// 05.03.07 V.Ivanchenko - add IfZAApplicable
31// 11.06.10 V. Grichine - update for antiprotons
32
33#include "G4GlauberGribovCrossSection.hh"
34
35#include "G4ParticleTable.hh"
36#include "G4IonTable.hh"
37#include "G4ParticleDefinition.hh"
38
39//////////////////////////////////////////////////////////////////////////////////////
40//
41//
42
43const G4double G4GlauberGribovCrossSection::fNeutronBarCorrectionTot[93] = {
44
451.0, 1.0, 1.118517e+00, 1.082002e+00, 1.116171e+00, 1.078747e+00, 1.061315e+00,
461.058205e+00, 1.082663e+00, 1.068500e+00, 1.076912e+00, 1.083475e+00, 1.079117e+00,
471.071856e+00, 1.071990e+00, 1.073774e+00, 1.079356e+00, 1.081314e+00, 1.082056e+00,
481.090772e+00, 1.096776e+00, 1.095828e+00, 1.097678e+00, 1.099157e+00, 1.103677e+00,
491.105132e+00, 1.109806e+00, 1.110816e+00, 1.117378e+00, 1.115165e+00, 1.115710e+00,
501.111855e+00, 1.110482e+00, 1.110112e+00, 1.106676e+00, 1.108706e+00, 1.105549e+00,
511.106318e+00, 1.106242e+00, 1.107672e+00, 1.107342e+00, 1.108119e+00, 1.106655e+00,
521.102588e+00, 1.096657e+00, 1.092920e+00, 1.086629e+00, 1.083592e+00, 1.076030e+00,
531.083777e+00, 1.089460e+00, 1.086545e+00, 1.079924e+00, 1.082218e+00, 1.077798e+00,
541.077062e+00, 1.072825e+00, 1.072241e+00, 1.072104e+00, 1.072490e+00, 1.069829e+00,
551.070398e+00, 1.065458e+00, 1.064968e+00, 1.060524e+00, 1.060048e+00, 1.057620e+00,
561.056428e+00, 1.055366e+00, 1.055017e+00, 1.052304e+00, 1.051767e+00, 1.049728e+00,
571.048745e+00, 1.047399e+00, 1.045876e+00, 1.042972e+00, 1.041824e+00, 1.039993e+00,
581.039021e+00, 1.036627e+00, 1.034176e+00, 1.032526e+00, 1.033633e+00, 1.036107e+00,
591.037803e+00, 1.031266e+00, 1.032991e+00, 1.033284e+00, 1.035015e+00, 1.033945e+00,
601.037075e+00, 1.034721e+00
61
62};
63
64const G4double G4GlauberGribovCrossSection::fNeutronBarCorrectionIn[93] = {
65
661.0, 1.0, 1.167421e+00, 1.156250e+00, 1.205364e+00, 1.154225e+00, 1.120391e+00,
671.124632e+00, 1.129460e+00, 1.107863e+00, 1.102152e+00, 1.104593e+00, 1.100285e+00,
681.098450e+00, 1.092677e+00, 1.101124e+00, 1.106461e+00, 1.115049e+00, 1.123903e+00,
691.126661e+00, 1.131259e+00, 1.133949e+00, 1.134185e+00, 1.133767e+00, 1.132813e+00,
701.131515e+00, 1.130338e+00, 1.134171e+00, 1.139206e+00, 1.141474e+00, 1.142189e+00,
711.140725e+00, 1.140100e+00, 1.139848e+00, 1.137674e+00, 1.138645e+00, 1.136339e+00,
721.136439e+00, 1.135946e+00, 1.136431e+00, 1.135702e+00, 1.135703e+00, 1.134113e+00,
731.131935e+00, 1.128381e+00, 1.126373e+00, 1.122453e+00, 1.120908e+00, 1.115953e+00,
741.115947e+00, 1.114426e+00, 1.111749e+00, 1.106207e+00, 1.107494e+00, 1.103622e+00,
751.102576e+00, 1.098816e+00, 1.097889e+00, 1.097306e+00, 1.097130e+00, 1.094578e+00,
761.094552e+00, 1.090222e+00, 1.089358e+00, 1.085409e+00, 1.084560e+00, 1.082182e+00,
771.080773e+00, 1.079464e+00, 1.078724e+00, 1.076121e+00, 1.075235e+00, 1.073159e+00,
781.071920e+00, 1.070395e+00, 1.069503e+00, 1.067525e+00, 1.066919e+00, 1.065779e+00,
791.065319e+00, 1.063730e+00, 1.062092e+00, 1.061085e+00, 1.059908e+00, 1.059815e+00,
801.059109e+00, 1.051920e+00, 1.051258e+00, 1.049473e+00, 1.048823e+00, 1.045984e+00,
811.046435e+00, 1.042614e+00
82
83};
84
85const G4double G4GlauberGribovCrossSection::fProtonBarCorrectionTot[93] = {
86
871.0, 1.0,
881.118515e+00, 1.082000e+00, 1.116169e+00, 1.078745e+00, 1.061313e+00, 1.058203e+00,
891.082661e+00, 1.068498e+00, 1.076910e+00, 1.083474e+00, 1.079115e+00, 1.071854e+00,
901.071988e+00, 1.073772e+00, 1.079355e+00, 1.081312e+00, 1.082054e+00, 1.090770e+00,
911.096774e+00, 1.095827e+00, 1.097677e+00, 1.099156e+00, 1.103676e+00, 1.105130e+00,
921.109805e+00, 1.110814e+00, 1.117377e+00, 1.115163e+00, 1.115708e+00, 1.111853e+00,
931.110480e+00, 1.110111e+00, 1.106674e+00, 1.108705e+00, 1.105548e+00, 1.106317e+00,
941.106241e+00, 1.107671e+00, 1.107341e+00, 1.108118e+00, 1.106654e+00, 1.102586e+00,
951.096655e+00, 1.092918e+00, 1.086628e+00, 1.083590e+00, 1.076028e+00, 1.083776e+00,
961.089458e+00, 1.086543e+00, 1.079923e+00, 1.082216e+00, 1.077797e+00, 1.077061e+00,
971.072824e+00, 1.072239e+00, 1.072103e+00, 1.072488e+00, 1.069828e+00, 1.070396e+00,
981.065456e+00, 1.064966e+00, 1.060523e+00, 1.060047e+00, 1.057618e+00, 1.056427e+00,
991.055365e+00, 1.055016e+00, 1.052303e+00, 1.051766e+00, 1.049727e+00, 1.048743e+00,
1001.047397e+00, 1.045875e+00, 1.042971e+00, 1.041823e+00, 1.039992e+00, 1.039019e+00,
1011.036626e+00, 1.034175e+00, 1.032525e+00, 1.033632e+00, 1.036106e+00, 1.037802e+00,
1021.031265e+00, 1.032990e+00, 1.033283e+00, 1.035014e+00, 1.033944e+00, 1.037074e+00,
1031.034720e+00
104
105};
106
107const G4double G4GlauberGribovCrossSection::fProtonBarCorrectionIn[93] = {
108
1091.0, 1.0,
1101.167419e+00, 1.156248e+00, 1.205362e+00, 1.154224e+00, 1.120390e+00, 1.124630e+00,
1111.129459e+00, 1.107861e+00, 1.102151e+00, 1.104591e+00, 1.100284e+00, 1.098449e+00,
1121.092675e+00, 1.101122e+00, 1.106460e+00, 1.115048e+00, 1.123902e+00, 1.126659e+00,
1131.131258e+00, 1.133948e+00, 1.134183e+00, 1.133766e+00, 1.132812e+00, 1.131514e+00,
1141.130337e+00, 1.134170e+00, 1.139205e+00, 1.141472e+00, 1.142188e+00, 1.140724e+00,
1151.140099e+00, 1.139847e+00, 1.137672e+00, 1.138644e+00, 1.136338e+00, 1.136438e+00,
1161.135945e+00, 1.136429e+00, 1.135701e+00, 1.135702e+00, 1.134112e+00, 1.131934e+00,
1171.128380e+00, 1.126371e+00, 1.122452e+00, 1.120907e+00, 1.115952e+00, 1.115946e+00,
1181.114425e+00, 1.111748e+00, 1.106205e+00, 1.107493e+00, 1.103621e+00, 1.102575e+00,
1191.098815e+00, 1.097888e+00, 1.097305e+00, 1.097129e+00, 1.094577e+00, 1.094551e+00,
1201.090221e+00, 1.089357e+00, 1.085408e+00, 1.084559e+00, 1.082181e+00, 1.080772e+00,
1211.079463e+00, 1.078723e+00, 1.076120e+00, 1.075234e+00, 1.073158e+00, 1.071919e+00,
1221.070394e+00, 1.069502e+00, 1.067524e+00, 1.066918e+00, 1.065778e+00, 1.065318e+00,
1231.063729e+00, 1.062091e+00, 1.061084e+00, 1.059907e+00, 1.059814e+00, 1.059108e+00,
1241.051919e+00, 1.051257e+00, 1.049472e+00, 1.048822e+00, 1.045983e+00, 1.046434e+00,
1251.042613e+00
126
127};
128
129
130const G4double G4GlauberGribovCrossSection::fPionPlusBarCorrectionTot[93] = {
131
1321.0, 1.0,
1331.075927e+00, 1.074407e+00, 1.126098e+00, 1.100127e+00, 1.089742e+00, 1.083536e+00,
1341.089988e+00, 1.103566e+00, 1.096922e+00, 1.126573e+00, 1.132734e+00, 1.136512e+00,
1351.136629e+00, 1.133086e+00, 1.132428e+00, 1.129299e+00, 1.125622e+00, 1.126992e+00,
1361.127840e+00, 1.162670e+00, 1.160392e+00, 1.157864e+00, 1.157227e+00, 1.154627e+00,
1371.192555e+00, 1.197243e+00, 1.197911e+00, 1.200326e+00, 1.220053e+00, 1.215019e+00,
1381.211703e+00, 1.209080e+00, 1.204248e+00, 1.203328e+00, 1.198671e+00, 1.196840e+00,
1391.194392e+00, 1.193037e+00, 1.190408e+00, 1.188583e+00, 1.206127e+00, 1.210028e+00,
1401.206434e+00, 1.204456e+00, 1.200547e+00, 1.199058e+00, 1.200174e+00, 1.200276e+00,
1411.198912e+00, 1.213048e+00, 1.207160e+00, 1.208020e+00, 1.203814e+00, 1.202380e+00,
1421.198306e+00, 1.197002e+00, 1.196027e+00, 1.195449e+00, 1.192563e+00, 1.192135e+00,
1431.187556e+00, 1.186308e+00, 1.182124e+00, 1.180900e+00, 1.178224e+00, 1.176471e+00,
1441.174811e+00, 1.173702e+00, 1.170827e+00, 1.169581e+00, 1.167205e+00, 1.165626e+00,
1451.180244e+00, 1.177626e+00, 1.175121e+00, 1.173903e+00, 1.172192e+00, 1.171128e+00,
1461.168997e+00, 1.166826e+00, 1.164130e+00, 1.165412e+00, 1.165504e+00, 1.165020e+00,
1471.158462e+00, 1.158014e+00, 1.156519e+00, 1.156081e+00, 1.153602e+00, 1.154190e+00,
1481.152974e+00
149
150};
151
152const G4double G4GlauberGribovCrossSection::fPionPlusBarCorrectionIn[93] = {
153
1541.0, 1.0,
1551.140246e+00, 1.097872e+00, 1.104301e+00, 1.068722e+00, 1.044495e+00, 1.062622e+00,
1561.047987e+00, 1.037032e+00, 1.035686e+00, 1.042870e+00, 1.052222e+00, 1.065100e+00,
1571.070480e+00, 1.078286e+00, 1.081488e+00, 1.089713e+00, 1.099105e+00, 1.098003e+00,
1581.102175e+00, 1.117707e+00, 1.121734e+00, 1.125229e+00, 1.126457e+00, 1.128905e+00,
1591.137312e+00, 1.126263e+00, 1.126459e+00, 1.115191e+00, 1.116986e+00, 1.117184e+00,
1601.117037e+00, 1.116777e+00, 1.115858e+00, 1.115745e+00, 1.114489e+00, 1.113993e+00,
1611.113226e+00, 1.112818e+00, 1.111890e+00, 1.111238e+00, 1.111209e+00, 1.111775e+00,
1621.110256e+00, 1.109414e+00, 1.107647e+00, 1.106980e+00, 1.106096e+00, 1.107331e+00,
1631.107849e+00, 1.106407e+00, 1.103426e+00, 1.103896e+00, 1.101756e+00, 1.101031e+00,
1641.098915e+00, 1.098260e+00, 1.097768e+00, 1.097487e+00, 1.095964e+00, 1.095773e+00,
1651.093348e+00, 1.092687e+00, 1.090465e+00, 1.089821e+00, 1.088394e+00, 1.087462e+00,
1661.086571e+00, 1.085997e+00, 1.084451e+00, 1.083798e+00, 1.082513e+00, 1.081670e+00,
1671.080735e+00, 1.075659e+00, 1.074341e+00, 1.073689e+00, 1.072787e+00, 1.072237e+00,
1681.071107e+00, 1.069955e+00, 1.064856e+00, 1.065873e+00, 1.065938e+00, 1.065694e+00,
1691.062192e+00, 1.061967e+00, 1.061180e+00, 1.060960e+00, 1.059646e+00, 1.059975e+00,
1701.059658e+00
171
172};
173
174
175const G4double G4GlauberGribovCrossSection::fPionMinusBarCorrectionTot[93] = {
176
1771.0, 1.0,
1781.075927e+00, 1.077959e+00, 1.129145e+00, 1.102088e+00, 1.089765e+00, 1.083542e+00,
1791.089995e+00, 1.104895e+00, 1.097154e+00, 1.127663e+00, 1.133063e+00, 1.137425e+00,
1801.136724e+00, 1.133859e+00, 1.132498e+00, 1.130276e+00, 1.127896e+00, 1.127656e+00,
1811.127905e+00, 1.164210e+00, 1.162259e+00, 1.160075e+00, 1.158978e+00, 1.156649e+00,
1821.194157e+00, 1.199177e+00, 1.198983e+00, 1.202325e+00, 1.221967e+00, 1.217548e+00,
1831.214389e+00, 1.211760e+00, 1.207335e+00, 1.206081e+00, 1.201766e+00, 1.199779e+00,
1841.197283e+00, 1.195706e+00, 1.193071e+00, 1.191115e+00, 1.208838e+00, 1.212681e+00,
1851.209235e+00, 1.207163e+00, 1.203451e+00, 1.201807e+00, 1.203283e+00, 1.203388e+00,
1861.202244e+00, 1.216509e+00, 1.211066e+00, 1.211504e+00, 1.207539e+00, 1.205991e+00,
1871.202143e+00, 1.200724e+00, 1.199595e+00, 1.198815e+00, 1.196025e+00, 1.195390e+00,
1881.191137e+00, 1.189791e+00, 1.185888e+00, 1.184575e+00, 1.181996e+00, 1.180229e+00,
1891.178545e+00, 1.177355e+00, 1.174616e+00, 1.173312e+00, 1.171016e+00, 1.169424e+00,
1901.184120e+00, 1.181478e+00, 1.179085e+00, 1.177817e+00, 1.176124e+00, 1.175003e+00,
1911.172947e+00, 1.170858e+00, 1.168170e+00, 1.169397e+00, 1.169304e+00, 1.168706e+00,
1921.162774e+00, 1.162217e+00, 1.160740e+00, 1.160196e+00, 1.157857e+00, 1.158220e+00,
1931.157267e+00
194};
195
196
197const G4double G4GlauberGribovCrossSection::fPionMinusBarCorrectionIn[93] = {
198
1991.0, 1.0,
2001.140246e+00, 1.100898e+00, 1.106773e+00, 1.070289e+00, 1.044514e+00, 1.062628e+00,
2011.047992e+00, 1.038041e+00, 1.035862e+00, 1.043679e+00, 1.052466e+00, 1.065780e+00,
2021.070551e+00, 1.078869e+00, 1.081541e+00, 1.090455e+00, 1.100847e+00, 1.098511e+00,
2031.102226e+00, 1.118865e+00, 1.123143e+00, 1.126904e+00, 1.127785e+00, 1.130444e+00,
2041.138502e+00, 1.127678e+00, 1.127244e+00, 1.116634e+00, 1.118347e+00, 1.118988e+00,
2051.118957e+00, 1.118696e+00, 1.118074e+00, 1.117722e+00, 1.116717e+00, 1.116111e+00,
2061.115311e+00, 1.114745e+00, 1.113814e+00, 1.113069e+00, 1.113141e+00, 1.113660e+00,
2071.112249e+00, 1.111343e+00, 1.109718e+00, 1.108942e+00, 1.108310e+00, 1.109549e+00,
2081.110227e+00, 1.108846e+00, 1.106183e+00, 1.106354e+00, 1.104388e+00, 1.103583e+00,
2091.101632e+00, 1.100896e+00, 1.100296e+00, 1.099873e+00, 1.098420e+00, 1.098082e+00,
2101.095892e+00, 1.095162e+00, 1.093144e+00, 1.092438e+00, 1.091083e+00, 1.090142e+00,
2111.089236e+00, 1.088604e+00, 1.087159e+00, 1.086465e+00, 1.085239e+00, 1.084388e+00,
2121.083473e+00, 1.078373e+00, 1.077136e+00, 1.076450e+00, 1.075561e+00, 1.074973e+00,
2131.073898e+00, 1.072806e+00, 1.067706e+00, 1.068684e+00, 1.068618e+00, 1.068294e+00,
2141.065241e+00, 1.064939e+00, 1.064166e+00, 1.063872e+00, 1.062659e+00, 1.062828e+00,
2151.062699e+00
216
217};
218
219
220
221
222////////////////////////////////////////////////////////////////////////////////
223//
224//
225
226G4GlauberGribovCrossSection::G4GlauberGribovCrossSection()
227: fUpperLimit( 100000 * GeV ),
228 fLowerLimit( 3 * GeV ),
229 fRadiusConst( 1.08*fermi ) // 1.1, 1.3 ?
230{
231 theGamma = G4Gamma::Gamma();
232 theProton = G4Proton::Proton();
233 theNeutron = G4Neutron::Neutron();
234 theAProton = G4AntiProton::AntiProton();
235 theANeutron = G4AntiNeutron::AntiNeutron();
236 thePiPlus = G4PionPlus::PionPlus();
237 thePiMinus = G4PionMinus::PionMinus();
238 thePiZero = G4PionZero::PionZero();
239 theKPlus = G4KaonPlus::KaonPlus();
240 theKMinus = G4KaonMinus::KaonMinus();
241 theK0S = G4KaonZeroShort::KaonZeroShort();
242 theK0L = G4KaonZeroLong::KaonZeroLong();
243 theL = G4Lambda::Lambda();
244 theAntiL = G4AntiLambda::AntiLambda();
245 theSPlus = G4SigmaPlus::SigmaPlus();
246 theASPlus = G4AntiSigmaPlus::AntiSigmaPlus();
247 theSMinus = G4SigmaMinus::SigmaMinus();
248 theASMinus = G4AntiSigmaMinus::AntiSigmaMinus();
249 theS0 = G4SigmaZero::SigmaZero();
250 theAS0 = G4AntiSigmaZero::AntiSigmaZero();
251 theXiMinus = G4XiMinus::XiMinus();
252 theXi0 = G4XiZero::XiZero();
253 theAXiMinus = G4AntiXiMinus::AntiXiMinus();
254 theAXi0 = G4AntiXiZero::AntiXiZero();
255 theOmega = G4OmegaMinus::OmegaMinus();
256 theAOmega = G4AntiOmegaMinus::AntiOmegaMinus();
257 theD = G4Deuteron::Deuteron();
258 theT = G4Triton::Triton();
259 theA = G4Alpha::Alpha();
260 theHe3 = G4He3::He3();
261}
262
263///////////////////////////////////////////////////////////////////////////////////////
264//
265//
266
267G4GlauberGribovCrossSection::~G4GlauberGribovCrossSection()
268{
269}
270
271
272////////////////////////////////////////////////////////////////////////////////////////
273//
274//
275
276
277G4bool
278G4GlauberGribovCrossSection::IsApplicable(const G4DynamicParticle* aDP,
279 const G4Element* anElement)
280{
281 return IsIsoApplicable(aDP, G4lrint(anElement->GetZ()),
282 G4lrint(anElement->GetN()));
283}
284
285////////////////////////////////////////////////////////////////////////////////////////
286//
287//
288
289G4bool
290G4GlauberGribovCrossSection::IsIsoApplicable(const G4DynamicParticle* aDP,
291 G4int Z, G4int)
292{
293 G4bool applicable = false;
294 // G4int baryonNumber = aDP->GetDefinition()->GetBaryonNumber();
295 G4double kineticEnergy = aDP->GetKineticEnergy();
296
297 const G4ParticleDefinition* theParticle = aDP->GetDefinition();
298
299 if ( ( kineticEnergy >= fLowerLimit &&
300 Z > 1 && // >= He
301 ( theParticle == theAProton ||
302 theParticle == theGamma ||
303 theParticle == theKPlus ||
304 theParticle == theKMinus ||
305 theParticle == theSMinus) ) ||
306
307 ( kineticEnergy >= fLowerLimit &&
308 Z > 1 && // >= He
309 ( theParticle == theProton ||
310 theParticle == theNeutron ||
311 theParticle == thePiPlus ||
312 theParticle == thePiMinus ) ) ) applicable = true;
313
314 return applicable;
315}
316
317////////////////////////////////////////////////////////////////////////////////////////
318//
319// Calculates total and inelastic Xsc, derives elastic as total - inelastic accordong to
320// Glauber model with Gribov correction calculated in the dipole approximation on
321// light cone. Gaussian density helps to calculate rest integrals of the model.
322// [1] B.Z. Kopeliovich, nucl-th/0306044
323
324
325G4double G4GlauberGribovCrossSection::
326GetCrossSection(const G4DynamicParticle* aParticle, const G4Element* anElement, G4double T)
327{
328 return GetZandACrossSection(aParticle, G4lrint(anElement->GetZ()),
329 G4lrint(anElement->GetN()), T);
330}
331
332////////////////////////////////////////////////////////////////////////////////////////
333//
334// Calculates total and inelastic Xsc, derives elastic as total - inelastic accordong to
335// Glauber model with Gribov correction calculated in the dipole approximation on
336// light cone. Gaussian density of point-like nucleons helps to calculate rest integrals of the model.
337// [1] B.Z. Kopeliovich, nucl-th/0306044 + simplification above
338
339
340
341G4double G4GlauberGribovCrossSection::
342GetZandACrossSection(const G4DynamicParticle* aParticle, G4int Z, G4int A, G4double)
343{
344 G4double xsection, sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
345 G4double R = GetNucleusRadius(A);
346
347 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
348
349 if( theParticle == theProton ||
350 theParticle == theNeutron ||
351 theParticle == thePiPlus ||
352 theParticle == thePiMinus )
353 {
354 sigma = GetHadronNucleonXscNS(aParticle, A, Z);
355 cofInelastic = 2.4;
356 cofTotal = 2.0;
357 }
358 else
359 {
360 sigma = GetHadronNucleonXscNS(aParticle, A, Z);
361 cofInelastic = 2.2;
362 cofTotal = 2.0;
363 }
364 // cofInelastic = 2.0;
365
366 if( A > 1 )
367 {
368 nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
369 ratio = sigma/nucleusSquare;
370
371 xsection = nucleusSquare*std::log( 1. + ratio );
372
373 xsection *= GetParticleBarCorTot(theParticle, Z);
374
375 fTotalXsc = xsection;
376
377
378
379 fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
380
381 fInelasticXsc *= GetParticleBarCorIn(theParticle, Z);
382
383 fElasticXsc = fTotalXsc - fInelasticXsc;
384
385
386 G4double difratio = ratio/(1.+ratio);
387
388 fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
389
390
391 sigma = GetHNinelasticXsc(aParticle, A, Z);
392 ratio = sigma/nucleusSquare;
393
394 fProductionXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
395
396 if (fElasticXsc < 0.) fElasticXsc = 0.;
397 }
398 else // H
399 {
400 fTotalXsc = sigma;
401 xsection = sigma;
402
403 if ( theParticle != theAProton )
404 {
405 sigma = GetHNinelasticXsc(aParticle, A, Z);
406 fInelasticXsc = sigma;
407 fElasticXsc = fTotalXsc - fInelasticXsc;
408 }
409 else
410 {
411 fElasticXsc = fTotalXsc - fInelasticXsc;
412 }
413 if (fElasticXsc < 0.) fElasticXsc = 0.;
414
415 }
416 return xsection;
417}
418
419//////////////////////////////////////////////////////////////////////////
420//
421// Return single-diffraction/inelastic cross-section ratio
422
423G4double G4GlauberGribovCrossSection::
424GetRatioSD(const G4DynamicParticle* aParticle, G4int A, G4int Z)
425{
426 G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
427 G4double R = GetNucleusRadius(A);
428
429 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
430
431 if( theParticle == theProton ||
432 theParticle == theNeutron ||
433 theParticle == thePiPlus ||
434 theParticle == thePiMinus )
435 {
436 sigma = GetHadronNucleonXscNS(aParticle, A, Z);
437 cofInelastic = 2.4;
438 cofTotal = 2.0;
439 }
440 else
441 {
442 sigma = GetHadronNucleonXscNS(aParticle, A, Z);
443 cofInelastic = 2.2;
444 cofTotal = 2.0;
445 }
446 nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
447 ratio = sigma/nucleusSquare;
448
449 fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
450
451 G4double difratio = ratio/(1.+ratio);
452
453 fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
454
455 if (fInelasticXsc > 0.) ratio = fDiffractionXsc/fInelasticXsc;
456 else ratio = 0.;
457
458 return ratio;
459}
460
461//////////////////////////////////////////////////////////////////////////
462//
463// Return suasi-elastic/inelastic cross-section ratio
464
465G4double G4GlauberGribovCrossSection::
466GetRatioQE(const G4DynamicParticle* aParticle, G4int A, G4int Z)
467{
468 G4double sigma, cofInelastic, cofTotal, nucleusSquare, ratio;
469 G4double R = GetNucleusRadius(A);
470
471 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
472
473 if( theParticle == theProton ||
474 theParticle == theNeutron ||
475 theParticle == thePiPlus ||
476 theParticle == thePiMinus )
477 {
478 sigma = GetHadronNucleonXscNS(aParticle, A, Z);
479 cofInelastic = 2.4;
480 cofTotal = 2.0;
481 }
482 else
483 {
484 sigma = GetHadronNucleonXscNS(aParticle, A, Z);
485 cofInelastic = 2.2;
486 cofTotal = 2.0;
487 }
488 nucleusSquare = cofTotal*pi*R*R; // basically 2piRR
489 ratio = sigma/nucleusSquare;
490
491 fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
492
493 sigma = GetHNinelasticXsc(aParticle, A, Z);
494 ratio = sigma/nucleusSquare;
495
496 fProductionXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
497
498 if (fInelasticXsc > fProductionXsc) ratio = (fInelasticXsc-fProductionXsc)/fInelasticXsc;
499 else ratio = 0.;
500 if ( ratio < 0. ) ratio = 0.;
501
502 return ratio;
503}
504
505/////////////////////////////////////////////////////////////////////////////////////
506//
507// Returns hadron-nucleon Xsc according to differnt parametrisations:
508// [2] E. Levin, hep-ph/9710546
509// [3] U. Dersch, et al, hep-ex/9910052
510// [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
511
512G4double
513G4GlauberGribovCrossSection::GetHadronNucleonXsc(const G4DynamicParticle* aParticle,
514 const G4Element* anElement)
515{
516 G4int At = G4lrint(anElement->GetN()); // number of nucleons
517 G4int Zt = G4lrint(anElement->GetZ()); // number of protons
518
519 return GetHadronNucleonXsc(aParticle, At, Zt);
520}
521
522
523
524
525/////////////////////////////////////////////////////////////////////////////////////
526//
527// Returns hadron-nucleon Xsc according to differnt parametrisations:
528// [2] E. Levin, hep-ph/9710546
529// [3] U. Dersch, et al, hep-ex/9910052
530// [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
531
532G4double
533G4GlauberGribovCrossSection::GetHadronNucleonXsc(const G4DynamicParticle* aParticle,
534 G4int At, G4int Zt)
535{
536 G4double xsection;
537
538 G4double targ_mass = G4ParticleTable::GetParticleTable()->
539 GetIonTable()->GetIonMass(Zt, At);
540// GetIonTable()->GetIonMass( G4int(Zt+0.5) , G4int(At+0.5) );
541
542 targ_mass = 0.939*GeV; // ~mean neutron and proton ???
543
544 G4double proj_mass = aParticle->GetMass();
545 G4double proj_momentum = aParticle->GetMomentum().mag();
546 G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
547
548 sMand /= GeV*GeV; // in GeV for parametrisation
549 proj_momentum /= GeV;
550
551 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
552
553 G4double aa = At;
554
555 if(theParticle == theGamma)
556 {
557 xsection = aa*(0.0677*std::pow(sMand,0.0808) + 0.129*std::pow(sMand,-0.4525));
558 }
559 else if(theParticle == theNeutron) // as proton ???
560 {
561 xsection = aa*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
562 }
563 else if(theParticle == theProton)
564 {
565 xsection = aa*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
566 // xsection = At*( 49.51*std::pow(sMand,-0.097) + 0.314*std::log(sMand)*std::log(sMand) );
567 // xsection = At*( 38.4 + 0.85*std::abs(std::pow(log(sMand),1.47)) );
568 }
569 else if(theParticle == theAProton)
570 {
571 xsection = aa*( 21.70*std::pow(sMand,0.0808) + 98.39*std::pow(sMand,-0.4525));
572 }
573 else if(theParticle == thePiPlus)
574 {
575 xsection = aa*(13.63*std::pow(sMand,0.0808) + 27.56*std::pow(sMand,-0.4525));
576 }
577 else if(theParticle == thePiMinus)
578 {
579 // xsection = At*( 55.2*std::pow(sMand,-0.255) + 0.346*std::log(sMand)*std::log(sMand) );
580 xsection = aa*(13.63*std::pow(sMand,0.0808) + 36.02*std::pow(sMand,-0.4525));
581 }
582 else if(theParticle == theKPlus)
583 {
584 xsection = aa*(11.82*std::pow(sMand,0.0808) + 8.15*std::pow(sMand,-0.4525));
585 }
586 else if(theParticle == theKMinus)
587 {
588 xsection = aa*(11.82*std::pow(sMand,0.0808) + 26.36*std::pow(sMand,-0.4525));
589 }
590 else // as proton ???
591 {
592 xsection = aa*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
593 }
594 xsection *= millibarn;
595 return xsection;
596}
597
598
599/////////////////////////////////////////////////////////////////////////////////////
600//
601// Returns hadron-nucleon Xsc according to PDG parametrisation (2005):
602// http://pdg.lbl.gov/2006/reviews/hadronicrpp.pdf
603
604G4double
605G4GlauberGribovCrossSection::GetHadronNucleonXscPDG(const G4DynamicParticle* aParticle,
606 const G4Element* anElement)
607{
608 G4int At = G4lrint(anElement->GetN()); // number of nucleons
609 G4int Zt = G4lrint(anElement->GetZ()); // number of protons
610
611 return GetHadronNucleonXscPDG(aParticle, At, Zt);
612}
613
614
615
616
617/////////////////////////////////////////////////////////////////////////////////////
618//
619// Returns hadron-nucleon Xsc according to PDG parametrisation (2005):
620// http://pdg.lbl.gov/2006/reviews/hadronicrpp.pdf
621// At = number of nucleons, Zt = number of protons
622
623G4double
624G4GlauberGribovCrossSection::GetHadronNucleonXscPDG(const G4DynamicParticle* aParticle,
625 G4int At, G4int Zt)
626{
627 G4double xsection;
628
629 G4int Nt = At-Zt; // number of neutrons
630 if (Nt < 0) Nt = 0;
631
632 G4double zz = Zt;
633 G4double aa = At;
634 G4double nn = Nt;
635
636 G4double targ_mass = G4ParticleTable::GetParticleTable()->
637 GetIonTable()->GetIonMass(Zt, At);
638
639 targ_mass = 0.939*GeV; // ~mean neutron and proton ???
640
641 G4double proj_mass = aParticle->GetMass();
642 G4double proj_momentum = aParticle->GetMomentum().mag();
643
644 G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
645
646 sMand /= GeV*GeV; // in GeV for parametrisation
647
648 // General PDG fit constants
649
650 G4double s0 = 5.38*5.38; // in Gev^2
651 G4double eta1 = 0.458;
652 G4double eta2 = 0.458;
653 G4double B = 0.308;
654
655
656 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
657
658
659 if(theParticle == theNeutron) // proton-neutron fit
660 {
661 xsection = zz*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
662 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
663 xsection += nn*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
664 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2)); // pp for nn
665 }
666 else if(theParticle == theProton)
667 {
668
669 xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
670 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
671
672 xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
673 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
674 }
675 else if(theParticle == theAProton)
676 {
677 xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
678 + 42.53*std::pow(sMand,-eta1) + 33.34*std::pow(sMand,-eta2));
679
680 xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
681 + 40.15*std::pow(sMand,-eta1) + 30.*std::pow(sMand,-eta2));
682 }
683 else if(theParticle == thePiPlus)
684 {
685 xsection = aa*( 20.86 + B*std::pow(std::log(sMand/s0),2.)
686 + 19.24*std::pow(sMand,-eta1) - 6.03*std::pow(sMand,-eta2));
687 }
688 else if(theParticle == thePiMinus)
689 {
690 xsection = aa*( 20.86 + B*std::pow(std::log(sMand/s0),2.)
691 + 19.24*std::pow(sMand,-eta1) + 6.03*std::pow(sMand,-eta2));
692 }
693 else if(theParticle == theKPlus)
694 {
695 xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
696 + 7.14*std::pow(sMand,-eta1) - 13.45*std::pow(sMand,-eta2));
697
698 xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
699 + 5.17*std::pow(sMand,-eta1) - 7.23*std::pow(sMand,-eta2));
700 }
701 else if(theParticle == theKMinus)
702 {
703 xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
704 + 7.14*std::pow(sMand,-eta1) + 13.45*std::pow(sMand,-eta2));
705
706 xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
707 + 5.17*std::pow(sMand,-eta1) + 7.23*std::pow(sMand,-eta2));
708 }
709 else if(theParticle == theSMinus)
710 {
711 xsection = aa*( 35.20 + B*std::pow(std::log(sMand/s0),2.)
712 - 199.*std::pow(sMand,-eta1) + 264.*std::pow(sMand,-eta2));
713 }
714 else if(theParticle == theGamma) // modify later on
715 {
716 xsection = aa*( 0.0 + B*std::pow(std::log(sMand/s0),2.)
717 + 0.032*std::pow(sMand,-eta1) - 0.0*std::pow(sMand,-eta2));
718
719 }
720 else // as proton ???
721 {
722 xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
723 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
724
725 xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
726 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
727 }
728 xsection *= millibarn; // parametrised in mb
729 return xsection;
730}
731
732
733/////////////////////////////////////////////////////////////////////////////////////
734//
735// Returns hadron-nucleon cross-section based on N. Starkov parametrisation of
736// data from mainly http://wwwppds.ihep.su:8001/c5-6A.html database
737
738G4double
739G4GlauberGribovCrossSection::GetHadronNucleonXscNS(const G4DynamicParticle* aParticle,
740 const G4Element* anElement)
741{
742 G4int At = G4lrint(anElement->GetN()); // number of nucleons
743 G4int Zt = G4lrint(anElement->GetZ()); // number of protons
744
745 return GetHadronNucleonXscNS(aParticle, At, Zt);
746}
747
748
749
750
751/////////////////////////////////////////////////////////////////////////////////////
752//
753// Returns hadron-nucleon cross-section based on N. Starkov parametrisation of
754// data from mainly http://wwwppds.ihep.su:8001/c5-6A.html database
755
756G4double
757G4GlauberGribovCrossSection::GetHadronNucleonXscNS(const G4DynamicParticle* aParticle,
758 G4int At, G4int Zt)
759{
760 G4double xsection(0), Delta, A0, B0;
761 G4double hpXsc(0);
762 G4double hnXsc(0);
763
764 G4int Nt = At-Zt; // number of neutrons
765 if (Nt < 0) Nt = 0;
766
767 G4double aa = At;
768 G4double zz = Zt;
769 G4double nn = Nt;
770
771 G4double targ_mass = G4ParticleTable::GetParticleTable()->
772 GetIonTable()->GetIonMass(Zt, At);
773
774 targ_mass = 0.939*GeV; // ~mean neutron and proton ???
775
776 G4double proj_mass = aParticle->GetMass();
777 G4double proj_energy = aParticle->GetTotalEnergy();
778 G4double proj_momentum = aParticle->GetMomentum().mag();
779
780 G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
781
782 sMand /= GeV*GeV; // in GeV for parametrisation
783 proj_momentum /= GeV;
784 proj_energy /= GeV;
785 proj_mass /= GeV;
786
787 // General PDG fit constants
788
789 G4double s0 = 5.38*5.38; // in Gev^2
790 G4double eta1 = 0.458;
791 G4double eta2 = 0.458;
792 G4double B = 0.308;
793
794
795 const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
796
797
798 if(theParticle == theNeutron)
799 {
800 if( proj_momentum >= 10.)
801 // if( proj_momentum >= 2.)
802 {
803 Delta = 1.;
804
805 if( proj_energy < 40. ) Delta = 0.916+0.0021*proj_energy;
806
807 if(proj_momentum >= 10.)
808 {
809 B0 = 7.5;
810 A0 = 100. - B0*std::log(3.0e7);
811
812 xsection = A0 + B0*std::log(proj_energy) - 11
813 + 103*std::pow(2*0.93827*proj_energy + proj_mass*proj_mass+
814 0.93827*0.93827,-0.165); // mb
815 }
816 xsection *= zz + nn;
817 }
818 else
819 {
820 // nn to be pp
821
822 if( proj_momentum < 0.73 )
823 {
824 hnXsc = 23 + 50*( std::pow( std::log(0.73/proj_momentum), 3.5 ) );
825 }
826 else if( proj_momentum < 1.05 )
827 {
828 hnXsc = 23 + 40*(std::log(proj_momentum/0.73))*
829 (std::log(proj_momentum/0.73));
830 }
831 else // if( proj_momentum < 10. )
832 {
833 hnXsc = 39.0+
834 75*(proj_momentum - 1.2)/(std::pow(proj_momentum,3.0) + 0.15);
835 }
836 // pn to be np
837
838 if( proj_momentum < 0.8 )
839 {
840 hpXsc = 33+30*std::pow(std::log(proj_momentum/1.3),4.0);
841 }
842 else if( proj_momentum < 1.4 )
843 {
844 hpXsc = 33+30*std::pow(std::log(proj_momentum/0.95),2.0);
845 }
846 else // if( proj_momentum < 10. )
847 {
848 hpXsc = 33.3+
849 20.8*(std::pow(proj_momentum,2.0)-1.35)/
850 (std::pow(proj_momentum,2.50)+0.95);
851 }
852 xsection = hpXsc*zz + hnXsc*nn;
853 }
854 }
855 else if(theParticle == theProton)
856 {
857 if( proj_momentum >= 10.)
858 // if( proj_momentum >= 2.)
859 {
860 Delta = 1.;
861
862 if( proj_energy < 40. ) Delta = 0.916+0.0021*proj_energy;
863
864 if(proj_momentum >= 10.)
865 {
866 B0 = 7.5;
867 A0 = 100. - B0*std::log(3.0e7);
868
869 xsection = A0 + B0*std::log(proj_energy) - 11
870 + 103*std::pow(2*0.93827*proj_energy + proj_mass*proj_mass+
871 0.93827*0.93827,-0.165); // mb
872 }
873 xsection *= zz + nn;
874 }
875 else
876 {
877 // pp
878
879 if( proj_momentum < 0.73 )
880 {
881 hpXsc = 23 + 50*( std::pow( std::log(0.73/proj_momentum), 3.5 ) );
882 }
883 else if( proj_momentum < 1.05 )
884 {
885 hpXsc = 23 + 40*(std::log(proj_momentum/0.73))*
886 (std::log(proj_momentum/0.73));
887 }
888 else // if( proj_momentum < 10. )
889 {
890 hpXsc = 39.0+
891 75*(proj_momentum - 1.2)/(std::pow(proj_momentum,3.0) + 0.15);
892 }
893 // pn to be np
894
895 if( proj_momentum < 0.8 )
896 {
897 hnXsc = 33+30*std::pow(std::log(proj_momentum/1.3),4.0);
898 }
899 else if( proj_momentum < 1.4 )
900 {
901 hnXsc = 33+30*std::pow(std::log(proj_momentum/0.95),2.0);
902 }
903 else // if( proj_momentum < 10. )
904 {
905 hnXsc = 33.3+
906 20.8*(std::pow(proj_momentum,2.0)-1.35)/
907 (std::pow(proj_momentum,2.50)+0.95);
908 }
909 xsection = hpXsc*zz + hnXsc*nn;
910 // xsection = hpXsc*(Zt + Nt);
911 // xsection = hnXsc*(Zt + Nt);
912 }
913 // xsection *= 0.95;
914 }
915 else if( theParticle == theAProton )
916 {
917 // xsection = Zt*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
918 // + 42.53*std::pow(sMand,-eta1) + 33.34*std::pow(sMand,-eta2));
919
920 // xsection += Nt*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
921 // + 40.15*std::pow(sMand,-eta1) + 30.*std::pow(sMand,-eta2));
922
923 G4double logP = std::log(proj_momentum);
924
925 if( proj_momentum <= 1.0 )
926 {
927 xsection = zz*(65.55 + 53.84/(proj_momentum+1.e-6) );
928 }
929 else
930 {
931 xsection = zz*( 41.1 + 77.2*std::pow( proj_momentum, -0.68)
932 + 0.293*logP*logP - 1.82*logP );
933 }
934 if ( nn > 0.)
935 {
936 xsection += nn*( 41.9 + 96.2*std::pow( proj_momentum, -0.99) - 0.154*logP);
937 }
938 else // H
939 {
940 fInelasticXsc = 38.0 + 38.0*std::pow( proj_momentum, -0.96)
941 - 0.169*logP*logP;
942 fInelasticXsc *= millibarn;
943 }
944 }
945 else if( theParticle == thePiPlus )
946 {
947 if(proj_momentum < 0.4)
948 {
949 G4double Ex3 = 180*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.085/0.085);
950 hpXsc = Ex3+20.0;
951 }
952 else if( proj_momentum < 1.15 )
953 {
954 G4double Ex4 = 88*(std::log(proj_momentum/0.75))*(std::log(proj_momentum/0.75));
955 hpXsc = Ex4+14.0;
956 }
957 else if(proj_momentum < 3.5)
958 {
959 G4double Ex1 = 3.2*std::exp(-(proj_momentum-2.55)*(proj_momentum-2.55)/0.55/0.55);
960 G4double Ex2 = 12*std::exp(-(proj_momentum-1.47)*(proj_momentum-1.47)/0.225/0.225);
961 hpXsc = Ex1+Ex2+27.5;
962 }
963 else // if(proj_momentum > 3.5) // mb
964 {
965 hpXsc = 10.6+2.*std::log(proj_energy)+25*std::pow(proj_energy,-0.43);
966 }
967 // pi+n = pi-p??
968
969 if(proj_momentum < 0.37)
970 {
971 hnXsc = 28.0 + 40*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.07/0.07);
972 }
973 else if(proj_momentum<0.65)
974 {
975 hnXsc = 26+110*(std::log(proj_momentum/0.48))*(std::log(proj_momentum/0.48));
976 }
977 else if(proj_momentum<1.3)
978 {
979 hnXsc = 36.1+
980 10*std::exp(-(proj_momentum-0.72)*(proj_momentum-0.72)/0.06/0.06)+
981 24*std::exp(-(proj_momentum-1.015)*(proj_momentum-1.015)/0.075/0.075);
982 }
983 else if(proj_momentum<3.0)
984 {
985 hnXsc = 36.1+0.079-4.313*std::log(proj_momentum)+
986 3*std::exp(-(proj_momentum-2.1)*(proj_momentum-2.1)/0.4/0.4)+
987 1.5*std::exp(-(proj_momentum-1.4)*(proj_momentum-1.4)/0.12/0.12);
988 }
989 else // mb
990 {
991 hnXsc = 10.6+2*std::log(proj_energy)+30*std::pow(proj_energy,-0.43);
992 }
993 xsection = hpXsc*zz + hnXsc*nn;
994 }
995 else if(theParticle == thePiMinus)
996 {
997 // pi-n = pi+p??
998
999 if(proj_momentum < 0.4)
1000 {
1001 G4double Ex3 = 180*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.085/0.085);
1002 hnXsc = Ex3+20.0;
1003 }
1004 else if(proj_momentum < 1.15)
1005 {
1006 G4double Ex4 = 88*(std::log(proj_momentum/0.75))*(std::log(proj_momentum/0.75));
1007 hnXsc = Ex4+14.0;
1008 }
1009 else if(proj_momentum < 3.5)
1010 {
1011 G4double Ex1 = 3.2*std::exp(-(proj_momentum-2.55)*(proj_momentum-2.55)/0.55/0.55);
1012 G4double Ex2 = 12*std::exp(-(proj_momentum-1.47)*(proj_momentum-1.47)/0.225/0.225);
1013 hnXsc = Ex1+Ex2+27.5;
1014 }
1015 else // if(proj_momentum > 3.5) // mb
1016 {
1017 hnXsc = 10.6+2.*std::log(proj_energy)+25*std::pow(proj_energy,-0.43);
1018 }
1019 // pi-p
1020
1021 if(proj_momentum < 0.37)
1022 {
1023 hpXsc = 28.0 + 40*std::exp(-(proj_momentum-0.29)*(proj_momentum-0.29)/0.07/0.07);
1024 }
1025 else if(proj_momentum<0.65)
1026 {
1027 hpXsc = 26+110*(std::log(proj_momentum/0.48))*(std::log(proj_momentum/0.48));
1028 }
1029 else if(proj_momentum<1.3)
1030 {
1031 hpXsc = 36.1+
1032 10*std::exp(-(proj_momentum-0.72)*(proj_momentum-0.72)/0.06/0.06)+
1033 24*std::exp(-(proj_momentum-1.015)*(proj_momentum-1.015)/0.075/0.075);
1034 }
1035 else if(proj_momentum<3.0)
1036 {
1037 hpXsc = 36.1+0.079-4.313*std::log(proj_momentum)+
1038 3*std::exp(-(proj_momentum-2.1)*(proj_momentum-2.1)/0.4/0.4)+
1039 1.5*std::exp(-(proj_momentum-1.4)*(proj_momentum-1.4)/0.12/0.12);
1040 }
1041 else // mb
1042 {
1043 hpXsc = 10.6+2*std::log(proj_energy)+30*std::pow(proj_energy,-0.43);
1044 }
1045 xsection = hpXsc*zz + hnXsc*nn;
1046 }
1047 else if(theParticle == theKPlus)
1048 {
1049 xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
1050 + 7.14*std::pow(sMand,-eta1) - 13.45*std::pow(sMand,-eta2));
1051
1052 xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
1053 + 5.17*std::pow(sMand,-eta1) - 7.23*std::pow(sMand,-eta2));
1054 }
1055 else if(theParticle == theKMinus)
1056 {
1057 xsection = zz*( 17.91 + B*std::pow(std::log(sMand/s0),2.)
1058 + 7.14*std::pow(sMand,-eta1) + 13.45*std::pow(sMand,-eta2));
1059
1060 xsection += nn*( 17.87 + B*std::pow(std::log(sMand/s0),2.)
1061 + 5.17*std::pow(sMand,-eta1) + 7.23*std::pow(sMand,-eta2));
1062 }
1063 else if(theParticle == theSMinus)
1064 {
1065 xsection = aa*( 35.20 + B*std::pow(std::log(sMand/s0),2.)
1066 - 199.*std::pow(sMand,-eta1) + 264.*std::pow(sMand,-eta2));
1067 }
1068 else if(theParticle == theGamma) // modify later on
1069 {
1070 xsection = aa*( 0.0 + B*std::pow(std::log(sMand/s0),2.)
1071 + 0.032*std::pow(sMand,-eta1) - 0.0*std::pow(sMand,-eta2));
1072
1073 }
1074 else // as proton ???
1075 {
1076 xsection = zz*( 35.45 + B*std::pow(std::log(sMand/s0),2.)
1077 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
1078
1079 xsection += nn*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
1080 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
1081 }
1082 xsection *= millibarn; // parametrised in mb
1083 return xsection;
1084}
1085
1086
1087/////////////////////////////////////////////////////////////////////////////////////
1088//
1089// Returns hadron-nucleon inelastic cross-section based on proper parametrisation
1090
1091G4double
1092G4GlauberGribovCrossSection::GetHNinelasticXsc(const G4DynamicParticle* aParticle,
1093 const G4Element* anElement)
1094{
1095 G4int At = G4lrint(anElement->GetN()); // number of nucleons
1096 G4int Zt = G4lrint(anElement->GetZ()); // number of protons
1097
1098 return GetHNinelasticXsc(aParticle, At, Zt);
1099}
1100
1101/////////////////////////////////////////////////////////////////////////////////////
1102//
1103// Returns hadron-nucleon inelastic cross-section based on FTF-parametrisation
1104
1105G4double
1106G4GlauberGribovCrossSection::GetHNinelasticXsc(const G4DynamicParticle* aParticle,
1107 G4int At, G4int Zt)
1108{
1109 G4ParticleDefinition* hadron = aParticle->GetDefinition();
1110 G4double sumInelastic;
1111 G4int Nt = At - Zt;
1112 if(Nt < 0) Nt = 0;
1113
1114 if( hadron == theKPlus )
1115 {
1116 sumInelastic = GetHNinelasticXscVU(aParticle, At, Zt);
1117 }
1118 else
1119 {
1120 //sumInelastic = Zt*GetHadronNucleonXscMK(aParticle, theProton);
1121 // sumInelastic += Nt*GetHadronNucleonXscMK(aParticle, theNeutron);
1122 sumInelastic = G4double(Zt)*GetHadronNucleonXscNS(aParticle, 1, 1);
1123 sumInelastic += G4double(Nt)*GetHadronNucleonXscNS(aParticle, 1, 0);
1124 }
1125 return sumInelastic;
1126}
1127
1128
1129/////////////////////////////////////////////////////////////////////////////////////
1130//
1131// Returns hadron-nucleon inelastic cross-section based on FTF-parametrisation
1132
1133G4double
1134G4GlauberGribovCrossSection::GetHNinelasticXscVU(const G4DynamicParticle* aParticle,
1135 G4int At, G4int Zt)
1136{
1137 G4int PDGcode = aParticle->GetDefinition()->GetPDGEncoding();
1138 G4int absPDGcode = std::abs(PDGcode);
1139
1140 G4double Elab = aParticle->GetTotalEnergy();
1141 // (s - 2*0.88*GeV*GeV)/(2*0.939*GeV)/GeV;
1142 G4double Plab = aParticle->GetMomentum().mag();
1143 // std::sqrt(Elab * Elab - 0.88);
1144
1145 Elab /= GeV;
1146 Plab /= GeV;
1147
1148 G4double LogPlab = std::log( Plab );
1149 G4double sqrLogPlab = LogPlab * LogPlab;
1150
1151 //G4cout<<"Plab = "<<Plab<<G4endl;
1152
1153 G4double NumberOfTargetProtons = G4double(Zt);
1154 G4double NumberOfTargetNucleons = G4double(At);
1155 G4double NumberOfTargetNeutrons = NumberOfTargetNucleons - NumberOfTargetProtons;
1156
1157 if(NumberOfTargetNeutrons < 0.0) NumberOfTargetNeutrons = 0.0;
1158
1159 G4double Xtotal, Xelastic, Xinelastic;
1160
1161 if( absPDGcode > 1000 ) //------Projectile is baryon --------
1162 {
1163 G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) +
1164 0.522*sqrLogPlab - 4.51*LogPlab;
1165
1166 G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) +
1167 0.513*sqrLogPlab - 4.27*LogPlab;
1168
1169 G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) +
1170 0.169*sqrLogPlab - 1.85*LogPlab;
1171
1172 G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) +
1173 0.169*sqrLogPlab - 1.85*LogPlab;
1174
1175 Xtotal = (NumberOfTargetProtons * XtotPP +
1176 NumberOfTargetNeutrons * XtotPN);
1177
1178 Xelastic = (NumberOfTargetProtons * XelPP +
1179 NumberOfTargetNeutrons * XelPN);
1180 }
1181 else if( PDGcode == 211 ) //------Projectile is PionPlus -------
1182 {
1183 G4double XtotPiP = 16.4 + 19.3 *std::pow(Plab,-0.42) +
1184 0.19 *sqrLogPlab - 0.0 *LogPlab;
1185
1186 G4double XtotPiN = 33.0 + 14.0 *std::pow(Plab,-1.36) +
1187 0.456*sqrLogPlab - 4.03*LogPlab;
1188
1189 G4double XelPiP = 0.0 + 11.4*std::pow(Plab,-0.40) +
1190 0.079*sqrLogPlab - 0.0 *LogPlab;
1191
1192 G4double XelPiN = 1.76 + 11.2*std::pow(Plab,-0.64) +
1193 0.043*sqrLogPlab - 0.0 *LogPlab;
1194
1195 Xtotal = ( NumberOfTargetProtons * XtotPiP +
1196 NumberOfTargetNeutrons * XtotPiN );
1197
1198 Xelastic = ( NumberOfTargetProtons * XelPiP +
1199 NumberOfTargetNeutrons * XelPiN );
1200 }
1201 else if( PDGcode == -211 ) //------Projectile is PionMinus -------
1202 {
1203 G4double XtotPiP = 33.0 + 14.0 *std::pow(Plab,-1.36) +
1204 0.456*sqrLogPlab - 4.03*LogPlab;
1205
1206 G4double XtotPiN = 16.4 + 19.3 *std::pow(Plab,-0.42) +
1207 0.19 *sqrLogPlab - 0.0 *LogPlab;
1208
1209 G4double XelPiP = 1.76 + 11.2*std::pow(Plab,-0.64) +
1210 0.043*sqrLogPlab - 0.0 *LogPlab;
1211
1212 G4double XelPiN = 0.0 + 11.4*std::pow(Plab,-0.40) +
1213 0.079*sqrLogPlab - 0.0 *LogPlab;
1214
1215 Xtotal = ( NumberOfTargetProtons * XtotPiP +
1216 NumberOfTargetNeutrons * XtotPiN );
1217
1218 Xelastic = ( NumberOfTargetProtons * XelPiP +
1219 NumberOfTargetNeutrons * XelPiN );
1220 }
1221 else if( PDGcode == 111 ) //------Projectile is PionZero -------
1222 {
1223 G4double XtotPiP =(16.4 + 19.3 *std::pow(Plab,-0.42) +
1224 0.19 *sqrLogPlab - 0.0 *LogPlab + //Pi+
1225 33.0 + 14.0 *std::pow(Plab,-1.36) +
1226 0.456*sqrLogPlab - 4.03*LogPlab)/2; //Pi-
1227
1228 G4double XtotPiN =(33.0 + 14.0 *std::pow(Plab,-1.36) +
1229 0.456*sqrLogPlab - 4.03*LogPlab + //Pi+
1230 16.4 + 19.3 *std::pow(Plab,-0.42) +
1231 0.19 *sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1232
1233 G4double XelPiP =( 0.0 + 11.4*std::pow(Plab,-0.40) +
1234 0.079*sqrLogPlab - 0.0 *LogPlab + //Pi+
1235 1.76 + 11.2*std::pow(Plab,-0.64) +
1236 0.043*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1237
1238 G4double XelPiN =( 1.76 + 11.2*std::pow(Plab,-0.64) +
1239 0.043*sqrLogPlab - 0.0 *LogPlab + //Pi+
1240 0.0 + 11.4*std::pow(Plab,-0.40) +
1241 0.079*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
1242
1243 Xtotal = ( NumberOfTargetProtons * XtotPiP +
1244 NumberOfTargetNeutrons * XtotPiN );
1245
1246 Xelastic = ( NumberOfTargetProtons * XelPiP +
1247 NumberOfTargetNeutrons * XelPiN );
1248 }
1249 else if( PDGcode == 321 ) //------Projectile is KaonPlus -------
1250 {
1251 G4double XtotKP = 18.1 + 0. *std::pow(Plab, 0. ) +
1252 0.26 *sqrLogPlab - 1.0 *LogPlab;
1253 G4double XtotKN = 18.7 + 0. *std::pow(Plab, 0. ) +
1254 0.21 *sqrLogPlab - 0.89*LogPlab;
1255
1256 G4double XelKP = 5.0 + 8.1*std::pow(Plab,-1.8 ) +
1257 0.16 *sqrLogPlab - 1.3 *LogPlab;
1258
1259 G4double XelKN = 7.3 + 0. *std::pow(Plab,-0. ) +
1260 0.29 *sqrLogPlab - 2.4 *LogPlab;
1261
1262 Xtotal = ( NumberOfTargetProtons * XtotKP +
1263 NumberOfTargetNeutrons * XtotKN );
1264
1265 Xelastic = ( NumberOfTargetProtons * XelKP +
1266 NumberOfTargetNeutrons * XelKN );
1267 }
1268 else if( PDGcode ==-321 ) //------Projectile is KaonMinus ------
1269 {
1270 G4double XtotKP = 32.1 + 0. *std::pow(Plab, 0. ) +
1271 0.66 *sqrLogPlab - 5.6 *LogPlab;
1272 G4double XtotKN = 25.2 + 0. *std::pow(Plab, 0. ) +
1273 0.38 *sqrLogPlab - 2.9 *LogPlab;
1274
1275 G4double XelKP = 7.3 + 0. *std::pow(Plab,-0. ) +
1276 0.29 *sqrLogPlab - 2.4 *LogPlab;
1277
1278 G4double XelKN = 5.0 + 8.1*std::pow(Plab,-1.8 ) +
1279 0.16 *sqrLogPlab - 1.3 *LogPlab;
1280
1281 Xtotal = ( NumberOfTargetProtons * XtotKP +
1282 NumberOfTargetNeutrons * XtotKN );
1283
1284 Xelastic = ( NumberOfTargetProtons * XelKP +
1285 NumberOfTargetNeutrons * XelKN );
1286 }
1287 else if( PDGcode == 311 ) //------Projectile is KaonZero ------
1288 {
1289 G4double XtotKP = ( 18.1 + 0. *std::pow(Plab, 0. ) +
1290 0.26 *sqrLogPlab - 1.0 *LogPlab + //K+
1291 32.1 + 0. *std::pow(Plab, 0. ) +
1292 0.66 *sqrLogPlab - 5.6 *LogPlab)/2; //K-
1293
1294 G4double XtotKN = ( 18.7 + 0. *std::pow(Plab, 0. ) +
1295 0.21 *sqrLogPlab - 0.89*LogPlab + //K+
1296 25.2 + 0. *std::pow(Plab, 0. ) +
1297 0.38 *sqrLogPlab - 2.9 *LogPlab)/2; //K-
1298
1299 G4double XelKP = ( 5.0 + 8.1*std::pow(Plab,-1.8 )
1300 + 0.16 *sqrLogPlab - 1.3 *LogPlab + //K+
1301 7.3 + 0. *std::pow(Plab,-0. ) +
1302 0.29 *sqrLogPlab - 2.4 *LogPlab)/2; //K-
1303
1304 G4double XelKN = ( 7.3 + 0. *std::pow(Plab,-0. ) +
1305 0.29 *sqrLogPlab - 2.4 *LogPlab + //K+
1306 5.0 + 8.1*std::pow(Plab,-1.8 ) +
1307 0.16 *sqrLogPlab - 1.3 *LogPlab)/2; //K-
1308
1309 Xtotal = ( NumberOfTargetProtons * XtotKP +
1310 NumberOfTargetNeutrons * XtotKN );
1311
1312 Xelastic = ( NumberOfTargetProtons * XelKP +
1313 NumberOfTargetNeutrons * XelKN );
1314 }
1315 else //------Projectile is undefined, Nucleon assumed
1316 {
1317 G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) +
1318 0.522*sqrLogPlab - 4.51*LogPlab;
1319
1320 G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) +
1321 0.513*sqrLogPlab - 4.27*LogPlab;
1322
1323 G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) +
1324 0.169*sqrLogPlab - 1.85*LogPlab;
1325 G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) +
1326 0.169*sqrLogPlab - 1.85*LogPlab;
1327
1328 Xtotal = ( NumberOfTargetProtons * XtotPP +
1329 NumberOfTargetNeutrons * XtotPN );
1330
1331 Xelastic = ( NumberOfTargetProtons * XelPP +
1332 NumberOfTargetNeutrons * XelPN );
1333 }
1334 Xinelastic = Xtotal - Xelastic;
1335
1336 if( Xinelastic < 0.) Xinelastic = 0.;
1337
1338 return Xinelastic*= millibarn;
1339}
1340
1341////////////////////////////////////////////////////////////////////////////////////
1342//
1343//
1344
1345G4double
1346G4GlauberGribovCrossSection::GetNucleusRadius(const G4DynamicParticle* ,
1347 const G4Element* anElement)
1348{
1349 G4int At = G4lrint(anElement->GetN());
1350 G4double oneThird = 1.0/3.0;
1351 G4double cubicrAt = std::pow(G4double(At), oneThird);
1352
1353 G4double R; // = fRadiusConst*cubicrAt;
1354 /*
1355 G4double tmp = std::pow( cubicrAt-1., 3.);
1356 tmp += At;
1357 tmp *= 0.5;
1358
1359 if (At > 20.) // 20.
1360 {
1361 R = fRadiusConst*std::pow (tmp, oneThird);
1362 }
1363 else
1364 {
1365 R = fRadiusConst*cubicrAt;
1366 }
1367 */
1368
1369 R = fRadiusConst*cubicrAt;
1370
1371 G4double meanA = 21.;
1372
1373 G4double tauA1 = 40.;
1374 G4double tauA2 = 10.;
1375 G4double tauA3 = 5.;
1376
1377 G4double a1 = 0.85;
1378 G4double b1 = 1. - a1;
1379
1380 G4double b2 = 0.3;
1381 G4double b3 = 4.;
1382
1383 if (At > 20) // 20.
1384 {
1385 R *= ( a1 + b1*std::exp( -(At - meanA)/tauA1) );
1386 }
1387 else if (At > 3)
1388 {
1389 R *= ( 1.0 + b2*( 1. - std::exp( (At - meanA)/tauA2) ) );
1390 }
1391 else
1392 {
1393 R *= ( 1.0 + b3*( 1. - std::exp( (At - meanA)/tauA3) ) );
1394 }
1395 return R;
1396
1397}
1398////////////////////////////////////////////////////////////////////////////////////
1399//
1400//
1401
1402G4double
1403G4GlauberGribovCrossSection::GetNucleusRadius(G4int At)
1404{
1405 G4double oneThird = 1.0/3.0;
1406 G4double cubicrAt = std::pow(G4double(At), oneThird);
1407
1408 G4double R; // = fRadiusConst*cubicrAt;
1409
1410 /*
1411 G4double tmp = std::pow( cubicrAt-1., 3.);
1412 tmp += At;
1413 tmp *= 0.5;
1414
1415 if (At > 20.)
1416 {
1417 R = fRadiusConst*std::pow (tmp, oneThird);
1418 }
1419 else
1420 {
1421 R = fRadiusConst*cubicrAt;
1422 }
1423 */
1424
1425 R = fRadiusConst*cubicrAt;
1426
1427 G4double meanA = 20.;
1428 G4double tauA = 20.;
1429
1430 if (At > 20) // 20.
1431 {
1432 R *= ( 0.8 + 0.2*std::exp( -(G4double(At) - meanA)/tauA) );
1433 }
1434 else
1435 {
1436 R *= ( 1.0 + 0.1*( 1. - std::exp( (G4double(At) - meanA)/tauA) ) );
1437 }
1438
1439 return R;
1440}
1441
1442////////////////////////////////////////////////////////////////////////////////////
1443//
1444//
1445
1446G4double G4GlauberGribovCrossSection::CalculateEcmValue( const G4double mp ,
1447 const G4double mt ,
1448 const G4double Plab )
1449{
1450 G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
1451 G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
1452 // G4double Pcm = Plab * mt / Ecm;
1453 // G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
1454
1455 return Ecm ; // KEcm;
1456}
1457
1458
1459////////////////////////////////////////////////////////////////////////////////////
1460//
1461//
1462
1463G4double G4GlauberGribovCrossSection::CalcMandelstamS( const G4double mp ,
1464 const G4double mt ,
1465 const G4double Plab )
1466{
1467 G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
1468 G4double sMand = mp*mp + mt*mt + 2*Elab*mt ;
1469
1470 return sMand;
1471}
1472
1473
1474//
1475//
1476///////////////////////////////////////////////////////////////////////////////////////
Note: See TracBrowser for help on using the repository browser.