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

#include <G4VPhononProcess.hh>

+ Inheritance diagram for G4VPhononProcess:

Public Member Functions

 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 &)
 
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
 
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 const G4VProcessGetCreatorProcess () const
 
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 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
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

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
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 41 of file G4VPhononProcess.hh.

Constructor & Destructor Documentation

◆ G4VPhononProcess()

G4VPhononProcess::G4VPhononProcess ( const G4String processName)

Definition at line 53 of file G4VPhononProcess.cc.

54 : G4VDiscreteProcess(processName, fPhonon),
56 currentTrack(0) {}
@ fPhonon
static G4PhononTrackMap * GetInstance()
const G4LatticePhysical * theLattice
G4PhononTrackMap * trackKmap

◆ ~G4VPhononProcess()

G4VPhononProcess::~G4VPhononProcess ( )
virtual

Definition at line 58 of file G4VPhononProcess.cc.

58{;}

Member Function Documentation

◆ ChoosePolarization()

G4int G4VPhononProcess::ChoosePolarization ( G4double  Ldos,
G4double  STdos,
G4double  FTdos 
) const
protectedvirtual

Definition at line 103 of file G4VPhononProcess.cc.

104 {
105 G4double norm = Ldos + STdos + FTdos;
106 G4double cProbST = STdos/norm;
107 G4double cProbFT = FTdos/norm + cProbST;
108
109 // NOTE: Order of selection done to match previous random sequences
110 G4double modeMixer = G4UniformRand();
111 if (modeMixer<cProbST) return G4PhononPolarization::TransSlow;
112 if (modeMixer<cProbFT) return G4PhononPolarization::TransFast;
114}
double G4double
Definition: G4Types.hh:83
#define G4UniformRand()
Definition: Randomize.hh:52

Referenced by G4PhononScattering::PostStepDoIt().

◆ CreateSecondary()

G4Track * G4VPhononProcess::CreateSecondary ( G4int  polarization,
const G4ThreeVector K,
G4double  energy 
) const
protectedvirtual

Definition at line 119 of file G4VPhononProcess.cc.

121 {
122 if (verboseLevel>1) {
123 G4cout << GetProcessName() << " CreateSecondary pol " << polarization
124 << " K " << waveVec << " E " << energy << G4endl;
125 }
126
127 G4ThreeVector vgroup = theLattice->MapKtoVDir(polarization, waveVec);
128 if (verboseLevel>1) G4cout << " MapKtoVDir returned " << vgroup << G4endl;
129
130 vgroup = theLattice->RotateToGlobal(vgroup);
131 if (verboseLevel>1) G4cout << " RotateToGlobal returned " << vgroup << G4endl;
132
133 if (verboseLevel && std::fabs(vgroup.mag()-1.) > 0.01) {
134 G4cout << "WARNING: " << GetProcessName() << " vgroup not a unit vector: "
135 << vgroup << G4endl;
136 }
137
138 G4ParticleDefinition* thePhonon = G4PhononPolarization::Get(polarization);
139
140 // Secondaries are created at the current track coordinates
141 G4Track* sec = new G4Track(new G4DynamicParticle(thePhonon, vgroup, energy),
142 currentTrack->GetGlobalTime(),
143 currentTrack->GetPosition());
144
145 // Store wavevector in lookup table for future tracking
146 trackKmap->SetK(sec, theLattice->RotateToGlobal(waveVec));
147
148 if (verboseLevel>1) {
149 G4cout << GetProcessName() << " secondary K rotated to "
150 << trackKmap->GetK(sec) << G4endl;
151 }
152
153 sec->SetVelocity(theLattice->MapKtoV(polarization, waveVec));
154 sec->UseGivenVelocity(true);
155
156 return sec;
157}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double mag() const
G4ThreeVector RotateToGlobal(const G4ThreeVector &dir) const
G4double MapKtoV(G4int, G4ThreeVector) const
G4ThreeVector MapKtoVDir(G4int, G4ThreeVector) const
const G4ThreeVector & GetK(const G4Track *track) const
void SetK(const G4Track *track, const G4ThreeVector &K)
void SetVelocity(G4double val)
G4bool UseGivenVelocity() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4int verboseLevel
Definition: G4VProcess.hh:360
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
G4double energy(const ThreeVector &p, const G4double m)
G4int Get(const G4ParticleDefinition *aPD)

Referenced by G4PhononScattering::PostStepDoIt().

◆ EndTracking()

void G4VPhononProcess::EndTracking ( )
virtual

Reimplemented from G4VProcess.

Definition at line 86 of file G4VPhononProcess.cc.

86 {
87 G4VProcess::EndTracking(); // Apply base class actions
88 trackKmap->RemoveTrack(currentTrack);
89 currentTrack = 0;
90 theLattice = 0;
91}
void RemoveTrack(const G4Track *track)
virtual void EndTracking()
Definition: G4VProcess.cc:102

◆ GetPolarization() [1/2]

G4int G4VPhononProcess::GetPolarization ( const G4Track track) const
protectedvirtual

Definition at line 96 of file G4VPhononProcess.cc.

96 {
98}
const G4ParticleDefinition * GetParticleDefinition() const

Referenced by GetPolarization(), and G4PhononReflection::PostStepDoIt().

◆ GetPolarization() [2/2]

virtual G4int G4VPhononProcess::GetPolarization ( const G4Track track) const
inlineprotectedvirtual

Definition at line 56 of file G4VPhononProcess.hh.

56 {
57 return GetPolarization(*track);
58 }
virtual G4int GetPolarization(const G4Track &track) const

◆ IsApplicable()

G4bool G4VPhononProcess::IsApplicable ( const G4ParticleDefinition aPD)
virtual

Reimplemented from G4VProcess.

Reimplemented in G4PhononDownconversion.

Definition at line 63 of file G4VPhononProcess.cc.

63 {
64 return (&aPD==G4PhononLong::Definition() ||
67}
static G4PhononLong * Definition()
Definition: G4PhononLong.cc:39
static G4PhononTransFast * Definition()
static G4PhononTransSlow * Definition()

◆ StartTracking()

void G4VPhononProcess::StartTracking ( G4Track track)
virtual

Reimplemented from G4VProcess.

Definition at line 72 of file G4VPhononProcess.cc.

72 {
73 G4VProcess::StartTracking(track); // Apply base class actions
74
75 // FIXME: THE WAVEVECTOR SHOULD BE COMPUTED BY INVERTING THE K/V MAP
76 if (!trackKmap->Find(track))
77 trackKmap->SetK(track, track->GetMomentumDirection());
78
79 currentTrack = track; // Save for use by EndTracking
80
81 // Fetch lattice for current track once, use in subsequent steps
83 theLattice = LM->GetLattice(track->GetVolume());
84}
G4LatticeLogical * GetLattice(G4Material *) const
static G4LatticeManager * GetLatticeManager()
G4bool Find(const G4Track *track) const
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetMomentumDirection() const
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:87

Member Data Documentation

◆ theLattice

◆ trackKmap

G4PhononTrackMap* G4VPhononProcess::trackKmap
protected

Definition at line 70 of file G4VPhononProcess.hh.

Referenced by CreateSecondary(), EndTracking(), and StartTracking().


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