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

#include <G4PhononDownconversion.hh>

+ Inheritance diagram for G4PhononDownconversion:

Public Member Functions

 G4PhononDownconversion (const G4String &processName="phononDownconversion")
 
virtual ~G4PhononDownconversion ()
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
- 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 &)
 
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 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 G4PhononDownconversion.hh.

Constructor & Destructor Documentation

◆ G4PhononDownconversion()

G4PhononDownconversion::G4PhononDownconversion ( const G4String processName = "phononDownconversion")

Definition at line 50 of file G4PhononDownconversion.cc.

51 : G4VPhononProcess(aName), fBeta(0.), fGamma(0.), fLambda(0.), fMu(0.) {;}

◆ ~G4PhononDownconversion()

G4PhononDownconversion::~G4PhononDownconversion ( )
virtual

Definition at line 53 of file G4PhononDownconversion.cc.

53{;}

Member Function Documentation

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 57 of file G4PhononDownconversion.cc.

59 {
60 //Determines mean free path for longitudinal phonons to split
62 G4double Eoverh = aTrack.GetKineticEnergy()/h_Planck;
63
64 //Calculate mean free path for anh. decay
65 G4double mfp = aTrack.GetVelocity()/(Eoverh*Eoverh*Eoverh*Eoverh*Eoverh*A);
66
67 if (verboseLevel > 1)
68 G4cout << "G4PhononDownconversion::GetMeanFreePath = " << mfp << G4endl;
69
71 return mfp;
72}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
double G4double
Definition: G4Types.hh:83
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetAnhDecConstant() const
G4double GetVelocity() const
G4double GetKineticEnergy() const
const G4LatticePhysical * theLattice
G4int verboseLevel
Definition: G4VProcess.hh:360

◆ IsApplicable()

G4bool G4PhononDownconversion::IsApplicable ( const G4ParticleDefinition aPD)
virtual

Reimplemented from G4VPhononProcess.

Definition at line 101 of file G4PhononDownconversion.cc.

101 {
102 //Only L-phonons decay
103 return (&aPD==G4PhononLong::PhononDefinition());
104}
static G4PhononLong * PhononDefinition()
Definition: G4PhononLong.cc:73

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 77 of file G4PhononDownconversion.cc.

78 {
80
81 //Obtain dynamical constants from this volume's lattice
82 fBeta=theLattice->GetBeta();
83 fGamma=theLattice->GetGamma();
84 fLambda=theLattice->GetLambda();
85 fMu=theLattice->GetMu();
86
87 //Destroy the parent phonon and create the daughter phonons.
88 //74% chance that daughter phonons are both transverse
89 //26% Transverse and Longitudinal
90 if (G4UniformRand()>0.740) MakeLTSecondaries(aTrack);
91 else MakeTTSecondaries(aTrack);
92
95
96 return &aParticleChange;
97}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
G4double GetGamma() const
G4double GetMu() const
G4double GetBeta() const
G4double GetLambda() const
void Initialize(const G4Track &) override
void ProposeEnergy(G4double finalEnergy)
void ProposeTrackStatus(G4TrackStatus status)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:331

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