Geant4  10.01.p02
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 //
28 //
29 // $Id: G4PhononScattering.cc 76693 2013-11-14 08:47:37Z gcosmo $
30 //
31 // 20131111 Add verbose output for MFP calculation
32 
33 #include "G4PhononScattering.hh"
34 #include "G4LatticePhysical.hh"
35 #include "G4PhononPolarization.hh"
36 #include "G4PhononLong.hh"
37 #include "G4PhononTrackMap.hh"
38 #include "G4PhononTransFast.hh"
39 #include "G4PhononTransSlow.hh"
40 #include "G4PhysicalConstants.hh"
41 #include "G4RandomDirection.hh"
42 #include "G4Step.hh"
43 #include "G4SystemOfUnits.hh"
44 #include "G4VParticleChange.hh"
45 #include "Randomize.hh"
46 
47 
49  : G4VPhononProcess(aName) {;}
50 
52 
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54 
56  G4double /*previousStepSize*/,
58  //Dynamical constants retrieved from PhysicalLattice
60  G4double Eoverh = aTrack.GetKineticEnergy()/h_Planck;
61 
62  //Calculate mean free path
63  G4double mfp = aTrack.GetVelocity()/(Eoverh*Eoverh*Eoverh*Eoverh*B);
64 
65  if (verboseLevel > 1)
66  G4cout << "G4PhononScattering::GetMeanFreePath = " << mfp << G4endl;
67 
68  *condition = NotForced;
69 
70  return mfp;
71 }
72 
73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74 
76  const G4Step& aStep) {
77  G4StepPoint* postStepPoint = aStep.GetPostStepPoint();
78  if (postStepPoint->GetStepStatus()==fGeomBoundary) {
79  return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
80  }
81 
82  //Initialize particle change
84 
85  //randomly generate a new direction and polarization state
87  G4int polarization = ChoosePolarization(theLattice->GetLDOS(),
89  theLattice->GetFTDOS());
90 
91  // Generate the new track after scattering
92  // FIXME: If polarization state is the same, just step the track!
93  G4Track* sec =
94  CreateSecondary(polarization, newDir, aTrack.GetKineticEnergy());
97 
98  // Scattered phonon replaces current track
101 
102  return &aParticleChange;
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double condition(const G4ErrorSymMatrix &m)
const G4LatticePhysical * theLattice
G4int verboseLevel
Definition: G4VProcess.hh:368
CLHEP::Hep3Vector G4ThreeVector
G4double GetVelocity() const
G4double GetSTDOS() const
G4StepStatus GetStepStatus() const
G4ThreeVector G4RandomDirection()
Definition of the G4PhononPolarization enum.
Definition of the G4PhononScattering class.
int G4int
Definition: G4Types.hh:78
G4double GetFTDOS() const
G4double GetKineticEnergy() const
Definition of the G4PhononTrackMap base class.
G4GLOB_DLL std::ostream G4cout
Definition: G4Step.hh:76
G4double GetScatteringConstant() const
G4PhononScattering(const G4String &processName="phononScattering")
virtual void Initialize(const G4Track &)
virtual G4int ChoosePolarization(G4double Ldos, G4double STdos, G4double FTdos) const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void SetNumberOfSecondaries(G4int totSecondaries)
G4StepPoint * GetPostStepPoint() const
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
Definition of the G4LatticePhysical class.
void AddSecondary(G4Track *aSecondary)
virtual G4Track * CreateSecondary(G4int polarization, const G4ThreeVector &K, G4double energy) const
virtual G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *)
#define G4endl
Definition: G4ios.hh:61
G4double GetLDOS() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)