Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EmCorrections Class Reference

#include <G4EmCorrections.hh>

Public Member Functions

 G4EmCorrections (G4int verb)
 
virtual ~G4EmCorrections ()
 
G4double HighOrderCorrections (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4double cutEnergy)
 
G4double IonHighOrderCorrections (const G4ParticleDefinition *, const G4MaterialCutsCouple *, G4double kineticEnergy)
 
G4double ComputeIonCorrections (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double IonBarkasCorrection (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double Bethe (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double SpinCorrection (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double KShellCorrection (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double LShellCorrection (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double ShellCorrection (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double ShellCorrectionSTD (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double DensityCorrection (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double BarkasCorrection (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double BlochCorrection (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double MottCorrection (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double NuclearDEDX (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy, G4bool fluct=true)
 
void AddStoppingData (G4int Z, G4int A, const G4String &materialName, G4PhysicsVector *dVector)
 
void InitialiseForNewRun ()
 
G4double EffectiveChargeCorrection (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
G4double EffectiveChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
void SetIonisationModels (G4VEmModel *m1=nullptr, G4VEmModel *m2=nullptr)
 
G4int GetNumberOfStoppingVectors () const
 
void SetVerbose (G4int verb)
 

Detailed Description

Definition at line 71 of file G4EmCorrections.hh.

Constructor & Destructor Documentation

G4EmCorrections::G4EmCorrections ( G4int  verb)
explicit

Definition at line 111 of file G4EmCorrections.cc.

112 {
113  particle = nullptr;
114  curParticle= nullptr;
115  material = nullptr;
116  curMaterial= nullptr;
117  theElementVector = nullptr;
118  atomDensity= nullptr;
119  curVector = nullptr;
120  ionLEModel = nullptr;
121  ionHEModel = nullptr;
122 
123  kinEnergy = 0.0;
124  verbose = verb;
125  massFactor = 1.0;
126  eth = 2.0*CLHEP::MeV;
127  nbinCorr = 20;
128  eCorrMin = 25.*CLHEP::keV;
129  eCorrMax = 250.*CLHEP::MeV;
130 
132  g4calc = G4Pow::GetInstance();
133 
134  nIons = ncouples = numberOfElements = idx = currentZ = 0;
135  mass = tau = gamma = bg2 = beta2 = beta = ba2 = tmax = charge = q2 = 0.0;
136 
137  // Constants
138  theZieglerFactor = CLHEP::eV*CLHEP::cm2*1.0e-15;
140  lossFlucFlag = true;
141 
142  // G.S. Khandelwal Nucl. Phys. A116(1968)97 - 111.
143  // "Shell corrections for K- and L- electrons
144 
145  nK = 20;
146  nL = 26;
147  nEtaK = 29;
148  nEtaL = 28;
149 
150  isMaster = false;
151 
152  // fill vectors
153  if(BarkasCorr == nullptr) { Initialise(); }
154 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
static constexpr double cm2
static constexpr double keV
G4IonTable * GetIonTable() const
static constexpr double MeV
static constexpr double eV
static G4ParticleTable * GetParticleTable()
static constexpr double fine_structure_const

Here is the call graph for this function:

G4EmCorrections::~G4EmCorrections ( )
virtual

Definition at line 158 of file G4EmCorrections.cc.

159 {
160  for(G4int i=0; i<nIons; ++i) {delete stopData[i];}
161  if(isMaster) {
162  delete BarkasCorr;
163  delete ThetaK;
164  delete ThetaL;
165  BarkasCorr = ThetaK = ThetaL = nullptr;
166  }
167 }
int G4int
Definition: G4Types.hh:78

Member Function Documentation

void G4EmCorrections::AddStoppingData ( G4int  Z,
G4int  A,
const G4String materialName,
G4PhysicsVector dVector 
)

Definition at line 991 of file G4EmCorrections.cc.

994 {
995  G4int i = 0;
996  for(; i<nIons; ++i) {
997  if(Z == Zion[i] && A == Aion[i] && mname == materialName[i]) break;
998  }
999  if(i == nIons) {
1000  Zion.push_back(Z);
1001  Aion.push_back(A);
1002  materialName.push_back(mname);
1003  materialList.push_back(nullptr);
1004  ionList.push_back(nullptr);
1005  stopData.push_back(dVector);
1006  nIons++;
1007  if(verbose > 1) {
1008  G4cout << "AddStoppingData Z= " << Z << " A= " << A << " " << mname
1009  << " idx= " << i << G4endl;
1010  }
1011  }
1012 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
#define G4endl
Definition: G4ios.hh:61

Here is the caller graph for this function:

G4double G4EmCorrections::BarkasCorrection ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)

Definition at line 659 of file G4EmCorrections.cc.

662 {
663  // . Z^3 Barkas effect in the stopping power of matter for charged particles
664  // J.C Ashley and R.H.Ritchie
665  // Physical review B Vol.5 No.7 1 April 1972 pp. 2393-2397
666  // valid for kineticEnergy > 0.5 MeV
667 
668  SetupKinematics(p, mat, e);
669  G4double BarkasTerm = 0.0;
670 
671  for (G4int i = 0; i<numberOfElements; ++i) {
672 
673  G4double Z = (*theElementVector)[i]->GetZ();
674  G4int iz = G4lrint(Z);
675  if(iz == 47) {
676  BarkasTerm += atomDensity[i]*0.006812*G4Exp(-G4Log(beta)*0.9);
677  } else if(iz >= 64) {
678  BarkasTerm += atomDensity[i]*0.002833*G4Exp(-G4Log(beta)*1.2);
679  } else {
680 
681  G4double X = ba2 / Z;
682  G4double b = 1.3;
683  if(1 == iz) {
684  if(material->GetName() == "G4_lH2") { b = 0.6; }
685  else { b = 1.8; }
686  }
687  else if(2 == iz) { b = 0.6; }
688  else if(10 >= iz) { b = 1.8; }
689  else if(17 >= iz) { b = 1.4; }
690  else if(18 == iz) { b = 1.8; }
691  else if(25 >= iz) { b = 1.4; }
692  else if(50 >= iz) { b = 1.35;}
693 
694  G4double W = b/std::sqrt(X);
695 
696  G4double val = BarkasCorr->Value(W);
697  if(W > BarkasCorr->Energy(46)) {
698  val *= BarkasCorr->Energy(46)/W;
699  }
700  // G4cout << "i= " << i << " b= " << b << " W= " << W
701  // << " Z= " << Z << " X= " << X << " val= " << val<< G4endl;
702  BarkasTerm += val*atomDensity[i] / (std::sqrt(Z*X)*X);
703  }
704  }
705 
706  BarkasTerm *= 1.29*charge/material->GetTotNbOfAtomsPerVolume();
707 
708  return BarkasTerm;
709 }
const G4String & GetName() const
Definition: G4Material.hh:178
int G4int
Definition: G4Types.hh:78
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:209
int G4lrint(double ad)
Definition: templates.hh:163
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmCorrections::Bethe ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)

Definition at line 308 of file G4EmCorrections.cc.

311 {
312  SetupKinematics(p, mat, e);
313  G4double eexc = material->GetIonisation()->GetMeanExcitationEnergy();
314  G4double eexc2 = eexc*eexc;
315  G4double dedx = 0.5*G4Log(2.0*electron_mass_c2*bg2*tmax/eexc2)-beta2;
316  return dedx;
317 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
static constexpr double electron_mass_c2
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double GetMeanExcitationEnergy() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4EmCorrections::BlochCorrection ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)

Definition at line 713 of file G4EmCorrections.cc.

716 {
717  SetupKinematics(p, mat, e);
718 
719  G4double y2 = q2/ba2;
720 
721  G4double term = 1.0/(1.0 + y2);
722  G4double del;
723  G4double j = 1.0;
724  do {
725  j += 1.0;
726  del = 1.0/(j* (j*j + y2));
727  term += del;
728  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
729  } while (del > 0.01*term);
730 
731  G4double res = -y2*term;
732  return res;
733 }
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

G4double G4EmCorrections::ComputeIonCorrections ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)

Definition at line 226 of file G4EmCorrections.cc.

229 {
230  // . Z^3 Barkas effect in the stopping power of matter for charged particles
231  // J.C Ashley and R.H.Ritchie
232  // Physical review B Vol.5 No.7 1 April 1972 pagg. 2393-2397
233  // and ICRU49 report
234  // valid for kineticEnergy < 0.5 MeV
235  // Other corrections from S.P.Ahlen Rev. Mod. Phys., Vol 52, No1, 1980
236  SetupKinematics(p, mat, e);
237  if(tau <= 0.0) { return 0.0; }
238 
239  G4double Barkas = BarkasCorrection (p, mat, e);
240  G4double Bloch = BlochCorrection (p, mat, e);
241  G4double Mott = MottCorrection (p, mat, e);
242 
243  G4double sum = 2.0*(Barkas*(charge - 1.0)/charge + Bloch) + Mott;
244 
245  if(verbose > 1) {
246  G4cout << "EmCorrections: E(MeV)= " << e/MeV << " Barkas= " << Barkas
247  << " Bloch= " << Bloch << " Mott= " << Mott
248  << " Sum= " << sum << G4endl;
249  }
250  sum *= material->GetElectronDensity() * q2 * twopi_mc2_rcl2 /beta2;
251 
252  if(verbose > 1) { G4cout << " Sum= " << sum << G4endl; }
253  return sum;
254 }
static constexpr double twopi_mc2_rcl2
G4double MottCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double BlochCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4GLOB_DLL std::ostream G4cout
G4double GetElectronDensity() const
Definition: G4Material.hh:217
G4double BarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmCorrections::DensityCorrection ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)

Definition at line 632 of file G4EmCorrections.cc.

635 {
636  SetupKinematics(p, mat, e);
637 
638  G4double cden = material->GetIonisation()->GetCdensity();
639  G4double mden = material->GetIonisation()->GetMdensity();
640  G4double aden = material->GetIonisation()->GetAdensity();
641  G4double x0den = material->GetIonisation()->GetX0density();
642  G4double x1den = material->GetIonisation()->GetX1density();
643 
644  G4double dedx = 0.0;
645 
646  // density correction
647  static const G4double twoln10 = 2.0*G4Log(10.0);
648  G4double x = G4Log(bg2)/twoln10;
649  if ( x >= x0den ) {
650  dedx = twoln10*x - cden ;
651  if ( x < x1den ) dedx += aden*G4Exp(G4Log(x1den-x)*mden) ;
652  }
653 
654  return dedx;
655 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4double GetAdensity() const
G4double GetX1density() const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double GetX0density() const
G4double GetCdensity() const
double G4double
Definition: G4Types.hh:76
G4double GetMdensity() const

Here is the call graph for this function:

G4double G4EmCorrections::EffectiveChargeCorrection ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)

Definition at line 942 of file G4EmCorrections.cc.

945 {
946  G4double factor = 1.0;
947  if(p->GetPDGCharge() <= 2.5*CLHEP::eplus || nIons <= 0) { return factor; }
948 
949  if(verbose > 1) {
950  G4cout << "EffectiveChargeCorrection: " << p->GetParticleName()
951  << " in " << mat->GetName()
952  << " ekin(MeV)= " << ekin/MeV << G4endl;
953  }
954 
955  if(p != curParticle || mat != curMaterial) {
956  curParticle = p;
957  curMaterial = mat;
958  curVector = nullptr;
959  currentZ = p->GetAtomicNumber();
960  if(verbose > 1) {
961  G4cout << "G4EmCorrections::EffectiveChargeCorrection: Zion= "
962  << currentZ << " Aion= " << p->GetPDGMass()/amu_c2 << G4endl;
963  }
964  massFactor = proton_mass_c2/p->GetPDGMass();
965  idx = -1;
966 
967  for(G4int i=0; i<nIons; ++i) {
968  if(materialList[i] == mat && currentZ == Zion[i]) {
969  idx = i;
970  break;
971  }
972  }
973  //G4cout << " idx= " << idx << " dz= " << G4endl;
974  if(idx >= 0) {
975  if(!ionList[idx]) { BuildCorrectionVector(); }
976  if(ionList[idx]) { curVector = stopData[idx]; }
977  } else { return factor; }
978  }
979  if(curVector) {
980  factor = curVector->Value(ekin*massFactor);
981  if(verbose > 1) {
982  G4cout << "E= " << ekin << " factor= " << factor << " massfactor= "
983  << massFactor << G4endl;
984  }
985  }
986  return factor;
987 }
static constexpr double proton_mass_c2
const char * p
Definition: xmltok.h:285
const G4String & GetName() const
Definition: G4Material.hh:178
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4int GetAtomicNumber() const
G4GLOB_DLL std::ostream G4cout
static constexpr double amu_c2
static constexpr double eplus
G4double Value(G4double theEnergy, size_t &lastidx) const
G4double GetPDGMass() const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmCorrections::EffectiveChargeSquareRatio ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)
inline

Definition at line 329 of file G4EmCorrections.hh.

332 {
333  return effCharge.EffectiveChargeSquareRatio(p,mat,kineticEnergy);
334 }
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *p, const G4Material *material, G4double kineticEnergy)

Here is the call graph for this function:

Here is the caller graph for this function:

G4int G4EmCorrections::GetNumberOfStoppingVectors ( ) const
inline

Definition at line 315 of file G4EmCorrections.hh.

316 {
317  return nIons;
318 }

Here is the caller graph for this function:

G4double G4EmCorrections::GetParticleCharge ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)
inline

Definition at line 321 of file G4EmCorrections.hh.

324 {
325  return effCharge.EffectiveCharge(p,mat,kineticEnergy);
326 }
G4double EffectiveCharge(const G4ParticleDefinition *p, const G4Material *material, G4double kineticEnergy)

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmCorrections::HighOrderCorrections ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy,
G4double  cutEnergy 
)

Definition at line 171 of file G4EmCorrections.cc.

174 {
175  // . Z^3 Barkas effect in the stopping power of matter for charged particles
176  // J.C Ashley and R.H.Ritchie
177  // Physical review B Vol.5 No.7 1 April 1972 pagg. 2393-2397
178  // and ICRU49 report
179  // valid for kineticEnergy < 0.5 MeV
180  // Other corrections from S.P.Ahlen Rev. Mod. Phys., Vol 52, No1, 1980
181 
182  SetupKinematics(p, mat, e);
183  if(tau <= 0.0) { return 0.0; }
184 
185  G4double Barkas = BarkasCorrection (p, mat, e);
186  G4double Bloch = BlochCorrection (p, mat, e);
187  G4double Mott = MottCorrection (p, mat, e);
188 
189  G4double sum = (2.0*(Barkas + Bloch) + Mott);
190 
191  if(verbose > 1) {
192  G4cout << "EmCorrections: E(MeV)= " << e/MeV << " Barkas= " << Barkas
193  << " Bloch= " << Bloch << " Mott= " << Mott
194  << " Sum= " << sum << " q2= " << q2 << G4endl;
195  G4cout << " ShellCorrection: " << ShellCorrection(p, mat, e)
196  << " Kshell= " << KShellCorrection(p, mat, e)
197  << " Lshell= " << LShellCorrection(p, mat, e)
198  << " " << mat->GetName() << G4endl;
199  }
200  sum *= material->GetElectronDensity() * q2 * twopi_mc2_rcl2 /beta2;
201  return sum;
202 }
static constexpr double twopi_mc2_rcl2
const G4String & GetName() const
Definition: G4Material.hh:178
G4double MottCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double BlochCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4GLOB_DLL std::ostream G4cout
G4double GetElectronDensity() const
Definition: G4Material.hh:217
G4double BarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double ShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double KShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double LShellCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmCorrections::InitialiseForNewRun ( )

Definition at line 1091 of file G4EmCorrections.cc.

1092 {
1094  ncouples = tb->GetTableSize();
1095  if(currmat.size() != ncouples) {
1096  currmat.resize(ncouples);
1097  for(std::map< G4int, std::vector<G4double> >::iterator it =
1098  thcorr.begin(); it != thcorr.end(); ++it){
1099  (it->second).clear();
1100  }
1101  thcorr.clear();
1102  for(size_t i=0; i<ncouples; ++i) {
1103  currmat[i] = tb->GetMaterialCutsCouple(i)->GetMaterial();
1104  G4String nam = currmat[i]->GetName();
1105  for(G4int j=0; j<nIons; ++j) {
1106  if(nam == materialName[j]) { materialList[j] = currmat[i]; }
1107  }
1108  }
1109  }
1110 }
int G4int
Definition: G4Types.hh:78
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
const G4Material * GetMaterial() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmCorrections::IonBarkasCorrection ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)

Definition at line 206 of file G4EmCorrections.cc.

209 {
210  // . Z^3 Barkas effect in the stopping power of matter for charged particles
211  // J.C Ashley and R.H.Ritchie
212  // Physical review B Vol.5 No.7 1 April 1972 pagg. 2393-2397
213  // and ICRU49 report
214  // valid for kineticEnergy < 0.5 MeV
215 
216  SetupKinematics(p, mat, e);
217  G4double res = 0.0;
218  if(tau > 0.0)
219  res = 2.0*BarkasCorrection(p, mat, e)*
220  material->GetElectronDensity() * q2 * twopi_mc2_rcl2 /beta2;
221  return res;
222 }
static constexpr double twopi_mc2_rcl2
G4double GetElectronDensity() const
Definition: G4Material.hh:217
G4double BarkasCorrection(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmCorrections::IonHighOrderCorrections ( const G4ParticleDefinition p,
const G4MaterialCutsCouple couple,
G4double  kineticEnergy 
)

Definition at line 258 of file G4EmCorrections.cc.

261 {
262  // . Z^3 Barkas effect in the stopping power of matter for charged particles
263  // J.C Ashley and R.H.Ritchie
264  // Physical review B Vol.5 No.7 1 April 1972 pagg. 2393-2397
265  // and ICRU49 report
266  // valid for kineticEnergy < 0.5 MeV
267  // Other corrections from S.P.Ahlen Rev. Mod. Phys., Vol 52, No1, 1980
268 
269  G4double sum = 0.0;
270 
271  if(ionHEModel) {
272  G4int Z = G4lrint(p->GetPDGCharge()*inveplus);
273  if(Z >= 100) Z = 99;
274  else if(Z < 1) Z = 1;
275 
276  G4double ethscaled = eth*p->GetPDGMass()/proton_mass_c2;
277  G4int ionPDG = p->GetPDGEncoding();
278  if(thcorr.find(ionPDG)==thcorr.end()) { // Not found: fill the map
279  std::vector<G4double> v;
280  for(size_t i=0; i<ncouples; ++i){
281  v.push_back(ethscaled*ComputeIonCorrections(p,currmat[i],ethscaled));
282  }
283  thcorr.insert(std::pair< G4int, std::vector<G4double> >(ionPDG,v));
284  }
285 
286  //G4cout << " map size=" << thcorr.size() << G4endl;
287  //for(std::map< G4int, std::vector<G4double> >::iterator
288  // it = thcorr.begin(); it != thcorr.end(); ++it){
289  // G4cout << "\t map element: first (key)=" << it->first
290  // << "\t second (vector): vec size=" << (it->second).size() << G4endl;
291  // for(size_t i=0; i<(it->second).size(); ++i){
292  // G4cout << "\t \t vec element: [" << i << "]=" << (it->second)[i]
293  //<< G4endl; } }
294 
295  G4double rest = (thcorr.find(ionPDG)->second)[couple->GetIndex()];
296 
297  sum = ComputeIonCorrections(p,couple->GetMaterial(),e) - rest/e;
298 
299  if(verbose > 1) {
300  G4cout << " Sum= " << sum << " dSum= " << rest/e << G4endl;
301  }
302  }
303  return sum;
304 }
static constexpr double proton_mass_c2
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4double GetPDGMass() const
int G4lrint(double ad)
Definition: templates.hh:163
#define G4endl
Definition: G4ios.hh:61
G4double ComputeIonCorrections(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
const G4Material * GetMaterial() const

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmCorrections::KShellCorrection ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)

Definition at line 332 of file G4EmCorrections.cc.

335 {
336  SetupKinematics(p, mat, e);
337  G4double term = 0.0;
338  for (G4int i = 0; i<numberOfElements; ++i) {
339 
340  G4double Z = (*theElementVector)[i]->GetZ();
341  G4int iz = G4lrint(Z);
342  G4double f = 1.0;
343  G4double Z2= (Z-0.3)*(Z-0.3);
344  if(1 == iz) {
345  f = 0.5;
346  Z2 = 1.0;
347  }
348  G4double eta = ba2/Z2;
349  G4double tet = Z2*(1. + Z2*0.25*alpha2);
350  if(11 < iz) { tet = ThetaK->Value(Z); }
351  term += f*atomDensity[i]*KShell(tet,eta)/Z;
352  }
353 
354  term /= material->GetTotNbOfAtomsPerVolume();
355 
356  return term;
357 }
int G4int
Definition: G4Types.hh:78
static G4double tet[DIM]
G4double Value(G4double theEnergy, size_t &lastidx) const
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:209
int G4lrint(double ad)
Definition: templates.hh:163
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmCorrections::LShellCorrection ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)

Definition at line 361 of file G4EmCorrections.cc.

364 {
365  SetupKinematics(p, mat, e);
366  G4double term = 0.0;
367  for (G4int i = 0; i<numberOfElements; ++i) {
368 
369  G4double Z = (*theElementVector)[i]->GetZ();
370  G4int iz = G4lrint(Z);
371  if(2 < iz) {
372  G4double Zeff = Z - ZD[10];
373  if(iz < 10) { Zeff = Z - ZD[iz]; }
374  G4double Z2= Zeff*Zeff;
375  G4double f = 0.125;
376  G4double eta = ba2/Z2;
377  G4double tet = ThetaL->Value(Z);
379  for(G4int j=1; j<nmax; ++j) {
381  if(15 >= iz) {
382  if(3 > j) { tet = 0.25*Z2*(1.0 + 5*Z2*alpha2/16.); }
383  else { tet = 0.25*Z2*(1.0 + Z2*alpha2/16.); }
384  }
385  //G4cout << " LShell: j= " << j << " ne= " << ne << " e(eV)= " << e/eV
386  // << " ThetaL= " << tet << G4endl;
387  term += f*ne*atomDensity[i]*LShell(tet,eta)/Z;
388  }
389  }
390  }
391 
392  term /= material->GetTotNbOfAtomsPerVolume();
393 
394  return term;
395 }
int G4int
Definition: G4Types.hh:78
static G4double tet[DIM]
G4double Value(G4double theEnergy, size_t &lastidx) const
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:209
int G4lrint(double ad)
Definition: templates.hh:163
static G4int GetNumberOfElectrons(G4int Z, G4int SubshellNb)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76
static G4int GetNumberOfShells(G4int Z)

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmCorrections::MottCorrection ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)

Definition at line 737 of file G4EmCorrections.cc.

740 {
741  SetupKinematics(p, mat, e);
742  G4double mterm = CLHEP::pi*fine_structure_const*beta*charge;
743  return mterm;
744 }
double G4double
Definition: G4Types.hh:76
static constexpr double fine_structure_const
static constexpr double pi
Definition: SystemOfUnits.h:54

Here is the caller graph for this function:

G4double G4EmCorrections::NuclearDEDX ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy,
G4bool  fluct = true 
)

Definition at line 748 of file G4EmCorrections.cc.

752 {
753  G4double nloss = 0.0;
754  if(e <= 0.0) return nloss;
755  SetupKinematics(p, mat, e);
756 
757  lossFlucFlag = fluct;
758 
759  // Projectile nucleus
760  G4double z1 = std::abs(particle->GetPDGCharge()*inveplus);
761  G4double mass1 = mass/amu_c2;
762 
763  // loop for the elements in the material
764  for (G4int iel=0; iel<numberOfElements; iel++) {
765  const G4Element* element = (*theElementVector)[iel] ;
766  G4double z2 = element->GetZ();
767  G4double mass2 = element->GetN();
768  nloss += (NuclearStoppingPower(kinEnergy, z1, z2, mass1, mass2))
769  * atomDensity[iel] ;
770  }
771  nloss *= theZieglerFactor;
772  return nloss;
773 }
G4double GetN() const
Definition: G4Element.hh:135
G4double GetZ() const
Definition: G4Element.hh:131
int G4int
Definition: G4Types.hh:78
static constexpr double amu_c2
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmCorrections::SetIonisationModels ( G4VEmModel m1 = nullptr,
G4VEmModel m2 = nullptr 
)
inline

Definition at line 309 of file G4EmCorrections.hh.

310 {
311  if(mod1) { ionLEModel = mod1; }
312  if(mod2) { ionHEModel = mod2; }
313 }

Here is the caller graph for this function:

void G4EmCorrections::SetVerbose ( G4int  verb)
inline

Definition at line 365 of file G4EmCorrections.hh.

366 {
367  verbose = verb;
368 }
G4double G4EmCorrections::ShellCorrection ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)

Definition at line 529 of file G4EmCorrections.cc.

532 {
533  SetupKinematics(p, mat, ekin);
534 
535  G4double term = 0.0;
536  //G4cout << "### G4EmCorrections::ShellCorrection " << mat->GetName()
537  // << " " << ekin/MeV << " MeV " << G4endl;
538  for (G4int i = 0; i<numberOfElements; ++i) {
539 
540  G4double res = 0.0;
541  G4double res0 = 0.0;
542  G4double Z = (*theElementVector)[i]->GetZ();
543  G4int iz = G4lrint(Z);
544  G4double Z2= (Z-0.3)*(Z-0.3);
545  G4double f = 1.0;
546  if(1 == iz) {
547  f = 0.5;
548  Z2 = 1.0;
549  }
550  G4double eta = ba2/Z2;
551  G4double tet = Z2*(1. + Z2*0.25*alpha2);
552  if(11 < iz) { tet = ThetaK->Value(Z); }
553  res0 = f*KShell(tet,eta);
554  res += res0;
555  //G4cout << " Z= " << iz << " Shell 0" << " tet= " << tet
556  // << " eta= " << eta << " resK= " << res0 << G4endl;
557  if(2 < iz) {
558  G4double Zeff = Z - ZD[10];
559  if(iz < 10) { Zeff = Z - ZD[iz]; }
560  Z2= Zeff*Zeff;
561  eta = ba2/Z2;
562  f = 0.125;
563  tet = ThetaL->Value(Z);
565  G4int nmax = std::min(4, ntot);
566  G4double norm = 0.0;
567  G4double eshell = 0.0;
568  for(G4int j=1; j<nmax; ++j) {
570  if(15 >= iz) {
571  if(3 > j) { tet = 0.25*Z2*(1.0 + 5*Z2*alpha2/16.); }
572  else { tet = 0.25*Z2*(1.0 + Z2*alpha2/16.); }
573  }
574  norm += ne;
575  eshell += tet*ne;
576  res0 = f*ne*LShell(tet,eta);
577  res += res0;
578  //G4cout << " Z= " << iz << " Shell " << j << " Ne= " << ne
579  // << " tet= " << tet << " eta= " << eta
580  // << " resL= " << res0 << G4endl;
581  }
582  if(ntot > nmax) {
583  eshell /= norm;
584 
585  static const G4double HM[53] = {
586  12.0, 12.0, 12.0, 12.0, 11.9, 11.7, 11.5, 11.2, 10.8, 10.4,
587  10.0, 9.51, 8.97, 8.52, 8.03, 7.46, 6.95, 6.53, 6.18, 5.87,
588  5.61, 5.39, 5.19, 5.01, 4.86, 4.72, 4.62, 4.53, 4.44, 4.38,
589  4.32, 4.26, 4.20, 4.15, 4.1, 4.04, 4.00, 3.95, 3.93, 3.91,
590  3.90, 3.89, 3.89, 3.88, 3.88, 3.88, 3.88, 3.88, 3.89, 3.89,
591  3.90, 3.92, 3.93 };
592  static const G4double HN[31] = {
593  75.5, 61.9, 52.2, 45.1, 39.6, 35.4, 31.9, 29.1, 27.2, 25.8,
594  24.5, 23.6, 22.7, 22.0, 21.4, 20.9, 20.5, 20.2, 19.9, 19.7,
595  19.5, 19.3, 19.2, 19.1, 18.4, 18.8, 18.7, 18.6, 18.5, 18.4,
596  18.2};
597 
598  // Add M-shell
599  if(28 > iz) {
600  res += f*(iz - 10)*LShell(eshell,HM[iz-11]*eta);
601  } else if(63 > iz) {
602  res += f*18*LShell(eshell,HM[iz-11]*eta);
603  } else {
604  res += f*18*LShell(eshell,HM[52]*eta);
605  }
606  // Add N-shell
607  if(32 < iz) {
608  if(60 > iz) {
609  res += f*(iz - 28)*LShell(eshell,HN[iz-33]*eta);
610  } else if(63 > iz) {
611  res += 4*LShell(eshell,HN[iz-33]*eta);
612  } else {
613  res += 4*LShell(eshell,HN[30]*eta);
614  }
615  // Add O-P-shells
616  if(60 < iz) {
617  res += f*(iz - 60)*LShell(eshell,150*eta);
618  }
619  }
620  }
621  }
622  term += res*atomDensity[i]/Z;
623  }
624 
625  term /= material->GetTotNbOfAtomsPerVolume();
626  //G4cout << "# Shell Correction= " << term << G4endl;
627  return term;
628 }
int G4int
Definition: G4Types.hh:78
static G4double tet[DIM]
G4double Value(G4double theEnergy, size_t &lastidx) const
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:209
int G4lrint(double ad)
Definition: templates.hh:163
static G4int GetNumberOfElectrons(G4int Z, G4int SubshellNb)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
double G4double
Definition: G4Types.hh:76
static G4int GetNumberOfShells(G4int Z)

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmCorrections::ShellCorrectionSTD ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)

Definition at line 495 of file G4EmCorrections.cc.

498 {
499  SetupKinematics(p, mat, e);
500  G4double taulim= 8.0*MeV/mass;
501  G4double bg2lim= taulim * (taulim+2.0);
502 
503  G4double* shellCorrectionVector =
505  G4double sh = 0.0;
506  G4double x = 1.0;
507  G4double taul = material->GetIonisation()->GetTaul();
508 
509  if ( bg2 >= bg2lim ) {
510  for (G4int k=0; k<3; ++k) {
511  x *= bg2 ;
512  sh += shellCorrectionVector[k]/x;
513  }
514 
515  } else {
516  for (G4int k=0; k<3; ++k) {
517  x *= bg2lim ;
518  sh += shellCorrectionVector[k]/x;
519  }
520  sh *= G4Log(tau/taul)/G4Log(taulim/taul);
521  }
522  sh *= 0.5;
523  return sh;
524 }
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:226
G4double GetTaul() const
int G4int
Definition: G4Types.hh:78
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double * GetShellCorrectionVector() const
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4EmCorrections::SpinCorrection ( const G4ParticleDefinition p,
const G4Material mat,
G4double  kineticEnergy 
)

Definition at line 321 of file G4EmCorrections.cc.

324 {
325  SetupKinematics(p, mat, e);
326  G4double dedx = 0.5*tmax/(kinEnergy + mass);
327  return 0.5*dedx*dedx;
328 }
double G4double
Definition: G4Types.hh:76

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