Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PhononScattering Class Reference

#include <G4PhononScattering.hh>

+ Inheritance diagram for G4PhononScattering:

Public Member Functions

 G4PhononScattering (const G4String &processName="phononScattering")
 
virtual ~G4PhononScattering ()
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VPhononProcess
 G4VPhononProcess (const G4String &processName)
 
virtual ~G4VPhononProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &aPD)
 
virtual void StartTracking (G4Track *track)
 
virtual void EndTracking ()
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

virtual G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *)
 
- Protected Member Functions inherited from G4VPhononProcess
virtual G4int GetPolarization (const G4Track &track) const
 
virtual G4int GetPolarization (const G4Track *track) const
 
virtual G4int ChoosePolarization (G4double Ldos, G4double STdos, G4double FTdos) const
 
virtual G4TrackCreateSecondary (G4int polarization, const G4ThreeVector &K, G4double energy) const
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VPhononProcess
G4PhononTrackMaptrackKmap
 
const G4LatticePhysicaltheLattice
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 35 of file G4PhononScattering.hh.

Constructor & Destructor Documentation

◆ G4PhononScattering()

G4PhononScattering::G4PhononScattering ( const G4String processName = "phononScattering")

Definition at line 47 of file G4PhononScattering.cc.

◆ ~G4PhononScattering()

G4PhononScattering::~G4PhononScattering ( )
virtual

Definition at line 50 of file G4PhononScattering.cc.

50{;}

Member Function Documentation

◆ GetMeanFreePath()

G4double G4PhononScattering::GetMeanFreePath ( const G4Track aTrack,
G4double  ,
G4ForceCondition condition 
)
protectedvirtual

Implements G4VDiscreteProcess.

Definition at line 54 of file G4PhononScattering.cc.

56 {
57 //Dynamical constants retrieved from PhysicalLattice
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}
double B(double temperature)
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetScatteringConstant() const
G4double GetVelocity() const
G4double GetKineticEnergy() const
const G4LatticePhysical * theLattice
G4int verboseLevel
Definition: G4VProcess.hh:356

◆ PostStepDoIt()

G4VParticleChange * G4PhononScattering::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 74 of file G4PhononScattering.cc.

75 {
76 G4StepPoint* postStepPoint = aStep.GetPostStepPoint();
77 if (postStepPoint->GetStepStatus()==fGeomBoundary) {
78 return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
79 }
80
81 //Initialize particle change
83
84 //randomly generate a new direction and polarization state
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());
96
97 // Scattered phonon replaces current track
100
101 return &aParticleChange;
102}
G4ThreeVector G4RandomDirection()
@ fGeomBoundary
Definition: G4StepStatus.hh:43
@ fStopAndKill
int G4int
Definition: G4Types.hh:85
G4double GetFTDOS() const
G4double GetLDOS() const
G4double GetSTDOS() const
void AddSecondary(G4Track *aSecondary)
void ProposeEnergy(G4double finalEnergy)
virtual void Initialize(const G4Track &)
G4StepStatus GetStepStatus() const
G4StepPoint * GetPostStepPoint() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void SetNumberOfSecondaries(G4int totSecondaries)
virtual G4int ChoosePolarization(G4double Ldos, G4double STdos, G4double FTdos) const
virtual G4Track * CreateSecondary(G4int polarization, const G4ThreeVector &K, G4double energy) const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327

The documentation for this class was generated from the following files: