Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhononScattering.cc
Go to the documentation of this file.
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/// \file processes/phonon/src/G4PhononScattering.cc
27/// \brief Implementation of the G4PhononScattering class
28//
29//
30// 20131111 Add verbose output for MFP calculation
31
32#include "G4PhononScattering.hh"
33#include "G4LatticePhysical.hh"
35#include "G4PhononLong.hh"
36#include "G4PhononTrackMap.hh"
37#include "G4PhononTransFast.hh"
38#include "G4PhononTransSlow.hh"
40#include "G4RandomDirection.hh"
41#include "G4Step.hh"
42#include "G4SystemOfUnits.hh"
43#include "G4VParticleChange.hh"
44#include "Randomize.hh"
45
46
49
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53
55 G4double /*previousStepSize*/,
57 //Dynamical constants retrieved from PhysicalLattice
58 G4double B = theLattice->GetScatteringConstant();
59 G4double Eoverh = aTrack.GetKineticEnergy()/h_Planck;
60
61 //Calculate mean free path
62 G4double mfp = aTrack.GetVelocity()/(Eoverh*Eoverh*Eoverh*Eoverh*B);
63
64 if (verboseLevel > 1)
65 G4cout << "G4PhononScattering::GetMeanFreePath = " << mfp << G4endl;
66
68
69 return mfp;
70}
71
72//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
73
75 const G4Step& aStep) {
76 G4StepPoint* postStepPoint = aStep.GetPostStepPoint();
77 if (postStepPoint->GetStepStatus()==fGeomBoundary) {
78 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
79 }
80
81 //Initialize particle change
82 aParticleChange.Initialize(aTrack);
83
84 //randomly generate a new direction and polarization state
86 G4int polarization = ChoosePolarization(theLattice->GetLDOS(),
87 theLattice->GetSTDOS(),
88 theLattice->GetFTDOS());
89
90 // Generate the new track after scattering
91 // FIXME: If polarization state is the same, just step the track!
92 G4Track* sec =
93 CreateSecondary(polarization, newDir, aTrack.GetKineticEnergy());
94 aParticleChange.SetNumberOfSecondaries(1);
95 aParticleChange.AddSecondary(sec);
96
97 // Scattered phonon replaces current track
98 aParticleChange.ProposeEnergy(0.);
99 aParticleChange.ProposeTrackStatus(fStopAndKill);
100
101 return &aParticleChange;
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
Definition of the G4LatticePhysical class.
G4ThreeVector G4RandomDirection()
@ fGeomBoundary
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4PhononScattering(const G4String &processName="phononScattering")
G4StepStatus GetStepStatus() const
G4StepPoint * GetPostStepPoint() const
G4double GetVelocity() const
G4double GetKineticEnergy() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4VPhononProcess(const G4String &processName)
virtual G4int ChoosePolarization(G4double Ldos, G4double STdos, G4double FTdos) const
virtual G4Track * CreateSecondary(G4int polarization, const G4ThreeVector &K, G4double energy) const
const G4LatticePhysical * theLattice
G4ParticleChange aParticleChange
G4int verboseLevel