Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/src/G4eCoulombScatteringModel.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

Diff markup

Differences between /processes/electromagnetic/standard/src/G4eCoulombScatteringModel.cc (Version 11.3.0) and /processes/electromagnetic/standard/src/G4eCoulombScatteringModel.cc (Version 10.0.p3)


  1 //                                                  1 //
  2 // *******************************************      2 // ********************************************************************
  3 // * License and Disclaimer                         3 // * License and Disclaimer                                           *
  4 // *                                                4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of th      5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided      6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License      7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/      8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.           9 // * include a list of copyright holders.                             *
 10 // *                                               10 // *                                                                  *
 11 // * Neither the authors of this software syst     11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing fin     12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warran     13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assum     14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file      15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitatio     16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                               17 // *                                                                  *
 18 // * This  code  implementation is the result      18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboratio     19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distri     20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  ag     21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publicati     22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Sof     23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // *******************************************     24 // ********************************************************************
 25 //                                                 25 //
                                                   >>  26 // $Id: G4eCoulombScatteringModel.cc 79188 2014-02-20 09:22:48Z gcosmo $
 26 //                                                 27 //
 27 // -------------------------------------------     28 // -------------------------------------------------------------------
 28 //                                                 29 //
 29 // GEANT4 Class file                               30 // GEANT4 Class file
 30 //                                                 31 //
 31 //                                                 32 //
 32 // File name:     G4eCoulombScatteringModel        33 // File name:     G4eCoulombScatteringModel
 33 //                                                 34 //
 34 // Author:        Vladimir Ivanchenko              35 // Author:        Vladimir Ivanchenko 
 35 //                                                 36 //
 36 // Creation date: 22.08.2005                       37 // Creation date: 22.08.2005
 37 //                                                 38 //
 38 // Modifications: V.Ivanchenko                 <<  39 // Modifications:
 39 //                                                 40 //
                                                   >>  41 // 01.08.06 V.Ivanchenko extend upper limit of table to TeV and review the
                                                   >>  42 //          logic of building - only elements from G4ElementTable
                                                   >>  43 // 08.08.06 V.Ivanchenko build internal table in ekin scale, introduce faclim
                                                   >>  44 // 19.08.06 V.Ivanchenko add inline function ScreeningParameter 
                                                   >>  45 // 09.10.07 V.Ivanchenko reorganized methods, add cut dependence in scattering off e- 
                                                   >>  46 // 09.06.08 V.Ivanchenko add SelectIsotope and sampling of the recoil ion 
                                                   >>  47 // 16.06.09 C.Consolandi fixed computation of effective mass
                                                   >>  48 // 27.05.10 V.Ivanchenko added G4WentzelOKandVIxSection class to
                                                   >>  49 //              compute cross sections and sample scattering angle
 40 //                                                 50 //
 41 //                                                 51 //
 42 // Class Description:                              52 // Class Description:
 43 //                                                 53 //
 44 // -------------------------------------------     54 // -------------------------------------------------------------------
 45 //                                                 55 //
 46 //....oooOO0OOooo........oooOO0OOooo........oo     56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 47 //....oooOO0OOooo........oooOO0OOooo........oo     57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 48                                                    58 
 49 #include "G4eCoulombScatteringModel.hh"            59 #include "G4eCoulombScatteringModel.hh"
 50 #include "G4PhysicalConstants.hh"                  60 #include "G4PhysicalConstants.hh"
 51 #include "G4SystemOfUnits.hh"                      61 #include "G4SystemOfUnits.hh"
 52 #include "Randomize.hh"                            62 #include "Randomize.hh"
 53 #include "G4DataVector.hh"                         63 #include "G4DataVector.hh"
 54 #include "G4ElementTable.hh"                       64 #include "G4ElementTable.hh"
 55 #include "G4ParticleChangeForGamma.hh"             65 #include "G4ParticleChangeForGamma.hh"
 56 #include "G4Proton.hh"                             66 #include "G4Proton.hh"
 57 #include "G4ParticleTable.hh"                      67 #include "G4ParticleTable.hh"
 58 #include "G4IonTable.hh"                           68 #include "G4IonTable.hh"
 59 #include "G4ProductionCutsTable.hh"                69 #include "G4ProductionCutsTable.hh"
 60 #include "G4NucleiProperties.hh"                   70 #include "G4NucleiProperties.hh"
 61 #include "G4Pow.hh"                                71 #include "G4Pow.hh"
                                                   >>  72 #include "G4LossTableManager.hh"
                                                   >>  73 #include "G4LossTableBuilder.hh"
 62 #include "G4NistManager.hh"                        74 #include "G4NistManager.hh"
 63                                                    75 
 64 //....oooOO0OOooo........oooOO0OOooo........oo     76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 65                                                    77 
 66 using namespace std;                               78 using namespace std;
 67                                                    79 
 68 G4eCoulombScatteringModel::G4eCoulombScatterin <<  80 G4eCoulombScatteringModel::G4eCoulombScatteringModel(const G4String& nam)
 69   : G4VEmModel("eCoulombScattering"), isCombin <<  81   : G4VEmModel(nam),
                                                   >>  82     cosThetaMin(1.0),
                                                   >>  83     cosThetaMax(-1.0),
                                                   >>  84     isInitialised(false)
 70 {                                                  85 {
                                                   >>  86   fParticleChange = 0;
 71   fNistManager = G4NistManager::Instance();        87   fNistManager = G4NistManager::Instance();
 72   theIonTable  = G4ParticleTable::GetParticleT     88   theIonTable  = G4ParticleTable::GetParticleTable()->GetIonTable();
 73   theProton    = G4Proton::Proton();               89   theProton    = G4Proton::Proton();
                                                   >>  90   currentMaterial = 0; 
                                                   >>  91   fixedCut = -1.0;
 74                                                    92 
 75   wokvi = new G4WentzelOKandVIxSection(isCombi <<  93   pCuts = 0;
 76                                                    94 
 77   mass = CLHEP::proton_mass_c2;                <<  95   lowEnergyThreshold = 1*keV;  // particle will be killed for lower energy
                                                   >>  96   recoilThreshold = 0.*keV; // by default does not work
                                                   >>  97 
                                                   >>  98   particle = 0;
                                                   >>  99   currentCouple = 0;
                                                   >> 100   wokvi = new G4WentzelOKandVIxSection();
                                                   >> 101 
                                                   >> 102   currentMaterialIndex = 0;
                                                   >> 103 
                                                   >> 104   cosTetMinNuc = 1.0;
                                                   >> 105   cosTetMaxNuc = -1.0;
                                                   >> 106   elecRatio = 0.0;
                                                   >> 107   mass = proton_mass_c2;
 78 }                                                 108 }
 79                                                   109 
 80 //....oooOO0OOooo........oooOO0OOooo........oo    110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 81                                                   111 
 82 G4eCoulombScatteringModel::~G4eCoulombScatteri    112 G4eCoulombScatteringModel::~G4eCoulombScatteringModel()
 83 {                                                 113 {
 84   delete wokvi;                                   114   delete wokvi;
 85 }                                                 115 }
 86                                                   116 
 87 //....oooOO0OOooo........oooOO0OOooo........oo    117 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 88                                                   118 
 89 void G4eCoulombScatteringModel::Initialise(con    119 void G4eCoulombScatteringModel::Initialise(const G4ParticleDefinition* part,
 90              const G4DataVector& cuts)            120              const G4DataVector& cuts)
 91 {                                                 121 {
 92   SetupParticle(part);                            122   SetupParticle(part);
 93   currentCouple = nullptr;                     << 123   currentCouple = 0;
 94                                                << 124   cosThetaMin = cos(PolarAngleLimit());
 95   G4double tet = PolarAngleLimit();            << 
 96                                                << 
 97   // defined theta limit between single and mu << 
 98   if(isCombined) {                             << 
 99     if(tet >= CLHEP::pi) { cosThetaMin = -1.0; << 
100     else if(tet > 0.0) { cosThetaMin = std::co << 
101                                                << 
102     // single scattering without multiple      << 
103   } else if(tet > 0.0) {                       << 
104     cosThetaMin = std::cos(std::min(tet, CLHEP << 
105   }                                            << 
106                                                << 
107   wokvi->Initialise(part, cosThetaMin);           125   wokvi->Initialise(part, cosThetaMin);
                                                   >> 126   /*      
                                                   >> 127   G4cout << "G4eCoulombScatteringModel: " << particle->GetParticleName()
                                                   >> 128          << "  1-cos(ThetaLimit)= " << 1 - cosThetaMin
                                                   >> 129    << "  cos(thetaMax)= " <<  cosThetaMax
                                                   >> 130    << G4endl;
                                                   >> 131   */
108   pCuts = &cuts;                                  132   pCuts = &cuts;
                                                   >> 133   //  G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3);
109   /*                                              134   /*
110   G4cout << "G4eCoulombScatteringModel::Initia << 135   G4cout << "!!! G4eCoulombScatteringModel::Initialise for " 
111      << part->GetParticleName() << " 1-cos(Tet << 136      << part->GetParticleName() << "  cos(TetMin)= " << cosThetaMin 
112      << " 1-cos(TetMax)= " << 1. - cosThetaMax << 137      << "  cos(TetMax)= " << cosThetaMax <<G4endl;
113   G4cout << "cut[0]= " << (*pCuts)[0] << G4end << 138   G4cout << "cut= " << pCuts[0] << "  cut1= " << pCuts[1] << G4endl;
114   */                                              139   */
115   if(nullptr == fParticleChange) {             << 140   if(!isInitialised) {
                                                   >> 141     isInitialised = true;
116     fParticleChange = GetParticleChangeForGamm    142     fParticleChange = GetParticleChangeForGamma();
117   }                                               143   }
118   if(IsMaster() && mass < GeV && part->GetPart    144   if(IsMaster() && mass < GeV && part->GetParticleName() != "GenericIon") {
119     InitialiseElementSelectors(part, cuts);    << 145     InitialiseElementSelectors(part,cuts);
120   }                                               146   }
121 }                                                 147 }
122                                                   148 
123 //....oooOO0OOooo........oooOO0OOooo........oo    149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
124                                                   150 
125 void G4eCoulombScatteringModel::InitialiseLoca    151 void G4eCoulombScatteringModel::InitialiseLocal(const G4ParticleDefinition*, 
126             G4VEmModel* masterModel)              152             G4VEmModel* masterModel)
127 {                                                 153 {
128   SetElementSelectors(masterModel->GetElementS    154   SetElementSelectors(masterModel->GetElementSelectors());
129 }                                                 155 }
130                                                   156 
131 //....oooOO0OOooo........oooOO0OOooo........oo    157 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132                                                   158 
133 G4double                                          159 G4double 
134 G4eCoulombScatteringModel::MinPrimaryEnergy(co    160 G4eCoulombScatteringModel::MinPrimaryEnergy(const G4Material* material,
135               const G4ParticleDefinition* part    161               const G4ParticleDefinition* part,
136               G4double)                           162               G4double)
137 {                                                 163 {
138   SetupParticle(part);                            164   SetupParticle(part);
139                                                   165 
140   // define cut using cuts for proton             166   // define cut using cuts for proton
141   G4double cut =                                  167   G4double cut = 
142     std::max(recoilThreshold, (*pCuts)[Current    168     std::max(recoilThreshold, (*pCuts)[CurrentCouple()->GetIndex()]);
143                                                   169 
144   // find out lightest element                    170   // find out lightest element
145   const G4ElementVector* theElementVector = ma    171   const G4ElementVector* theElementVector = material->GetElementVector();
146   std::size_t nelm = material->GetNumberOfElem << 172   G4int nelm = material->GetNumberOfElements();
147                                                << 
148   // select lightest element                   << 
149   G4int Z = 300;                                  173   G4int Z = 300;
150   for (std::size_t j=0; j<nelm; ++j) {         << 174   for (G4int j=0; j<nelm; ++j) {        
151     Z = std::min(Z,(*theElementVector)[j]->Get << 175     G4int iz = (G4int)(*theElementVector)[j]->GetZ();
                                                   >> 176     if(iz < Z) { Z = iz; }
152   }                                               177   }
153   G4int A = G4lrint(fNistManager->GetAtomicMas    178   G4int A = G4lrint(fNistManager->GetAtomicMassAmu(Z));
154   G4double targetMass = G4NucleiProperties::Ge    179   G4double targetMass = G4NucleiProperties::GetNuclearMass(A, Z);
155   G4double t = std::max(cut, 0.5*(cut + sqrt(2    180   G4double t = std::max(cut, 0.5*(cut + sqrt(2*cut*targetMass)));
156                                                   181 
157   return t;                                    << 182   return std::max(lowEnergyThreshold, t);
158 }                                                 183 }
159                                                   184 
160 //....oooOO0OOooo........oooOO0OOooo........oo    185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
161                                                   186 
162 G4double G4eCoulombScatteringModel::ComputeCro    187 G4double G4eCoulombScatteringModel::ComputeCrossSectionPerAtom(
163                 const G4ParticleDefinition* p,    188                 const G4ParticleDefinition* p,
164     G4double kinEnergy,                           189     G4double kinEnergy,
165     G4double Z, G4double,                         190     G4double Z, G4double,
166     G4double cutEnergy, G4double)                 191     G4double cutEnergy, G4double)
167 {                                                 192 {
168   /*                                           << 193   //G4cout << "### G4eCoulombScatteringModel::ComputeCrossSectionPerAtom  for " 
169   G4cout << "### G4eCoulombScatteringModel::Co << 194   //<< p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl; 
170    << p->GetParticleName()<<" Z= "<<Z<<" e(MeV << 
171    << G4endl;                                  << 
172   */                                           << 
173   G4double cross = 0.0;                           195   G4double cross = 0.0;
174   elecRatio = 0.0;                             << 
175   if(p != particle) { SetupParticle(p); }         196   if(p != particle) { SetupParticle(p); }
176                                                   197 
177   // cross section is set to zero to avoid pro    198   // cross section is set to zero to avoid problems in sample secondary
178   if(kinEnergy <= 0.0) { return cross; }          199   if(kinEnergy <= 0.0) { return cross; }
179   DefineMaterial(CurrentCouple());                200   DefineMaterial(CurrentCouple());
180   G4double costmin = wokvi->SetupKinematic(kin << 201   cosTetMinNuc = wokvi->SetupKinematic(kinEnergy, currentMaterial);
181                                                << 202   if(cosThetaMax < cosTetMinNuc) {
182   //G4cout << "cosThetaMax= "<<cosThetaMax<<"  << 203     G4int iz = G4int(Z);
183                                                << 204     G4double cut = cutEnergy;
184   if(cosThetaMax < costmin) {                  << 205     if(fixedCut > 0.0) { cut = fixedCut; }
185     G4int iz = G4lrint(Z);                     << 206     cosTetMinNuc = wokvi->SetupTarget(iz, cut);
186     G4double cut = (0.0 < fixedCut) ? fixedCut << 207     cosTetMaxNuc = cosThetaMax; 
187     costmin = wokvi->SetupTarget(iz, cut);     << 208     if(iz == 1 && cosTetMaxNuc < 0.0 && particle == theProton) { 
188     //G4cout << "SetupTarget: Z= " << iz << "  << 209       cosTetMaxNuc = 0.0; 
189     //     << costmin << G4endl;               << 
190     G4double costmax = (1 == iz && particle == << 
191       ? 0.0 : cosThetaMax;                     << 
192     if(costmin > costmax) {                    << 
193       cross = wokvi->ComputeNuclearCrossSectio << 
194         + wokvi->ComputeElectronCrossSection(c << 
195     }                                             210     }
196     /*                                         << 211     cross =  wokvi->ComputeNuclearCrossSection(cosTetMinNuc, cosTetMaxNuc);
197     if(p->GetParticleName() == "e-")           << 212     elecRatio = wokvi->ComputeElectronCrossSection(cosTetMinNuc, cosThetaMax);
198     G4cout << "Z= " << Z << " e(MeV)= " << kin << 213     cross += elecRatio;
199      << " cross(b)= " << cross/barn << " 1-cos << 214     if(cross > 0.0) { elecRatio /= cross; }  
200      << " 1-costmax= " << 1-costmax            << 
201      << " 1-cosThetaMax= " << 1-cosThetaMax    << 
202      << "  " << currentMaterial->GetName()     << 
203      << G4endl;                                << 
204     */                                         << 
205   }                                               215   }
206   //G4cout << "====== cross= " << cross << G4e << 216   /*
                                                   >> 217   if(p->GetParticleName() == "e-") 
                                                   >> 218   G4cout << "e(MeV)= " << kinEnergy/MeV << " cross(b)= " << cross/barn  
                                                   >> 219    << " 1-cosTetMinNuc= " << 1-cosTetMinNuc
                                                   >> 220    << " 1-cosTetMaxNuc= " << 1-cosTetMaxNuc
                                                   >> 221    << G4endl;
                                                   >> 222   */
207   return cross;                                   223   return cross;  
208 }                                                 224 }
209                                                   225 
210 //....oooOO0OOooo........oooOO0OOooo........oo    226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
211                                                   227 
212 void G4eCoulombScatteringModel::SampleSecondar    228 void G4eCoulombScatteringModel::SampleSecondaries(
213                 std::vector<G4DynamicParticle*    229                 std::vector<G4DynamicParticle*>* fvect,
214     const G4MaterialCutsCouple* couple,           230     const G4MaterialCutsCouple* couple,
215     const G4DynamicParticle* dp,                  231     const G4DynamicParticle* dp,
216     G4double cutEnergy,                           232     G4double cutEnergy,
217     G4double)                                     233     G4double)
218 {                                                 234 {
219   G4double kinEnergy = dp->GetKineticEnergy();    235   G4double kinEnergy = dp->GetKineticEnergy();
                                                   >> 236 
                                                   >> 237   // absorb particle below low-energy limit to avoid situation
                                                   >> 238   // when a particle has no energy loss
                                                   >> 239   if(kinEnergy < lowEnergyThreshold) { 
                                                   >> 240     fParticleChange->SetProposedKineticEnergy(0.0);
                                                   >> 241     fParticleChange->ProposeLocalEnergyDeposit(kinEnergy);
                                                   >> 242     fParticleChange->ProposeNonIonizingEnergyDeposit(kinEnergy);
                                                   >> 243     return; 
                                                   >> 244   }
220   SetupParticle(dp->GetDefinition());             245   SetupParticle(dp->GetDefinition());
221   DefineMaterial(couple);                         246   DefineMaterial(couple);
222   /*                                              247   /*
223   G4cout << "G4eCoulombScatteringModel::Sample    248   G4cout << "G4eCoulombScatteringModel::SampleSecondaries e(MeV)= " 
224      << kinEnergy << "  " << particle->GetPart    249      << kinEnergy << "  " << particle->GetParticleName() 
225      << " cut= " << cutEnergy<< G4endl;           250      << " cut= " << cutEnergy<< G4endl;
226   */                                              251   */
227   // Choose nucleus                               252   // Choose nucleus
228   G4double cut = (0.0 < fixedCut) ? fixedCut : << 253   G4double cut = cutEnergy;
                                                   >> 254   if(fixedCut > 0.0) { cut = fixedCut; }
229                                                   255 
230   wokvi->SetupKinematic(kinEnergy, currentMate << 256   const G4Element* currentElement = 
                                                   >> 257     SelectRandomAtom(couple,particle,kinEnergy,cut,kinEnergy);
231                                                   258 
232   const G4Element* currentElement = SelectTarg << 259   G4double Z = currentElement->GetZ();
233                                        dp->Get << 
234   G4int iz = currentElement->GetZasInt();      << 
235                                                << 
236   G4double costmin = wokvi->SetupTarget(iz, cu << 
237   G4double costmax = (1 == iz && particle == t << 
238     ? 0.0 :  cosThetaMax;                      << 
239   if(costmin <= costmax) { return; }           << 
240                                                << 
241   G4double cross = wokvi->ComputeNuclearCrossS << 
242   G4double ecross = wokvi->ComputeElectronCros << 
243   G4double ratio = ecross/(cross + ecross);    << 
244                                                   260 
                                                   >> 261   if(ComputeCrossSectionPerAtom(particle,kinEnergy, Z,
                                                   >> 262           kinEnergy, cut, kinEnergy) == 0.0) 
                                                   >> 263     { return; }
                                                   >> 264 
                                                   >> 265   G4int iz = G4int(Z);
245   G4int ia = SelectIsotopeNumber(currentElemen    266   G4int ia = SelectIsotopeNumber(currentElement);
246   G4double targetMass = G4NucleiProperties::Ge    267   G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz);
247   wokvi->SetTargetMass(targetMass);               268   wokvi->SetTargetMass(targetMass);
248                                                   269 
249   G4ThreeVector newDirection =                    270   G4ThreeVector newDirection = 
250     wokvi->SampleSingleScattering(costmin, cos << 271     wokvi->SampleSingleScattering(cosTetMinNuc, cosThetaMax, elecRatio);
251   G4double cost = newDirection.z();               272   G4double cost = newDirection.z();
252     /*                                         << 273 
253       G4cout << "SampleSec: e(MeV)= " << kinEn << 
254              << " 1-costmin= " << 1-costmin    << 
255              << " 1-costmax= " << 1-costmax    << 
256              << " 1-cost= " << 1-cost          << 
257              << " ratio= " << ratio            << 
258              << G4endl;                        << 
259     */                                         << 
260   G4ThreeVector direction = dp->GetMomentumDir    274   G4ThreeVector direction = dp->GetMomentumDirection(); 
261   newDirection.rotateUz(direction);               275   newDirection.rotateUz(direction);   
262                                                   276 
263   fParticleChange->ProposeMomentumDirection(ne    277   fParticleChange->ProposeMomentumDirection(newDirection);   
264                                                   278 
265   // recoil sampling assuming a small recoil      279   // recoil sampling assuming a small recoil
266   // and first order correction to primary 4-m    280   // and first order correction to primary 4-momentum
267   G4double mom2 = wokvi->GetMomentumSquare();     281   G4double mom2 = wokvi->GetMomentumSquare();
268   G4double trec = mom2*(1.0 - cost)               282   G4double trec = mom2*(1.0 - cost)
269     /(targetMass + (mass + kinEnergy)*(1.0 - c    283     /(targetMass + (mass + kinEnergy)*(1.0 - cost));
270                                                   284 
271   // the check likely not needed                  285   // the check likely not needed
272   trec = std::min(trec, kinEnergy);            << 286   if(trec > kinEnergy) { trec = kinEnergy; }
273   G4double finalT = kinEnergy - trec;             287   G4double finalT = kinEnergy - trec; 
274   G4double edep = 0.0;                            288   G4double edep = 0.0;
275     /*                                         << 289   //G4cout<<"G4eCoulombScatteringModel: finalT= "<<finalT<<" Trec= "
276     G4cout<<"G4eCoulombScatteringModel: finalT << 290   //  <<trec << " Z= " << iz << " A= " << ia<<G4endl;
277     <<trec << " Z= " << iz << " A= " << ia     << 291 
278     << " tcut(keV)= " << (*pCuts)[currentMater << 
279     */                                         << 
280   G4double tcut = recoilThreshold;                292   G4double tcut = recoilThreshold;
281   if(pCuts) { tcut= std::max(tcut,(*pCuts)[cur    293   if(pCuts) { tcut= std::max(tcut,(*pCuts)[currentMaterialIndex]); }
282                                                   294 
283   if(trec > tcut) {                               295   if(trec > tcut) {
284     G4ParticleDefinition* ion = theIonTable->G    296     G4ParticleDefinition* ion = theIonTable->GetIon(iz, ia, 0);
285     G4ThreeVector dir = (direction*sqrt(mom2)     297     G4ThreeVector dir = (direction*sqrt(mom2) - 
286        newDirection*sqrt(finalT*(2*mass + fina    298        newDirection*sqrt(finalT*(2*mass + finalT))).unit();
287     auto newdp = new G4DynamicParticle(ion, di << 299     G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, trec);
288     fvect->push_back(newdp);                      300     fvect->push_back(newdp);
289   } else {                                        301   } else {
290     edep = trec;                                  302     edep = trec;
291     fParticleChange->ProposeNonIonizingEnergyD    303     fParticleChange->ProposeNonIonizingEnergyDeposit(edep);
292   }                                               304   }
293                                                   305 
294     // finelize primary energy and energy bala << 306   // finelize primary energy and energy balance
295     // this threshold may be applied only beca << 307   // this threshold may be applied only because for low-enegry
296     // e+e- msc model is applied               << 308   // e+e- msc model is applied
297   if(finalT < 0.0) {                           << 309   if(finalT <= lowEnergyThreshold) { 
298     edep += finalT;                               310     edep += finalT;  
299     finalT = 0.0;                                 311     finalT = 0.0;
300   }                                               312   } 
301   edep = std::max(edep, 0.0);                  << 
302   fParticleChange->SetProposedKineticEnergy(fi    313   fParticleChange->SetProposedKineticEnergy(finalT);
303   fParticleChange->ProposeLocalEnergyDeposit(e    314   fParticleChange->ProposeLocalEnergyDeposit(edep);
304 }                                                 315 }
305                                                   316 
306 //....oooOO0OOooo........oooOO0OOooo........oo    317 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
                                                   >> 318 
                                                   >> 319 
307                                                   320