Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/src/G4eDPWACoulombScatteringModel.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 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 //
 27 // -------------------------------------------------------------------
 28 //
 29 // GEANT4 Class header file
 30 //
 31 //
 32 // File name:     G4eDPWACoulombScatteringModel
 33 //
 34 // Author:        Mihaly Novak
 35 //
 36 // Creation date: 02.07.2020
 37 //
 38 // Modifications:
 39 //
 40 // -------------------------------------------------------------------
 41 
 42 #include "G4eDPWACoulombScatteringModel.hh"
 43 
 44 #include "G4eDPWAElasticDCS.hh"
 45 #include "G4ParticleChangeForGamma.hh"
 46 #include "G4ParticleDefinition.hh"
 47 #include "G4DataVector.hh"
 48 
 49 #include "G4ProductionCutsTable.hh"
 50 #include "G4Material.hh"
 51 #include "G4Element.hh"
 52 #include "G4ElementVector.hh"
 53 
 54 #include "G4Electron.hh"
 55 
 56 #include "G4PhysicalConstants.hh"
 57 #include "G4SystemOfUnits.hh"
 58 #include "Randomize.hh"
 59 #include "G4ThreeVector.hh"
 60 
 61 
 62 G4eDPWACoulombScatteringModel::G4eDPWACoulombScatteringModel(G4bool ismixed, G4bool isscpcor, G4double mumin)
 63 : G4VEmModel("eDPWACoulombScattering"),
 64   fIsMixedModel(ismixed),
 65   fIsScpCorrection(isscpcor),
 66   fMuMin(mumin),
 67   fTheDCS(nullptr),
 68   fParticleChange(nullptr)
 69 {
 70   SetLowEnergyLimit (  0.0*CLHEP::eV);  // ekin = 10 eV   is used if (E< 10  eV)
 71   SetHighEnergyLimit(100.0*CLHEP::MeV); // ekin = 100 MeV is used if (E>100 MeV)
 72 }
 73 
 74 
 75 G4eDPWACoulombScatteringModel::~G4eDPWACoulombScatteringModel()
 76 {
 77   if (IsMaster()) {
 78     delete fTheDCS;
 79   }
 80 }
 81 
 82 
 83 void G4eDPWACoulombScatteringModel::Initialise(const G4ParticleDefinition* pdef,
 84                                                const G4DataVector& prodcuts)
 85 {
 86   if(!fParticleChange) {
 87     fParticleChange = GetParticleChangeForGamma();
 88   }
 89   fMuMin        = 0.5*(1.0-std::cos(PolarAngleLimit()));
 90   fIsMixedModel = (fMuMin > 0.0);
 91   if(IsMaster()) {
 92     // clean the G4eDPWAElasticDCS object if any
 93     delete fTheDCS;
 94     fTheDCS = new G4eDPWAElasticDCS(pdef==G4Electron::Electron(), fIsMixedModel);
 95     // init only for the elements that are used in the geometry
 96     G4ProductionCutsTable* theCpTable = G4ProductionCutsTable::GetProductionCutsTable();
 97     G4int numOfCouples = (G4int)theCpTable->GetTableSize();
 98     for(G4int j=0; j<numOfCouples; ++j) {
 99       const G4Material* mat = theCpTable->GetMaterialCutsCouple(j)->GetMaterial();
100       const G4ElementVector* elV = mat->GetElementVector();
101       std::size_t numOfElem = mat->GetNumberOfElements();
102       for (std::size_t ie = 0; ie < numOfElem; ++ie) {
103         fTheDCS->InitialiseForZ((*elV)[ie]->GetZasInt());
104       }
105     }
106     // init scattering power correction
107     if (fIsScpCorrection) {
108       fTheDCS->InitSCPCorrection(LowEnergyLimit(), HighEnergyLimit());
109     }
110     // will make use of the cross sections so the above needs to be done before
111     InitialiseElementSelectors(pdef, prodcuts);
112   }
113 }
114 
115 
116 void G4eDPWACoulombScatteringModel::InitialiseLocal(const G4ParticleDefinition*,
117                                                     G4VEmModel* masterModel)
118 {
119   SetElementSelectors(masterModel->GetElementSelectors());
120   SetTheDCS(static_cast<G4eDPWACoulombScatteringModel*>(masterModel)->GetTheDCS());
121 }
122 
123 
124 G4double
125 G4eDPWACoulombScatteringModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
126                                                           G4double ekin,
127                                                           G4double Z,
128                                                           G4double /*A*/,
129                                                           G4double /*prodcut*/,
130                                                           G4double /*emax*/)
131 {
132   // Cross sections are computed by numerical integration of the pre-computed
133   // DCS data between the muMin, muMax limits where mu(theta)=0.5[1-cos(theta)].
134   // In case of single scattering model (i.e. when fMuMin=0): [muMin=0, muMax=1]
135   // In case of mixed simulation model  (i.e. when fMuMin>0): [fMuMin , muMax=1]
136   // NOTE: cross sections will be zero if the kinetic enrgy is out of the
137   //       [10 eV-100 MeV] range for which DCS data has been computed.
138   //
139   G4double elCS  = 0.0;          // elastic cross section
140   G4double tr1CS = 0.0;          // first transport cross section
141   G4double tr2CS = 0.0;          // second transport cross section
142   const G4double muMin = fMuMin;
143   const G4double muMax = 1.0;
144   fTheDCS->ComputeCSPerAtom((G4int)Z, ekin, elCS, tr1CS, tr2CS, muMin, muMax);
145   // scattering power correction: should be only in condensed history ioni!
146   if (fIsScpCorrection && CurrentCouple()) {
147     const G4double theScpCor = fTheDCS->ComputeScatteringPowerCorrection(CurrentCouple(), ekin);
148     elCS *= (theScpCor*(1.0+1.0/Z));
149   }
150   return std::max(0.0, elCS);
151 }
152 
153 
154 void
155 G4eDPWACoulombScatteringModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
156                                                  const G4MaterialCutsCouple* cp,
157                                                  const G4DynamicParticle* dp,
158                                                  G4double, G4double)
159 {
160   const G4double    ekin   = dp->GetKineticEnergy();
161   const G4double    lekin  = dp->GetLogKineticEnergy();
162   const G4Element*  target = SelectTargetAtom(cp, dp->GetParticleDefinition(), ekin, lekin);
163   const G4int       izet   = target->GetZasInt();
164   // sample cosine of the polar scattering angle in (hard) elastic insteraction
165   CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
166   G4double cost = 1.0;
167   if (!fIsMixedModel) {
168     G4double rndm[3];
169     rndmEngine->flatArray(3, rndm);
170     cost = fTheDCS->SampleCosineTheta(izet, lekin, rndm[0], rndm[1], rndm[2]);
171   } else {
172     //sample cost between costMax,costMin where costMax = 1-2xfMuMin;
173     const G4double costMax = 1.0-2.0*fMuMin;
174     const G4double costMin = -1.0;
175     G4double rndm[2];
176     rndmEngine->flatArray(2, rndm);
177     cost = fTheDCS->SampleCosineThetaRestricted(izet, lekin, rndm[0], rndm[1], costMin, costMax);
178   }
179   // compute the new direction in the scattering frame
180   const G4double sint = std::sqrt((1.0-cost)*(1.0+cost));
181   const G4double phi  = CLHEP::twopi*rndmEngine->flat();
182   G4ThreeVector theNewDirection(sint*std::cos(phi), sint*std::sin(phi), cost);
183   // get original direction in lab frame and rotate new direction to lab frame
184   G4ThreeVector theOrgDirectionLab = dp->GetMomentumDirection();
185   theNewDirection.rotateUz(theOrgDirectionLab);
186   // set new direction
187   fParticleChange->ProposeMomentumDirection(theNewDirection);
188 }
189 
190 
191 
192 
193 
194