Geant4  10.02.p03
G4hhElastic Class Reference

#include <G4hhElastic.hh>

Inheritance diagram for G4hhElastic:
Collaboration diagram for G4hhElastic:

Public Member Functions

 G4hhElastic ()
 
 G4hhElastic (G4ParticleDefinition *target, G4ParticleDefinition *projectile, G4double plab)
 
 G4hhElastic (G4ParticleDefinition *target, G4ParticleDefinition *projectile)
 
virtual ~G4hhElastic ()
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
void Initialise ()
 
void BuildTableT (G4ParticleDefinition *target, G4ParticleDefinition *projectile)
 
void BuildTableTest (G4ParticleDefinition *target, G4ParticleDefinition *projectile, G4double plab)
 
G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int, G4int)
 
G4double SampleBisectionalT (const G4ParticleDefinition *p, G4double plab)
 
G4double SampleTest (G4double tMin)
 
G4double GetTransfer (G4int iMomentum, G4int iTransfer, G4double position)
 
void SetParameters ()
 
void SetSigmaTot (G4double stot)
 
void SetSpp (G4double spp)
 
G4double GetSpp ()
 
void SetParametersCMS (G4double plab)
 
G4double GetBq ()
 
G4double GetBQ ()
 
G4double GetBqQ ()
 
void SetBq (G4double b)
 
void SetBQ (G4double b)
 
void SetBqQ (G4double b)
 
G4double GetRhoReIm ()
 
void CalculateBQ (G4double b)
 
void CalculateBqQ13 (G4double b)
 
void CalculateBqQ12 (G4double b)
 
void CalculateBqQ123 (G4double b)
 
void SetRA (G4double rn, G4double pq, G4double pQ)
 
void SetRB (G4double rn, G4double pq, G4double pQ)
 
void SetAlphaP (G4double a)
 
void SetImCof (G4double a)
 
G4double GetImCof ()
 
void SetLambda (G4double L)
 
void SetEta (G4double E)
 
void SetCofF2 (G4double f)
 
void SetCofF3 (G4double f)
 
G4double GetCofF2 ()
 
G4double GetCofF3 ()
 
G4double GetRA ()
 
G4double GetRq ()
 
G4double GetRQ ()
 
G4double GetRB ()
 
G4double GetRg ()
 
G4double GetRG ()
 
G4complex Pomeron ()
 
G4complex Phi13 ()
 
G4complex Phi14 ()
 
G4complex Phi23 ()
 
G4complex Phi24 ()
 
G4complex GetF1qQgG (G4double qp)
 
G4double GetdsdtF1qQgG (G4double s, G4double q)
 
G4complex GetF2qQgG (G4double qp)
 
G4double GetdsdtF12qQgG (G4double s, G4double q)
 
G4complex GetF3qQgG (G4double qp)
 
G4double GetdsdtF123qQgG (G4double q)
 
G4double GetdsdtF13qQG (G4double s, G4double q)
 
G4complex GetAqq ()
 
G4complex GetAQQ ()
 
G4complex GetAqQ ()
 
G4double GetCofS1 ()
 
G4double GetCofS2 ()
 
G4double GetCofS3 ()
 
G4double GetOpticalRatio ()
 
G4complex GetF1 (G4double qp)
 
G4double GetdsdtF1 (G4double s, G4double q)
 
G4complex GetF2 (G4double qp)
 
G4double GetdsdtF12 (G4double s, G4double q)
 
G4complex GetF3 (G4double qp)
 
G4double GetdsdtF123 (G4double q)
 
G4double GetExpRatioF123 (G4double s, G4double q)
 
- Public Member Functions inherited from G4HadronElastic
 G4HadronElastic (const G4String &name="hElasticLHEP")
 
virtual ~G4HadronElastic ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
void SetLowestEnergyLimit (G4double value)
 
G4double LowestEnergyLimit () const
 
G4double ComputeMomentumCMS (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual void ModelDescription (std::ostream &) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
const G4HadronicInteractionGetMyPointer () const
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 

Private Attributes

G4ParticleDefinitionfTarget
 
G4ParticleDefinitionfProjectile
 
G4ParticleDefinitiontheProton
 
G4ParticleDefinitiontheNeutron
 
G4ParticleDefinitionthePionPlus
 
G4ParticleDefinitionthePionMinus
 
G4double lowEnergyRecoilLimit
 
G4double lowEnergyLimitHE
 
G4double lowEnergyLimitQ
 
G4double lowestEnergyLimit
 
G4double plabLowLimit
 
G4int fEnergyBin
 
G4int fBinT
 
G4PhysicsLogVectorfEnergyVector
 
G4PhysicsTablefTableT
 
std::vector< G4PhysicsTable * > fBankT
 
G4double fMff2
 
G4double fMQ
 
G4double fMq
 
G4double fMassTarg
 
G4double fMassProj
 
G4double fMassSum2
 
G4double fMassDif2
 
G4double fRA
 
G4double fRQ
 
G4double fRq
 
G4double fAlpha
 
G4double fBeta
 
G4double fRB
 
G4double fRG
 
G4double fRg
 
G4double fGamma
 
G4double fDelta
 
G4double fAlphaP
 
G4double fLambdaFF
 
G4double fLambda
 
G4double fEta
 
G4double fImCof
 
G4double fCofF2
 
G4double fCofF3
 
G4double fRhoReIm
 
G4double fExpSlope
 
G4double fSo
 
G4double fSigmaTot
 
G4double fBq
 
G4double fBQ
 
G4double fBqQ
 
G4double fOptRatio
 
G4double fSpp
 
G4double fPcms
 
G4double fQcof
 
G4int fInTkin
 
G4double fOldTkin
 
G4HadronNucleonXscfHadrNuclXsc
 

Static Private Attributes

static const G4double theNuclNuclData [18][6]
 
static const G4double thePiKaNuclData [8][6]
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 63 of file G4hhElastic.hh.

Constructor & Destructor Documentation

◆ G4hhElastic() [1/3]

G4hhElastic::G4hhElastic ( )

Definition at line 73 of file G4hhElastic.cc.

74  : G4HadronElastic("HadrHadrElastic")
75 {
76  SetMinEnergy( 1.*GeV );
77  SetMaxEnergy( 10000.*TeV );
78  verboseLevel = 0;
79  lowEnergyRecoilLimit = 100.*keV;
80  lowEnergyLimitQ = 0.0*GeV;
81  lowEnergyLimitHE = 0.0*GeV;
82  lowestEnergyLimit= 0.0*keV;
83  plabLowLimit = 20.0*MeV;
84 
86  fInTkin=0;
91 
93  fProjectile = 0;
95 
96  fEnergyBin = 200;
97  fBinT = 514; // 514; // 500; // 200;
98 
100 
101  fTableT = 0;
102  fOldTkin = 0.;
103  SetParameters();
104 
105  Initialise();
106 }
void SetParameters()
Definition: G4hhElastic.hh:273
static const double MeV
Definition: G4SIunits.hh:211
G4ParticleDefinition * thePionMinus
Definition: G4hhElastic.hh:101
G4PhysicsLogVector * fEnergyVector
Definition: G4hhElastic.hh:113
G4HadronElastic(const G4String &name="hElasticLHEP")
G4double fOldTkin
Definition: G4hhElastic.hh:246
G4double lowEnergyLimitQ
Definition: G4hhElastic.hh:106
G4double lowEnergyRecoilLimit
Definition: G4hhElastic.hh:104
G4double fSpp
Definition: G4hhElastic.hh:159
G4double lowestEnergyLimit
Definition: G4hhElastic.hh:107
G4ParticleDefinition * theProton
Definition: G4hhElastic.hh:98
void SetMinEnergy(G4double anEnergy)
G4double lowEnergyLimitHE
Definition: G4hhElastic.hh:105
G4double fOptRatio
Definition: G4hhElastic.hh:158
G4double fPcms
Definition: G4hhElastic.hh:160
G4HadronNucleonXsc * fHadrNuclXsc
Definition: G4hhElastic.hh:249
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static const double GeV
Definition: G4SIunits.hh:214
G4int fEnergyBin
Definition: G4hhElastic.hh:110
G4ParticleDefinition * fProjectile
Definition: G4hhElastic.hh:97
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4PhysicsTable * fTableT
Definition: G4hhElastic.hh:114
G4double fSigmaTot
Definition: G4hhElastic.hh:154
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
void SetMaxEnergy(const G4double anEnergy)
G4ParticleDefinition * thePionPlus
Definition: G4hhElastic.hh:100
static const double TeV
Definition: G4SIunits.hh:215
static const double keV
Definition: G4SIunits.hh:213
void Initialise()
Definition: G4hhElastic.cc:231
G4ParticleDefinition * fTarget
Definition: G4hhElastic.hh:96
G4ParticleDefinition * theNeutron
Definition: G4hhElastic.hh:99
G4double plabLowLimit
Definition: G4hhElastic.hh:108
G4double fRhoReIm
Definition: G4hhElastic.hh:150
Here is the call graph for this function:

◆ G4hhElastic() [2/3]

G4hhElastic::G4hhElastic ( G4ParticleDefinition target,
G4ParticleDefinition projectile,
G4double  plab 
)

Definition at line 114 of file G4hhElastic.cc.

115  : G4HadronElastic("HadrHadrElastic")
116 {
117  SetMinEnergy( 1.*GeV );
118  SetMaxEnergy( 10000.*TeV );
119  verboseLevel = 0;
120  lowEnergyRecoilLimit = 100.*keV;
121  lowEnergyLimitQ = 0.0*GeV;
122  lowEnergyLimitHE = 0.0*GeV;
123  lowestEnergyLimit = 0.0*keV;
124  plabLowLimit = 20.0*MeV;
125 
127  fInTkin=0;
132 
133  fTarget = target;
134  fProjectile = projectile;
140 
141  fEnergyBin = 200;
142  fBinT = 514; // 200;
143 
145  fTableT = 0;
146  fOldTkin = 0.;
147 
148 
149  SetParameters();
150  SetParametersCMS( plab);
151 }
void SetParameters()
Definition: G4hhElastic.hh:273
void SetParametersCMS(G4double plab)
Definition: G4hhElastic.hh:318
static const double MeV
Definition: G4SIunits.hh:211
G4ParticleDefinition * thePionMinus
Definition: G4hhElastic.hh:101
G4PhysicsLogVector * fEnergyVector
Definition: G4hhElastic.hh:113
G4HadronElastic(const G4String &name="hElasticLHEP")
G4double fOldTkin
Definition: G4hhElastic.hh:246
G4double lowEnergyLimitQ
Definition: G4hhElastic.hh:106
G4double lowEnergyRecoilLimit
Definition: G4hhElastic.hh:104
G4double fMassProj
Definition: G4hhElastic.hh:126
G4double fSpp
Definition: G4hhElastic.hh:159
G4double lowestEnergyLimit
Definition: G4hhElastic.hh:107
G4ParticleDefinition * theProton
Definition: G4hhElastic.hh:98
void SetMinEnergy(G4double anEnergy)
G4double lowEnergyLimitHE
Definition: G4hhElastic.hh:105
G4double fOptRatio
Definition: G4hhElastic.hh:158
G4double fPcms
Definition: G4hhElastic.hh:160
G4HadronNucleonXsc * fHadrNuclXsc
Definition: G4hhElastic.hh:249
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static const double GeV
Definition: G4SIunits.hh:214
G4int fEnergyBin
Definition: G4hhElastic.hh:110
G4ParticleDefinition * fProjectile
Definition: G4hhElastic.hh:97
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double fMassSum2
Definition: G4hhElastic.hh:127
G4PhysicsTable * fTableT
Definition: G4hhElastic.hh:114
G4double fSigmaTot
Definition: G4hhElastic.hh:154
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
cout<< "-> Edep in the target
Definition: analysis.C:54
G4double fMassDif2
Definition: G4hhElastic.hh:128
void SetMaxEnergy(const G4double anEnergy)
G4double fMassTarg
Definition: G4hhElastic.hh:125
G4ParticleDefinition * thePionPlus
Definition: G4hhElastic.hh:100
static const double TeV
Definition: G4SIunits.hh:215
static const double keV
Definition: G4SIunits.hh:213
G4ParticleDefinition * fTarget
Definition: G4hhElastic.hh:96
G4ParticleDefinition * theNeutron
Definition: G4hhElastic.hh:99
G4double plabLowLimit
Definition: G4hhElastic.hh:108
G4double fRhoReIm
Definition: G4hhElastic.hh:150
Here is the call graph for this function:

◆ G4hhElastic() [3/3]

G4hhElastic::G4hhElastic ( G4ParticleDefinition target,
G4ParticleDefinition projectile 
)

Definition at line 159 of file G4hhElastic.cc.

160  : G4HadronElastic("HadrHadrElastic")
161 {
162  SetMinEnergy( 1.*GeV );
163  SetMaxEnergy( 10000.*TeV );
164  verboseLevel = 0;
165  lowEnergyRecoilLimit = 100.*keV;
166  lowEnergyLimitQ = 0.0*GeV;
167  lowEnergyLimitHE = 0.0*GeV;
168  lowestEnergyLimit= 0.0*keV;
169  plabLowLimit = 20.0*MeV;
170 
172  fInTkin=0;
173 
174  fTarget = target; // later vmg
175  fProjectile = projectile;
180 
181  fTarget = G4Proton::Proton(); // later vmg
182  fProjectile = 0;
188 
189  fEnergyBin = 200;
190  fBinT = 514; // 514; // 500; // 200;
191 
193 
194  fTableT = 0;
195  fOldTkin = 0.;
196 
197  SetParameters();
198 }
void SetParameters()
Definition: G4hhElastic.hh:273
static const double MeV
Definition: G4SIunits.hh:211
G4ParticleDefinition * thePionMinus
Definition: G4hhElastic.hh:101
G4PhysicsLogVector * fEnergyVector
Definition: G4hhElastic.hh:113
G4HadronElastic(const G4String &name="hElasticLHEP")
G4double fOldTkin
Definition: G4hhElastic.hh:246
G4double lowEnergyLimitQ
Definition: G4hhElastic.hh:106
G4double lowEnergyRecoilLimit
Definition: G4hhElastic.hh:104
G4double fMassProj
Definition: G4hhElastic.hh:126
G4double fSpp
Definition: G4hhElastic.hh:159
G4double lowestEnergyLimit
Definition: G4hhElastic.hh:107
G4ParticleDefinition * theProton
Definition: G4hhElastic.hh:98
void SetMinEnergy(G4double anEnergy)
G4double lowEnergyLimitHE
Definition: G4hhElastic.hh:105
G4double fOptRatio
Definition: G4hhElastic.hh:158
G4double fPcms
Definition: G4hhElastic.hh:160
G4HadronNucleonXsc * fHadrNuclXsc
Definition: G4hhElastic.hh:249
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static const double GeV
Definition: G4SIunits.hh:214
G4int fEnergyBin
Definition: G4hhElastic.hh:110
G4ParticleDefinition * fProjectile
Definition: G4hhElastic.hh:97
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double fMassSum2
Definition: G4hhElastic.hh:127
G4PhysicsTable * fTableT
Definition: G4hhElastic.hh:114
G4double fSigmaTot
Definition: G4hhElastic.hh:154
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
cout<< "-> Edep in the target
Definition: analysis.C:54
G4double fMassDif2
Definition: G4hhElastic.hh:128
void SetMaxEnergy(const G4double anEnergy)
G4double fMassTarg
Definition: G4hhElastic.hh:125
G4ParticleDefinition * thePionPlus
Definition: G4hhElastic.hh:100
static const double TeV
Definition: G4SIunits.hh:215
static const double keV
Definition: G4SIunits.hh:213
G4ParticleDefinition * fTarget
Definition: G4hhElastic.hh:96
G4ParticleDefinition * theNeutron
Definition: G4hhElastic.hh:99
G4double plabLowLimit
Definition: G4hhElastic.hh:108
G4double fRhoReIm
Definition: G4hhElastic.hh:150
Here is the call graph for this function:

◆ ~G4hhElastic()

G4hhElastic::~G4hhElastic ( )
virtual

Definition at line 206 of file G4hhElastic.cc.

207 {
208  if ( fEnergyVector ) {
209  delete fEnergyVector;
210  fEnergyVector = 0;
211  }
212 
213  for ( std::vector<G4PhysicsTable*>::iterator it = fBankT.begin();
214  it != fBankT.end(); ++it ) {
215  if ( (*it) ) (*it)->clearAndDestroy();
216  delete *it;
217  *it = 0;
218  }
219  fTableT = 0;
220  if(fHadrNuclXsc) delete fHadrNuclXsc;
221 }
G4PhysicsLogVector * fEnergyVector
Definition: G4hhElastic.hh:113
std::vector< G4PhysicsTable * > fBankT
Definition: G4hhElastic.hh:115
G4HadronNucleonXsc * fHadrNuclXsc
Definition: G4hhElastic.hh:249
G4PhysicsTable * fTableT
Definition: G4hhElastic.hh:114

Member Function Documentation

◆ BuildTableT()

void G4hhElastic::BuildTableT ( G4ParticleDefinition target,
G4ParticleDefinition projectile 
)

Definition at line 255 of file G4hhElastic.cc.

256 {
257  G4int iTkin, jTransfer;
258  G4double plab, Tkin, tMax;
259  G4double t1, t2, dt, delta = 0., sum = 0.;
260 
261  fTarget = target;
262  fProjectile = projectile;
267 
269  // G4HadronNucleonXsc* hnXsc = new G4HadronNucleonXsc();
271 
272  for( iTkin = 0; iTkin < fEnergyBin; iTkin++)
273  {
274  Tkin = fEnergyVector->GetLowEdgeEnergy(iTkin);
275  plab = std::sqrt( Tkin*( Tkin + 2*fMassProj ) );
276  // G4DynamicParticle* theDynamicParticle = new G4DynamicParticle(projectile,
277  // G4ParticleMomentum(0.,0.,1.),
278  // Tkin);
279  // fSigmaTot = fHadrNuclXsc->GetHadronNucleonXscNS( theDynamicParticle, target );
280 
281  SetParametersCMS( plab );
282 
283  tMax = 4.*fPcms*fPcms;
284  if( tMax > 15.*GeV*GeV ) tMax = 15.*GeV*GeV; // Check vs. energy ???
285 
287  sum = 0.;
288  dt = tMax/fBinT;
289 
290  // for(j = 1; j < fBinT; j++)
291 
292  for( jTransfer = fBinT-1; jTransfer >= 1; jTransfer--)
293  {
294  t1 = dt*(jTransfer-1);
295  t2 = t1 + dt;
296 
297  if( fMassProj > 900.*MeV ) // pp, pn
298  {
299  delta = integral.Legendre10(this, &G4hhElastic::GetdsdtF123, t1, t2);
300  // delta = integral.Legendre96(this, &G4hhElastic::GetdsdtF123, t1, t2);
301  }
302  else // pi+-p, K+-p
303  {
304  delta = integral.Legendre10(this, &G4hhElastic::GetdsdtF123qQgG, t1, t2);
305  // delta = integral.Legendre96(this, &G4hhElastic::GetdsdtF123qQgG, t1, t2);
306  }
307  sum += delta;
308  vectorT->PutValue( jTransfer-1, t1, sum ); // t2
309  }
310  // vectorT->PutValue( fBinT-1, dt*(fBinT-1), 0. ); // t2
311  fTableT->insertAt( iTkin, vectorT );
312  // delete theDynamicParticle;
313  }
314  // delete hnXsc;
315 
316  return;
317 }
G4double Legendre10(T &typeT, F f, G4double a, G4double b)
void SetParametersCMS(G4double plab)
Definition: G4hhElastic.hh:318
TTree * t1
Definition: plottest35.C:26
static const double MeV
Definition: G4SIunits.hh:211
void PutValue(size_t binNumber, G4double binValue, G4double dataValue)
G4PhysicsLogVector * fEnergyVector
Definition: G4hhElastic.hh:113
G4double fMassProj
Definition: G4hhElastic.hh:126
int G4int
Definition: G4Types.hh:78
G4double GetLowEdgeEnergy(size_t binNumber) const
G4double fPcms
Definition: G4hhElastic.hh:160
G4double GetdsdtF123qQgG(G4double q)
Definition: G4hhElastic.hh:748
static const double GeV
Definition: G4SIunits.hh:214
G4int fEnergyBin
Definition: G4hhElastic.hh:110
G4double GetdsdtF123(G4double q)
G4ParticleDefinition * fProjectile
Definition: G4hhElastic.hh:97
G4double fMassSum2
Definition: G4hhElastic.hh:127
G4PhysicsTable * fTableT
Definition: G4hhElastic.hh:114
cout<< "-> Edep in the target
Definition: analysis.C:54
G4double fMassDif2
Definition: G4hhElastic.hh:128
void insertAt(size_t, G4PhysicsVector *)
G4double fMassTarg
Definition: G4hhElastic.hh:125
TTree * t2
Definition: plottest35.C:36
double G4double
Definition: G4Types.hh:76
G4ParticleDefinition * fTarget
Definition: G4hhElastic.hh:96
Here is the call graph for this function:
Here is the caller graph for this function:

◆ BuildTableTest()

void G4hhElastic::BuildTableTest ( G4ParticleDefinition target,
G4ParticleDefinition projectile,
G4double  plab 
)

Definition at line 521 of file G4hhElastic.cc.

522 {
523  G4int jTransfer;
524  G4double tMax; // , sQq, sQG;
525  G4double t1, t2, dt, delta = 0., sum = 0. ; // , threshold;
526 
527  fTarget = target;
528  fProjectile = projectile;
533  fSpp = fMassProj*fMassProj + fMassTarg*fMassTarg + 2.*fMassTarg*std::sqrt(plab*plab + fMassProj*fMassProj);
534  fPcms = std::sqrt( (fSpp - fMassSum2)*(fSpp - fMassDif2)/4./fSpp);
535 
536  G4cout<<"fMassTarg = "<<fMassTarg<<" MeV; fMassProj = "<<fMassProj<<" MeV"<<G4endl;
537  tMax = 4.*fPcms*fPcms;
538  if( tMax > 15.*GeV*GeV ) tMax = 15.*GeV*GeV; // Check vs. energy ???
539 
540 
542  fTableT = new G4PhysicsTable(1);
544 
545  sum = 0.;
546  dt = tMax/G4double(fBinT);
547  G4cout<<"s = "<<std::sqrt(fSpp)/GeV<<" GeV; fPcms = "<<fPcms/GeV
548  <<" GeV; qMax = "<<tMax/GeV/GeV<<" GeV2; dt = "<<dt/GeV/GeV<<" GeV2"<<G4endl;
549 
550  // G4cout<<"fRA = "<<fRA*GeV<<"; fRB = "<<fRB*GeV<<G4endl;
551 
552  // for(jTransfer = 1; jTransfer < fBinT; jTransfer++)
553  for( jTransfer = fBinT-1; jTransfer >= 1; jTransfer-- )
554  {
555  t1 = dt*(jTransfer-1);
556  t2 = t1 + dt;
557 
558  if( fMassProj > 900.*MeV ) // pp, pn
559  {
560  delta = integral.Legendre10(this, &G4hhElastic::GetdsdtF123, t1, t2);
561  // threshold = integral.Legendre96(this, &G4hhElastic::GetdsdtF123, t1, tMax);
562  }
563  else // pi+-p, K+-p
564  {
565  delta = integral.Legendre10(this, &G4hhElastic::GetdsdtF123qQgG, t1, t2);
566  // threshold = integral.Legendre96(this, &G4hhElastic::GetdsdtF123qQgG, t1, tMax);
567  // delta = integral.Legendre96(this, &G4hhElastic::GetdsdtF123, t1, t2);
568  }
569  sum += delta;
570  // G4cout<<delta<<"\t"<<sum<<"\t"<<threshold<<G4endl;
571 
572  // sQq = GetdsdtF123(q1);
573  // sQG = GetdsdtF123qQgG(q1);
574  // G4cout<<q1/GeV<<"\t"<<sQG*GeV*GeV/millibarn<<"\t"<<sQq*GeV*GeV/millibarn<<G4endl;
575  // G4cout<<"sum = "<<sum<<", ";
576 
577  vectorT->PutValue( jTransfer-1, t1, sum ); // t2
578  }
579  // vectorT->PutValue( fBinT-1, dt*(fBinT-1), 0. ); // t2
580  fTableT->insertAt( 0, vectorT );
581  fBankT.push_back( fTableT ); // 0
582 
583  // for(jTransfer = 0; jTransfer < fBinT-1; jTransfer++)
584  // G4cout<<(*(*fTableT)(0))(jTransfer)/sum<<"\t\t"<<G4Pow::GetInstance()->powN(2.,-jTransfer)<<G4endl;
585 
586  return;
587 }
G4double Legendre10(T &typeT, F f, G4double a, G4double b)
TTree * t1
Definition: plottest35.C:26
static const double MeV
Definition: G4SIunits.hh:211
void PutValue(size_t binNumber, G4double binValue, G4double dataValue)
G4double fMassProj
Definition: G4hhElastic.hh:126
G4double fSpp
Definition: G4hhElastic.hh:159
int G4int
Definition: G4Types.hh:78
std::vector< G4PhysicsTable * > fBankT
Definition: G4hhElastic.hh:115
G4GLOB_DLL std::ostream G4cout
G4double fPcms
Definition: G4hhElastic.hh:160
G4double GetdsdtF123qQgG(G4double q)
Definition: G4hhElastic.hh:748
static const double GeV
Definition: G4SIunits.hh:214
G4double GetdsdtF123(G4double q)
G4ParticleDefinition * fProjectile
Definition: G4hhElastic.hh:97
G4double fMassSum2
Definition: G4hhElastic.hh:127
G4PhysicsTable * fTableT
Definition: G4hhElastic.hh:114
cout<< "-> Edep in the target
Definition: analysis.C:54
G4double fMassDif2
Definition: G4hhElastic.hh:128
void insertAt(size_t, G4PhysicsVector *)
G4double fMassTarg
Definition: G4hhElastic.hh:125
TTree * t2
Definition: plottest35.C:36
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4ParticleDefinition * fTarget
Definition: G4hhElastic.hh:96
Here is the call graph for this function:

◆ CalculateBQ()

void G4hhElastic::CalculateBQ ( G4double  b)
inline

Definition at line 959 of file G4hhElastic.hh.

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 }
Double_t y2[nxs]
G4double GetCofS1()
Definition: G4hhElastic.hh:905
Double_t y1[nxs]
Double_t x2[nxs]
Float_t d
G4double GetCofS3()
Definition: G4hhElastic.hh:931
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4double fOptRatio
Definition: G4hhElastic.hh:158
Double_t x1[nxs]
G4double GetCofS2()
Definition: G4hhElastic.hh:918
G4double fBQ
Definition: G4hhElastic.hh:156
double D(double temp)
static const G4double b1
#define G4endl
Definition: G4ios.hh:61
G4double fBq
Definition: G4hhElastic.hh:155
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CalculateBqQ12()

void G4hhElastic::CalculateBqQ12 ( G4double  b)
inline

Definition at line 788 of file G4hhElastic.hh.

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 }
G4complex Phi14()
Definition: G4hhElastic.hh:532
G4complex Phi13()
Definition: G4hhElastic.hh:522
G4double fLambda
Definition: G4hhElastic.hh:145
static const double hbarc
G4double fEta
Definition: G4hhElastic.hh:146
G4complex Phi23()
Definition: G4hhElastic.hh:542
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4GLOB_DLL std::ostream G4cout
static const double pi
Definition: SystemOfUnits.h:53
G4double fSigmaTot
Definition: G4hhElastic.hh:154
G4double fBQ
Definition: G4hhElastic.hh:156
static const G4double b1
#define G4endl
Definition: G4ios.hh:61
G4complex Phi24()
Definition: G4hhElastic.hh:552
G4double fBq
Definition: G4hhElastic.hh:155
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CalculateBqQ123()

void G4hhElastic::CalculateBqQ123 ( G4double  b)
inline

Definition at line 814 of file G4hhElastic.hh.

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 }
G4complex Phi14()
Definition: G4hhElastic.hh:532
G4complex Phi13()
Definition: G4hhElastic.hh:522
Double_t x2[nxs]
G4double fLambda
Definition: G4hhElastic.hh:145
static const double hbarc
G4double fEta
Definition: G4hhElastic.hh:146
G4complex Phi23()
Definition: G4hhElastic.hh:542
double C(double temp)
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4GLOB_DLL std::ostream G4cout
G4double fCofF2
Definition: G4hhElastic.hh:148
double A(double temperature)
G4double fOptRatio
Definition: G4hhElastic.hh:158
static const double pi
Definition: SystemOfUnits.h:53
Double_t x1[nxs]
G4double fSigmaTot
Definition: G4hhElastic.hh:154
G4double fBQ
Definition: G4hhElastic.hh:156
static const G4double b1
#define G4endl
Definition: G4ios.hh:61
G4complex Phi24()
Definition: G4hhElastic.hh:552
G4double fBq
Definition: G4hhElastic.hh:155
double G4double
Definition: G4Types.hh:76
G4double fCofF3
Definition: G4hhElastic.hh:149
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CalculateBqQ13()

void G4hhElastic::CalculateBqQ13 ( G4double  b)
inline

Definition at line 765 of file G4hhElastic.hh.

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 }
G4complex Phi14()
Definition: G4hhElastic.hh:532
G4double fAlpha
Definition: G4hhElastic.hh:134
static const double hbarc
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4GLOB_DLL std::ostream G4cout
static const double pi
Definition: SystemOfUnits.h:53
static const G4double b2
G4double fBqQ
Definition: G4hhElastic.hh:157
G4double fSigmaTot
Definition: G4hhElastic.hh:154
G4double fBQ
Definition: G4hhElastic.hh:156
#define G4endl
Definition: G4ios.hh:61
G4complex Phi24()
Definition: G4hhElastic.hh:552
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetAqq()

G4complex G4hhElastic::GetAqq ( )
inline

Definition at line 869 of file G4hhElastic.hh.

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 }
G4double fRq
Definition: G4hhElastic.hh:133
G4double fImCof
Definition: G4hhElastic.hh:147
G4double fLambda
Definition: G4hhElastic.hh:145
G4double fSpp
Definition: G4hhElastic.hh:159
std::complex< G4double > G4complex
Definition: G4Types.hh:81
static const double pi
Definition: SystemOfUnits.h:53
G4double fSo
Definition: G4hhElastic.hh:152
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double fAlphaP
Definition: G4hhElastic.hh:143
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetAQQ()

G4complex G4hhElastic::GetAQQ ( )
inline

Definition at line 882 of file G4hhElastic.hh.

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 }
G4double fRQ
Definition: G4hhElastic.hh:132
G4double fImCof
Definition: G4hhElastic.hh:147
G4double fLambda
Definition: G4hhElastic.hh:145
G4double fSpp
Definition: G4hhElastic.hh:159
std::complex< G4double > G4complex
Definition: G4Types.hh:81
static const double pi
Definition: SystemOfUnits.h:53
G4double fSo
Definition: G4hhElastic.hh:152
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double fAlphaP
Definition: G4hhElastic.hh:143
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetAqQ()

G4complex G4hhElastic::GetAqQ ( )
inline

Definition at line 895 of file G4hhElastic.hh.

896 {
897  G4complex z = 0.5*( GetAqq() + GetAQQ() );
898  return z;
899 }
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4complex GetAQQ()
Definition: G4hhElastic.hh:882
G4complex GetAqq()
Definition: G4hhElastic.hh:869
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetBq()

G4double G4hhElastic::GetBq ( )
inline

Definition at line 171 of file G4hhElastic.hh.

171 { return fBq;};
G4double fBq
Definition: G4hhElastic.hh:155

◆ GetBQ()

G4double G4hhElastic::GetBQ ( )
inline

Definition at line 172 of file G4hhElastic.hh.

172 { return fBQ;};
G4double fBQ
Definition: G4hhElastic.hh:156

◆ GetBqQ()

G4double G4hhElastic::GetBqQ ( )
inline

Definition at line 173 of file G4hhElastic.hh.

173 { return fBqQ;};
G4double fBqQ
Definition: G4hhElastic.hh:157

◆ GetCofF2()

G4double G4hhElastic::GetCofF2 ( )
inline

Definition at line 193 of file G4hhElastic.hh.

193 {return fCofF2;};
G4double fCofF2
Definition: G4hhElastic.hh:148

◆ GetCofF3()

G4double G4hhElastic::GetCofF3 ( )
inline

Definition at line 194 of file G4hhElastic.hh.

194 {return fCofF3;};
G4double fCofF3
Definition: G4hhElastic.hh:149

◆ GetCofS1()

G4double G4hhElastic::GetCofS1 ( )
inline

Definition at line 905 of file G4hhElastic.hh.

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 }
G4double fLambda
Definition: G4hhElastic.hh:145
G4complex GetAqQ()
Definition: G4hhElastic.hh:895
static const double hbarc
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4double fCofF2
Definition: G4hhElastic.hh:148
static const double pi
Definition: SystemOfUnits.h:53
G4double fSigmaTot
Definition: G4hhElastic.hh:154
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetCofS2()

G4double G4hhElastic::GetCofS2 ( )
inline

Definition at line 918 of file G4hhElastic.hh.

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 }
G4double fLambda
Definition: G4hhElastic.hh:145
G4complex GetAqQ()
Definition: G4hhElastic.hh:895
static const double hbarc
std::complex< G4double > G4complex
Definition: G4Types.hh:81
static const double pi
Definition: SystemOfUnits.h:53
G4complex GetAqq()
Definition: G4hhElastic.hh:869
G4double fSigmaTot
Definition: G4hhElastic.hh:154
double G4double
Definition: G4Types.hh:76
G4double fCofF3
Definition: G4hhElastic.hh:149
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetCofS3()

G4double G4hhElastic::GetCofS3 ( )
inline

Definition at line 931 of file G4hhElastic.hh.

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 }
G4double fLambda
Definition: G4hhElastic.hh:145
G4complex GetAqQ()
Definition: G4hhElastic.hh:895
static const double hbarc
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4complex GetAQQ()
Definition: G4hhElastic.hh:882
static const double pi
Definition: SystemOfUnits.h:53
G4double fSigmaTot
Definition: G4hhElastic.hh:154
double G4double
Definition: G4Types.hh:76
G4double fCofF3
Definition: G4hhElastic.hh:149
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetdsdtF1()

G4double G4hhElastic::GetdsdtF1 ( G4double  s,
G4double  q 
)
inline

Definition at line 1044 of file G4hhElastic.hh.

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 }
G4double fSpp
Definition: G4hhElastic.hh:159
std::complex< G4double > G4complex
Definition: G4Types.hh:81
static const double proton_mass_c2
static const double pi
Definition: SystemOfUnits.h:53
double G4double
Definition: G4Types.hh:76
G4complex GetF1(G4double qp)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetdsdtF12()

G4double G4hhElastic::GetdsdtF12 ( G4double  s,
G4double  q 
)
inline

Definition at line 1087 of file G4hhElastic.hh.

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 }
G4double fSpp
Definition: G4hhElastic.hh:159
std::complex< G4double > G4complex
Definition: G4Types.hh:81
static const double proton_mass_c2
static const double pi
Definition: SystemOfUnits.h:53
G4complex GetF2(G4double qp)
double G4double
Definition: G4Types.hh:76
G4complex GetF1(G4double qp)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetdsdtF123()

G4double G4hhElastic::GetdsdtF123 ( G4double  q)
inline

Definition at line 1130 of file G4hhElastic.hh.

1131 {
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 }
G4double fSpp
Definition: G4hhElastic.hh:159
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4double fCofF2
Definition: G4hhElastic.hh:148
static const double proton_mass_c2
static const double pi
Definition: SystemOfUnits.h:53
G4complex GetF2(G4double qp)
G4complex GetF3(G4double qp)
double G4double
Definition: G4Types.hh:76
G4double fCofF3
Definition: G4hhElastic.hh:149
G4complex GetF1(G4double qp)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetdsdtF123qQgG()

G4double G4hhElastic::GetdsdtF123qQgG ( G4double  q)
inline

Definition at line 748 of file G4hhElastic.hh.

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 }
G4complex GetF3qQgG(G4double qp)
Definition: G4hhElastic.hh:694
G4complex GetF2qQgG(G4double qp)
Definition: G4hhElastic.hh:645
G4double fSpp
Definition: G4hhElastic.hh:159
G4complex GetF1qQgG(G4double qp)
Definition: G4hhElastic.hh:564
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4double fCofF2
Definition: G4hhElastic.hh:148
static const double pi
Definition: SystemOfUnits.h:53
G4double fMassSum2
Definition: G4hhElastic.hh:127
G4double fMassDif2
Definition: G4hhElastic.hh:128
double G4double
Definition: G4Types.hh:76
G4double fCofF3
Definition: G4hhElastic.hh:149
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetdsdtF12qQgG()

G4double G4hhElastic::GetdsdtF12qQgG ( G4double  s,
G4double  q 
)
inline

Definition at line 678 of file G4hhElastic.hh.

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 }
G4complex GetF2qQgG(G4double qp)
Definition: G4hhElastic.hh:645
G4double fSpp
Definition: G4hhElastic.hh:159
G4complex GetF1qQgG(G4double qp)
Definition: G4hhElastic.hh:564
std::complex< G4double > G4complex
Definition: G4Types.hh:81
static const double pi
Definition: SystemOfUnits.h:53
G4double fMassSum2
Definition: G4hhElastic.hh:127
G4double fMassDif2
Definition: G4hhElastic.hh:128
#define F12
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetdsdtF13qQG()

G4double G4hhElastic::GetdsdtF13qQG ( G4double  s,
G4double  q 
)
inline

Definition at line 586 of file G4hhElastic.hh.

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.);
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 }
G4complex Phi14()
Definition: G4hhElastic.hh:532
G4double fAlpha
Definition: G4hhElastic.hh:134
#define F13
G4double fLambda
Definition: G4hhElastic.hh:145
static const double hbarc
G4double fEta
Definition: G4hhElastic.hh:146
G4double fSpp
Definition: G4hhElastic.hh:159
G4double fGamma
Definition: G4hhElastic.hh:140
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4double fBeta
Definition: G4hhElastic.hh:135
static const double pi
Definition: SystemOfUnits.h:53
G4double fMassSum2
Definition: G4hhElastic.hh:127
G4double fBqQ
Definition: G4hhElastic.hh:157
G4double fSigmaTot
Definition: G4hhElastic.hh:154
G4double fBQ
Definition: G4hhElastic.hh:156
G4double fMassDif2
Definition: G4hhElastic.hh:128
G4complex Phi24()
Definition: G4hhElastic.hh:552
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetdsdtF1qQgG()

G4double G4hhElastic::GetdsdtF1qQgG ( G4double  s,
G4double  q 
)
inline

Definition at line 629 of file G4hhElastic.hh.

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 }
G4double fSpp
Definition: G4hhElastic.hh:159
G4complex GetF1qQgG(G4double qp)
Definition: G4hhElastic.hh:564
std::complex< G4double > G4complex
Definition: G4Types.hh:81
static const double pi
Definition: SystemOfUnits.h:53
G4double fMassSum2
Definition: G4hhElastic.hh:127
G4double fMassDif2
Definition: G4hhElastic.hh:128
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetExpRatioF123()

G4double G4hhElastic::GetExpRatioF123 ( G4double  s,
G4double  q 
)
inline

Definition at line 1145 of file G4hhElastic.hh.

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 }
G4double fSpp
Definition: G4hhElastic.hh:159
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4double fCofF2
Definition: G4hhElastic.hh:148
static const double proton_mass_c2
static const double pi
Definition: SystemOfUnits.h:53
#define F10
G4complex GetF2(G4double qp)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double fExpSlope
Definition: G4hhElastic.hh:151
G4complex GetF3(G4double qp)
double G4double
Definition: G4Types.hh:76
G4double fCofF3
Definition: G4hhElastic.hh:149
G4complex GetF1(G4double qp)
G4double fRhoReIm
Definition: G4hhElastic.hh:150
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetF1()

G4complex G4hhElastic::GetF1 ( G4double  qp)
inline

Definition at line 1026 of file G4hhElastic.hh.

1027 {
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 }
G4complex GetAqQ()
Definition: G4hhElastic.hh:895
static const double hbarc
G4double fSpp
Definition: G4hhElastic.hh:159
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4complex GetAQQ()
Definition: G4hhElastic.hh:882
static const double proton_mass_c2
static const double pi
Definition: SystemOfUnits.h:53
G4complex GetAqq()
Definition: G4hhElastic.hh:869
G4double fSigmaTot
Definition: G4hhElastic.hh:154
G4double fBQ
Definition: G4hhElastic.hh:156
G4double fBq
Definition: G4hhElastic.hh:155
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetF1qQgG()

G4complex G4hhElastic::GetF1qQgG ( G4double  qp)
inline

Definition at line 564 of file G4hhElastic.hh.

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 }
G4complex Phi14()
Definition: G4hhElastic.hh:532
G4double fAlpha
Definition: G4hhElastic.hh:134
G4complex Phi13()
Definition: G4hhElastic.hh:522
G4double fLambda
Definition: G4hhElastic.hh:145
static const double hbarc
G4double fEta
Definition: G4hhElastic.hh:146
G4complex Phi23()
Definition: G4hhElastic.hh:542
G4double fSpp
Definition: G4hhElastic.hh:159
G4double fGamma
Definition: G4hhElastic.hh:140
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4double fBeta
Definition: G4hhElastic.hh:135
static const double pi
Definition: SystemOfUnits.h:53
G4double fMassSum2
Definition: G4hhElastic.hh:127
G4double fDelta
Definition: G4hhElastic.hh:141
G4double fSigmaTot
Definition: G4hhElastic.hh:154
G4double fBQ
Definition: G4hhElastic.hh:156
G4double fMassDif2
Definition: G4hhElastic.hh:128
G4complex Phi24()
Definition: G4hhElastic.hh:552
G4double fBq
Definition: G4hhElastic.hh:155
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetF2()

G4complex G4hhElastic::GetF2 ( G4double  qp)
inline

Definition at line 1059 of file G4hhElastic.hh.

1060 {
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;
1079 
1080  return res;
1081 }
G4double fLambda
Definition: G4hhElastic.hh:145
G4complex GetAqQ()
Definition: G4hhElastic.hh:895
static const double hbarc
G4double fSpp
Definition: G4hhElastic.hh:159
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4complex GetAQQ()
Definition: G4hhElastic.hh:882
static const double proton_mass_c2
static const double pi
Definition: SystemOfUnits.h:53
G4complex GetAqq()
Definition: G4hhElastic.hh:869
G4double fSigmaTot
Definition: G4hhElastic.hh:154
G4double fBQ
Definition: G4hhElastic.hh:156
G4double fBq
Definition: G4hhElastic.hh:155
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetF2qQgG()

G4complex G4hhElastic::GetF2qQgG ( G4double  qp)
inline

Definition at line 645 of file G4hhElastic.hh.

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 }
G4complex Phi14()
Definition: G4hhElastic.hh:532
G4double fAlpha
Definition: G4hhElastic.hh:134
G4complex Phi13()
Definition: G4hhElastic.hh:522
G4double fLambda
Definition: G4hhElastic.hh:145
static const double hbarc
G4double fEta
Definition: G4hhElastic.hh:146
G4complex Phi23()
Definition: G4hhElastic.hh:542
G4double fSpp
Definition: G4hhElastic.hh:159
G4double fGamma
Definition: G4hhElastic.hh:140
std::complex< G4double > G4complex
Definition: G4Types.hh:81
static const double pi
Definition: SystemOfUnits.h:53
G4double fMassSum2
Definition: G4hhElastic.hh:127
G4double fDelta
Definition: G4hhElastic.hh:141
G4double fSigmaTot
Definition: G4hhElastic.hh:154
G4double fBQ
Definition: G4hhElastic.hh:156
G4double fMassDif2
Definition: G4hhElastic.hh:128
G4complex Phi24()
Definition: G4hhElastic.hh:552
G4double fBq
Definition: G4hhElastic.hh:155
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetF3()

G4complex G4hhElastic::GetF3 ( G4double  qp)
inline

Definition at line 1102 of file G4hhElastic.hh.

1103 {
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 }
G4double fLambda
Definition: G4hhElastic.hh:145
G4complex GetAqQ()
Definition: G4hhElastic.hh:895
static const double hbarc
G4double fSpp
Definition: G4hhElastic.hh:159
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4complex GetAQQ()
Definition: G4hhElastic.hh:882
static const double proton_mass_c2
static const double pi
Definition: SystemOfUnits.h:53
G4complex GetAqq()
Definition: G4hhElastic.hh:869
G4double fSigmaTot
Definition: G4hhElastic.hh:154
G4double fBQ
Definition: G4hhElastic.hh:156
G4double fBq
Definition: G4hhElastic.hh:155
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetF3qQgG()

G4complex G4hhElastic::GetF3qQgG ( G4double  qp)
inline

Definition at line 694 of file G4hhElastic.hh.

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.);
740 
741  return res;
742 }
G4complex Phi14()
Definition: G4hhElastic.hh:532
G4double fAlpha
Definition: G4hhElastic.hh:134
G4complex Phi13()
Definition: G4hhElastic.hh:522
G4double fLambda
Definition: G4hhElastic.hh:145
static const double hbarc
G4double fEta
Definition: G4hhElastic.hh:146
G4complex Phi23()
Definition: G4hhElastic.hh:542
G4double fSpp
Definition: G4hhElastic.hh:159
G4double fGamma
Definition: G4hhElastic.hh:140
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4double fBeta
Definition: G4hhElastic.hh:135
static const double pi
Definition: SystemOfUnits.h:53
G4double fMassSum2
Definition: G4hhElastic.hh:127
G4double fDelta
Definition: G4hhElastic.hh:141
G4double fSigmaTot
Definition: G4hhElastic.hh:154
G4double fBQ
Definition: G4hhElastic.hh:156
G4double fMassDif2
Definition: G4hhElastic.hh:128
G4complex Phi24()
Definition: G4hhElastic.hh:552
G4double fBq
Definition: G4hhElastic.hh:155
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetImCof()

G4double G4hhElastic::GetImCof ( )
inline

Definition at line 188 of file G4hhElastic.hh.

188 {return fImCof;};
G4double fImCof
Definition: G4hhElastic.hh:147

◆ GetOpticalRatio()

G4double G4hhElastic::GetOpticalRatio ( )
inline

Definition at line 944 of file G4hhElastic.hh.

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 }
G4double fOptRatio
Definition: G4hhElastic.hh:158
Here is the caller graph for this function:

◆ GetRA()

G4double G4hhElastic::GetRA ( )
inline

Definition at line 196 of file G4hhElastic.hh.

196 { return fRA;};
G4double fRA
Definition: G4hhElastic.hh:131
Here is the caller graph for this function:

◆ GetRB()

G4double G4hhElastic::GetRB ( )
inline

Definition at line 200 of file G4hhElastic.hh.

200 { return fRB;};
G4double fRB
Definition: G4hhElastic.hh:137
Here is the caller graph for this function:

◆ GetRg()

G4double G4hhElastic::GetRg ( )
inline

Definition at line 201 of file G4hhElastic.hh.

201 { return fRg;};
G4double fRg
Definition: G4hhElastic.hh:139

◆ GetRG()

G4double G4hhElastic::GetRG ( )
inline

Definition at line 202 of file G4hhElastic.hh.

202 { return fRG;};
G4double fRG
Definition: G4hhElastic.hh:138
Here is the call graph for this function:

◆ GetRhoReIm()

G4double G4hhElastic::GetRhoReIm ( )
inline

Definition at line 177 of file G4hhElastic.hh.

177 { return fRhoReIm;};
G4double fRhoReIm
Definition: G4hhElastic.hh:150
Here is the call graph for this function:

◆ GetRq()

G4double G4hhElastic::GetRq ( )
inline

Definition at line 197 of file G4hhElastic.hh.

197 { return fRq;};
G4double fRq
Definition: G4hhElastic.hh:133

◆ GetRQ()

G4double G4hhElastic::GetRQ ( )
inline

Definition at line 198 of file G4hhElastic.hh.

198 { return fRQ;};
G4double fRQ
Definition: G4hhElastic.hh:132

◆ GetSpp()

G4double G4hhElastic::GetSpp ( )
inline

Definition at line 168 of file G4hhElastic.hh.

168 {return fSpp;};
G4double fSpp
Definition: G4hhElastic.hh:159
Here is the call graph for this function:

◆ GetTransfer()

G4double G4hhElastic::GetTransfer ( G4int  iMomentum,
G4int  iTransfer,
G4double  position 
)

Definition at line 630 of file G4hhElastic.cc.

631 {
632  G4double x1, x2, y1, y2, randTransfer, delta, mean, epsilon = 1.e-6;
633 
634  if( iTransfer == 0 )
635  {
636  randTransfer = (*fTableT)(iTkin)->GetLowEdgeEnergy(iTransfer);
637  // iTransfer++;
638  }
639  else
640  {
641  if ( iTransfer >= G4int((*fTableT)(iTkin)->GetVectorLength()) )
642  {
643  iTransfer = (*fTableT)(iTkin)->GetVectorLength() - 1;
644  }
645  y1 = (*(*fTableT)(iTkin))(iTransfer-1);
646  y2 = (*(*fTableT)(iTkin))(iTransfer);
647 
648  x1 = (*fTableT)(iTkin)->GetLowEdgeEnergy(iTransfer-1);
649  x2 = (*fTableT)(iTkin)->GetLowEdgeEnergy(iTransfer);
650 
651  delta = y2 - y1;
652  mean = y2 + y1;
653 
654  if ( x1 == x2 ) randTransfer = x2;
655  else
656  {
657  // if ( y1 == y2 )
658  if ( delta < epsilon*mean )
659  randTransfer = x1 + ( x2 - x1 )*G4UniformRand();
660  else randTransfer = x1 + ( position - y1 )*( x2 - x1 )/delta; // ( y2 - y1 );
661  }
662  }
663  return randTransfer;
664 }
Double_t y2[nxs]
Double_t y1[nxs]
Double_t x2[nxs]
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
Double_t x1[nxs]
G4PhysicsTable * fTableT
Definition: G4hhElastic.hh:114
double G4double
Definition: G4Types.hh:76
double epsilon(double density, double temperature)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Initialise()

void G4hhElastic::Initialise ( )

Definition at line 231 of file G4hhElastic.cc.

232 {
233  // pp,pn
234 
237  fBankT.push_back(fTableT); // 0
238 
239  // pi+-p
240 
243  fBankT.push_back(fTableT); // 1
244  //K+-p
247  fBankT.push_back(fTableT); // 2
248 
249 }
void BuildTableT(G4ParticleDefinition *target, G4ParticleDefinition *projectile)
Definition: G4hhElastic.cc:255
std::vector< G4PhysicsTable * > fBankT
Definition: G4hhElastic.hh:115
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
G4ParticleDefinition * fProjectile
Definition: G4hhElastic.hh:97
G4PhysicsTable * fTableT
Definition: G4hhElastic.hh:114
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
G4ParticleDefinition * fTarget
Definition: G4hhElastic.hh:96
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsApplicable()

G4bool G4hhElastic::IsApplicable ( const G4HadProjectile projectile,
G4Nucleus nucleus 
)
inlinevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 258 of file G4hhElastic.hh.

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 }
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
const G4ParticleDefinition * GetDefinition() const
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
Here is the call graph for this function:

◆ Phi13()

G4complex G4hhElastic::Phi13 ( )
inline

Definition at line 522 of file G4hhElastic.hh.

523 {
524  G4double re = (fRq*fRq + fRg*fRg)/16.;
525  G4complex result(re,0.);
526  result += Pomeron();
527  return result;
528 }
G4double fRq
Definition: G4hhElastic.hh:133
G4complex Pomeron()
Definition: G4hhElastic.hh:511
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4double fRg
Definition: G4hhElastic.hh:139
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Phi14()

G4complex G4hhElastic::Phi14 ( )
inline

Definition at line 532 of file G4hhElastic.hh.

533 {
534  G4double re = (fRq*fRq + fRG*fRG)/16.;
535  G4complex result(re,0.);
536  result += Pomeron();
537  return result;
538 }
G4double fRq
Definition: G4hhElastic.hh:133
G4complex Pomeron()
Definition: G4hhElastic.hh:511
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4double fRG
Definition: G4hhElastic.hh:138
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Phi23()

G4complex G4hhElastic::Phi23 ( )
inline

Definition at line 542 of file G4hhElastic.hh.

543 {
544  G4double re = (fRQ*fRQ + fRg*fRg)/16.;
545  G4complex result(re,0.);
546  result += Pomeron();
547  return result;
548 }
G4double fRQ
Definition: G4hhElastic.hh:132
G4complex Pomeron()
Definition: G4hhElastic.hh:511
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4double fRg
Definition: G4hhElastic.hh:139
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Phi24()

G4complex G4hhElastic::Phi24 ( )
inline

Definition at line 552 of file G4hhElastic.hh.

553 {
554  G4double re = (fRQ*fRQ + fRG*fRG)/16.;
555  G4complex result(re,0.);
556  result += Pomeron();
557  return result;
558 }
G4double fRQ
Definition: G4hhElastic.hh:132
G4complex Pomeron()
Definition: G4hhElastic.hh:511
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4double fRG
Definition: G4hhElastic.hh:138
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Pomeron()

G4complex G4hhElastic::Pomeron ( )
inline

Definition at line 511 of file G4hhElastic.hh.

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 }
G4double fImCof
Definition: G4hhElastic.hh:147
G4double fSpp
Definition: G4hhElastic.hh:159
std::complex< G4double > G4complex
Definition: G4Types.hh:81
static const double pi
Definition: SystemOfUnits.h:53
G4double fSo
Definition: G4hhElastic.hh:152
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double fAlphaP
Definition: G4hhElastic.hh:143
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SampleBisectionalT()

G4double G4hhElastic::SampleBisectionalT ( const G4ParticleDefinition p,
G4double  plab 
)

Definition at line 440 of file G4hhElastic.cc.

441 {
442  G4int iTkin, iTransfer;
443  G4double t, position, m1 = aParticle->GetPDGMass();
444  G4double Tkin = std::sqrt(m1*m1+p*p) - m1;
445 
446  if( aParticle == G4Proton::Proton() || aParticle == G4Neutron::Neutron() )
447  {
448  fTableT = fBankT[0];
449  }
450  if( aParticle == G4PionPlus::PionPlus() || aParticle == G4PionMinus::PionMinus() )
451  {
452  fTableT = fBankT[1];
453  }
454  if( aParticle == G4KaonPlus::KaonPlus() || aParticle == G4KaonMinus::KaonMinus() )
455  {
456  fTableT = fBankT[2];
457  }
458  G4double delta = std::abs(Tkin - fOldTkin)/(Tkin + fOldTkin);
459  G4double deltaMax = 1.e-2;
460 
461  if ( delta < deltaMax ) iTkin = fInTkin;
462  else
463  {
464  for( iTkin = 0; iTkin < fEnergyBin; iTkin++ )
465  {
466  if( Tkin < fEnergyVector->GetLowEdgeEnergy(iTkin) ) break;
467  }
468  }
469  if ( iTkin >= fEnergyBin ) iTkin = fEnergyBin-1; // Tkin is more then theMaxEnergy
470  if ( iTkin < 0 ) iTkin = 0; // against negative index, Tkin < theMinEnergy
471 
472  fOldTkin = Tkin;
473  fInTkin = iTkin;
474 
475  if (iTkin == fEnergyBin -1 || iTkin == 0 ) // the table edges
476  {
477  position = (*(*fTableT)(iTkin))(0)*G4UniformRand();
478 
479  for(iTransfer = 0; iTransfer < fBinT-1; iTransfer++)
480  {
481  if( position >= (*(*fTableT)(iTkin))(iTransfer) ) break;
482  }
483  if (iTransfer >= fBinT-1) iTransfer = fBinT-2;
484 
485  t = GetTransfer(iTkin, iTransfer, position);
486 
487 
488  }
489  else // Tkin inside between energy table edges
490  {
491  G4double rand = G4UniformRand();
492  position = (*(*fTableT)(iTkin))(0)*rand;
493 
494  //
495  // (*fTableT)(iTkin)->GetLowEdgeEnergy(fBinT-2);
496  G4int sTransfer = 0, fTransfer = fBinT - 2, dTransfer = fTransfer - sTransfer;
497  G4double y2;
498 
499  for( iTransfer = 0; iTransfer < fBinT - 1; iTransfer++ )
500  {
501  // dTransfer %= 2;
502  dTransfer /= 2;
503  // dTransfer *= 0.5;
504  y2 = (*(*fTableT)(iTkin))( sTransfer + dTransfer );
505 
506  if( y2 > position ) sTransfer += dTransfer;
507 
508  // if( dTransfer <= 1 ) break;
509  if( dTransfer < 1 ) break;
510  }
511  t = (*fTableT)(iTkin)->GetLowEdgeEnergy(sTransfer); // +(-0.5+rand)*(*fTableT)(iTkin)->GetLowEdgeEnergy(3);
512  }
513  return t;
514 }
Double_t y2[nxs]
G4double fOldTkin
Definition: G4hhElastic.hh:246
int G4int
Definition: G4Types.hh:78
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
std::vector< G4PhysicsTable * > fBankT
Definition: G4hhElastic.hh:115
#define position
Definition: xmlparse.cc:622
#define G4UniformRand()
Definition: Randomize.hh:97
G4double GetTransfer(G4int iMomentum, G4int iTransfer, G4double position)
Definition: G4hhElastic.cc:630
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
G4int fEnergyBin
Definition: G4hhElastic.hh:110
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4PhysicsTable * fTableT
Definition: G4hhElastic.hh:114
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
double G4double
Definition: G4Types.hh:76
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
Here is the call graph for this function:

◆ SampleInvariantT()

G4double G4hhElastic::SampleInvariantT ( const G4ParticleDefinition p,
G4double  plab,
G4int  ,
G4int   
)
virtual

Reimplemented from G4HadronElastic.

Definition at line 323 of file G4hhElastic.cc.

325 {
326  G4int iTkin, iTransfer;
327  G4double t, t2, position, m1 = aParticle->GetPDGMass();
328  G4double Tkin = std::sqrt(m1*m1+p*p) - m1;
329 
330  if( aParticle == G4Proton::Proton() || aParticle == G4Neutron::Neutron() )
331  {
332  fTableT = fBankT[0];
333  }
334  if( aParticle == G4PionPlus::PionPlus() || aParticle == G4PionMinus::PionMinus() )
335  {
336  fTableT = fBankT[1];
337  }
338  if( aParticle == G4KaonPlus::KaonPlus() || aParticle == G4KaonMinus::KaonMinus() )
339  {
340  fTableT = fBankT[2];
341  }
342 
343  G4double delta = std::abs(Tkin - fOldTkin)/(Tkin + fOldTkin);
344  G4double deltaMax = 1.e-2;
345 
346  if ( delta < deltaMax ) iTkin = fInTkin;
347  else
348  {
349  for( iTkin = 0; iTkin < fEnergyBin; iTkin++)
350  {
351  if( Tkin < fEnergyVector->GetLowEdgeEnergy(iTkin) ) break;
352  }
353  }
354  if ( iTkin >= fEnergyBin ) iTkin = fEnergyBin-1; // Tkin is more then theMaxEnergy
355  if ( iTkin < 0 ) iTkin = 0; // against negative index, Tkin < theMinEnergy
356 
357  fOldTkin = Tkin;
358  fInTkin = iTkin;
359 
360  if (iTkin == fEnergyBin -1 || iTkin == 0 ) // the table edges
361  {
362  position = (*(*fTableT)(iTkin))(0)*G4UniformRand();
363 
364  // G4cout<<"position = "<<position<<G4endl;
365 
366  for(iTransfer = 0; iTransfer < fBinT-1; iTransfer++)
367  {
368  if( position >= (*(*fTableT)(iTkin))(iTransfer) ) break;
369  }
370  if (iTransfer >= fBinT-1) iTransfer = fBinT-2;
371 
372  // G4cout<<"iTransfer = "<<iTransfer<<G4endl;
373 
374  t = GetTransfer(iTkin, iTransfer, position);
375 
376  // G4cout<<"t = "<<t<<G4endl;
377  }
378  else // Tkin inside between energy table edges
379  {
380  // position = (*(*fTableT)(iTkin))(fBinT-2)*G4UniformRand();
381  position = (*(*fTableT)(iTkin))(0)*G4UniformRand();
382 
383  // G4cout<<"position = "<<position<<G4endl;
384 
385  for(iTransfer = 0; iTransfer < fBinT-1; iTransfer++)
386  {
387  // if( position < (*(*fTableT)(iTkin))(iTransfer) ) break;
388  if( position >= (*(*fTableT)(iTkin))(iTransfer) ) break;
389  }
390  if (iTransfer >= fBinT-1) iTransfer = fBinT-2;
391 
392  // G4cout<<"iTransfer = "<<iTransfer<<G4endl;
393 
394  t2 = GetTransfer(iTkin, iTransfer, position);
395  return t2;
396  /*
397  G4double t1, E1, E2, W, W1, W2;
398  // G4cout<<"t2 = "<<t2<<G4endl;
399 
400  E2 = fEnergyVector->GetLowEdgeEnergy(iTkin);
401 
402  // G4cout<<"E2 = "<<E2<<G4endl;
403 
404  iTkin--;
405 
406  // position = (*(*fTableT)(iTkin))(fBinT-2)*G4UniformRand();
407 
408  // G4cout<<"position = "<<position<<G4endl;
409 
410  for(iTransfer = 0; iTransfer < fBinT-1; iTransfer++)
411  {
412  // if( position < (*(*fTableT)(iTkin))(iTransfer) ) break;
413  if( position >= (*(*fTableT)(iTkin))(iTransfer) ) break;
414  }
415  if (iTransfer >= fBinT-1) iTransfer = fBinT-2;
416 
417  t1 = GetTransfer(iTkin, iTransfer, position);
418 
419  // G4cout<<"t1 = "<<t1<<G4endl;
420 
421  E1 = fEnergyVector->GetLowEdgeEnergy(iTkin);
422 
423  // G4cout<<"E1 = "<<E1<<G4endl;
424 
425  W = 1.0/(E2 - E1);
426  W1 = (E2 - Tkin)*W;
427  W2 = (Tkin - E1)*W;
428 
429  t = W1*t1 + W2*t2;
430  */
431  }
432  return t;
433 }
G4double fOldTkin
Definition: G4hhElastic.hh:246
int G4int
Definition: G4Types.hh:78
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
std::vector< G4PhysicsTable * > fBankT
Definition: G4hhElastic.hh:115
#define position
Definition: xmlparse.cc:622
#define G4UniformRand()
Definition: Randomize.hh:97
G4double GetTransfer(G4int iMomentum, G4int iTransfer, G4double position)
Definition: G4hhElastic.cc:630
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
G4int fEnergyBin
Definition: G4hhElastic.hh:110
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4PhysicsTable * fTableT
Definition: G4hhElastic.hh:114
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
TTree * t2
Definition: plottest35.C:36
double G4double
Definition: G4Types.hh:76
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
Here is the call graph for this function:

◆ SampleTest()

G4double G4hhElastic::SampleTest ( G4double  tMin)

Definition at line 594 of file G4hhElastic.cc.

595 {
596  G4int iTkin, iTransfer, iTmin;
597  G4double t, position;
598  // G4double qMin = std::sqrt(tMin);
599 
600  fTableT = fBankT[0];
601  iTkin = 0;
602 
603  for(iTransfer = 0; iTransfer < fBinT-1; iTransfer++)
604  {
605  // if( qMin <= (*fTableT)(iTkin)->GetLowEdgeEnergy(iTransfer) ) break;
606  if( tMin <= (*fTableT)(iTkin)->GetLowEdgeEnergy(iTransfer) ) break;
607  }
608  iTmin = iTransfer-1;
609  if(iTmin < 0 ) iTmin = 0;
610 
611  position = (*(*fTableT)(iTkin))(iTmin)*G4UniformRand();
612 
613  for( iTmin = 0; iTransfer < fBinT-1; iTransfer++)
614  {
615  if( position > (*(*fTableT)(iTkin))(iTransfer) ) break;
616  }
617  if (iTransfer >= fBinT-1) iTransfer = fBinT-2;
618 
619  t = GetTransfer(iTkin, iTransfer, position);
620 
621  return t;
622 }
int G4int
Definition: G4Types.hh:78
std::vector< G4PhysicsTable * > fBankT
Definition: G4hhElastic.hh:115
#define position
Definition: xmlparse.cc:622
#define G4UniformRand()
Definition: Randomize.hh:97
G4double GetTransfer(G4int iMomentum, G4int iTransfer, G4double position)
Definition: G4hhElastic.cc:630
G4PhysicsTable * fTableT
Definition: G4hhElastic.hh:114
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ SetAlphaP()

void G4hhElastic::SetAlphaP ( G4double  a)
inline

Definition at line 186 of file G4hhElastic.hh.

186 {fAlphaP = a;};
G4double fAlphaP
Definition: G4hhElastic.hh:143

◆ SetBq()

void G4hhElastic::SetBq ( G4double  b)
inline

Definition at line 174 of file G4hhElastic.hh.

174 {fBq = b;};
G4double fBq
Definition: G4hhElastic.hh:155
Here is the caller graph for this function:

◆ SetBQ()

void G4hhElastic::SetBQ ( G4double  b)
inline

Definition at line 175 of file G4hhElastic.hh.

175 {fBQ = b;};
G4double fBQ
Definition: G4hhElastic.hh:156
Here is the caller graph for this function:

◆ SetBqQ()

void G4hhElastic::SetBqQ ( G4double  b)
inline

Definition at line 176 of file G4hhElastic.hh.

176 {fBqQ = b;};
G4double fBqQ
Definition: G4hhElastic.hh:157

◆ SetCofF2()

void G4hhElastic::SetCofF2 ( G4double  f)
inline

Definition at line 191 of file G4hhElastic.hh.

191 {fCofF2 = f;};
G4double fCofF2
Definition: G4hhElastic.hh:148
Here is the caller graph for this function:

◆ SetCofF3()

void G4hhElastic::SetCofF3 ( G4double  f)
inline

Definition at line 192 of file G4hhElastic.hh.

192 {fCofF3 = f;};
G4double fCofF3
Definition: G4hhElastic.hh:149
Here is the caller graph for this function:

◆ SetEta()

void G4hhElastic::SetEta ( G4double  E)
inline

Definition at line 190 of file G4hhElastic.hh.

190 {fEta = E;};
G4double fEta
Definition: G4hhElastic.hh:146
Here is the caller graph for this function:

◆ SetImCof()

void G4hhElastic::SetImCof ( G4double  a)
inline

Definition at line 187 of file G4hhElastic.hh.

187 {fImCof = a;};
G4double fImCof
Definition: G4hhElastic.hh:147
Here is the caller graph for this function:

◆ SetLambda()

void G4hhElastic::SetLambda ( G4double  L)
inline

Definition at line 189 of file G4hhElastic.hh.

189 {fLambda = L;};
G4double fLambda
Definition: G4hhElastic.hh:145
static const double L
Definition: G4SIunits.hh:123
Here is the caller graph for this function:

◆ SetParameters()

void G4hhElastic::SetParameters ( )
inline

Definition at line 273 of file G4hhElastic.hh.

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 }
G4double fRQ
Definition: G4hhElastic.hh:132
G4double fAlpha
Definition: G4hhElastic.hh:134
G4double fRq
Definition: G4hhElastic.hh:133
G4double fImCof
Definition: G4hhElastic.hh:147
G4double fLambda
Definition: G4hhElastic.hh:145
G4double fEta
Definition: G4hhElastic.hh:146
G4double fRA
Definition: G4hhElastic.hh:131
static const double GeV
G4double fMq
Definition: G4hhElastic.hh:123
G4double fGamma
Definition: G4hhElastic.hh:140
G4double fCofF2
Definition: G4hhElastic.hh:148
G4double fBeta
Definition: G4hhElastic.hh:135
G4double fRG
Definition: G4hhElastic.hh:138
G4double fQcof
Definition: G4hhElastic.hh:161
G4double fSo
Definition: G4hhElastic.hh:152
G4double fBqQ
Definition: G4hhElastic.hh:157
G4double fAlphaP
Definition: G4hhElastic.hh:143
G4double fRg
Definition: G4hhElastic.hh:139
G4double fDelta
Definition: G4hhElastic.hh:141
G4double fBQ
Definition: G4hhElastic.hh:156
G4double fExpSlope
Definition: G4hhElastic.hh:151
G4double fLambdaFF
Definition: G4hhElastic.hh:144
G4double fMff2
Definition: G4hhElastic.hh:121
G4double fBq
Definition: G4hhElastic.hh:155
G4double fMQ
Definition: G4hhElastic.hh:122
G4double fCofF3
Definition: G4hhElastic.hh:149
G4double fRB
Definition: G4hhElastic.hh:137
Here is the caller graph for this function:

◆ SetParametersCMS()

void G4hhElastic::SetParametersCMS ( G4double  plab)
inline

Definition at line 318 of file G4hhElastic.hh.

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 }
void SetRA(G4double rn, G4double pq, G4double pQ)
Definition: G4hhElastic.hh:489
void SetBQ(G4double b)
Definition: G4hhElastic.hh:175
static const G4double theNuclNuclData[18][6]
Definition: G4hhElastic.hh:247
G4double fMassProj
Definition: G4hhElastic.hh:126
void SetImCof(G4double a)
Definition: G4hhElastic.hh:187
G4double fSpp
Definition: G4hhElastic.hh:159
int G4int
Definition: G4Types.hh:78
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
static const double GeV
G4double fGamma
Definition: G4hhElastic.hh:140
void SetCofF3(G4double f)
Definition: G4hhElastic.hh:192
void SetEta(G4double E)
Definition: G4hhElastic.hh:190
G4double fPcms
Definition: G4hhElastic.hh:160
void SetCofF2(G4double f)
Definition: G4hhElastic.hh:191
G4HadronNucleonXsc * fHadrNuclXsc
Definition: G4hhElastic.hh:249
G4ParticleDefinition * fProjectile
Definition: G4hhElastic.hh:97
G4double fMassSum2
Definition: G4hhElastic.hh:127
void SetBq(G4double b)
Definition: G4hhElastic.hh:174
void SetRB(G4double rn, G4double pq, G4double pQ)
Definition: G4hhElastic.hh:500
G4double fDelta
Definition: G4hhElastic.hh:141
G4double fSigmaTot
Definition: G4hhElastic.hh:154
G4double fMassDif2
Definition: G4hhElastic.hh:128
G4double fMassTarg
Definition: G4hhElastic.hh:125
static const G4double thePiKaNuclData[8][6]
Definition: G4hhElastic.hh:248
double G4double
Definition: G4Types.hh:76
void SetLambda(G4double L)
Definition: G4hhElastic.hh:189
G4ThreeVector G4ParticleMomentum
G4ParticleDefinition * fTarget
Definition: G4hhElastic.hh:96
G4double GetRA()
Definition: G4hhElastic.hh:196
static const double MeV
G4double GetRB()
Definition: G4hhElastic.hh:200
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetRA()

void G4hhElastic::SetRA ( G4double  rn,
G4double  pq,
G4double  pQ 
)
inline

Definition at line 489 of file G4hhElastic.hh.

490 {
491  fRA = rA;
492  fRq = fRA*pq;
493  fRQ = fRA*pQ;
494 }
G4double fRQ
Definition: G4hhElastic.hh:132
G4double fRq
Definition: G4hhElastic.hh:133
G4double fRA
Definition: G4hhElastic.hh:131
Here is the caller graph for this function:

◆ SetRB()

void G4hhElastic::SetRB ( G4double  rn,
G4double  pq,
G4double  pQ 
)
inline

Definition at line 500 of file G4hhElastic.hh.

501 {
502  fRB = rB;
503  fRg = fRB*pg;
504  fRG = fRB*pG;
505 }
G4double fRG
Definition: G4hhElastic.hh:138
G4double fRg
Definition: G4hhElastic.hh:139
G4double fRB
Definition: G4hhElastic.hh:137
Here is the caller graph for this function:

◆ SetSigmaTot()

void G4hhElastic::SetSigmaTot ( G4double  stot)
inline

Definition at line 166 of file G4hhElastic.hh.

166 {fSigmaTot = stot;};
G4double fSigmaTot
Definition: G4hhElastic.hh:154

◆ SetSpp()

void G4hhElastic::SetSpp ( G4double  spp)
inline

Definition at line 167 of file G4hhElastic.hh.

167 {fSpp = spp;};
G4double fSpp
Definition: G4hhElastic.hh:159

Member Data Documentation

◆ fAlpha

G4double G4hhElastic::fAlpha
private

Definition at line 134 of file G4hhElastic.hh.

◆ fAlphaP

G4double G4hhElastic::fAlphaP
private

Definition at line 143 of file G4hhElastic.hh.

◆ fBankT

std::vector<G4PhysicsTable*> G4hhElastic::fBankT
private

Definition at line 115 of file G4hhElastic.hh.

◆ fBeta

G4double G4hhElastic::fBeta
private

Definition at line 135 of file G4hhElastic.hh.

◆ fBinT

G4int G4hhElastic::fBinT
private

Definition at line 111 of file G4hhElastic.hh.

◆ fBq

G4double G4hhElastic::fBq
private

Definition at line 155 of file G4hhElastic.hh.

◆ fBQ

G4double G4hhElastic::fBQ
private

Definition at line 156 of file G4hhElastic.hh.

◆ fBqQ

G4double G4hhElastic::fBqQ
private

Definition at line 157 of file G4hhElastic.hh.

◆ fCofF2

G4double G4hhElastic::fCofF2
private

Definition at line 148 of file G4hhElastic.hh.

◆ fCofF3

G4double G4hhElastic::fCofF3
private

Definition at line 149 of file G4hhElastic.hh.

◆ fDelta

G4double G4hhElastic::fDelta
private

Definition at line 141 of file G4hhElastic.hh.

◆ fEnergyBin

G4int G4hhElastic::fEnergyBin
private

Definition at line 110 of file G4hhElastic.hh.

◆ fEnergyVector

G4PhysicsLogVector* G4hhElastic::fEnergyVector
private

Definition at line 113 of file G4hhElastic.hh.

◆ fEta

G4double G4hhElastic::fEta
private

Definition at line 146 of file G4hhElastic.hh.

◆ fExpSlope

G4double G4hhElastic::fExpSlope
private

Definition at line 151 of file G4hhElastic.hh.

◆ fGamma

G4double G4hhElastic::fGamma
private

Definition at line 140 of file G4hhElastic.hh.

◆ fHadrNuclXsc

G4HadronNucleonXsc* G4hhElastic::fHadrNuclXsc
private

Definition at line 249 of file G4hhElastic.hh.

◆ fImCof

G4double G4hhElastic::fImCof
private

Definition at line 147 of file G4hhElastic.hh.

◆ fInTkin

G4int G4hhElastic::fInTkin
private

Definition at line 245 of file G4hhElastic.hh.

◆ fLambda

G4double G4hhElastic::fLambda
private

Definition at line 145 of file G4hhElastic.hh.

◆ fLambdaFF

G4double G4hhElastic::fLambdaFF
private

Definition at line 144 of file G4hhElastic.hh.

◆ fMassDif2

G4double G4hhElastic::fMassDif2
private

Definition at line 128 of file G4hhElastic.hh.

◆ fMassProj

G4double G4hhElastic::fMassProj
private

Definition at line 126 of file G4hhElastic.hh.

◆ fMassSum2

G4double G4hhElastic::fMassSum2
private

Definition at line 127 of file G4hhElastic.hh.

◆ fMassTarg

G4double G4hhElastic::fMassTarg
private

Definition at line 125 of file G4hhElastic.hh.

◆ fMff2

G4double G4hhElastic::fMff2
private

Definition at line 121 of file G4hhElastic.hh.

◆ fMQ

G4double G4hhElastic::fMQ
private

Definition at line 122 of file G4hhElastic.hh.

◆ fMq

G4double G4hhElastic::fMq
private

Definition at line 123 of file G4hhElastic.hh.

◆ fOldTkin

G4double G4hhElastic::fOldTkin
private

Definition at line 246 of file G4hhElastic.hh.

◆ fOptRatio

G4double G4hhElastic::fOptRatio
private

Definition at line 158 of file G4hhElastic.hh.

◆ fPcms

G4double G4hhElastic::fPcms
private

Definition at line 160 of file G4hhElastic.hh.

◆ fProjectile

G4ParticleDefinition* G4hhElastic::fProjectile
private

Definition at line 97 of file G4hhElastic.hh.

◆ fQcof

G4double G4hhElastic::fQcof
private

Definition at line 161 of file G4hhElastic.hh.

◆ fRA

G4double G4hhElastic::fRA
private

Definition at line 131 of file G4hhElastic.hh.

◆ fRB

G4double G4hhElastic::fRB
private

Definition at line 137 of file G4hhElastic.hh.

◆ fRG

G4double G4hhElastic::fRG
private

Definition at line 138 of file G4hhElastic.hh.

◆ fRg

G4double G4hhElastic::fRg
private

Definition at line 139 of file G4hhElastic.hh.

◆ fRhoReIm

G4double G4hhElastic::fRhoReIm
private

Definition at line 150 of file G4hhElastic.hh.

◆ fRQ

G4double G4hhElastic::fRQ
private

Definition at line 132 of file G4hhElastic.hh.

◆ fRq

G4double G4hhElastic::fRq
private

Definition at line 133 of file G4hhElastic.hh.

◆ fSigmaTot

G4double G4hhElastic::fSigmaTot
private

Definition at line 154 of file G4hhElastic.hh.

◆ fSo

G4double G4hhElastic::fSo
private

Definition at line 152 of file G4hhElastic.hh.

◆ fSpp

G4double G4hhElastic::fSpp
private

Definition at line 159 of file G4hhElastic.hh.

◆ fTableT

G4PhysicsTable* G4hhElastic::fTableT
private

Definition at line 114 of file G4hhElastic.hh.

◆ fTarget

G4ParticleDefinition* G4hhElastic::fTarget
private

Definition at line 96 of file G4hhElastic.hh.

◆ lowEnergyLimitHE

G4double G4hhElastic::lowEnergyLimitHE
private

Definition at line 105 of file G4hhElastic.hh.

◆ lowEnergyLimitQ

G4double G4hhElastic::lowEnergyLimitQ
private

Definition at line 106 of file G4hhElastic.hh.

◆ lowEnergyRecoilLimit

G4double G4hhElastic::lowEnergyRecoilLimit
private

Definition at line 104 of file G4hhElastic.hh.

◆ lowestEnergyLimit

G4double G4hhElastic::lowestEnergyLimit
private

Definition at line 107 of file G4hhElastic.hh.

◆ plabLowLimit

G4double G4hhElastic::plabLowLimit
private

Definition at line 108 of file G4hhElastic.hh.

◆ theNeutron

G4ParticleDefinition* G4hhElastic::theNeutron
private

Definition at line 99 of file G4hhElastic.hh.

◆ theNuclNuclData

const G4double G4hhElastic::theNuclNuclData
staticprivate
Initial value:
=
{
{ 2.76754, 4.8, 4.8, 0.05, 0.742441, 10.5 },
{ 3.07744, 5.4, 5.4, 0.02, 0.83818, 6.5 },
{ 3.36305, 5.2, 5.2, 0.02, 0.838893, 7.5 },
{ 4.32941, 6, 6, 0.03, 0.769389, 7.5 },
{ 4.62126, 6, 6, 0.03, 0.770111, 6.5 },
{ 5.47416, 4.5, 4.5, 0.03, 0.813185, 7.5 },
{ 6.15088, 6.5, 6.5, 0.02, 0.799539, 6.5 },
{ 6.77474, 5.2, 5.2, 0.03, 0.784901, 7.5 },
{ 9.77775, 7, 7, 0.03, 0.742531, 6.5 },
{ 10.4728, 5.2, 5.2, 0.03, 0.780439, 7.5 },
{ 13.7631, 7, 7, 0.008, 0.8664, 5.0 },
{ 19.4184, 6.8, 6.8, 0.009, 0.861337, 2.5 },
{ 23.5, 6.8, 6.8, 0.007, 0.878112, 1.5 },
{ 24.1362, 7.2, 7.2, 0.008, 0.864745, 5.5 },
{ 52.8, 6.8, 6.8, 0.008, 0.871929, 1.5 },
{ 546, 7.4, 7.4, 0.013, 0.845877, 5.5 },
{ 1960, 7.8, 7.8, 0.022, 0.809062, 7.5 },
{ 7000, 8, 8, 0.024, 0.820441, 5.5 }
}

Definition at line 247 of file G4hhElastic.hh.

◆ thePiKaNuclData

const G4double G4hhElastic::thePiKaNuclData
staticprivate
Initial value:
=
{
{ 2.5627, 3.8, 3.3, 0.22, 0.222, 1.5 },
{ 2.93928, 4.3, 3.8, 0.2, 0.250601, 1.3 },
{ 3.22326, 4.8, 4.3, 0.13, 0.32751, 2.5 },
{ 7.80704, 5.5, 5, 0.13, 0.340631, 2.5 },
{ 9.7328, 5, 4.5, 0.05, 0.416319, 5.5 },
{ 13.7315, 5.3, 4.8, 0.05, 0.418426, 5.5 },
{ 16.6359, 6.3, 5.8, 0.05, 0.423817, 5.5 },
{ 19.3961, 5, 4.5, 0.05, 0.413477, 3.5 }
}

Definition at line 248 of file G4hhElastic.hh.

◆ thePionMinus

G4ParticleDefinition* G4hhElastic::thePionMinus
private

Definition at line 101 of file G4hhElastic.hh.

◆ thePionPlus

G4ParticleDefinition* G4hhElastic::thePionPlus
private

Definition at line 100 of file G4hhElastic.hh.

◆ theProton

G4ParticleDefinition* G4hhElastic::theProton
private

Definition at line 98 of file G4hhElastic.hh.


The documentation for this class was generated from the following files: