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

Last change on this file since 1347 was 1347, checked in by garnier, 13 years ago

geant4 tag 9.4

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