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

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

maj sur la beta de geant 4.9.3

File size: 28.3 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//
27// 24.11.08 V. Grichine - first implementation
28//
29
30#include "G4GGNuclNuclCrossSection.hh"
31
32#include "G4ParticleTable.hh"
33#include "G4IonTable.hh"
34#include "G4ParticleDefinition.hh"
35
36
37
38////////////////////////////////////////////////////////////////////////////////
39//
40//
41
42G4GGNuclNuclCrossSection::G4GGNuclNuclCrossSection() 
43: fUpperLimit( 100000 * GeV ),
44  fLowerLimit( 0.1 * MeV ),
45  fRadiusConst( 1.08*fermi )  // 1.1, 1.3 ?
46{
47  theProton   = G4Proton::Proton();
48  theNeutron  = G4Neutron::Neutron();
49}
50
51///////////////////////////////////////////////////////////////////////////////////////
52//
53//
54
55G4GGNuclNuclCrossSection::~G4GGNuclNuclCrossSection()
56{
57}
58
59
60////////////////////////////////////////////////////////////////////////////////////////
61//
62//
63
64
65G4bool
66G4GGNuclNuclCrossSection::IsApplicable(const G4DynamicParticle* aDP, 
67                                          const G4Element*  anElement)
68{
69  return IsZAApplicable(aDP, anElement->GetZ(), anElement->GetN());
70} 
71
72////////////////////////////////////////////////////////////////////////////////////////
73//
74//
75
76G4bool
77G4GGNuclNuclCrossSection::IsZAApplicable(const G4DynamicParticle* aDP, 
78                                            G4double Z, G4double)
79{
80  G4bool applicable      = false;
81  // G4int baryonNumber     = aDP->GetDefinition()->GetBaryonNumber();
82  G4double kineticEnergy = aDP->GetKineticEnergy();
83
84  //  const G4ParticleDefinition* theParticle = aDP->GetDefinition();
85 
86  if ( kineticEnergy  >= fLowerLimit && Z > 1.5 ) applicable = true;
87
88  return applicable;
89}
90
91////////////////////////////////////////////////////////////////////////////////////////
92//
93// Calculates total and inelastic Xsc, derives elastic as total - inelastic accordong to
94// Glauber model with Gribov correction calculated in the dipole approximation on
95// light cone. Gaussian density helps to calculate rest integrals of the model.
96// [1] B.Z. Kopeliovich, nucl-th/0306044
97
98
99G4double G4GGNuclNuclCrossSection::
100GetCrossSection(const G4DynamicParticle* aParticle, const G4Element* anElement, G4double T)
101{
102  return GetIsoZACrossSection(aParticle, anElement->GetZ(), anElement->GetN(), T);
103}
104
105////////////////////////////////////////////////////////////////////////////////////////
106//
107// Calculates total and inelastic Xsc, derives elastic as total - inelastic accordong to
108// Glauber model with Gribov correction calculated in the dipole approximation on
109// light cone. Gaussian density of point-like nucleons helps to calculate rest integrals of the model.
110// [1] B.Z. Kopeliovich, nucl-th/0306044 + simplification above
111
112
113
114G4double G4GGNuclNuclCrossSection::
115GetIsoZACrossSection(const G4DynamicParticle* aParticle, G4double tZ, G4double tA, G4double)
116{
117  G4double xsection, sigma, cofInelastic = 2.4, cofTotal = 2.0, nucleusSquare, cB, ratio;
118
119  G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
120  G4double pA = aParticle->GetDefinition()->GetBaryonNumber();
121
122  G4double pTkin = aParticle->GetKineticEnergy(); 
123  pTkin /= pA;
124
125  G4double pN = pA - pZ;
126  if( pN < 0. ) pN = 0.;
127
128  G4double tN = tA - tZ;
129  if( tN < 0. ) tN = 0.;
130
131  G4double tR = GetNucleusRadius(tA); 
132  G4double pR = GetNucleusRadius(pA); 
133
134  cB = GetCoulombBarier(aParticle, tZ, tA, pR, tR);
135
136  if(cB > 0.)
137  {
138
139    sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscNS(theProton, pTkin, theProton) +
140          (pZ*tN+pN*tZ)*GetHadronNucleonXscNS(theProton, pTkin, theNeutron);
141
142    nucleusSquare = cofTotal*pi*( pR*pR + tR*tR );   // basically 2piRR
143
144    ratio = sigma/nucleusSquare;
145
146    xsection =  nucleusSquare*std::log( 1. + ratio );
147
148    fTotalXsc = xsection;
149
150    fTotalXsc *= cB;
151
152    fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
153
154    fInelasticXsc *= cB;
155
156    fElasticXsc   = fTotalXsc - fInelasticXsc;
157
158    // if (fElasticXsc < DBL_MIN) fElasticXsc = DBL_MIN;
159    /*   
160    G4double difratio = ratio/(1.+ratio);
161
162    fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
163    */
164    // production to be checked !!! edit MK xsc
165
166    sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscMK(theProton, pTkin, theProton) +
167          (pZ*tN+pN*tZ)*GetHadronNucleonXscMK(theProton, pTkin, theNeutron);
168 
169    ratio = sigma/nucleusSquare;
170
171    fProductionXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
172
173    if (fElasticXsc < 0.) fElasticXsc = 0.;
174  }
175  else
176  {
177    fInelasticXsc  = 0.;
178    fTotalXsc      = 0.;
179    fElasticXsc    = 0.;
180    fProductionXsc = 0.;
181  }
182  return fInelasticXsc;   // xsection;
183}
184
185////////////////////////////////////////////////////////////////////////////////////////
186//
187//
188
189G4double G4GGNuclNuclCrossSection::
190GetCoulombBarier(const G4DynamicParticle* aParticle, G4double tZ, G4double tA, G4double pR, G4double tR)
191{
192  G4double ratio;
193  G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
194
195
196  G4double pTkin = aParticle->GetKineticEnergy();
197  // G4double pPlab = aParticle->GetTotalMomentum();
198  G4double pM    = aParticle->GetDefinition()->GetPDGMass();
199  // G4double tM    = tZ*proton_mass_c2 + (tA-tZ)*neutron_mass_c2; // ~ 1% accuracy
200  G4double tM    = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass( G4int(tZ), G4int(tA) );
201  G4double pElab = pTkin + pM;
202  G4double totEcm  = std::sqrt(pM*pM + tM*tM + 2.*pElab*tM);
203  // G4double pPcm  = pPlab*tM/totEcm;
204  // G4double pTcm  = std::sqrt(pM*pM + pPcm*pPcm) - pM;
205  G4double totTcm  = totEcm - pM -tM;
206
207  G4double bC    = fine_structure_const*hbarc*pZ*tZ;
208           bC   /= pR + tR;
209           bC   /= 2.;  // 4., 2. parametrisation cof ??? vmg
210
211           // G4cout<<"pTkin = "<<pTkin/GeV<<"; pPlab = "
212           // <<pPlab/GeV<<"; bC = "<<bC/GeV<<"; pTcm = "<<pTcm/GeV<<G4endl;
213
214  if( totTcm <= bC ) ratio = 0.;
215  else             ratio = 1. - bC/totTcm;
216
217  // if(ratio < DBL_MIN) ratio = DBL_MIN;
218  if( ratio < 0.) ratio = 0.;
219
220  // G4cout <<"ratio = "<<ratio<<G4endl;
221  return ratio;
222}
223
224
225//////////////////////////////////////////////////////////////////////////
226//
227// Return single-diffraction/inelastic cross-section ratio
228
229G4double G4GGNuclNuclCrossSection::
230GetRatioSD(const G4DynamicParticle* aParticle, G4double tA, G4double tZ)
231{
232  G4double sigma, cofInelastic = 2.4, cofTotal = 2.0, nucleusSquare, ratio;
233
234  G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
235  G4double pA = aParticle->GetDefinition()->GetBaryonNumber();
236
237  G4double pTkin = aParticle->GetKineticEnergy(); 
238  pTkin /= pA;
239
240  G4double pN = pA - pZ;
241  if( pN < 0. ) pN = 0.;
242
243  G4double tN = tA - tZ;
244  if( tN < 0. ) tN = 0.;
245
246  G4double tR = GetNucleusRadius(tA); 
247  G4double pR = GetNucleusRadius(pA); 
248
249  sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscNS(theProton, pTkin, theProton) +
250          (pZ*tN+pN*tZ)*GetHadronNucleonXscNS(theProton, pTkin, theNeutron);
251
252  nucleusSquare = cofTotal*pi*( pR*pR + tR*tR );   // basically 2piRR
253
254  ratio = sigma/nucleusSquare;
255
256
257  fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
258   
259  G4double difratio = ratio/(1.+ratio);
260
261  fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
262
263  if (fInelasticXsc > 0.) ratio = fDiffractionXsc/fInelasticXsc;
264  else                    ratio = 0.;
265
266  return ratio; 
267}
268
269//////////////////////////////////////////////////////////////////////////
270//
271// Return suasi-elastic/inelastic cross-section ratio
272
273G4double G4GGNuclNuclCrossSection::
274GetRatioQE(const G4DynamicParticle* aParticle, G4double tA, G4double tZ)
275{
276  G4double sigma, cofInelastic = 2.4, cofTotal = 2.0, nucleusSquare, ratio;
277
278  G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
279  G4double pA = aParticle->GetDefinition()->GetBaryonNumber();
280
281  G4double pTkin = aParticle->GetKineticEnergy(); 
282  pTkin /= pA;
283
284  G4double pN = pA - pZ;
285  if( pN < 0. ) pN = 0.;
286
287  G4double tN = tA - tZ;
288  if( tN < 0. ) tN = 0.;
289
290  G4double tR = GetNucleusRadius(tA); 
291  G4double pR = GetNucleusRadius(pA); 
292
293  sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscNS(theProton, pTkin, theProton) +
294          (pZ*tN+pN*tZ)*GetHadronNucleonXscNS(theProton, pTkin, theNeutron);
295
296  nucleusSquare = cofTotal*pi*( pR*pR + tR*tR );   // basically 2piRR
297
298  ratio = sigma/nucleusSquare;
299
300  fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
301
302  //  sigma = GetHNinelasticXsc(aParticle, tA, tZ);
303  ratio = sigma/nucleusSquare;
304
305  fProductionXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
306
307  if (fInelasticXsc > fProductionXsc) ratio = (fInelasticXsc-fProductionXsc)/fInelasticXsc;
308  else                                ratio = 0.;
309  if ( ratio < 0. )                   ratio = 0.;
310
311  return ratio; 
312}
313
314/////////////////////////////////////////////////////////////////////////////////////
315//
316// Returns hadron-nucleon Xsc according to differnt parametrisations:
317// [2] E. Levin, hep-ph/9710546
318// [3] U. Dersch, et al, hep-ex/9910052
319// [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
320
321G4double
322G4GGNuclNuclCrossSection::GetHadronNucleonXsc(const G4DynamicParticle* aParticle, 
323                                                  const G4Element* anElement          )
324{
325  G4double At = anElement->GetN();  // number of nucleons
326  G4double Zt = anElement->GetZ();  // number of protons
327
328
329  return GetHadronNucleonXsc( aParticle, At, Zt );
330}
331
332
333
334
335/////////////////////////////////////////////////////////////////////////////////////
336//
337// Returns hadron-nucleon Xsc according to differnt parametrisations:
338// [2] E. Levin, hep-ph/9710546
339// [3] U. Dersch, et al, hep-ex/9910052
340// [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
341
342G4double
343G4GGNuclNuclCrossSection::GetHadronNucleonXsc(const G4DynamicParticle* aParticle, 
344                                                   G4double At,  G4double Zt       )
345{
346  G4double xsection = 0.;
347
348
349  G4double targ_mass = G4ParticleTable::GetParticleTable()->
350  GetIonTable()->GetIonMass( G4int(Zt+0.5) , G4int(At+0.5) );
351
352  targ_mass = 0.939*GeV;  // ~mean neutron and proton ???
353
354  G4double proj_mass     = aParticle->GetMass();
355  G4double proj_momentum = aParticle->GetMomentum().mag();
356  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
357
358  sMand /= GeV*GeV;  // in GeV for parametrisation
359  proj_momentum /= GeV;
360
361  const G4ParticleDefinition* pParticle = aParticle->GetDefinition();
362 
363
364  if(pParticle == theNeutron) // as proton ???
365  {
366    xsection = At*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
367  } 
368  else if(pParticle == theProton) 
369  {
370    xsection = At*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
371    // xsection = At*( 49.51*std::pow(sMand,-0.097) + 0.314*std::log(sMand)*std::log(sMand) );
372    // xsection = At*( 38.4 + 0.85*std::abs(std::pow(log(sMand),1.47)) );
373  } 
374 
375  xsection *= millibarn;
376  return xsection;
377}
378
379
380/////////////////////////////////////////////////////////////////////////////////////
381//
382// Returns hadron-nucleon Xsc according to PDG parametrisation (2005):
383// http://pdg.lbl.gov/2006/reviews/hadronicrpp.pdf
384
385G4double
386G4GGNuclNuclCrossSection::GetHadronNucleonXscPDG(const G4DynamicParticle* aParticle, 
387                                                  const G4Element* anElement          )
388{
389  G4double At = anElement->GetN();  // number of nucleons
390  G4double Zt = anElement->GetZ();  // number of protons
391
392
393  return GetHadronNucleonXscPDG( aParticle, At, Zt );
394}
395
396
397
398
399/////////////////////////////////////////////////////////////////////////////////////
400//
401// Returns hadron-nucleon Xsc according to PDG parametrisation (2005):
402// http://pdg.lbl.gov/2006/reviews/hadronicrpp.pdf
403//  At = number of nucleons,  Zt = number of protons
404
405G4double
406G4GGNuclNuclCrossSection::GetHadronNucleonXscPDG(const G4DynamicParticle* aParticle, 
407                                                     G4double At,  G4double Zt )
408{
409  G4double xsection = 0.;
410
411  G4double Nt = At-Zt;              // number of neutrons
412  if (Nt < 0.) Nt = 0.; 
413
414
415  G4double targ_mass = G4ParticleTable::GetParticleTable()->
416  GetIonTable()->GetIonMass( G4int(Zt+0.5) , G4int(At+0.5) );
417
418  targ_mass = 0.939*GeV;  // ~mean neutron and proton ???
419
420  G4double proj_mass     = aParticle->GetMass(); 
421  G4double proj_momentum = aParticle->GetMomentum().mag();
422
423  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
424
425  sMand         /= GeV*GeV;  // in GeV for parametrisation
426
427  // General PDG fit constants
428
429  G4double s0   = 5.38*5.38; // in Gev^2
430  G4double eta1 = 0.458;
431  G4double eta2 = 0.458;
432  G4double B    = 0.308;
433
434
435  const G4ParticleDefinition* pParticle = aParticle->GetDefinition();
436 
437
438  if(pParticle == theNeutron) // proton-neutron fit
439  {
440    xsection = Zt*( 35.80 + B*std::pow(std::log(sMand/s0),2.) 
441                          + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
442    xsection  += Nt*( 35.45 + B*std::pow(std::log(sMand/s0),2.) 
443                      + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2)); // pp for nn
444  } 
445  else if(pParticle == theProton) 
446  {
447     
448      xsection  = Zt*( 35.45 + B*std::pow(std::log(sMand/s0),2.) 
449                          + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
450
451      xsection += Nt*( 35.80 + B*std::pow(std::log(sMand/s0),2.) 
452                          + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
453  } 
454  xsection *= millibarn; // parametrised in mb
455  return xsection;
456}
457
458
459
460
461/////////////////////////////////////////////////////////////////////////////////////
462//
463// Returns nucleon-nucleon cross-section based on N. Starkov parametrisation of
464// data from mainly http://wwwppds.ihep.su:8001/c5-6A.html database
465// projectile nucleon is pParticle with pTkin shooting target nucleon tParticle
466
467G4double
468G4GGNuclNuclCrossSection::GetHadronNucleonXscNS( G4ParticleDefinition* pParticle, 
469                                                 G4double pTkin, 
470                                                 G4ParticleDefinition* tParticle)
471{
472  G4double xsection(0), Delta, A0, B0;
473  G4double hpXsc(0);
474  G4double hnXsc(0);
475
476
477  G4double targ_mass     = tParticle->GetPDGMass();
478  G4double proj_mass     = pParticle->GetPDGMass(); 
479
480  G4double proj_energy   = proj_mass + pTkin; 
481  G4double proj_momentum = std::sqrt(pTkin*(pTkin+2*proj_mass));
482
483  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
484
485  sMand         /= GeV*GeV;  // in GeV for parametrisation
486  proj_momentum /= GeV;
487  proj_energy   /= GeV;
488  proj_mass     /= GeV;
489
490  // General PDG fit constants
491
492  //  G4double s0   = 5.38*5.38; // in Gev^2
493  //  G4double eta1 = 0.458;
494  //  G4double eta2 = 0.458;
495  //  G4double B    = 0.308;
496
497 
498
499
500 
501  if( proj_momentum >= 10. ) // high energy: pp = nn = np
502    // if( proj_momentum >= 2.)
503  {
504    Delta = 1.;
505
506    if( proj_energy < 40. ) Delta = 0.916+0.0021*proj_energy;
507
508    if( proj_momentum >= 10.)
509    {
510        B0 = 7.5;
511        A0 = 100. - B0*std::log(3.0e7);
512
513        xsection = A0 + B0*std::log(proj_energy) - 11
514                  + 103*std::pow(2*0.93827*proj_energy + proj_mass*proj_mass+
515                     0.93827*0.93827,-0.165);        //  mb
516    }
517  }
518  else // low energy pp = nn != np
519  {
520      if(pParticle == tParticle) // pp or nn      // nn to be pp
521      {
522        if( proj_momentum < 0.73 )
523        {
524          hnXsc = 23 + 50*( std::pow( std::log(0.73/proj_momentum), 3.5 ) );
525        }
526        else if( proj_momentum < 1.05  )
527        {
528          hnXsc = 23 + 40*(std::log(proj_momentum/0.73))*
529                         (std::log(proj_momentum/0.73));
530        }
531        else  // if( proj_momentum < 10.  )
532        {
533          hnXsc = 39.0 + 
534              75*(proj_momentum - 1.2)/(std::pow(proj_momentum,3.0) + 0.15);
535        }
536        xsection = hnXsc;
537      }
538      else  // pn to be np
539      {
540        if( proj_momentum < 0.8 )
541        {
542          hpXsc = 33+30*std::pow(std::log(proj_momentum/1.3),4.0);
543        }     
544        else if( proj_momentum < 1.4 )
545        {
546          hpXsc = 33+30*std::pow(std::log(proj_momentum/0.95),2.0);
547        }
548        else    // if( proj_momentum < 10.  )
549        {
550          hpXsc = 33.3+
551              20.8*(std::pow(proj_momentum,2.0)-1.35)/
552                 (std::pow(proj_momentum,2.50)+0.95);
553        }
554        xsection = hpXsc;
555      }
556  }
557  xsection *= millibarn; // parametrised in mb
558  return xsection;
559}
560
561/*
562/////////////////////////////////////////////////////////////////////////////////////
563//
564// Returns hadron-nucleon inelastic cross-section based on proper parametrisation
565
566G4double
567G4GGNuclNuclCrossSection::GetHNinelasticXsc(const G4DynamicParticle* aParticle,
568                                                  const G4Element* anElement          )
569{
570  G4double At = anElement->GetN();  // number of nucleons
571  G4double Zt = anElement->GetZ();  // number of protons
572
573
574  return GetHNinelasticXsc( aParticle, At, Zt );
575}
576
577/////////////////////////////////////////////////////////////////////////////////////
578//
579// Returns hadron-nucleon inelastic cross-section based on FTF-parametrisation
580
581G4double
582G4GGNuclNuclCrossSection::GetHNinelasticXsc(const G4DynamicParticle* aParticle,
583                                                     G4double At,  G4double Zt )
584{
585  //  G4ParticleDefinition* hadron = aParticle->GetDefinition();
586  G4double sumInelastic, Nt = At - Zt;
587
588  if(Nt < 0.) Nt = 0.;
589 
590  sumInelastic  = Zt*GetHadronNucleonXscMK(aParticle, theProton);
591  sumInelastic += Nt*GetHadronNucleonXscMK(aParticle, theNeutron);   
592 
593  return sumInelastic;
594}
595*/
596
597/////////////////////////////////////////////////////////////////////////////////////
598//
599// Returns hadron-nucleon inelastic cross-section based on FTF-parametrisation
600
601G4double
602G4GGNuclNuclCrossSection::GetHNinelasticXscVU(const G4DynamicParticle* aParticle, 
603                                                     G4double At,  G4double Zt )
604{
605  G4int PDGcode    = aParticle->GetDefinition()->GetPDGEncoding();
606  G4int absPDGcode = std::abs(PDGcode);
607
608  G4double Elab = aParticle->GetTotalEnergy();             
609                          // (s - 2*0.88*GeV*GeV)/(2*0.939*GeV)/GeV;
610  G4double Plab = aParticle->GetMomentum().mag();           
611                          // std::sqrt(Elab * Elab - 0.88);
612
613  Elab /= GeV;
614  Plab /= GeV;
615
616  G4double LogPlab    = std::log( Plab );
617  G4double sqrLogPlab = LogPlab * LogPlab;
618
619  //G4cout<<"Plab = "<<Plab<<G4endl;
620
621  G4double NumberOfTargetProtons  = Zt; 
622  G4double NumberOfTargetNucleons = At;
623  G4double NumberOfTargetNeutrons = NumberOfTargetNucleons - NumberOfTargetProtons;
624
625  if(NumberOfTargetNeutrons < 0.) NumberOfTargetNeutrons = 0.;
626
627  G4double Xtotal = 0., Xelastic = 0., Xinelastic =0.;
628
629  if( absPDGcode > 1000 )  //------Projectile is baryon --------
630  {
631       G4double XtotPP = 48.0 +  0. *std::pow(Plab, 0.  ) +
632                         0.522*sqrLogPlab - 4.51*LogPlab;
633
634       G4double XtotPN = 47.3 +  0. *std::pow(Plab, 0.  ) +
635                         0.513*sqrLogPlab - 4.27*LogPlab;
636
637       G4double XelPP  = 11.9 + 26.9*std::pow(Plab,-1.21) +
638                         0.169*sqrLogPlab - 1.85*LogPlab;
639
640       G4double XelPN  = 11.9 + 26.9*std::pow(Plab,-1.21) +
641                         0.169*sqrLogPlab - 1.85*LogPlab;
642
643       Xtotal          = ( NumberOfTargetProtons  * XtotPP +
644                           NumberOfTargetNeutrons * XtotPN  );
645
646       Xelastic        = ( NumberOfTargetProtons  * XelPP  +
647                           NumberOfTargetNeutrons * XelPN   );
648  }
649  Xinelastic = Xtotal - Xelastic;
650
651  if(Xinelastic < 0.) Xinelastic = 0.;
652
653  return Xinelastic*= millibarn;
654}
655
656/////////////////////////////////////////////////////////////////////////////////////
657//
658// Returns hadron-nucleon cross-section based on Mikhail Kossov CHIPS parametrisation of
659// data from G4QuasiFreeRatios class
660
661G4double
662G4GGNuclNuclCrossSection::GetHadronNucleonXscMK(G4ParticleDefinition* pParticle, G4double pTkin,
663                                                G4ParticleDefinition* nucleon  )
664{
665  G4int I = -1;
666  G4int PDG = pParticle->GetPDGEncoding();
667  G4double totalXsc = 0;
668  G4double elasticXsc = 0;
669  G4double inelasticXsc;
670  // G4int absPDG = std::abs(PDG);
671
672  G4double pM = pParticle->GetPDGMass();
673  G4double p  = std::sqrt(pTkin*(pTkin+2*pM))/GeV;
674
675  G4bool F = false;           
676  if(nucleon == theProton)       F = true;
677  else if(nucleon == theNeutron) F = false;
678  else
679  {
680    G4cout << "nucleon is not proton or neutron, return xsc for proton" << G4endl;
681    F = true;
682  }
683
684  G4bool kfl = true;                             // Flag of K0/aK0 oscillation
685  G4bool kf  = false;
686
687  if( PDG == 130 || PDG == 310 )
688  {
689    kf = true;
690    if( G4UniformRand() > .5 ) kfl = false;
691  }
692  if     ( (PDG == 2212 && F) || (PDG == 2112 && !F) ) I = 0; // pp/nn
693  else if( (PDG == 2112 && F) || (PDG == 2212 && !F) ) I = 1; // np/pn
694  else
695  {
696    G4cout<<"MK PDG = "<<PDG
697          <<", while it is defined only for p,n,hyperons,anti-baryons,pi,K/antiK"<<G4endl;
698    G4Exception("G4QuasiFreeRatio::FetchElTot:","22",FatalException,"CHIPScrash");
699  }
700
701  // Each parameter set can have not more than nPoints = 128 parameters
702
703  static const G4double lmi = 3.5;       // min of (lnP-lmi)^2 parabola
704  static const G4double pbe = .0557;     // elastic (lnP-lmi)^2 parabola coefficient
705  static const G4double pbt = .3;        // total (lnP-lmi)^2 parabola coefficient
706  static const G4double pmi = .1;        // Below that fast LE calculation is made
707  static const G4double pma = 1000.;     // Above that fast HE calculation is made
708                 
709  if( p <= 0.)
710  {
711    G4cout<<" p = "<<p<<" is zero or negative"<<G4endl;
712
713    elasticXsc   = 0.;
714    inelasticXsc = 0.;
715    totalXsc     = 0.;
716
717    return totalXsc;
718  }
719  if (!I)                          // pp/nn
720  {
721    if( p < pmi )
722    {
723      G4double p2 = p*p;
724      elasticXsc          = 1./(.00012 + p2*.2);
725      totalXsc          = elasticXsc;
726    }
727    else if(p>pma)
728    {
729      G4double lp  = std::log(p)-lmi;
730      G4double lp2 = lp*lp;
731      elasticXsc  = pbe*lp2 + 6.72;
732      totalXsc    = pbt*lp2 + 38.2;
733    }
734    else
735    {
736      G4double p2  = p*p;
737      G4double LE  = 1./( .00012 + p2*.2);
738      G4double lp  = std::log(p) - lmi;
739      G4double lp2 = lp*lp;
740      G4double rp2 = 1./p2;
741      elasticXsc  = LE + ( pbe*lp2 + 6.72+32.6/p)/( 1. + rp2/p);
742      totalXsc    = LE + ( pbt*lp2 + 38.2+52.7*rp2)/( 1. + 2.72*rp2*rp2);
743    }
744  }
745  else if( I==1 )                        // np/pn
746  {
747    if( p < pmi )
748    {
749      G4double p2 = p*p;
750      elasticXsc = 1./( .00012 + p2*( .051 + .1*p2));
751      totalXsc   = elasticXsc;
752    }
753    else if( p > pma )
754    {
755      G4double lp  = std::log(p) - lmi;
756      G4double lp2 = lp*lp;
757      elasticXsc  = pbe*lp2 + 6.72;
758      totalXsc    = pbt*lp2 + 38.2;
759    }
760    else
761    {
762      G4double p2  = p*p;
763      G4double LE  = 1./( .00012 + p2*( .051 + .1*p2 ) );
764      G4double lp  = std::log(p) - lmi;
765      G4double lp2 = lp*lp;
766      G4double rp2 = 1./p2;
767      elasticXsc  = LE + (pbe*lp2 + 6.72 + 30./p)/( 1. + .49*rp2/p);
768      totalXsc    = LE + (pbt*lp2 + 38.2)/( 1. + .54*rp2*rp2);
769    }
770  }
771  else
772  {
773    G4cout<<"PDG incoding = "<<I<<" is not defined (0-1)"<<G4endl;
774 
775  }
776  if( elasticXsc > totalXsc ) elasticXsc = totalXsc;
777
778  totalXsc   *= millibarn;
779  elasticXsc *= millibarn;
780  inelasticXsc   = totalXsc - elasticXsc;
781  if (inelasticXsc < 0.) inelasticXsc = 0.;
782
783  return inelasticXsc;
784}
785
786////////////////////////////////////////////////////////////////////////////////////
787//
788//
789
790G4double
791G4GGNuclNuclCrossSection::GetNucleusRadius( const G4DynamicParticle* , 
792                                               const G4Element* anElement)
793{
794  G4double At       = anElement->GetN();
795  G4double oneThird = 1.0/3.0;
796  G4double cubicrAt = std::pow (At, oneThird); 
797
798
799  G4double R;  // = fRadiusConst*cubicrAt;
800  /* 
801  G4double tmp = std::pow( cubicrAt-1., 3.);
802  tmp         += At;
803  tmp         *= 0.5;
804
805  if (At > 20.)   // 20.
806  {
807    R = fRadiusConst*std::pow (tmp, oneThird);
808  }
809  else
810  {
811    R = fRadiusConst*cubicrAt;
812  }
813  */
814 
815  R = fRadiusConst*cubicrAt;
816
817  // return R;  // !!!!
818
819
820 
821  G4double meanA  = 21.;
822
823  G4double tauA1  = 40.; 
824  G4double tauA2  = 10.; 
825  G4double tauA3  = 5.; 
826
827  G4double a1 = 0.85;
828  G4double b1 = 1. - a1;
829
830  G4double b2 = 0.3;
831  G4double b3 = 4.;
832
833  if (At > 20.)   // 20.
834  {
835    R *= ( a1 + b1*std::exp( -(At - meanA)/tauA1) ); 
836  }
837  else if (At > 3.5)
838  {
839    R *= ( 1.0 + b2*( 1. - std::exp( (At - meanA)/tauA2) ) ); 
840  }
841  else 
842  {
843    R *= ( 1.0 + b3*( 1. - std::exp( (At - meanA)/tauA3) ) ); 
844  } 
845  return R;
846 
847}
848
849////////////////////////////////////////////////////////////////////////////////////
850//
851//
852
853G4double
854G4GGNuclNuclCrossSection::GetNucleusRadius(G4double At)
855{
856  G4double R;
857
858  // R = GetNucleusRadiusGG(At);
859
860  R = GetNucleusRadiusDE(At);
861
862  return R;
863}
864
865///////////////////////////////////////////////////////////////////
866
867G4double
868G4GGNuclNuclCrossSection::GetNucleusRadiusGG(G4double At)
869{
870 
871  G4double oneThird = 1.0/3.0;
872  G4double cubicrAt = std::pow (At, oneThird); 
873
874
875  G4double R;  // = fRadiusConst*cubicrAt;
876 
877  /*
878  G4double tmp = std::pow( cubicrAt-1., 3.);
879  tmp         += At;
880  tmp         *= 0.5;
881
882  if (At > 20.)
883  {
884    R = fRadiusConst*std::pow (tmp, oneThird);
885  }
886  else
887  {
888    R = fRadiusConst*cubicrAt;
889  }
890  */
891 
892  R = fRadiusConst*cubicrAt;
893
894  G4double meanA = 20.;
895  G4double tauA  = 20.; 
896
897  if ( At > 20.)   // 20.
898  {
899    R *= ( 0.8 + 0.2*std::exp( -(At - meanA)/tauA) ); 
900  }
901  else
902  {
903    R *= ( 1.0 + 0.1*( 1. - std::exp( (At - meanA)/tauA) ) ); 
904  }
905
906  return R;
907
908}
909
910
911G4double
912G4GGNuclNuclCrossSection::GetNucleusRadiusDE(G4double A)
913{
914
915  // algorithm from diffuse-elastic
916
917  G4double R, r0, a11, a12, a13, a2, a3;
918
919  a11 = 1.26;  // 1.08, 1.16
920  a12 = 1.;  // 1.08, 1.16
921  a13 = 1.12;  // 1.08, 1.16
922  a2 = 1.1;
923  a3 = 1.;
924
925
926  if( A < 50. )
927  {
928    if( 10  < A && A <= 15. )     r0  = a11*( 1 - std::pow(A, -2./3.) )*fermi;   // 1.08*fermi;
929    else if( 15  < A && A <= 20 ) r0  = a12*( 1 - std::pow(A, -2./3.) )*fermi;
930    else if( 20  < A && A <= 30 ) r0  = a13*( 1 - std::pow(A, -2./3.) )*fermi;
931    else                          r0  = a2*fermi;
932
933    R = r0*std::pow( A, 1./3. );
934  }
935  else
936  {
937    r0 = a3*fermi;
938
939    R  = r0*std::pow(A, 0.27);
940  }
941  return R;
942
943
944
945}
946
947
948
949
950
951
952////////////////////////////////////////////////////////////////////////////////////
953//
954//
955
956G4double G4GGNuclNuclCrossSection::CalculateEcmValue( const G4double mp , 
957                                                         const G4double mt , 
958                                                         const G4double Plab )
959{
960  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
961  G4double Ecm  = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
962  // G4double Pcm  = Plab * mt / Ecm;
963  // G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
964
965  return Ecm ; // KEcm;
966}
967
968
969////////////////////////////////////////////////////////////////////////////////////
970//
971//
972
973G4double G4GGNuclNuclCrossSection::CalcMandelstamS( const G4double mp , 
974                                                       const G4double mt , 
975                                                       const G4double Plab )
976{
977  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
978  G4double sMand  = mp*mp + mt*mt + 2*Elab*mt ;
979
980  return sMand;
981}
982
983
984//
985//
986///////////////////////////////////////////////////////////////////////////////////////
Note: See TracBrowser for help on using the repository browser.