Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4ComponentGGNuclNuclXsc.cc
Go to the documentation of this file.
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
30
32#include "G4SystemOfUnits.hh"
33#include "G4ParticleTable.hh"
34#include "G4IonTable.hh"
36#include "G4HadTmpUtil.hh"
37#include "G4HadronNucleonXsc.hh"
38
39
41 : G4VComponentCrossSection("Glauber-Gribov nucleus nucleus"),
42 fUpperLimit(100000*GeV), fLowerLimit(0.1*MeV),
43 fRadiusConst(1.08*fermi), // 1.1, 1.3 ?
44 fTotalXsc(0.0), fElasticXsc(0.0), fInelasticXsc(0.0), fProductionXsc(0.0),
45 fDiffractionXsc(0.0), fHadronNucleonXsc(0.0)
46{
47 theProton = G4Proton::Proton();
48 theNeutron = G4Neutron::Neutron();
49 hnXsc = new G4HadronNucleonXsc();
50}
51
52
54{}
55
56////////////////////////////////////////////////////////////////////
57
59 G4double kinEnergy,
60 G4int Z, G4int A)
61{
62 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
63 kinEnergy);
64 fInelasticXsc = GetZandACrossSection(aDP, Z, A);
65 delete aDP;
66
67 return fTotalXsc;
68}
69
70//////////////////////////////////////////////////////////////////////
71
73 G4double kinEnergy,
74 G4int Z, G4double A)
75{
76 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
77 kinEnergy);
78 fInelasticXsc = GetZandACrossSection(aDP, Z, G4int(A));
79 delete aDP;
80
81 return fTotalXsc;
82}
83
84////////////////////////////////////////////////////////////////////
85
87 G4double kinEnergy,
88 G4int Z, G4int A)
89{
90 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
91 kinEnergy);
92 fInelasticXsc = GetZandACrossSection(aDP, Z, A);
93 delete aDP;
94
95 return fInelasticXsc;
96}
97
98/////////////////////////////////////////////////////////////////////
99
101 G4double kinEnergy,
102 G4int Z, G4double A)
103{
104 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
105 kinEnergy);
106 fInelasticXsc = GetZandACrossSection(aDP, Z, G4int(A));
107 delete aDP;
108
109 return fInelasticXsc;
110}
111
112//////////////////////////////////////////////////////////////////
113
115 G4double kinEnergy,
116 G4int Z, G4double A)
117{
118 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
119 kinEnergy);
120 fInelasticXsc = GetZandACrossSection(aDP, Z, G4int(A));
121 delete aDP;
122
123 return fElasticXsc;
124}
125
126///////////////////////////////////////////////////////////////////
127
129 G4double kinEnergy,
130 G4int Z, G4int A)
131{
132 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
133 kinEnergy);
134 fInelasticXsc = GetZandACrossSection(aDP, Z, A);
135 delete aDP;
136
137 return fElasticXsc;
138}
139
140////////////////////////////////////////////////////////////////
141
143 G4double kinEnergy,
144 G4int Z, G4int A)
145{
146 G4DynamicParticle* aDP = new G4DynamicParticle(aParticle,G4ParticleMomentum(1.,0.,0.),
147 kinEnergy);
148 fInelasticXsc = GetZandACrossSection(aDP, Z, A);
149 delete aDP;
150 G4double ratio = 0.;
151
152 if(fInelasticXsc > 0.)
153 {
154 ratio = (fInelasticXsc - fProductionXsc)/fInelasticXsc;
155 if(ratio < 0.) ratio = 0.;
156 }
157 return ratio;
158}
159
160//////////////////////////////////////////////////////////////////////
161
162void
164{
165 outFile << "G4ComponentGGNuclNuclXsc calculates total, inelastic and\n"
166 << "elastic cross sections for nucleus-nucleus collisions using\n"
167 << "the Glauber model with Gribov corrections. It is valid for\n"
168 << "all incident energies above 100 keV./n";
169}
170
171/////////////////////////////////////////////////////////////////////
172
173G4bool
175 G4int Z, const G4Material*)
176{
177 G4bool applicable = false;
178 G4double kineticEnergy = aDP->GetKineticEnergy();
179
180 if (kineticEnergy >= fLowerLimit && Z > 1) applicable = true;
181 return applicable;
182}
183
184///////////////////////////////////////////////////////////////////////////////
185//
186// Calculates total and inelastic Xsc, derives elastic as total - inelastic
187// accordong to Glauber model with Gribov correction calculated in the dipole
188// approximation on light cone. Gaussian density helps to calculate rest
189// integrals of the model. [1] B.Z. Kopeliovich, nucl-th/0306044
190
191
194 const G4Material*)
195{
196 G4int A = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
197 return GetZandACrossSection(aParticle, Z, A);
198}
199
200///////////////////////////////////////////////////////////////////////////////
201//
202// Calculates total and inelastic Xsc, derives elastic as total - inelastic
203// accordong to Glauber model with Gribov correction calculated in the dipole
204// approximation on light cone. Gaussian density of point-like nucleons helps
205// to calculate rest integrals of the model. [1] B.Z. Kopeliovich,
206// nucl-th/0306044 + simplification above
207
208
211 G4int tZ, G4int tA)
212{
213 G4double xsection;
214 G4double sigma;
215 G4double cofInelastic = 2.4;
216 G4double cofTotal = 2.0;
217 G4double nucleusSquare;
218 G4double cB;
219 G4double ratio;
220
221 G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
222 G4double pA = aParticle->GetDefinition()->GetBaryonNumber();
223
224 G4double pTkin = aParticle->GetKineticEnergy();
225 pTkin /= pA;
226
227 G4double pN = pA - pZ;
228 if( pN < 0. ) pN = 0.;
229
230 G4double tN = tA - tZ;
231 if( tN < 0. ) tN = 0.;
232
234 G4double pR = GetNucleusRadius(pZ,pA);
235
236 cB = GetCoulombBarier(aParticle, G4double(tZ), G4double(tA), pR, tR);
237
238 if ( cB > 0. )
239 {
240 G4DynamicParticle* dProton = new G4DynamicParticle(theProton,
241 G4ParticleMomentum(1.,0.,0.),
242 pTkin);
243
244 G4DynamicParticle* dNeutron = new G4DynamicParticle(theNeutron,
245 G4ParticleMomentum(1.,0.,0.),
246 pTkin);
247
248 sigma = (pZ*tZ+pN*tN)*hnXsc->GetHadronNucleonXscNS(dProton, theProton);
249
250 G4double ppInXsc = hnXsc->GetInelasticHadronNucleonXsc();
251
252 sigma += (pZ*tN+pN*tZ)*hnXsc->GetHadronNucleonXscNS(dNeutron, theProton);
253
254 G4double npInXsc = hnXsc->GetInelasticHadronNucleonXsc();
255
256 // G4cout<<"ppInXsc = "<<ppInXsc/millibarn<<"; npInXsc = "<<npInXsc/millibarn<<G4endl;
257 // G4cout<<"npTotXsc = "<<hnXsc->GetTotalHadronNucleonXsc()/millibarn<<"; npElXsc = "
258 // <<hnXsc->GetElasticHadronNucleonXsc()/millibarn<<G4endl;
259
260 nucleusSquare = cofTotal*pi*( pR*pR + tR*tR ); // basically 2piRR
261
262 ratio = sigma/nucleusSquare;
263 xsection = nucleusSquare*std::log( 1. + ratio );
264 fTotalXsc = xsection;
265 fTotalXsc *= cB;
266
267 fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
268
269 fInelasticXsc *= cB;
270 fElasticXsc = fTotalXsc - fInelasticXsc;
271
272 // if (fElasticXsc < DBL_MIN) fElasticXsc = DBL_MIN;
273 /*
274 G4double difratio = ratio/(1.+ratio);
275
276 fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
277 */
278 // production to be checked !!! edit MK xsc
279
280 //sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscMK(theProton, pTkin, theProton) +
281 // (pZ*tN+pN*tZ)*GetHadronNucleonXscMK(theProton, pTkin, theNeutron);
282
283 sigma = (pZ*tZ+pN*tN)*ppInXsc + (pZ*tN+pN*tZ)*npInXsc;
284
285 ratio = sigma/nucleusSquare;
286 fProductionXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
287
288 if (fElasticXsc < 0.) fElasticXsc = 0.;
289 }
290 else
291 {
292 fInelasticXsc = 0.;
293 fTotalXsc = 0.;
294 fElasticXsc = 0.;
295 fProductionXsc = 0.;
296 }
297 return fInelasticXsc; // xsection;
298}
299
300///////////////////////////////////////////////////////////////////////////////
301//
302//
303
305GetCoulombBarier(const G4DynamicParticle* aParticle, G4double tZ, G4double tA,
306 G4double pR, G4double tR)
307{
308 G4double ratio;
309 G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
310
311 G4double pTkin = aParticle->GetKineticEnergy();
312 // G4double pPlab = aParticle->GetTotalMomentum();
313 G4double pM = aParticle->GetDefinition()->GetPDGMass();
314 // G4double tM = tZ*proton_mass_c2 + (tA-tZ)*neutron_mass_c2; // ~ 1% accuracy
316 G4double pElab = pTkin + pM;
317 G4double totEcm = std::sqrt(pM*pM + tM*tM + 2.*pElab*tM);
318 // G4double pPcm = pPlab*tM/totEcm;
319 // G4double pTcm = std::sqrt(pM*pM + pPcm*pPcm) - pM;
320 G4double totTcm = totEcm - pM -tM;
321
322 G4double bC = fine_structure_const*hbarc*pZ*tZ;
323 bC /= pR + tR;
324 bC /= 2.; // 4., 2. parametrisation cof ??? vmg
325
326 // G4cout<<"pTkin = "<<pTkin/GeV<<"; pPlab = "
327 // <<pPlab/GeV<<"; bC = "<<bC/GeV<<"; pTcm = "<<pTcm/GeV<<G4endl;
328
329 if( totTcm <= bC ) ratio = 0.;
330 else ratio = 1. - bC/totTcm;
331
332 // if(ratio < DBL_MIN) ratio = DBL_MIN;
333 if( ratio < 0.) ratio = 0.;
334
335 // G4cout <<"ratio = "<<ratio<<G4endl;
336 return ratio;
337}
338
339
340//////////////////////////////////////////////////////////////////////////
341//
342// Return single-diffraction/inelastic cross-section ratio
343
345GetRatioSD(const G4DynamicParticle* aParticle, G4double tA, G4double tZ)
346{
347 G4double sigma, cofInelastic = 2.4, cofTotal = 2.0, nucleusSquare, ratio;
348
349 G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
350 G4double pA = aParticle->GetDefinition()->GetBaryonNumber();
351
352 G4double pTkin = aParticle->GetKineticEnergy();
353 pTkin /= pA;
354
355 G4double pN = pA - pZ;
356 if( pN < 0. ) pN = 0.;
357
358 G4double tN = tA - tZ;
359 if( tN < 0. ) tN = 0.;
360
361 G4double tR = GetNucleusRadius(tZ,tA);
362 G4double pR = GetNucleusRadius(pZ,pA);
363
364 sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscNS(theProton, pTkin, theProton) +
365 (pZ*tN+pN*tZ)*GetHadronNucleonXscNS(theProton, pTkin, theNeutron);
366
367 nucleusSquare = cofTotal*pi*( pR*pR + tR*tR ); // basically 2piRR
368 ratio = sigma/nucleusSquare;
369 fInelasticXsc = nucleusSquare*std::log(1. + cofInelastic*ratio)/cofInelastic;
370 G4double difratio = ratio/(1.+ratio);
371
372 fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
373
374 if (fInelasticXsc > 0.) ratio = fDiffractionXsc/fInelasticXsc;
375 else ratio = 0.;
376
377 return ratio;
378}
379
380//////////////////////////////////////////////////////////////////////////
381//
382// Return quasi-elastic/inelastic cross-section ratio
383
385GetRatioQE(const G4DynamicParticle* aParticle, G4double tA, G4double tZ)
386{
387 G4double sigma, cofInelastic = 2.4, cofTotal = 2.0, nucleusSquare, ratio;
388
389 G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
390 G4double pA = aParticle->GetDefinition()->GetBaryonNumber();
391
392 G4double pTkin = aParticle->GetKineticEnergy();
393 pTkin /= pA;
394
395 G4double pN = pA - pZ;
396 if( pN < 0. ) pN = 0.;
397
398 G4double tN = tA - tZ;
399 if( tN < 0. ) tN = 0.;
400
401 G4double tR = GetNucleusRadius(tZ,tA);
402 G4double pR = GetNucleusRadius(pZ,pA);
403
404 sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscNS(theProton, pTkin, theProton) +
405 (pZ*tN+pN*tZ)*GetHadronNucleonXscNS(theProton, pTkin, theNeutron);
406
407 nucleusSquare = cofTotal*pi*( pR*pR + tR*tR ); // basically 2piRR
408 ratio = sigma/nucleusSquare;
409 fInelasticXsc = nucleusSquare*std::log(1. + cofInelastic*ratio)/cofInelastic;
410
411 // sigma = GetHNinelasticXsc(aParticle, tA, tZ);
412 ratio = sigma/nucleusSquare;
413 fProductionXsc = nucleusSquare*std::log(1. + cofInelastic*ratio)/cofInelastic;
414
415 if (fInelasticXsc > fProductionXsc) ratio = (fInelasticXsc-fProductionXsc)/fInelasticXsc;
416 else ratio = 0.;
417 if ( ratio < 0. ) ratio = 0.;
418
419 return ratio;
420}
421
422///////////////////////////////////////////////////////////////////////////////
423//
424// Returns hadron-nucleon Xsc according to differnt parametrisations:
425// [2] E. Levin, hep-ph/9710546
426// [3] U. Dersch, et al, hep-ex/9910052
427// [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
428
431 const G4Element* anElement)
432{
433 G4int At = G4lrint(anElement->GetN()); // number of nucleons
434 G4int Zt = G4lrint(anElement->GetZ()); // number of protons
435 return GetHadronNucleonXsc(aParticle, At, Zt);
436}
437
438
439
440
441///////////////////////////////////////////////////////////////////////////////
442//
443// Returns hadron-nucleon Xsc according to differnt parametrisations:
444// [2] E. Levin, hep-ph/9710546
445// [3] U. Dersch, et al, hep-ex/9910052
446// [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
447
450 G4int At, G4int Zt)
451{
452 G4double xsection = 0.;
453
455 GetIonTable()->GetIonMass(Zt, At);
456 targ_mass = 0.939*GeV; // ~mean neutron and proton ???
457
458 G4double proj_mass = aParticle->GetMass();
459 G4double proj_momentum = aParticle->GetMomentum().mag();
460 G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
461
462 sMand /= GeV*GeV; // in GeV for parametrisation
463 proj_momentum /= GeV;
464 const G4ParticleDefinition* pParticle = aParticle->GetDefinition();
465
466 if(pParticle == theNeutron) // as proton ???
467 {
468 xsection = G4double(At)*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
469 }
470 else if(pParticle == theProton)
471 {
472 xsection = G4double(At)*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
473 }
474
475 xsection *= millibarn;
476 return xsection;
477}
478
479
480
481///////////////////////////////////////////////////////////////////////////////
482//
483// Returns hadron-nucleon Xsc according to PDG parametrisation (2005):
484// http://pdg.lbl.gov/2006/reviews/hadronicrpp.pdf
485// At = number of nucleons, Zt = number of protons
486
489 G4double sMand,
490 G4ParticleDefinition* tParticle)
491{
492 G4double xsection = 0.;
493 // G4bool pORn = (tParticle == theProton || nucleon == theNeutron );
494 G4bool proton = (tParticle == theProton);
495 G4bool neutron = (tParticle == theNeutron);
496
497 // General PDG fit constants
498
499 G4double s0 = 5.38*5.38; // in Gev^2
500 G4double eta1 = 0.458;
501 G4double eta2 = 0.458;
502 G4double B = 0.308;
503
504 // const G4ParticleDefinition* pParticle = aParticle->GetDefinition();
505
506 if(pParticle == theNeutron) // proton-neutron fit
507 {
508 if ( proton )
509 {
510 xsection = ( 35.80 + B*std::pow(std::log(sMand/s0),2.)
511 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
512 }
513 if ( neutron )
514 {
515 xsection = (35.45 + B*std::pow(std::log(sMand/s0),2.)
516 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2)); // pp for nn
517 }
518 }
519 else if(pParticle == theProton)
520 {
521 if ( proton )
522 {
523 xsection = (35.45 + B*std::pow(std::log(sMand/s0),2.)
524 + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
525
526 }
527 if ( neutron )
528 {
529 xsection = (35.80 + B*std::pow(std::log(sMand/s0),2.)
530 + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
531 }
532 }
533 xsection *= millibarn; // parametrised in mb
534 return xsection;
535}
536
537
538///////////////////////////////////////////////////////////////////////////////
539//
540// Returns nucleon-nucleon cross-section based on N. Starkov parametrisation of
541// data from mainly http://wwwppds.ihep.su:8001/c5-6A.html database
542// projectile nucleon is pParticle with pTkin shooting target nucleon tParticle
543
546 G4double pTkin,
547 G4ParticleDefinition* tParticle)
548{
549 G4double xsection(0);
550 // G4double Delta; DHW 19 May 2011: variable set but not used
551 G4double A0, B0;
552 G4double hpXscv(0);
553 G4double hnXscv(0);
554
555 G4double targ_mass = tParticle->GetPDGMass();
556 G4double proj_mass = pParticle->GetPDGMass();
557
558 G4double proj_energy = proj_mass + pTkin;
559 G4double proj_momentum = std::sqrt(pTkin*(pTkin+2*proj_mass));
560
561 G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
562
563 sMand /= GeV*GeV; // in GeV for parametrisation
564 proj_momentum /= GeV;
565 proj_energy /= GeV;
566 proj_mass /= GeV;
567
568 // General PDG fit constants
569
570 // G4double s0 = 5.38*5.38; // in Gev^2
571 // G4double eta1 = 0.458;
572 // G4double eta2 = 0.458;
573 // G4double B = 0.308;
574
575 if( proj_momentum >= 373.)
576 {
577 return GetHadronNucleonXscPDG(pParticle,sMand,tParticle);
578 }
579 else if( proj_momentum >= 10. ) // high energy: pp = nn = np
580 // if( proj_momentum >= 2.)
581 {
582 // Delta = 1.; // DHW 19 May 2011: variable set but not used
583 // if (proj_energy < 40.) Delta = 0.916+0.0021*proj_energy;
584
585 if (proj_momentum >= 10.) {
586 B0 = 7.5;
587 A0 = 100. - B0*std::log(3.0e7);
588
589 xsection = A0 + B0*std::log(proj_energy) - 11
590 + 103*std::pow(2*0.93827*proj_energy + proj_mass*proj_mass+
591 0.93827*0.93827,-0.165); // mb
592 }
593 }
594 else // low energy pp = nn != np
595 {
596 if(pParticle == tParticle) // pp or nn // nn to be pp
597 {
598 if( proj_momentum < 0.73 )
599 {
600 hnXscv = 23 + 50*( std::pow( std::log(0.73/proj_momentum), 3.5 ) );
601 }
602 else if( proj_momentum < 1.05 )
603 {
604 hnXscv = 23 + 40*(std::log(proj_momentum/0.73))*
605 (std::log(proj_momentum/0.73));
606 }
607 else // if( proj_momentum < 10. )
608 {
609 hnXscv = 39.0 +
610 75*(proj_momentum - 1.2)/(std::pow(proj_momentum,3.0) + 0.15);
611 }
612 xsection = hnXscv;
613 }
614 else // pn to be np
615 {
616 if( proj_momentum < 0.8 )
617 {
618 hpXscv = 33+30*std::pow(std::log(proj_momentum/1.3),4.0);
619 }
620 else if( proj_momentum < 1.4 )
621 {
622 hpXscv = 33+30*std::pow(std::log(proj_momentum/0.95),2.0);
623 }
624 else // if( proj_momentum < 10. )
625 {
626 hpXscv = 33.3+
627 20.8*(std::pow(proj_momentum,2.0)-1.35)/
628 (std::pow(proj_momentum,2.50)+0.95);
629 }
630 xsection = hpXscv;
631 }
632 }
633 xsection *= millibarn; // parametrised in mb
634 return xsection;
635}
636
637/////////////////////////////////////////////////////////////////////////////////
638//
639// Returns hadron-nucleon inelastic cross-section based on FTF-parametrisation
640
643 G4int At, G4int Zt)
644{
645 G4int PDGcode = aParticle->GetDefinition()->GetPDGEncoding();
646 G4int absPDGcode = std::abs(PDGcode);
647 G4double Elab = aParticle->GetTotalEnergy();
648 // (s - 2*0.88*GeV*GeV)/(2*0.939*GeV)/GeV;
649 G4double Plab = aParticle->GetMomentum().mag();
650 // std::sqrt(Elab * Elab - 0.88);
651
652 Elab /= GeV;
653 Plab /= GeV;
654
655 G4double LogPlab = std::log( Plab );
656 G4double sqrLogPlab = LogPlab * LogPlab;
657
658 //G4cout<<"Plab = "<<Plab<<G4endl;
659
660 G4double NumberOfTargetProtons = Zt;
661 G4double NumberOfTargetNucleons = At;
662 G4double NumberOfTargetNeutrons = NumberOfTargetNucleons - NumberOfTargetProtons;
663
664 if(NumberOfTargetNeutrons < 0.) NumberOfTargetNeutrons = 0.;
665
666 G4double Xtotal = 0., Xelastic = 0., Xinelastic =0.;
667
668 if( absPDGcode > 1000 ) //------Projectile is baryon --------
669 {
670 G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) +
671 0.522*sqrLogPlab - 4.51*LogPlab;
672
673 G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) +
674 0.513*sqrLogPlab - 4.27*LogPlab;
675
676 G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) +
677 0.169*sqrLogPlab - 1.85*LogPlab;
678
679 G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) +
680 0.169*sqrLogPlab - 1.85*LogPlab;
681
682 Xtotal = ( NumberOfTargetProtons * XtotPP +
683 NumberOfTargetNeutrons * XtotPN );
684
685 Xelastic = ( NumberOfTargetProtons * XelPP +
686 NumberOfTargetNeutrons * XelPN );
687 }
688
689 Xinelastic = Xtotal - Xelastic;
690 if(Xinelastic < 0.) Xinelastic = 0.;
691
692 return Xinelastic*= millibarn;
693}
694
695///////////////////////////////////////////////////////////////////////////////
696//
697//
698
701 const G4Element* anElement)
702{
703 G4double At = anElement->GetN();
704 G4double oneThird = 1.0/3.0;
705 G4double cubicrAt = std::pow (At, oneThird);
706
707 G4double R; // = fRadiusConst*cubicrAt;
708 R = fRadiusConst*cubicrAt;
709
710 G4double meanA = 21.;
711 G4double tauA1 = 40.;
712 G4double tauA2 = 10.;
713 G4double tauA3 = 5.;
714
715 G4double a1 = 0.85;
716 G4double b1 = 1. - a1;
717
718 G4double b2 = 0.3;
719 G4double b3 = 4.;
720
721 if (At > 20.) // 20.
722 {
723 R *= ( a1 + b1*std::exp( -(At - meanA)/tauA1) );
724 }
725 else if (At > 3.5)
726 {
727 R *= ( 1.0 + b2*( 1. - std::exp( (At - meanA)/tauA2) ) );
728 }
729 else
730 {
731 R *= ( 1.0 + b3*( 1. - std::exp( (At - meanA)/tauA3) ) );
732 }
733
734 return R;
735}
736
737///////////////////////////////////////////////////////////////////////////////
738//
739//
740
743{
744 G4double R;
745 R = GetNucleusRadiusDE(Zt,At);
746 // R = GetNucleusRadiusRMS(Zt,At);
747
748 return R;
749}
750
751///////////////////////////////////////////////////////////////////
752
755{
756 G4double oneThird = 1.0/3.0;
757 G4double cubicrAt = std::pow (At, oneThird);
758
759 G4double R; // = fRadiusConst*cubicrAt;
760 R = fRadiusConst*cubicrAt;
761
762 G4double meanA = 20.;
763 G4double tauA = 20.;
764
765 if ( At > 20.) // 20.
766 {
767 R *= ( 0.8 + 0.2*std::exp( -(At - meanA)/tauA) );
768 }
769 else
770 {
771 R *= ( 1.0 + 0.1*( 1. - std::exp( (At - meanA)/tauA) ) );
772 }
773
774 return R;
775}
776
777/////////////////////////////////////////////////////////////////////////////
778//
779//
780
783{
784 // algorithm from diffuse-elastic
785
786 G4double R, r0, a11, a12, a13, a2, a3;
787
788 a11 = 1.26; // 1.08, 1.16
789 a12 = 1.; // 1.08, 1.16
790 a13 = 1.12; // 1.08, 1.16
791 a2 = 1.1;
792 a3 = 1.;
793
794 // Special rms radii for light nucleii
795
796 if (A < 50.)
797 {
798 if (std::abs(A-1.) < 0.5) return 0.89*fermi; // p
799 else if(std::abs(A-2.) < 0.5) return 2.13*fermi; // d
800 else if(std::abs(Z-1.) < 0.5 && std::abs(A-3.) < 0.5) return 1.80*fermi; // t
801
802 else if(std::abs(Z-2.) < 0.5 && std::abs(A-3.) < 0.5) return 1.96*fermi; // He3
803 else if(std::abs(Z-2.) < 0.5 && std::abs(A-4.) < 0.5) return 1.68*fermi; // He4
804
805 else if(std::abs(Z-3.) < 0.5) return 2.40*fermi; // Li7
806 else if(std::abs(Z-4.) < 0.5) return 2.51*fermi; // Be9
807
808 else if( 10. < A && A <= 16. ) r0 = a11*( 1 - std::pow(A, -2./3.) )*fermi; // 1.08*fermi;
809 else if( 15. < A && A <= 20. ) r0 = a12*( 1 - std::pow(A, -2./3.) )*fermi;
810 else if( 20. < A && A <= 30. ) r0 = a13*( 1 - std::pow(A, -2./3.) )*fermi;
811 else r0 = a2*fermi;
812
813 R = r0*std::pow( A, 1./3. );
814 }
815 else
816 {
817 r0 = a3*fermi;
818
819 R = r0*std::pow(A, 0.27);
820 }
821 return R;
822}
823
824
825/////////////////////////////////////////////////////////////////////////////
826//
827// RMS radii from e-A scattering data
828
831{
832
833 if (std::abs(A-1.) < 0.5) return 0.89*fermi; // p
834 else if(std::abs(A-2.) < 0.5) return 2.13*fermi; // d
835 else if(std::abs(Z-1.) < 0.5 && std::abs(A-3.) < 0.5) return 1.80*fermi; // t
836
837 else if(std::abs(Z-2.) < 0.5 && std::abs(A-3.) < 0.5) return 1.96*fermi; // He3
838 else if(std::abs(Z-2.) < 0.5 && std::abs(A-4.) < 0.5) return 1.68*fermi; // He4
839
840 else if(std::abs(Z-3.) < 0.5) return 2.40*fermi; // Li7
841 else if(std::abs(Z-4.) < 0.5) return 2.51*fermi; // Be9
842
843 else return 1.24*std::pow(A, 0.28 )*fermi; // A > 9
844}
845
846
847///////////////////////////////////////////////////////////////////////////////
848//
849//
850
852 const G4double mt,
853 const G4double Plab)
854{
855 G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
856 G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
857 // G4double Pcm = Plab * mt / Ecm;
858 // G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
859
860 return Ecm ; // KEcm;
861}
862
863
864///////////////////////////////////////////////////////////////////////////////
865//
866//
867
869 const G4double mt,
870 const G4double Plab)
871{
872 G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
873 G4double sMand = mp*mp + mt*mt + 2*Elab*mt ;
874
875 return sMand;
876}
877
878//
879//
880///////////////////////////////////////////////////////////////////////////////
@ neutron
G4ThreeVector G4ParticleMomentum
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
double mag() const
virtual G4double GetInelasticElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
G4double GetHadronNucleonXscPDG(G4ParticleDefinition *, G4double sMand, G4ParticleDefinition *)
G4double CalculateEcmValue(const G4double, const G4double, const G4double)
virtual G4double ComputeQuasiElasticRatio(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
G4double CalcMandelstamS(const G4double, const G4double, const G4double)
G4double GetNucleusRadiusGG(G4double At)
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *)
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
G4double GetHadronNucleonXsc(const G4DynamicParticle *, const G4Element *)
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *)
G4double GetRatioQE(const G4DynamicParticle *, G4double At, G4double Zt)
G4double GetRatioSD(const G4DynamicParticle *, G4double At, G4double Zt)
virtual G4double GetTotalElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
G4double GetNucleusRadiusRMS(G4double Zt, G4double At)
virtual G4double GetElasticElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double A)
virtual G4double GetTotalIsotopeCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
G4double GetNucleusRadiusDE(G4double Zt, G4double At)
G4double GetZandACrossSection(const G4DynamicParticle *, G4int Z, G4int A)
G4double GetHadronNucleonXscNS(G4ParticleDefinition *, G4double pTkin, G4ParticleDefinition *)
G4double GetHNinelasticXscVU(const G4DynamicParticle *, G4int At, G4int Zt)
G4double GetCoulombBarier(const G4DynamicParticle *, G4double Z, G4double A, G4double pR, G4double tR)
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4int A)
virtual void CrossSectionDescription(std::ostream &) const
G4double GetMass() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
G4double GetZ() const
Definition: G4Element.hh:131
G4double GetN() const
Definition: G4Element.hh:134
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetInelasticHadronNucleonXsc()
G4double GetIonMass(G4int Z, G4int A, G4int L=0) const
!! Only ground states are supported now
Definition: G4IonTable.cc:774
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4NistManager * Instance()
G4double GetPDGCharge() const
static G4ParticleTable * GetParticleTable()
G4IonTable * GetIonTable()
static G4Proton * Proton()
Definition: G4Proton.cc:93
int G4lrint(double ad)
Definition: templates.hh:163