source: trunk/source/processes/hadronic/cross_sections/src/G4GGNuclNuclCrossSection.cc@ 1350

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

geant4 tag 9.4

File size: 22.5 KB
RevLine 
[968]1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// 24.11.08 V. Grichine - first implementation
27//
28
29#include "G4GGNuclNuclCrossSection.hh"
30
31#include "G4ParticleTable.hh"
32#include "G4IonTable.hh"
33#include "G4ParticleDefinition.hh"
[1340]34#include "G4HadTmpUtil.hh"
[968]35
36
37G4GGNuclNuclCrossSection::G4GGNuclNuclCrossSection()
38: fUpperLimit( 100000 * GeV ),
[1055]39 fLowerLimit( 0.1 * MeV ),
[1347]40 fRadiusConst( 1.08*fermi ), // 1.1, 1.3 ?
41 fTotalXsc(0.0), fElasticXsc(0.0), fInelasticXsc(0.0), fProductionXsc(0.0),
42 fDiffractionXsc(0.0), fHadronNucleonXsc(0.0)
[968]43{
44 theProton = G4Proton::Proton();
45 theNeutron = G4Neutron::Neutron();
46}
47
48
49G4GGNuclNuclCrossSection::~G4GGNuclNuclCrossSection()
[1340]50{}
[968]51
52
53G4bool
54G4GGNuclNuclCrossSection::IsApplicable(const G4DynamicParticle* aDP,
[1340]55 const G4Element* anElement)
[968]56{
[1340]57 G4int Z = G4lrint(anElement->GetZ());
58 G4int N = G4lrint(anElement->GetN());
59 return IsIsoApplicable(aDP, Z, N);
[968]60}
61
[1340]62///////////////////////////////////////////////////////////////////////////////
[968]63//
64//
65
66G4bool
[1340]67G4GGNuclNuclCrossSection::IsIsoApplicable(const G4DynamicParticle* aDP,
68 G4int Z, G4int)
[968]69{
[1340]70 G4bool applicable = false;
[968]71 G4double kineticEnergy = aDP->GetKineticEnergy();
72
[1340]73 if (kineticEnergy >= fLowerLimit && Z > 1) applicable = true;
[968]74 return applicable;
75}
76
[1340]77///////////////////////////////////////////////////////////////////////////////
[968]78//
[1340]79// Calculates total and inelastic Xsc, derives elastic as total - inelastic
80// accordong to Glauber model with Gribov correction calculated in the dipole
81// approximation on light cone. Gaussian density helps to calculate rest
82// integrals of the model. [1] B.Z. Kopeliovich, nucl-th/0306044
[968]83
84
85G4double G4GGNuclNuclCrossSection::
[1340]86GetCrossSection(const G4DynamicParticle* aParticle, const G4Element* anElement,
87 G4double T)
[968]88{
[1340]89 G4int Z = G4lrint(anElement->GetZ());
90 G4int N = G4lrint(anElement->GetN());
91 return GetZandACrossSection(aParticle, Z, N, T);
[968]92}
93
[1340]94///////////////////////////////////////////////////////////////////////////////
[968]95//
[1340]96// Calculates total and inelastic Xsc, derives elastic as total - inelastic
97// accordong to Glauber model with Gribov correction calculated in the dipole
98// approximation on light cone. Gaussian density of point-like nucleons helps
99// to calculate rest integrals of the model. [1] B.Z. Kopeliovich,
100// nucl-th/0306044 + simplification above
[968]101
102
103G4double G4GGNuclNuclCrossSection::
[1340]104GetZandACrossSection(const G4DynamicParticle* aParticle,
105 G4int tZ, G4int tA, G4double)
[968]106{
[1340]107 G4double xsection;
108 G4double sigma;
109 G4double cofInelastic = 2.4;
110 G4double cofTotal = 2.0;
111 G4double nucleusSquare;
112 G4double cB;
113 G4double ratio;
[968]114
115 G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
116 G4double pA = aParticle->GetDefinition()->GetBaryonNumber();
117
118 G4double pTkin = aParticle->GetKineticEnergy();
119 pTkin /= pA;
120
121 G4double pN = pA - pZ;
122 if( pN < 0. ) pN = 0.;
123
124 G4double tN = tA - tZ;
125 if( tN < 0. ) tN = 0.;
126
127 G4double tR = GetNucleusRadius(tA);
128 G4double pR = GetNucleusRadius(pA);
129
[1340]130 cB = GetCoulombBarier(aParticle, G4double(tZ), G4double(tA), pR, tR);
131 if (cB > 0.) {
[1055]132
133 sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscNS(theProton, pTkin, theProton) +
[968]134 (pZ*tN+pN*tZ)*GetHadronNucleonXscNS(theProton, pTkin, theNeutron);
135
[1055]136 nucleusSquare = cofTotal*pi*( pR*pR + tR*tR ); // basically 2piRR
[968]137
[1055]138 ratio = sigma/nucleusSquare;
139 xsection = nucleusSquare*std::log( 1. + ratio );
140 fTotalXsc = xsection;
141 fTotalXsc *= cB;
[968]142
[1055]143 fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
[968]144
[1055]145 fInelasticXsc *= cB;
146 fElasticXsc = fTotalXsc - fInelasticXsc;
[968]147
[1055]148 // if (fElasticXsc < DBL_MIN) fElasticXsc = DBL_MIN;
149 /*
150 G4double difratio = ratio/(1.+ratio);
151
152 fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
153 */
154 // production to be checked !!! edit MK xsc
155
[1228]156 //sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscMK(theProton, pTkin, theProton) +
157 // (pZ*tN+pN*tZ)*GetHadronNucleonXscMK(theProton, pTkin, theNeutron);
158
159 sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscNS(theProton, pTkin, theProton) +
160 (pZ*tN+pN*tZ)*GetHadronNucleonXscNS(theProton, pTkin, theNeutron);
[968]161
[1055]162 ratio = sigma/nucleusSquare;
163 fProductionXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
[968]164
[1055]165 if (fElasticXsc < 0.) fElasticXsc = 0.;
166 }
167 else
168 {
169 fInelasticXsc = 0.;
170 fTotalXsc = 0.;
171 fElasticXsc = 0.;
172 fProductionXsc = 0.;
173 }
174 return fInelasticXsc; // xsection;
175}
[968]176
[1340]177///////////////////////////////////////////////////////////////////////////////
[1055]178//
179//
180
181G4double G4GGNuclNuclCrossSection::
[1340]182GetCoulombBarier(const G4DynamicParticle* aParticle, G4double tZ, G4double tA,
183 G4double pR, G4double tR)
[1055]184{
185 G4double ratio;
186 G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
187
188 G4double pTkin = aParticle->GetKineticEnergy();
189 // G4double pPlab = aParticle->GetTotalMomentum();
190 G4double pM = aParticle->GetDefinition()->GetPDGMass();
191 // G4double tM = tZ*proton_mass_c2 + (tA-tZ)*neutron_mass_c2; // ~ 1% accuracy
192 G4double tM = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass( G4int(tZ), G4int(tA) );
193 G4double pElab = pTkin + pM;
194 G4double totEcm = std::sqrt(pM*pM + tM*tM + 2.*pElab*tM);
195 // G4double pPcm = pPlab*tM/totEcm;
196 // G4double pTcm = std::sqrt(pM*pM + pPcm*pPcm) - pM;
197 G4double totTcm = totEcm - pM -tM;
198
199 G4double bC = fine_structure_const*hbarc*pZ*tZ;
200 bC /= pR + tR;
201 bC /= 2.; // 4., 2. parametrisation cof ??? vmg
202
203 // G4cout<<"pTkin = "<<pTkin/GeV<<"; pPlab = "
204 // <<pPlab/GeV<<"; bC = "<<bC/GeV<<"; pTcm = "<<pTcm/GeV<<G4endl;
205
206 if( totTcm <= bC ) ratio = 0.;
207 else ratio = 1. - bC/totTcm;
208
209 // if(ratio < DBL_MIN) ratio = DBL_MIN;
210 if( ratio < 0.) ratio = 0.;
211
212 // G4cout <<"ratio = "<<ratio<<G4endl;
213 return ratio;
[968]214}
215
[1055]216
[968]217//////////////////////////////////////////////////////////////////////////
218//
219// Return single-diffraction/inelastic cross-section ratio
220
221G4double G4GGNuclNuclCrossSection::
222GetRatioSD(const G4DynamicParticle* aParticle, G4double tA, G4double tZ)
223{
224 G4double sigma, cofInelastic = 2.4, cofTotal = 2.0, nucleusSquare, ratio;
225
226 G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
227 G4double pA = aParticle->GetDefinition()->GetBaryonNumber();
228
229 G4double pTkin = aParticle->GetKineticEnergy();
230 pTkin /= pA;
231
232 G4double pN = pA - pZ;
233 if( pN < 0. ) pN = 0.;
234
235 G4double tN = tA - tZ;
236 if( tN < 0. ) tN = 0.;
237
238 G4double tR = GetNucleusRadius(tA);
239 G4double pR = GetNucleusRadius(pA);
240
241 sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscNS(theProton, pTkin, theProton) +
242 (pZ*tN+pN*tZ)*GetHadronNucleonXscNS(theProton, pTkin, theNeutron);
243
244 nucleusSquare = cofTotal*pi*( pR*pR + tR*tR ); // basically 2piRR
245 ratio = sigma/nucleusSquare;
[1340]246 fInelasticXsc = nucleusSquare*std::log(1. + cofInelastic*ratio)/cofInelastic;
[968]247 G4double difratio = ratio/(1.+ratio);
248
249 fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
250
251 if (fInelasticXsc > 0.) ratio = fDiffractionXsc/fInelasticXsc;
252 else ratio = 0.;
253
254 return ratio;
255}
256
257//////////////////////////////////////////////////////////////////////////
258//
[1340]259// Return quasi-elastic/inelastic cross-section ratio
[968]260
261G4double G4GGNuclNuclCrossSection::
262GetRatioQE(const G4DynamicParticle* aParticle, G4double tA, G4double tZ)
263{
264 G4double sigma, cofInelastic = 2.4, cofTotal = 2.0, nucleusSquare, ratio;
265
266 G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
267 G4double pA = aParticle->GetDefinition()->GetBaryonNumber();
268
269 G4double pTkin = aParticle->GetKineticEnergy();
270 pTkin /= pA;
271
272 G4double pN = pA - pZ;
273 if( pN < 0. ) pN = 0.;
274
275 G4double tN = tA - tZ;
276 if( tN < 0. ) tN = 0.;
277
278 G4double tR = GetNucleusRadius(tA);
279 G4double pR = GetNucleusRadius(pA);
280
281 sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscNS(theProton, pTkin, theProton) +
282 (pZ*tN+pN*tZ)*GetHadronNucleonXscNS(theProton, pTkin, theNeutron);
283
284 nucleusSquare = cofTotal*pi*( pR*pR + tR*tR ); // basically 2piRR
285 ratio = sigma/nucleusSquare;
[1340]286 fInelasticXsc = nucleusSquare*std::log(1. + cofInelastic*ratio)/cofInelastic;
[968]287
288 // sigma = GetHNinelasticXsc(aParticle, tA, tZ);
289 ratio = sigma/nucleusSquare;
[1340]290 fProductionXsc = nucleusSquare*std::log(1. + cofInelastic*ratio)/cofInelastic;
[968]291
292 if (fInelasticXsc > fProductionXsc) ratio = (fInelasticXsc-fProductionXsc)/fInelasticXsc;
293 else ratio = 0.;
294 if ( ratio < 0. ) ratio = 0.;
295
296 return ratio;
297}
298
[1340]299///////////////////////////////////////////////////////////////////////////////
[968]300//
301// Returns hadron-nucleon Xsc according to differnt parametrisations:
302// [2] E. Levin, hep-ph/9710546
303// [3] U. Dersch, et al, hep-ex/9910052
304// [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
305
306G4double
307G4GGNuclNuclCrossSection::GetHadronNucleonXsc(const G4DynamicParticle* aParticle,
[1340]308 const G4Element* anElement)
[968]309{
[1340]310 G4int At = G4lrint(anElement->GetN()); // number of nucleons
311 G4int Zt = G4lrint(anElement->GetZ()); // number of protons
312 return GetHadronNucleonXsc(aParticle, At, Zt);
[968]313}
314
315
316
317
[1340]318///////////////////////////////////////////////////////////////////////////////
[968]319//
320// Returns hadron-nucleon Xsc according to differnt parametrisations:
321// [2] E. Levin, hep-ph/9710546
322// [3] U. Dersch, et al, hep-ex/9910052
323// [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
324
325G4double
326G4GGNuclNuclCrossSection::GetHadronNucleonXsc(const G4DynamicParticle* aParticle,
[1340]327 G4int At, G4int Zt)
[968]328{
329 G4double xsection = 0.;
330
331 G4double targ_mass = G4ParticleTable::GetParticleTable()->
[1340]332 GetIonTable()->GetIonMass(Zt, At);
[968]333 targ_mass = 0.939*GeV; // ~mean neutron and proton ???
334
[1340]335 G4double proj_mass = aParticle->GetMass();
[968]336 G4double proj_momentum = aParticle->GetMomentum().mag();
337 G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
338
339 sMand /= GeV*GeV; // in GeV for parametrisation
340 proj_momentum /= GeV;
341 const G4ParticleDefinition* pParticle = aParticle->GetDefinition();
342
343 if(pParticle == theNeutron) // as proton ???
344 {
[1340]345 xsection = G4double(At)*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
[968]346 }
347 else if(pParticle == theProton)
348 {
[1340]349 xsection = G4double(At)*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
[968]350 }
351
352 xsection *= millibarn;
353 return xsection;
354}
355
356
[1340]357///////////////////////////////////////////////////////////////////////////////
[968]358//
359// Returns hadron-nucleon Xsc according to PDG parametrisation (2005):
360// http://pdg.lbl.gov/2006/reviews/hadronicrpp.pdf
361
362G4double
363G4GGNuclNuclCrossSection::GetHadronNucleonXscPDG(const G4DynamicParticle* aParticle,
[1340]364 const G4Element* anElement)
[968]365{
[1340]366 G4int At = G4lrint(anElement->GetN()); // number of nucleons
367 G4int Zt = G4lrint(anElement->GetZ()); // number of protons
[968]368 return GetHadronNucleonXscPDG( aParticle, At, Zt );
369}
370
371
[1340]372///////////////////////////////////////////////////////////////////////////////
[968]373//
374// Returns hadron-nucleon Xsc according to PDG parametrisation (2005):
375// http://pdg.lbl.gov/2006/reviews/hadronicrpp.pdf
376// At = number of nucleons, Zt = number of protons
377
378G4double
379G4GGNuclNuclCrossSection::GetHadronNucleonXscPDG(const G4DynamicParticle* aParticle,
[1340]380 G4int At, G4int Zt)
[968]381{
382 G4double xsection = 0.;
383
384 G4double Nt = At-Zt; // number of neutrons
385 if (Nt < 0.) Nt = 0.;
386
387 G4double targ_mass = G4ParticleTable::GetParticleTable()->
[1340]388 GetIonTable()->GetIonMass(Zt, At);
[968]389 targ_mass = 0.939*GeV; // ~mean neutron and proton ???
390
391 G4double proj_mass = aParticle->GetMass();
392 G4double proj_momentum = aParticle->GetMomentum().mag();
393 G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
394 sMand /= GeV*GeV; // in GeV for parametrisation
395
396 // General PDG fit constants
397
398 G4double s0 = 5.38*5.38; // in Gev^2
399 G4double eta1 = 0.458;
400 G4double eta2 = 0.458;
401 G4double B = 0.308;
402
403 const G4ParticleDefinition* pParticle = aParticle->GetDefinition();
404
405 if(pParticle == theNeutron) // proton-neutron fit
406 {
[1340]407 xsection = G4double(Zt)*( 35.80 + B*std::pow(std::log(sMand/s0),2.)
408 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
409 xsection += Nt*(35.45 + B*std::pow(std::log(sMand/s0),2.)
410 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2)); // pp for nn
[968]411 }
412 else if(pParticle == theProton)
413 {
[1340]414 xsection = G4double(Zt)*(35.45 + B*std::pow(std::log(sMand/s0),2.)
415 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
[968]416
[1340]417 xsection += Nt*(35.80 + B*std::pow(std::log(sMand/s0),2.)
418 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
[968]419 }
420 xsection *= millibarn; // parametrised in mb
421 return xsection;
422}
423
424
[1340]425///////////////////////////////////////////////////////////////////////////////
[968]426//
427// Returns nucleon-nucleon cross-section based on N. Starkov parametrisation of
428// data from mainly http://wwwppds.ihep.su:8001/c5-6A.html database
429// projectile nucleon is pParticle with pTkin shooting target nucleon tParticle
430
431G4double
[1340]432G4GGNuclNuclCrossSection::GetHadronNucleonXscNS(G4ParticleDefinition* pParticle,
[968]433 G4double pTkin,
434 G4ParticleDefinition* tParticle)
435{
436 G4double xsection(0), Delta, A0, B0;
437 G4double hpXsc(0);
438 G4double hnXsc(0);
439
[1340]440 G4double targ_mass = tParticle->GetPDGMass();
441 G4double proj_mass = pParticle->GetPDGMass();
[968]442
443 G4double proj_energy = proj_mass + pTkin;
444 G4double proj_momentum = std::sqrt(pTkin*(pTkin+2*proj_mass));
445
446 G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
447
448 sMand /= GeV*GeV; // in GeV for parametrisation
449 proj_momentum /= GeV;
450 proj_energy /= GeV;
451 proj_mass /= GeV;
452
453 // General PDG fit constants
454
455 // G4double s0 = 5.38*5.38; // in Gev^2
456 // G4double eta1 = 0.458;
457 // G4double eta2 = 0.458;
458 // G4double B = 0.308;
459
460 if( proj_momentum >= 10. ) // high energy: pp = nn = np
461 // if( proj_momentum >= 2.)
462 {
463 Delta = 1.;
[1340]464 if (proj_energy < 40.) Delta = 0.916+0.0021*proj_energy;
[968]465
[1340]466 if (proj_momentum >= 10.) {
467 B0 = 7.5;
468 A0 = 100. - B0*std::log(3.0e7);
[968]469
[1340]470 xsection = A0 + B0*std::log(proj_energy) - 11
471 + 103*std::pow(2*0.93827*proj_energy + proj_mass*proj_mass+
472 0.93827*0.93827,-0.165); // mb
[968]473 }
474 }
475 else // low energy pp = nn != np
476 {
477 if(pParticle == tParticle) // pp or nn // nn to be pp
478 {
479 if( proj_momentum < 0.73 )
480 {
481 hnXsc = 23 + 50*( std::pow( std::log(0.73/proj_momentum), 3.5 ) );
482 }
483 else if( proj_momentum < 1.05 )
484 {
485 hnXsc = 23 + 40*(std::log(proj_momentum/0.73))*
486 (std::log(proj_momentum/0.73));
487 }
488 else // if( proj_momentum < 10. )
489 {
490 hnXsc = 39.0 +
491 75*(proj_momentum - 1.2)/(std::pow(proj_momentum,3.0) + 0.15);
492 }
493 xsection = hnXsc;
494 }
495 else // pn to be np
496 {
497 if( proj_momentum < 0.8 )
498 {
499 hpXsc = 33+30*std::pow(std::log(proj_momentum/1.3),4.0);
500 }
501 else if( proj_momentum < 1.4 )
502 {
503 hpXsc = 33+30*std::pow(std::log(proj_momentum/0.95),2.0);
504 }
505 else // if( proj_momentum < 10. )
506 {
507 hpXsc = 33.3+
508 20.8*(std::pow(proj_momentum,2.0)-1.35)/
509 (std::pow(proj_momentum,2.50)+0.95);
510 }
511 xsection = hpXsc;
512 }
513 }
514 xsection *= millibarn; // parametrised in mb
515 return xsection;
516}
517
[1340]518/////////////////////////////////////////////////////////////////////////////////
[968]519//
520// Returns hadron-nucleon inelastic cross-section based on FTF-parametrisation
521
522G4double
523G4GGNuclNuclCrossSection::GetHNinelasticXscVU(const G4DynamicParticle* aParticle,
[1340]524 G4int At, G4int Zt)
[968]525{
[1340]526 G4int PDGcode = aParticle->GetDefinition()->GetPDGEncoding();
[968]527 G4int absPDGcode = std::abs(PDGcode);
528 G4double Elab = aParticle->GetTotalEnergy();
529 // (s - 2*0.88*GeV*GeV)/(2*0.939*GeV)/GeV;
530 G4double Plab = aParticle->GetMomentum().mag();
531 // std::sqrt(Elab * Elab - 0.88);
532
533 Elab /= GeV;
534 Plab /= GeV;
535
536 G4double LogPlab = std::log( Plab );
537 G4double sqrLogPlab = LogPlab * LogPlab;
538
539 //G4cout<<"Plab = "<<Plab<<G4endl;
540
541 G4double NumberOfTargetProtons = Zt;
542 G4double NumberOfTargetNucleons = At;
543 G4double NumberOfTargetNeutrons = NumberOfTargetNucleons - NumberOfTargetProtons;
544
545 if(NumberOfTargetNeutrons < 0.) NumberOfTargetNeutrons = 0.;
546
547 G4double Xtotal = 0., Xelastic = 0., Xinelastic =0.;
548
549 if( absPDGcode > 1000 ) //------Projectile is baryon --------
550 {
551 G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) +
552 0.522*sqrLogPlab - 4.51*LogPlab;
553
554 G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) +
555 0.513*sqrLogPlab - 4.27*LogPlab;
556
557 G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) +
558 0.169*sqrLogPlab - 1.85*LogPlab;
559
560 G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) +
561 0.169*sqrLogPlab - 1.85*LogPlab;
562
563 Xtotal = ( NumberOfTargetProtons * XtotPP +
564 NumberOfTargetNeutrons * XtotPN );
565
566 Xelastic = ( NumberOfTargetProtons * XelPP +
567 NumberOfTargetNeutrons * XelPN );
568 }
[1340]569
[968]570 Xinelastic = Xtotal - Xelastic;
571 if(Xinelastic < 0.) Xinelastic = 0.;
572
573 return Xinelastic*= millibarn;
574}
575
[1340]576///////////////////////////////////////////////////////////////////////////////
[968]577//
578//
579
580G4double
[1340]581G4GGNuclNuclCrossSection::GetNucleusRadius(const G4DynamicParticle* ,
582 const G4Element* anElement)
[968]583{
[1340]584 G4double At = anElement->GetN();
[968]585 G4double oneThird = 1.0/3.0;
586 G4double cubicrAt = std::pow (At, oneThird);
587
588 G4double R; // = fRadiusConst*cubicrAt;
589 R = fRadiusConst*cubicrAt;
590
591 G4double meanA = 21.;
592 G4double tauA1 = 40.;
593 G4double tauA2 = 10.;
594 G4double tauA3 = 5.;
595
596 G4double a1 = 0.85;
597 G4double b1 = 1. - a1;
598
599 G4double b2 = 0.3;
600 G4double b3 = 4.;
601
602 if (At > 20.) // 20.
603 {
604 R *= ( a1 + b1*std::exp( -(At - meanA)/tauA1) );
605 }
606 else if (At > 3.5)
607 {
608 R *= ( 1.0 + b2*( 1. - std::exp( (At - meanA)/tauA2) ) );
609 }
610 else
611 {
612 R *= ( 1.0 + b3*( 1. - std::exp( (At - meanA)/tauA3) ) );
[1340]613 }
614
[968]615 return R;
616}
617
[1340]618///////////////////////////////////////////////////////////////////////////////
[968]619//
620//
621
622G4double
623G4GGNuclNuclCrossSection::GetNucleusRadius(G4double At)
624{
625 G4double R;
626 R = GetNucleusRadiusDE(At);
627
628 return R;
629}
630
631///////////////////////////////////////////////////////////////////
632
633G4double
634G4GGNuclNuclCrossSection::GetNucleusRadiusGG(G4double At)
635{
636 G4double oneThird = 1.0/3.0;
637 G4double cubicrAt = std::pow (At, oneThird);
638
[1340]639 G4double R; // = fRadiusConst*cubicrAt;
[968]640 R = fRadiusConst*cubicrAt;
641
642 G4double meanA = 20.;
[1340]643 G4double tauA = 20.;
[968]644
645 if ( At > 20.) // 20.
646 {
647 R *= ( 0.8 + 0.2*std::exp( -(At - meanA)/tauA) );
648 }
649 else
650 {
651 R *= ( 1.0 + 0.1*( 1. - std::exp( (At - meanA)/tauA) ) );
652 }
653
654 return R;
655}
656
657
658G4double
659G4GGNuclNuclCrossSection::GetNucleusRadiusDE(G4double A)
660{
661 // algorithm from diffuse-elastic
662
663 G4double R, r0, a11, a12, a13, a2, a3;
664
665 a11 = 1.26; // 1.08, 1.16
666 a12 = 1.; // 1.08, 1.16
667 a13 = 1.12; // 1.08, 1.16
668 a2 = 1.1;
669 a3 = 1.;
670
671
[1340]672 if (A < 50.)
[968]673 {
674 if( 10 < A && A <= 15. ) r0 = a11*( 1 - std::pow(A, -2./3.) )*fermi; // 1.08*fermi;
675 else if( 15 < A && A <= 20 ) r0 = a12*( 1 - std::pow(A, -2./3.) )*fermi;
676 else if( 20 < A && A <= 30 ) r0 = a13*( 1 - std::pow(A, -2./3.) )*fermi;
677 else r0 = a2*fermi;
678
679 R = r0*std::pow( A, 1./3. );
680 }
681 else
682 {
683 r0 = a3*fermi;
684
685 R = r0*std::pow(A, 0.27);
686 }
687 return R;
688}
689
690
[1340]691///////////////////////////////////////////////////////////////////////////////
[968]692//
693//
694
[1340]695G4double G4GGNuclNuclCrossSection::CalculateEcmValue(const G4double mp,
696 const G4double mt,
697 const G4double Plab)
[968]698{
699 G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
700 G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
701 // G4double Pcm = Plab * mt / Ecm;
702 // G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
703
704 return Ecm ; // KEcm;
705}
706
707
[1340]708///////////////////////////////////////////////////////////////////////////////
[968]709//
710//
711
[1340]712G4double G4GGNuclNuclCrossSection::CalcMandelstamS(const G4double mp,
713 const G4double mt,
714 const G4double Plab)
[968]715{
716 G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
717 G4double sMand = mp*mp + mt*mt + 2*Elab*mt ;
718
719 return sMand;
720}
721
722//
723//
[1340]724///////////////////////////////////////////////////////////////////////////////
Note: See TracBrowser for help on using the repository browser.