Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/highenergy/src/G4DynamicParticleMSC.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 file
 30 //
 31 //
 32 // File name:     G4DynamicParticleMSC
 33 //
 34 // Author:        Vladimir Ivanchenko
 35 //
 36 // Creation date: 17.08.2024
 37 //
 38 // -------------------------------------------------------------------
 39 //
 40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 41 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 42 
 43 #include "G4DynamicParticleMSC.hh"
 44 #include "G4PhysicalConstants.hh"
 45 #include "G4SystemOfUnits.hh"
 46 #include "G4EmProcessSubType.hh"
 47 #include "G4LossTableManager.hh"
 48 #include "G4Step.hh"
 49 #include "G4Track.hh"
 50 #include "G4Log.hh"
 51 #include "G4Exp.hh"
 52 
 53 namespace
 54 {
 55   constexpr G4double c_highland = 13.6*CLHEP::MeV;
 56 }
 57 
 58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 59 
 60 G4DynamicParticleMSC::G4DynamicParticleMSC()
 61   : G4VContinuousDiscreteProcess("dynPartMSC")
 62 {
 63   SetVerboseLevel(1);
 64   SetProcessSubType(fDynamicMultipleScattering);
 65   lManager = G4LossTableManager::Instance();  
 66   lManager->Register(this);
 67 }
 68 
 69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 70 
 71 G4DynamicParticleMSC::~G4DynamicParticleMSC()
 72 {
 73   lManager->DeRegister(this);
 74 }
 75 
 76 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 77 
 78 void G4DynamicParticleMSC::PreStepInitialisation(const G4Track& track)
 79 {
 80   fMaterial = track.GetMaterial();
 81   fZeff = fMaterial->GetIonisation()->GetZeffective();
 82   auto dpart = track.GetDynamicParticle();
 83   fEkinPreStep = dpart->GetKineticEnergy();
 84   fBeta = dpart->GetBeta();
 85   fCharge = dpart->GetCharge()/CLHEP::eplus;
 86   fMass = std::max(dpart->GetMass(), CLHEP::electron_mass_c2);
 87 }
 88 
 89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 90 
 91 G4double G4DynamicParticleMSC::AlongStepGetPhysicalInteractionLength(
 92                                   const G4Track& track, G4double, G4double, G4double&,
 93                                   G4GPILSelection* selection)
 94 {
 95   *selection = CandidateForSelection;
 96   PreStepInitialisation(track);
 97 
 98   // no step limit for the time being
 99   return DBL_MAX;
100 }
101 
102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
103 
104 G4double G4DynamicParticleMSC::PostStepGetPhysicalInteractionLength(
105                                   const G4Track&, G4double,
106                                   G4ForceCondition* condition)
107 {
108   *condition = NotForced;
109   return DBL_MAX;
110 }
111 
112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
113 
114 G4VParticleChange* G4DynamicParticleMSC::AlongStepDoIt(const G4Track& track,
115                           const G4Step& step)
116 {
117   fParticleChange.InitialiseMSC(track, step);
118 
119   // no energy loss
120   if (fCharge == 0.0) { return &fParticleChange; }
121 
122   G4double geomLength = step.GetStepLength();
123   G4double y = geomLength/fMaterial->GetRadlen();
124   G4double theta0 = c_highland*std::abs(fCharge)*std::sqrt(y)*
125     (1.0 + 0.038*G4Log(y*fCharge*fCharge/(fBeta*fBeta)))/fBeta;
126 
127   if (theta0 < 0.001) { return &fParticleChange; }
128   G4double cost = 1.0;
129   G4double r = G4UniformRand();
130   if (theta0 < 1.0) {
131     G4double theta2 = theta0*theta0;
132     cost -= theta2*G4Log(1.0 + r*(G4Exp(2.0/theta2) - 1.0));
133   } else {
134     cost -= 2.0*r;
135   }
136   G4double phi = CLHEP::twopi*G4UniformRand();
137   G4double sint = std::sqrt((1.0 - cost)*(1.0 + cost));
138   fNewDir.set(sint*std::cos(phi), sint*std::sin(phi), cost);
139   fNewDir.rotateUz(step.GetPostStepPoint()->GetMomentumDirection());
140   
141   fParticleChange.ProposeMomentumDirection(fNewDir);
142   fParticleChange.ProposeTrueStepLength(geomLength);
143   return &fParticleChange;
144 }
145 
146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
147 
148 G4double G4DynamicParticleMSC::GetMeanFreePath(const G4Track&, G4double,
149                  G4ForceCondition* condition)
150 {
151   *condition = Forced;
152   return DBL_MAX;
153 }
154 
155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
156 
157 G4double G4DynamicParticleMSC::GetContinuousStepLimit(const G4Track& track,
158                   G4double previousStepSize,
159                                                       G4double currentMinimalStep,
160                   G4double& currentSafety)
161 {
162   G4GPILSelection selection = NotCandidateForSelection;
163   G4double x = AlongStepGetPhysicalInteractionLength(track, previousStepSize,
164                  currentMinimalStep,
165                        currentSafety, &selection);
166   return x;
167 }
168 
169 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
170 
171 void G4DynamicParticleMSC::ProcessDescription(std::ostream& out) const
172 {
173   out << "G4DynamicParticleMSC: no delta rays" << G4endl;
174 }
175 
176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
177