Geant4  10.02.p02
G4hhElastic.hh
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 //
27 // $Id: G4hhElastic.hh,v 1.21 2010-11-09 09:04:29 grichine Exp $
28 // GEANT4 tag $Name: not supported by cvs2svn $
29 //
30 //
31 // G4 Model: hadron diffraction elastic scattering with 4-momentum balance
32 //
33 // Class Description
34 // Final state production model for hadron-hadron elastic scattering
35 // in the framework of quark-diquark model with springy Pomeron.
36 // Projectiles are proton, neutron, pions, kaons.
37 // Targets are proton (and neutron).
38 // Class Description - End
39 //
40 // 02.05.14 V. Grichine - 1-st implementation
41 // 10.10.14 V. Grichine - change to combine with low mass diffraction
42 
43 #ifndef G4hhElastic_h
44 #define G4hhElastic_h 1
45 
46 #include "globals.hh"
47 #include <complex>
48 #include "G4Integrator.hh"
49 
50 #include "G4HadronElastic.hh"
51 #include "G4HadProjectile.hh"
52 #include "G4Nucleus.hh"
53 #include "G4HadronNucleonXsc.hh"
54 
55 #include "G4Exp.hh"
56 #include "G4Log.hh"
57 
58 
60 class G4PhysicsTable;
61 class G4PhysicsLogVector;
62 
64 {
65 public:
66  // PL constructor
67  G4hhElastic();
68  // test constructor
70  G4double plab );
71  // constructor used for low mass diffraction
73 
74  virtual ~G4hhElastic();
75 
76  virtual G4bool IsApplicable(const G4HadProjectile &/*aTrack*/,
77  G4Nucleus & /*targetNucleus*/);
78 
79 
80  void Initialise();
81 
82  void BuildTableT( G4ParticleDefinition* target, G4ParticleDefinition* projectile); // , G4double plab );
83  void BuildTableTest( G4ParticleDefinition* target, G4ParticleDefinition* projectile,
84  G4double plab );
85 
88 
89  G4double SampleTest(G4double tMin ); // const G4ParticleDefinition* p, );
90 
91  G4double GetTransfer( G4int iMomentum, G4int iTransfer, G4double position );
92 
93 private:
94 
95 
102 
103 
109 
112 
115  std::vector<G4PhysicsTable*> fBankT;
116 
117 
118  // Gauss model parameters
119 
120 
124 
129 
130 
131  G4double fRA; // hadron A
136 
137  G4double fRB; // hadron B
142 
153 
161  G4double fQcof; // q prime when integrate
162 
163 public: // Gauss model methods
164 
165  void SetParameters();
166  void SetSigmaTot(G4double stot){fSigmaTot = stot;};
167  void SetSpp(G4double spp){fSpp = spp;};
168  G4double GetSpp(){return fSpp;};
169  void SetParametersCMS(G4double plab);
170 
171  G4double GetBq(){ return fBq;};
172  G4double GetBQ(){ return fBQ;};
173  G4double GetBqQ(){ return fBqQ;};
174  void SetBq(G4double b){fBq = b;};
175  void SetBQ(G4double b){fBQ = b;};
176  void SetBqQ(G4double b){fBqQ = b;};
178 
179  void CalculateBQ(G4double b);
180  void CalculateBqQ13(G4double b);
181  void CalculateBqQ12(G4double b);
182  void CalculateBqQ123(G4double b);
183  void SetRA(G4double rn, G4double pq, G4double pQ);
184  void SetRB(G4double rn, G4double pq, G4double pQ);
185 
186  void SetAlphaP(G4double a){fAlphaP = a;};
187  void SetImCof(G4double a){fImCof = a;};
189  void SetLambda(G4double L){fLambda = L;};
190  void SetEta(G4double E){fEta = E;};
191  void SetCofF2(G4double f){fCofF2 = f;};
192  void SetCofF3(G4double f){fCofF3 = f;};
195 
196  G4double GetRA(){ return fRA;};
197  G4double GetRq(){ return fRq;};
198  G4double GetRQ(){ return fRQ;};
199 
200  G4double GetRB(){ return fRB;};
201  G4double GetRg(){ return fRg;};
202  G4double GetRG(){ return fRG;};
203 
204  // FqQgG stuff
205 
206  G4complex Pomeron();
207 
208  G4complex Phi13();
209  G4complex Phi14();
210  G4complex Phi23();
211  G4complex Phi24();
212 
218  G4double GetdsdtF123qQgG(G4double q); // sampling ds/dt
220 
221 
222  // F123 stuff
223 
224  G4complex GetAqq();
225  G4complex GetAQQ();
226  G4complex GetAqQ();
227 
228  G4double GetCofS1();
229  G4double GetCofS2();
230  G4double GetCofS3();
232 
238  G4double GetdsdtF123(G4double q); // sampling ds/dt
240 
241  // parameter arrays
242 
243 private:
244 
247  static const G4double theNuclNuclData[18][6];
248  static const G4double thePiKaNuclData[8][6];
250 };
251 
255 
256 
257 
259  G4Nucleus & nucleus)
260 {
261  if( ( projectile.GetDefinition() == G4Proton::Proton() ||
262  projectile.GetDefinition() == G4Neutron::Neutron() ||
263  projectile.GetDefinition() == G4PionPlus::PionPlus() ||
264  projectile.GetDefinition() == G4PionMinus::PionMinus() ||
265  projectile.GetDefinition() == G4KaonPlus::KaonPlus() ||
266  projectile.GetDefinition() == G4KaonMinus::KaonMinus() ) &&
267 
268  nucleus.GetZ_asInt() < 2 ) return true;
269  else return false;
270 }
271 
272 
274 {
275  // masses
276 
277  fMq = 0.36*CLHEP::GeV; // 0.441*GeV; // 0.36*GeV;
278  fMQ = 0.441*CLHEP::GeV;
279  fMff2 = 0.26*CLHEP::GeV*CLHEP::GeV; // 0.25*GeV*GeV; // 0.5*GeV*GeV;
280 
281  fAlpha = 1./3.;
282  fBeta = 1. - fAlpha;
283 
284  fGamma = 1./2.; // 1./3.; //
285  fDelta = 1. - fGamma; // 1./2.;
286 
287  // radii and exp cof
288 
289  fRA = 6.5/CLHEP::GeV; // 7.3/GeV; // 3.25/GeV; // 7./GeV; // 2./GeV; // 1./GeV;
290  fRq = 0.173*fRA; // 2.4/GeV;
291  fRQ = 0.316*fRA; // 1./GeV; // 2./GeV; // 1./GeV;
292  fRB = 6.5/CLHEP::GeV; // 7.3/GeV; // 3.25/GeV; // 7./GeV; // 2./GeV; // 1./GeV;
293  fRg = 0.173*fRA; // 2.4/GeV;
294  fRG = 0.173*fRA; // 1./GeV; // 2./GeV; // 1./GeV;
295 
296  fAlphaP = 0.15/CLHEP::GeV/CLHEP::GeV; // 0.15/GeV/GeV;
297  fLambda = 0.25*fRA*fRA; // 0.25
298  fEta = 0.25*fRB*fRB; // 0.25
299  fImCof = 6.5;
300  fCofF2 = 1.;
301  fCofF3 = 1.;
302 
303  fBq = 0.02; // 0.21; // 1./3.;
304  fBQ = 1. + fBq - 2*std::sqrt(fBq); // 1 - fBq; // 2./3.;
305  fBqQ = std::sqrt(fBq*fBQ);
306 
307  fLambdaFF = 1.5/CLHEP::GeV/CLHEP::GeV; // 0.15/GeV/GeV;
308  fSo = 1.*CLHEP::GeV*CLHEP::GeV;
309  fQcof = 0.009*CLHEP::GeV;
310  fExpSlope = 19.9/CLHEP::GeV/CLHEP::GeV;
311 }
312 
313 
315 //
316 // Set target and projectile masses and calculate mass sum and difference squared for Pcms
317 
319 {
320  G4int i;
321  G4double trMass = 900.*CLHEP::MeV, Tkin;
322  G4double sl, sh, ds, rAl, rAh, drA, rBl, rBh, drB, bql, bqh, dbq, bQl, bQh, dbQ, cIl, cIh, dcI;
323 
324  Tkin = std::sqrt(fMassProj*fMassProj + plab*plab) - fMassProj;
325 
326  G4DynamicParticle* theDynamicParticle = new G4DynamicParticle(fProjectile,
327  G4ParticleMomentum(0.,0.,1.),
328  Tkin);
329  fSigmaTot = fHadrNuclXsc->GetHadronNucleonXscNS( theDynamicParticle, fTarget );
330 
331  delete theDynamicParticle;
332 
333  fSpp = fMassProj*fMassProj + fMassTarg*fMassTarg + 2.*fMassTarg*std::sqrt(plab*plab + fMassProj*fMassProj);
334  fPcms = std::sqrt( (fSpp - fMassSum2)*(fSpp - fMassDif2)/4./fSpp);
335 
336  G4double sCMS = std::sqrt(fSpp);
337 
338  if( fMassProj > trMass ) // p,n,pb on p
339  {
340  this->SetCofF2(1.);
341  this->SetCofF3(1.);
342  fGamma = 1./3.; // 1./3.; //
343  fDelta = 1. - fGamma; // 1./2.;
344 
345  if( sCMS <= theNuclNuclData[0][0]*CLHEP::GeV ) // low edge, as s=2.76754
346  {
347  this->SetRA(theNuclNuclData[0][1]/CLHEP::GeV,0.173,0.316);
348  this->SetRB(theNuclNuclData[0][2]/CLHEP::GeV,0.173,0.316);
349 
350  this->SetBq(theNuclNuclData[0][3]);
351  this->SetBQ(theNuclNuclData[0][4]);
352  this->SetImCof(theNuclNuclData[0][5]);
353 
354  this->SetLambda(0.25*this->GetRA()*this->GetRA());
355  this->SetEta(0.25*this->GetRB()*this->GetRB());
356  }
357  else if( sCMS >= theNuclNuclData[17][0]*CLHEP::GeV ) // high edge, as s=7000 ???
358  {
359  this->SetRA(theNuclNuclData[17][1]/CLHEP::GeV,0.173,0.316);
360  this->SetRB(theNuclNuclData[17][2]/CLHEP::GeV,0.173,0.316);
361 
362  this->SetBq(theNuclNuclData[17][3]);
363  this->SetBQ(theNuclNuclData[17][4]);
364  this->SetImCof(theNuclNuclData[17][5]);
365 
366  this->SetLambda(0.25*this->GetRA()*this->GetRA());
367  this->SetEta(0.25*this->GetRB()*this->GetRB());
368  }
369  else // in approximation between array points
370  {
371  for( i = 0; i < 18; i++ ) if( sCMS <= theNuclNuclData[i][0]*CLHEP::GeV ) break;
372  if( i == 0 ) i++;
373  if( i == 18 ) i--;
374 
375  sl = theNuclNuclData[i-1][0]*CLHEP::GeV;
376  sh = theNuclNuclData[i][0]*CLHEP::GeV;
377  ds = (sCMS - sl)/(sh - sl);
378 
379  rAl = theNuclNuclData[i-1][1]/CLHEP::GeV;
380  rAh = theNuclNuclData[i][1]/CLHEP::GeV;
381  drA = rAh - rAl;
382 
383  rBl = theNuclNuclData[i-1][2]/CLHEP::GeV;
384  rBh = theNuclNuclData[i][2]/CLHEP::GeV;
385  drB = rBh - rBl;
386 
387  bql = theNuclNuclData[i-1][3];
388  bqh = theNuclNuclData[i][3];
389  dbq = bqh - bql;
390 
391  bQl = theNuclNuclData[i-1][4];
392  bQh = theNuclNuclData[i][4];
393  dbQ = bQh - bQl;
394 
395  cIl = theNuclNuclData[i-1][5];
396  cIh = theNuclNuclData[i][5];
397  dcI = cIh - cIl;
398 
399  this->SetRA(rAl+drA*ds,0.173,0.316);
400  this->SetRB(rBl+drB*ds,0.173,0.316);
401 
402  this->SetBq(bql+dbq*ds);
403  this->SetBQ(bQl+dbQ*ds);
404  this->SetImCof(cIl+dcI*ds);
405 
406  this->SetLambda(0.25*this->GetRA()*this->GetRA());
407  this->SetEta(0.25*this->GetRB()*this->GetRB());
408  }
409  }
410  else // pi, K
411  {
412  this->SetCofF2(1.);
413  this->SetCofF3(-1.);
414  fGamma = 1./2.; // 1./3.; //
415  fDelta = 1. - fGamma; // 1./2.;
416 
417  if( sCMS <= thePiKaNuclData[0][0]*CLHEP::GeV ) // low edge, as s=2.76754
418  {
419  this->SetRA(thePiKaNuclData[0][1]/CLHEP::GeV,0.173,0.316);
420  this->SetRB(thePiKaNuclData[0][2]/CLHEP::GeV,0.173,0.173);
421 
422  this->SetBq(thePiKaNuclData[0][3]);
423  this->SetBQ(thePiKaNuclData[0][4]);
424  this->SetImCof(thePiKaNuclData[0][5]);
425 
426  this->SetLambda(0.25*this->GetRA()*this->GetRA());
427  this->SetEta(this->GetRB()*this->GetRB()/6.);
428  }
429  else if( sCMS >= thePiKaNuclData[7][0]*CLHEP::GeV ) // high edge, as s=7000 ???
430  {
431  this->SetRA(thePiKaNuclData[7][1]/CLHEP::GeV,0.173,0.316);
432  this->SetRB(thePiKaNuclData[7][2]/CLHEP::GeV,0.173,0.173);
433 
434  this->SetBq(thePiKaNuclData[7][3]);
435  this->SetBQ(thePiKaNuclData[7][4]);
436  this->SetImCof(thePiKaNuclData[7][5]);
437 
438  this->SetLambda(0.25*this->GetRA()*this->GetRA());
439  this->SetEta(this->GetRB()*this->GetRB()/6.);
440  }
441  else // in approximation between array points
442  {
443  for( i = 0; i < 8; i++ ) if( sCMS <= thePiKaNuclData[i][0]*CLHEP::GeV ) break;
444  if( i == 0 ) i++;
445  if( i == 8 ) i--;
446 
447  sl = thePiKaNuclData[i-1][0]*CLHEP::GeV;
448  sh = thePiKaNuclData[i][0]*CLHEP::GeV;
449  ds = (sCMS - sl)/(sh - sl);
450 
451  rAl = thePiKaNuclData[i-1][1]/CLHEP::GeV;
452  rAh = thePiKaNuclData[i][1]/CLHEP::GeV;
453  drA = rAh - rAl;
454 
455  rBl = thePiKaNuclData[i-1][2]/CLHEP::GeV;
456  rBh = thePiKaNuclData[i][2]/CLHEP::GeV;
457  drB = rBh - rBl;
458 
459  bql = thePiKaNuclData[i-1][3];
460  bqh = thePiKaNuclData[i][3];
461  dbq = bqh - bql;
462 
463  bQl = thePiKaNuclData[i-1][4];
464  bQh = thePiKaNuclData[i][4];
465  dbQ = bQh - bQl;
466 
467  cIl = thePiKaNuclData[i-1][5];
468  cIh = thePiKaNuclData[i][5];
469  dcI = cIh - cIl;
470 
471  this->SetRA(rAl+drA*ds,0.173,0.316);
472  this->SetRB(rBl+drB*ds,0.173,0.173);
473 
474  this->SetBq(bql+dbq*ds);
475  this->SetBQ(bQl+dbQ*ds);
476  this->SetImCof(cIl+dcI*ds);
477 
478  this->SetLambda(0.25*this->GetRA()*this->GetRA());
479  this->SetEta(this->GetRB()*this->GetRB()/6.);
480  }
481  }
482  return;
483 }
484 
486 //
487 // RA for qQ
488 
490 {
491  fRA = rA;
492  fRq = fRA*pq;
493  fRQ = fRA*pQ;
494 }
495 
497 //
498 // RB for gG
499 
501 {
502  fRB = rB;
503  fRg = fRB*pg;
504  fRG = fRB*pG;
505 }
506 
508 //
509 // Returns Pomeron parametrization with Im part modified, *= fImCof
510 
512 {
513  G4double re, im;
514 
515  re = fAlphaP*G4Log(fSpp/fSo);
516  im = -0.5*fAlphaP*fImCof*CLHEP::pi;
517  return G4complex(re,im);
518 }
519 
521 
523 {
524  G4double re = (fRq*fRq + fRg*fRg)/16.;
525  G4complex result(re,0.);
526  result += Pomeron();
527  return result;
528 }
529 
531 
533 {
534  G4double re = (fRq*fRq + fRG*fRG)/16.;
535  G4complex result(re,0.);
536  result += Pomeron();
537  return result;
538 }
539 
541 
543 {
544  G4double re = (fRQ*fRQ + fRg*fRg)/16.;
545  G4complex result(re,0.);
546  result += Pomeron();
547  return result;
548 }
549 
551 
553 {
554  G4double re = (fRQ*fRQ + fRG*fRG)/16.;
555  G4complex result(re,0.);
556  result += Pomeron();
557  return result;
558 }
559 
561 //
562 // F1, case qQ-gG
563 
565 {
566  G4double p = std::sqrt((fSpp - fMassSum2)*(fSpp - fMassDif2)/4./fSpp);
567  G4double k = p/CLHEP::hbarc;
568 
569  G4complex exp13 = fBq*std::exp(-(Phi13() + fBeta*fBeta*fLambda + fDelta*fDelta*fEta)*t);
570  G4complex exp14 = fBq*std::exp(-(Phi14() + fBeta*fBeta*fLambda + fGamma*fGamma*fEta)*t);
571  G4complex exp23 = fBQ*std::exp(-(Phi23() + fAlpha*fAlpha*fLambda + fDelta*fDelta*fEta)*t);
572  G4complex exp24 = fBQ*std::exp(-(Phi24() + fAlpha*fAlpha*fLambda + fGamma*fGamma*fEta)*t);
573 
574  G4complex res = exp13 + exp14 + exp23 + exp24;
575 
576  res *= 0.25*k*fSigmaTot/CLHEP::pi;
577  res *= G4complex(0.,1.);
578 
579  return res;
580 }
581 
583 //
584 //
585 
587 {
588  fSpp = spp;
589  G4double p = std::sqrt((fSpp - fMassSum2)*(fSpp - fMassDif2)/4./fSpp);
590  G4double k = p/CLHEP::hbarc;
591 
592  G4complex exp14 = fBqQ*std::exp(-(Phi14() + fBeta*fBeta*fLambda + fGamma*fGamma*fEta)*t);
593  G4complex exp24 = fBQ*std::exp(-(Phi24() + fAlpha*fAlpha*fLambda + fGamma*fGamma*fEta)*t);
594 
595  G4complex F1 = exp14 + exp24;
596 
597  F1 *= 0.25*k*fSigmaTot/CLHEP::pi;
598  F1 *= G4complex(0.,1.);
599 
600  // 1424
601 
602  G4complex z1424 = -(Phi24() + fAlpha*fLambda)*(Phi24() + fAlpha*fLambda);
603  z1424 /= Phi14() + Phi24() + fLambda;
604  z1424 += Phi24() + fAlpha*fAlpha*fLambda + fGamma*fGamma*fEta;
605 
606 
607  G4complex exp1424 = std::exp(-z1424*t);
608  exp1424 /= Phi14() + Phi24() + fLambda;
609 
610  G4complex F3 = fBqQ*fBQ*exp1424;
611 
612 
613  F3 *= 0.25*k/CLHEP::pi;
614  F3 *= G4complex(0.,1.);
615  F3 *= fSigmaTot*fSigmaTot/(8.*CLHEP::pi*CLHEP::hbarc*CLHEP::hbarc);
616 
617  G4complex F13 = F1 - F3;
618 
619  G4double dsdt = CLHEP::pi/p/p;
620  dsdt *= real(F13)*real(F13) + imag(F13)*imag(F13);
621 
622  return dsdt;
623 }
624 
626 //
627 // dsigma/dt(s,t) F1qQgG
628 
630 {
631  fSpp = spp;
632  G4double p = std::sqrt((fSpp - fMassSum2)*(fSpp - fMassDif2)/4./fSpp);
633 
634  G4complex F1 = GetF1qQgG(t);
635 
636  G4double dsdt = CLHEP::pi/p/p;
637  dsdt *= real(F1)*real(F1) + imag(F1)*imag(F1);
638  return dsdt;
639 }
640 
642 //
643 //
644 
646 {
647  G4double p = std::sqrt((fSpp - fMassSum2)*(fSpp - fMassDif2)/4./fSpp);
648  G4double k = p/CLHEP::hbarc;
649 
651  z1324 /= Phi13() + Phi24() + fLambda + fEta;
652  z1324 += Phi24() + fAlpha*fAlpha*fLambda + fGamma*fGamma*fEta;
653 
654  G4complex exp1324 = std::exp(-z1324*t);
655  exp1324 /= Phi13() + Phi24() + fLambda + fEta;
656 
657  G4complex z1423 = -(Phi23() + fAlpha*fLambda + fDelta*fEta)*(Phi24() + fAlpha*fLambda + fDelta*fEta);;
658  z1423 /= Phi14() + Phi23() + fLambda + fEta;
659  z1423 += Phi23() + fAlpha*fAlpha*fLambda + fDelta*fDelta*fEta;
660 
661  G4complex exp1423 = std::exp(-z1423*t);
662  exp1423 /= Phi14() + Phi23() + fLambda + fEta;
663 
664  G4complex res = exp1324 + exp1423;
665 
666 
667  res *= 0.25*k/CLHEP::pi;
668  res *= G4complex(0.,1.);
669  res *= fBq*fBQ*fSigmaTot*fSigmaTot/(8.*CLHEP::pi*CLHEP::hbarc*CLHEP::hbarc); // or 4. ???
670 
671  return res;
672 }
673 
675 //
676 // dsigma/dt(s,t) F12
677 
679 {
680  fSpp = spp;
681  G4double p = std::sqrt((fSpp - fMassSum2)*(fSpp - fMassDif2)/4./fSpp);
682 
683  G4complex F12 = GetF1qQgG(t) - GetF2qQgG(t);
684 
685  G4double dsdt = CLHEP::pi/p/p;
686  dsdt *= real(F12)*real(F12) + imag(F12)*imag(F12);
687  return dsdt;
688 }
689 
691 //
692 //
693 
695 {
696  G4double p = std::sqrt( (fSpp - fMassSum2)*(fSpp - fMassDif2)/4./fSpp);
697  G4double k = p/CLHEP::hbarc;
698 
699  // 1314
700 
701  G4complex z1314 = -(Phi14() + fGamma*fEta)*(Phi14() + fGamma*fEta);
702  z1314 /= Phi13() + Phi14() + fEta;
703  z1314 += Phi14() + fBeta*fBeta*fLambda + fGamma*fGamma*fEta;
704 
705  G4complex exp1314 = std::exp(-z1314*t);
706  exp1314 /= Phi13() + Phi14() + fEta;
707 
708  // 2324
709 
710  G4complex z2324 = -(Phi24() + fGamma*fEta)*(Phi24() + fGamma*fEta);;
711  z2324 /= Phi24() + Phi23() + fEta;
712  z2324 += Phi24() + fAlpha*fAlpha*fLambda + fGamma*fGamma*fEta;
713 
714  G4complex exp2324 = std::exp(-z2324*t);
715  exp2324 /= Phi24() + Phi23() + fEta;
716 
717  // 1323
718 
719  G4complex z1323 = -(Phi23() + fAlpha*fLambda)*(Phi23() + fAlpha*fLambda);
720  z1323 /= Phi13() + Phi23() + fLambda;
721  z1323 += Phi23() + fAlpha*fAlpha*fLambda + fDelta*fDelta*fEta;
722 
723  G4complex exp1323 = std::exp(-z1323*t);
724  exp1323 /= Phi13() + Phi23() + fLambda;
725 
726  // 1424
727 
728  G4complex z1424 = -(Phi24() + fAlpha*fLambda)*(Phi24() + fAlpha*fLambda);
729  z1424 /= Phi14() + Phi24() + fLambda;
730  z1424 += Phi24() + fAlpha*fAlpha*fLambda + fGamma*fGamma*fEta;
731 
732  G4complex exp1424 = std::exp(-z1424*t);
733  exp1424 /= Phi14() + Phi24() + fLambda;
734 
735  G4complex res = fBq*fBq*exp1314 + fBQ*fBQ*exp2324 + fBq*fBQ*exp1323 + fBq*fBQ*exp1424;
736 
737  res *= 0.25*k/CLHEP::pi;
738  res *= G4complex(0.,1.);
739  res *= fSigmaTot*fSigmaTot/(8.*CLHEP::pi*CLHEP::hbarc*CLHEP::hbarc);
740 
741  return res;
742 }
743 
745 //
746 // dsigma/dt(s,t) F123 sampling ds/dt
747 
749 {
750  G4double p = std::sqrt( (fSpp - fMassSum2)*(fSpp - fMassDif2)/4./fSpp );
751 
752  G4complex F123 = GetF1qQgG(t); // - fCofF2*GetF2qQgG(t) - fCofF3*GetF3qQgG(t);
753  F123 -= fCofF2*GetF2qQgG(t);
754  F123 -= fCofF3*GetF3qQgG(t);
755 
756  G4double dsdt = CLHEP::pi/p/p;
757  dsdt *= real(F123)*real(F123) + imag(F123)*imag(F123);
758  return dsdt;
759 }
760 
762 //
763 // Set fBqQ at a given fBQ=b2 according to the optical theorem,qQ-G
764 
766 {
767  fBQ = b2;
768 
769  G4complex z1424 = G4complex(1./8./CLHEP::pi,0.);
770  z1424 /= Phi14() + Phi24() + fAlpha;
771  G4double c1424 = real(z1424)/(CLHEP::hbarc*CLHEP::hbarc);
772 
773  fBqQ = 1. - fBQ;
774  fBQ /= 1. - fSigmaTot*fBQ*c1424;
775 
776  G4cout<<"fSigmaTot*fBQ*c1424 = "<<fSigmaTot*fBQ*c1424<<G4endl;
777 
778  G4double ratio = fBqQ + fBQ - fSigmaTot*fBqQ*fBQ*c1424;
779  G4cout<<"ratio = "<<ratio<<G4endl;
780 
781  return ;
782 }
783 
785 //
786 // Set fBQ at a given fBq=b according to the optical theorem, F1-F2
787 
789 {
790  fBq = b1;
791 
792  G4complex z1324 = G4complex(1./8./CLHEP::pi,0.);
793  z1324 /= Phi13() + Phi24() + fLambda + fEta;
794  G4double c1324 = real(z1324)/(CLHEP::hbarc*CLHEP::hbarc);
795 
796  G4complex z1423 = G4complex(1./8./CLHEP::pi,0.);
797  z1423 /= Phi14() + Phi23() + fLambda + fEta;
798  G4double c1423 = real(z1423)/(CLHEP::hbarc*CLHEP::hbarc);
799 
800  fBQ = 1. - 2.*fBq;
801  fBQ /= 2. - fSigmaTot*fBq*(c1324+1423);
802 
803  G4double ratio = 2.*(fBq + fBQ) - fSigmaTot*fBq*fBQ*(c1324 + c1423);
804  G4cout<<"ratio = "<<ratio<<G4endl;
805 
806  return ;
807 }
808 
810 //
811 // Set fBQ at a given fBq=b according to the optical theorem, F1-F2-F3,
812 // simplified meson-barion case g=G=q
813 
815 {
816  fBq = b1;
817 
818  G4complex z1324 = fCofF2*G4complex(1./8./CLHEP::pi,0.);
819  z1324 /= Phi13() + Phi24() + fLambda + fEta;
820  G4double c1324 = real(z1324)/(CLHEP::hbarc*CLHEP::hbarc);
821 
822  G4complex z1423 = fCofF2*G4complex(1./8./CLHEP::pi,0.);
823  z1423 /= Phi14() + Phi23() + fLambda + fEta;
824  G4double c1423 = real(z1423)/(CLHEP::hbarc*CLHEP::hbarc);
825 
826  G4complex z1314 = fCofF3*G4complex(1./8./CLHEP::pi,0.);
827  z1314 /= Phi13() + Phi14() + fEta;
828  G4double c1314 = real(z1314)/(CLHEP::hbarc*CLHEP::hbarc);
829 
830  G4complex z2324 = fCofF3*G4complex(1./8./CLHEP::pi,0.);
831  z2324 /= Phi23() + Phi24() + fEta;
832  G4double c2324 = real(z2324)/(CLHEP::hbarc*CLHEP::hbarc);
833 
834  G4complex z1323 = fCofF3*G4complex(1./8./CLHEP::pi,0.);
835  z1323 /= Phi13() + Phi23() + fLambda;
836  G4double c1323 = real(z1323)/(CLHEP::hbarc*CLHEP::hbarc);
837 
838  G4complex z1424 = fCofF3*G4complex(1./8./CLHEP::pi,0.);
839  z1424 /= Phi14() + Phi24() + fLambda;
840  G4double c1424 = real(z1424)/(CLHEP::hbarc*CLHEP::hbarc);
841 
842  G4double A = fSigmaTot*c2324;
843  G4double B = fSigmaTot*fBq*(c1324 + c1423 + c1323 + c1424) - 2.;
844  G4double C = 1. + fSigmaTot*fBq*fBq*c1314 - 2*fBq;
845  G4cout<<"A = "<<A<<"; B = "<<B<<"; C = "<<C<<G4endl;
846  G4cout<<"determinant = "<<B*B-4.*A*C<<G4endl;
847 
848  G4double x1 = ( -B - std::sqrt(B*B-4.*A*C) )/2./A;
849  G4double x2 = ( -B + std::sqrt(B*B-4.*A*C) )/2./A;
850  G4cout<<"x1 = "<<x1<<"; x2 = "<<x2<<G4endl;
851 
852  if( B*B-4.*A*C < 1.e-6 ) fBQ = std::abs(-B/2./A);
853  else if ( B < 0.) fBQ = std::abs( ( -B - std::sqrt(B*B-4.*A*C) )/2./A);
854  else fBQ = std::abs( ( -B + std::sqrt(B*B-4.*A*C) )/2./A);
855 
856  fOptRatio = 2*(fBq+fBQ) - fSigmaTot*fBq*fBQ*(c1324 + c1423 + c1323 + c1424);
857  fOptRatio -= fSigmaTot*fBq*fBq*c1314 + fSigmaTot*c2324*fBQ*fBQ;
858  G4cout<<"BqQ123, fOptRatio = "<<fOptRatio<<G4endl;
859 
860  return ;
861 }
862 
863 
864 
868 
870 {
871  G4double re, im;
872 
873  re = fRq*fRq/8. + fAlphaP*G4Log(fSpp/fSo) + 8.*fLambda/9.;
874  im = -0.5*fAlphaP*fImCof*CLHEP::pi;
875  return G4complex(re,im);
876 }
877 
879 //
880 //
881 
883 {
884  G4double re, im;
885 
886  re = fRQ*fRQ/8. + fAlphaP*G4Log(fSpp/fSo) + 2.*fLambda/9.;
887  im = -0.5*fAlphaP*fImCof*CLHEP::pi;
888  return G4complex(re,im);
889 }
890 
892 //
893 //
894 
896 {
897  G4complex z = 0.5*( GetAqq() + GetAQQ() );
898  return z;
899 }
900 
902 //
903 //
904 
906 {
907  G4complex z = 1./( GetAqQ() + 4.*fLambda/9. );
908  G4double result = real(z);
909  result /= 4.*CLHEP::pi*CLHEP::hbarc*CLHEP::hbarc;
910  result *= fSigmaTot*fCofF2;
911  return result;
912 }
913 
915 //
916 //
917 
919 {
920  G4complex z = 1./( GetAqq() + GetAqQ() - 4.*fLambda/9. );
921  G4double result = real(z);
922  result /= 4.*CLHEP::pi*CLHEP::hbarc*CLHEP::hbarc;
923  result *= fSigmaTot*fCofF3;
924  return result;
925 }
926 
928 //
929 //
930 
932 {
933  G4complex z = 1./( GetAQQ() + GetAqQ() + 2.*fLambda/9. );
934  G4double result = real(z);
935  result /= 4.*CLHEP::pi*CLHEP::hbarc*CLHEP::hbarc;
936  result *= fSigmaTot*fCofF3;
937  return result;
938 }
939 
941 //
942 //
943 
945 {
946  return fOptRatio;
947  // G4double sqrtBqBQ = std::sqrt(fBq*fBQ);
948  // G4double result = fBq + fBQ + 2.*sqrtBqBQ - 1.;
949  // result /= sqrtBqBQ*( GetCofS1()*sqrtBqBQ + GetCofS2()*fBq + GetCofS3()*fBQ );
950  // return result;
951 }
952 
953 
954 
956 //
957 // Set fBQ at a given fBq=b according to the optical theorem
958 
960 {
961  fBq = b1;
962  G4double s1 = GetCofS1();
963  G4double s2 = GetCofS2();
964  G4double s3 = GetCofS3();
965  G4double sqrtBq = std::sqrt(fBq);
966 
967  // cofs of the fBQ 3rd equation
968 
969  G4double a = s3*sqrtBq;
970  G4double b = s1*fBq - 1.;
971  G4double c = (s2*fBq - 2.)*sqrtBq;
972  G4double d = 1. - fBq;
973 
974  // cofs of the incomplete 3rd equation
975 
976  G4double p = c/a;
977  p -= b*b/a/a/3.;
978  G4double q = d/a;
979  q -= b*c/a/a/3.;
980  q += 2*b*b*b/a/a/a/27.;
981 
982  // cofs for the incomplete colutions
983 
984  G4double D = p*p*p/3./3./3.;
985  D += q*q/2./2.;
986  G4complex A1 = G4complex(- q/2., std::sqrt(-D) );
987  G4complex A = std::pow(A1,1./3.);
988 
989  G4complex B1 = G4complex(- q/2., -std::sqrt(-D) );
990  G4complex B = std::pow(B1,1./3.);
991 
992  // roots of the incomplete 3rd equation
993 
994  G4complex y1 = A + B;
995  G4complex y2 = -0.5*(A + B) + 0.5*std::sqrt(3.)*(A - B)*G4complex(0.,1.);
996  G4complex y3 = -0.5*(A + B) - 0.5*std::sqrt(3.)*(A - B)*G4complex(0.,1.);
997 
998  G4complex x1 = y1 - b/a/3.;
999  G4complex x2 = y2 - b/a/3.;
1000  G4complex x3 = y3 - b/a/3.;
1001 
1002  G4cout<<"re_x1 = "<<real(x1)<<"; re_x2 = "<<real(x2)<<"; re_x3 = "<<real(x3)<<G4endl;
1003  G4cout<<"im_x1 = "<<imag(x1)<<"; im_x2 = "<<imag(x2)<<"; im_x3 = "<<imag(x3)<<G4endl;
1004 
1005  G4double r1 = real(x1)*real(x1);
1006  G4double r2 = real(x2)*real(x2);
1007  G4double r3 = real(x3)*real(x3);
1008 
1009  if( r1 <= 1. && r1 >= 0. ) fBQ = r1;
1010  else if( r2 <= 1. && r2 >= 0. ) fBQ = r2;
1011  else if( r3 <= 1. && r3 >= 0. ) fBQ = r3;
1012  else fBQ = 1.;
1013  // fBQ = real(x3)*real(x3);
1014  G4double sqrtBqBQ = std::sqrt(fBq*fBQ);
1015  fOptRatio = fBq + fBQ + 2.*sqrtBqBQ - 1.;
1016  fOptRatio /= sqrtBqBQ*( GetCofS1()*sqrtBqBQ + GetCofS2()*fBq + GetCofS3()*fBQ );
1017  G4cout<<"F123, fOptRatio = "<<fOptRatio<<G4endl;
1018 
1019  return ;
1020 }
1021 
1023 //
1024 //
1025 
1027 {
1028  G4double p = std::sqrt(0.25*fSpp - CLHEP::proton_mass_c2*CLHEP::proton_mass_c2);
1029  G4double k = p/CLHEP::hbarc;
1030  G4complex exp1 = fBq*std::exp(-GetAqq()*t);
1031  G4complex exp2 = fBQ*std::exp(-GetAQQ()*t);
1032  G4complex exp3 = 2.*std::sqrt(fBq*fBQ)*std::exp(-GetAqQ()*t);
1033 
1034  G4complex res = exp1 + exp2 + exp3;
1035  res *= 0.25*k*fSigmaTot/CLHEP::pi;
1036  res *= G4complex(0.,1.);
1037  return res;
1038 }
1039 
1041 //
1042 // dsigma/dt(s,t) F1
1043 
1045 {
1046  fSpp = spp;
1047  G4double p = std::sqrt(0.25*spp - CLHEP::proton_mass_c2*CLHEP::proton_mass_c2);
1048  G4complex F1 = GetF1(t);
1049 
1050  G4double dsdt = CLHEP::pi/p/p;
1051  dsdt *= real(F1)*real(F1) + imag(F1)*imag(F1);
1052  return dsdt;
1053 }
1054 
1056 //
1057 //
1058 
1060 {
1061  G4double p = std::sqrt(0.25*fSpp - CLHEP::proton_mass_c2*CLHEP::proton_mass_c2);
1062  G4double k = p/CLHEP::hbarc;
1063  G4complex z1 = GetAqq()*GetAQQ() - 16.*fLambda*fLambda/81.;
1064  z1 /= 2.*(GetAqQ() + 4.*fLambda/9.);
1065  G4complex exp1 = std::exp(-z1*t);
1066 
1067  G4complex z2 = 0.5*( GetAqQ() - 4.*fLambda/9.);
1068 
1069  G4complex exp2 = std::exp(-z2*t);
1070 
1071  G4complex res = exp1 + exp2;
1072 
1073  G4complex z3 = GetAqQ() + 4.*fLambda/9.;
1074 
1075  res *= 0.25*k/CLHEP::pi;
1076  res *= G4complex(0.,1.);
1077  res /= z3;
1078  res *= fBq*fBQ*fSigmaTot*fSigmaTot/(8.*CLHEP::pi*CLHEP::hbarc*CLHEP::hbarc);
1079 
1080  return res;
1081 }
1082 
1084 //
1085 // dsigma/dt(s,t) F12
1086 
1088 {
1089  fSpp = spp;
1090  G4double p = std::sqrt(0.25*spp - CLHEP::proton_mass_c2*CLHEP::proton_mass_c2);
1091  G4complex F1 = GetF1(t) - GetF2(t);
1092 
1093  G4double dsdt = CLHEP::pi/p/p;
1094  dsdt *= real(F1)*real(F1) + imag(F1)*imag(F1);
1095  return dsdt;
1096 }
1097 
1099 //
1100 //
1101 
1103 {
1104  G4double p = std::sqrt(0.25*fSpp - CLHEP::proton_mass_c2*CLHEP::proton_mass_c2);
1105  G4double k = p/CLHEP::hbarc;
1106  G4complex z1 = GetAqq()*GetAqQ() - 4.*fLambda*fLambda/81.;
1107  z1 /= GetAqq() + GetAqQ() - 4.*fLambda/9.;
1108 
1109  G4complex exp1 = std::exp(-z1*t)*fBq/(GetAqq() + GetAqQ() - 4.*fLambda/9.);
1110 
1111  G4complex z2 = GetAqQ()*GetAQQ() - 1.*fLambda*fLambda/81.;
1112  z2 /= GetAQQ() + GetAqQ() + 2.*fLambda/9.;
1113 
1114  G4complex exp2 = std::exp(-z2*t)*fBQ/(GetAQQ() + GetAqQ() + 2.*fLambda/9.);
1115 
1116  G4complex res = exp1 + exp2;
1117 
1118 
1119  res *= 0.25*k/CLHEP::pi;
1120  res *= G4complex(0.,1.);
1121  res *= std::sqrt(fBq*fBQ)*fSigmaTot*fSigmaTot/(4.*CLHEP::pi*CLHEP::hbarc*CLHEP::hbarc);
1122 
1123  return res;
1124 }
1125 
1127 //
1128 // dsigma/dt(s,t) F123, sampling ds/dt
1129 
1131 {
1132  G4double p = std::sqrt(0.25*fSpp - CLHEP::proton_mass_c2*CLHEP::proton_mass_c2);
1133  G4complex F1 = GetF1(t);
1134  F1 -= fCofF2*GetF2(t);
1135  F1 -= fCofF3*GetF3(t);
1136  G4double dsdt = CLHEP::pi/p/p;
1137  dsdt *= real(F1)*real(F1) + imag(F1)*imag(F1);
1138  return dsdt;
1139 }
1140 
1142 //
1143 // dsigma/dt(s,t) F123
1144 
1146 {
1147  fSpp = spp;
1148  G4double p = std::sqrt(0.25*spp - CLHEP::proton_mass_c2*CLHEP::proton_mass_c2);
1149 
1150  // qQ-ds/dt
1151 
1152  G4complex F1 = GetF1(t) - fCofF2*GetF2(t) - fCofF3*GetF3(t);
1153 
1154  G4double dsdt = CLHEP::pi/p/p;
1155  dsdt *= real(F1)*real(F1) + imag(F1)*imag(F1);
1156 
1157  // exponent ds/dt
1158 
1159  G4complex F10 = GetF1(0.) - fCofF2*GetF2(0.) - fCofF3*GetF3(0.);
1160 
1161  fRhoReIm = real(F10)/imag(F10);
1162 
1163  G4double dsdt0 = CLHEP::pi/p/p;
1164  dsdt0 *= real(F10)*real(F10) + imag(F10)*imag(F10);
1165 
1166  dsdt0 *= G4Exp(-fExpSlope*t);
1167 
1168  G4double ratio = dsdt/dsdt0;
1169 
1170  return ratio;
1171 }
1172 
1173 
1174 //
1175 //
1177 
1178 #endif
1179 
1180 
1181 
1182 
1183 
1184 
G4double fRQ
Definition: G4hhElastic.hh:132
G4complex Phi14()
Definition: G4hhElastic.hh:532
void SetParameters()
Definition: G4hhElastic.hh:273
G4double fAlpha
Definition: G4hhElastic.hh:134
void SetParametersCMS(G4double plab)
Definition: G4hhElastic.hh:318
static const double MeV
Definition: G4SIunits.hh:211
G4double GetCofS1()
Definition: G4hhElastic.hh:905
G4double fRq
Definition: G4hhElastic.hh:133
void SetRA(G4double rn, G4double pq, G4double pQ)
Definition: G4hhElastic.hh:489
void SetBqQ(G4double b)
Definition: G4hhElastic.hh:176
G4ParticleDefinition * thePionMinus
Definition: G4hhElastic.hh:101
void CalculateBQ(G4double b)
Definition: G4hhElastic.hh:959
void SetSigmaTot(G4double stot)
Definition: G4hhElastic.hh:166
virtual ~G4hhElastic()
Definition: G4hhElastic.cc:206
G4double GetRg()
Definition: G4hhElastic.hh:201
G4complex Phi13()
Definition: G4hhElastic.hh:522
G4complex Pomeron()
Definition: G4hhElastic.hh:511
G4double fImCof
Definition: G4hhElastic.hh:147
#define F13
G4PhysicsLogVector * fEnergyVector
Definition: G4hhElastic.hh:113
G4double fLambda
Definition: G4hhElastic.hh:145
void CalculateBqQ13(G4double b)
Definition: G4hhElastic.hh:765
G4double GetRq()
Definition: G4hhElastic.hh:197
G4complex GetAqQ()
Definition: G4hhElastic.hh:895
G4double GetCofS3()
Definition: G4hhElastic.hh:931
G4complex GetF3qQgG(G4double qp)
Definition: G4hhElastic.hh:694
G4double GetBqQ()
Definition: G4hhElastic.hh:173
G4double fOldTkin
Definition: G4hhElastic.hh:246
G4double z
Definition: TRTMaterials.hh:39
G4complex GetF2qQgG(G4double qp)
Definition: G4hhElastic.hh:645
G4double fEta
Definition: G4hhElastic.hh:146
G4complex Phi23()
Definition: G4hhElastic.hh:542
void SetBQ(G4double b)
Definition: G4hhElastic.hh:175
static const G4double theNuclNuclData[18][6]
Definition: G4hhElastic.hh:247
G4double lowEnergyLimitQ
Definition: G4hhElastic.hh:106
G4double GetExpRatioF123(G4double s, G4double q)
G4double lowEnergyRecoilLimit
Definition: G4hhElastic.hh:104
G4double GetBq()
Definition: G4hhElastic.hh:171
G4double a
Definition: TRTMaterials.hh:39
double C(double temp)
G4double fRA
Definition: G4hhElastic.hh:131
G4double fMassProj
Definition: G4hhElastic.hh:126
double B(double temperature)
void SetSpp(G4double spp)
Definition: G4hhElastic.hh:167
void SetImCof(G4double a)
Definition: G4hhElastic.hh:187
G4double GetRG()
Definition: G4hhElastic.hh:202
G4double fSpp
Definition: G4hhElastic.hh:159
int G4int
Definition: G4Types.hh:78
static const double L
Definition: G4SIunits.hh:123
G4double lowestEnergyLimit
Definition: G4hhElastic.hh:107
void CalculateBqQ123(G4double b)
Definition: G4hhElastic.hh:814
G4ParticleDefinition * theProton
Definition: G4hhElastic.hh:98
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double fMq
Definition: G4hhElastic.hh:123
G4complex GetF1qQgG(G4double qp)
Definition: G4hhElastic.hh:564
G4double GetdsdtF12(G4double s, G4double q)
static const double s
Definition: G4SIunits.hh:168
G4double fGamma
Definition: G4hhElastic.hh:140
void BuildTableT(G4ParticleDefinition *target, G4ParticleDefinition *projectile)
Definition: G4hhElastic.cc:255
virtual G4bool IsApplicable(const G4HadProjectile &, G4Nucleus &)
Definition: G4hhElastic.hh:258
G4double GetBQ()
Definition: G4hhElastic.hh:172
std::vector< G4PhysicsTable * > fBankT
Definition: G4hhElastic.hh:115
void BuildTableTest(G4ParticleDefinition *target, G4ParticleDefinition *projectile, G4double plab)
Definition: G4hhElastic.cc:521
std::complex< G4double > G4complex
Definition: G4Types.hh:81
#define position
Definition: xmlparse.cc:622
G4complex GetAQQ()
Definition: G4hhElastic.hh:882
G4GLOB_DLL std::ostream G4cout
void SetCofF3(G4double f)
Definition: G4hhElastic.hh:192
G4double fCofF2
Definition: G4hhElastic.hh:148
double A(double temperature)
const G4ParticleDefinition * GetDefinition() const
G4double fBeta
Definition: G4hhElastic.hh:135
G4double GetTransfer(G4int iMomentum, G4int iTransfer, G4double position)
Definition: G4hhElastic.cc:630
G4double lowEnergyLimitHE
Definition: G4hhElastic.hh:105
void SetAlphaP(G4double a)
Definition: G4hhElastic.hh:186
bool G4bool
Definition: G4Types.hh:79
G4double GetCofF2()
Definition: G4hhElastic.hh:193
G4double fOptRatio
Definition: G4hhElastic.hh:158
void SetEta(G4double E)
Definition: G4hhElastic.hh:190
G4double fRG
Definition: G4hhElastic.hh:138
G4double fPcms
Definition: G4hhElastic.hh:160
void SetCofF2(G4double f)
Definition: G4hhElastic.hh:191
G4double GetdsdtF1(G4double s, G4double q)
G4HadronNucleonXsc * fHadrNuclXsc
Definition: G4hhElastic.hh:249
G4double fQcof
Definition: G4hhElastic.hh:161
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetCofS2()
Definition: G4hhElastic.hh:918
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
G4double GetdsdtF123qQgG(G4double q)
Definition: G4hhElastic.hh:748
static const double GeV
Definition: G4SIunits.hh:214
static const G4double b2
G4int fEnergyBin
Definition: G4hhElastic.hh:110
G4double GetdsdtF123(G4double q)
G4ParticleDefinition * fProjectile
Definition: G4hhElastic.hh:97
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
#define F10
G4complex GetF2(G4double qp)
G4double fMassSum2
Definition: G4hhElastic.hh:127
G4double GetSpp()
Definition: G4hhElastic.hh:168
void SetBq(G4double b)
Definition: G4hhElastic.hh:174
G4double fSo
Definition: G4hhElastic.hh:152
G4complex GetAqq()
Definition: G4hhElastic.hh:869
void SetRB(G4double rn, G4double pq, G4double pQ)
Definition: G4hhElastic.hh:500
G4double GetRQ()
Definition: G4hhElastic.hh:198
G4double SampleInvariantT(const G4ParticleDefinition *p, G4double plab, G4int, G4int)
Definition: G4hhElastic.cc:323
G4PhysicsTable * fTableT
Definition: G4hhElastic.hh:114
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double fBqQ
Definition: G4hhElastic.hh:157
G4double fAlphaP
Definition: G4hhElastic.hh:143
G4double fRg
Definition: G4hhElastic.hh:139
static const double pi
Definition: G4SIunits.hh:74
G4double fDelta
Definition: G4hhElastic.hh:141
G4double fSigmaTot
Definition: G4hhElastic.hh:154
G4double GetdsdtF1qQgG(G4double s, G4double q)
Definition: G4hhElastic.hh:629
G4double SampleBisectionalT(const G4ParticleDefinition *p, G4double plab)
Definition: G4hhElastic.cc:440
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
G4double fBQ
Definition: G4hhElastic.hh:156
G4double fExpSlope
Definition: G4hhElastic.hh:151
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double GetOpticalRatio()
Definition: G4hhElastic.hh:944
G4double fMassDif2
Definition: G4hhElastic.hh:128
double D(double temp)
static const G4double b1
G4complex GetF3(G4double qp)
G4double fLambdaFF
Definition: G4hhElastic.hh:144
void CalculateBqQ12(G4double b)
Definition: G4hhElastic.hh:788
G4double fMassTarg
Definition: G4hhElastic.hh:125
G4double GetdsdtF13qQG(G4double s, G4double q)
Definition: G4hhElastic.hh:586
#define G4endl
Definition: G4ios.hh:61
G4ParticleDefinition * thePionPlus
Definition: G4hhElastic.hh:100
G4complex Phi24()
Definition: G4hhElastic.hh:552
G4double fMff2
Definition: G4hhElastic.hh:121
#define F12
G4double GetdsdtF12qQgG(G4double s, G4double q)
Definition: G4hhElastic.hh:678
static const G4double thePiKaNuclData[8][6]
Definition: G4hhElastic.hh:248
void Initialise()
Definition: G4hhElastic.cc:231
G4double GetCofF3()
Definition: G4hhElastic.hh:194
G4double fBq
Definition: G4hhElastic.hh:155
G4double fMQ
Definition: G4hhElastic.hh:122
double G4double
Definition: G4Types.hh:76
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
G4double GetImCof()
Definition: G4hhElastic.hh:188
void SetLambda(G4double L)
Definition: G4hhElastic.hh:189
G4double fCofF3
Definition: G4hhElastic.hh:149
G4ThreeVector G4ParticleMomentum
G4ParticleDefinition * fTarget
Definition: G4hhElastic.hh:96
G4double GetRhoReIm()
Definition: G4hhElastic.hh:177
G4ParticleDefinition * theNeutron
Definition: G4hhElastic.hh:99
G4double plabLowLimit
Definition: G4hhElastic.hh:108
G4double fRB
Definition: G4hhElastic.hh:137
G4double GetRA()
Definition: G4hhElastic.hh:196
G4double SampleTest(G4double tMin)
Definition: G4hhElastic.cc:594
G4complex GetF1(G4double qp)
G4double GetRB()
Definition: G4hhElastic.hh:200
G4double fRhoReIm
Definition: G4hhElastic.hh:150