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

#include <G4PhononReflection.hh>

+ Inheritance diagram for G4PhononReflection:

Public Member Functions

 G4PhononReflection (const G4String &processName="phononReflection")
 
virtual ~G4PhononReflection ()
 
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 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 &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- 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
 
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 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 const G4VProcessGetCreatorProcess () const
 
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
 
- Protected Member Functions inherited from G4VDiscreteProcess
- 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 36 of file G4PhononReflection.hh.

Constructor & Destructor Documentation

◆ G4PhononReflection()

G4PhononReflection::G4PhononReflection ( const G4String & processName = "phononReflection")

Definition at line 52 of file G4PhononReflection.cc.

53 : G4VPhononProcess(aName),
54 kCarTolerance(G4GeometryTolerance::GetInstance()->GetSurfaceTolerance()) {;}
static G4GeometryTolerance * GetInstance()
G4VPhononProcess(const G4String &processName)

◆ ~G4PhononReflection()

G4PhononReflection::~G4PhononReflection ( )
virtual

Definition at line 56 of file G4PhononReflection.cc.

56{;}

Member Function Documentation

◆ GetMeanFreePath()

G4double G4PhononReflection::GetMeanFreePath ( const G4Track & ,
G4double ,
G4ForceCondition * condition )
protectedvirtual

Implements G4VDiscreteProcess.

Definition at line 65 of file G4PhononReflection.cc.

66 {
68 return DBL_MAX;
69}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced
#define DBL_MAX
Definition templates.hh:62

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 79 of file G4PhononReflection.cc.

80 {
82
83 //Check if current step is limited by a volume boundary
84 G4StepPoint* postStepPoint = aStep.GetPostStepPoint();
85 if (postStepPoint->GetStepStatus()!=fGeomBoundary) {
86 //make sure that correct phonon velocity is used after the step
87 int pol = GetPolarization(aTrack);
88 if (pol < 0 || pol > 2) {
89 G4Exception("G4PhononReflection::PostStepDoIt","Phonon001",
90 EventMustBeAborted, "Track is not a phonon");
91 return &aParticleChange; // NOTE: Will never get here
92 }
93
94 // FIXME: This should be using wave-vector, shouldn't it?
96
97 //Since step was not a volume boundary, just set correct phonon velocity and return
99 return &aParticleChange;
100 }
101
102 // do nothing but return is the step is too short
103 // This is to allow actual reflection where after
104 // the first boundary crossing a second, infinitesimal
105 // step occurs crossing back into the original volume
106 if (aTrack.GetStepLength()<=kCarTolerance/2) {
107 return &aParticleChange;
108 }
109
110 G4double eKin = aTrack.GetKineticEnergy();
113
114 return &aParticleChange;
115}
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
@ fGeomBoundary
@ fStopAndKill
double G4double
Definition G4Types.hh:83
G4double MapKtoV(G4int, G4ThreeVector) const
void Initialize(const G4Track &) override
void ProposeVelocity(G4double finalVelocity)
G4StepStatus GetStepStatus() const
G4StepPoint * GetPostStepPoint() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetStepLength() const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
virtual G4int GetPolarization(const G4Track &track) const
const G4LatticePhysical * theLattice
G4ParticleChange aParticleChange

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