Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PhononReflection Class Reference

#include <G4PhononReflection.hh>

Inheritance diagram for G4PhononReflection:
Collaboration 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 &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
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 &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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 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
 
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 G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
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
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 37 of file G4PhononReflection.hh.

Constructor & Destructor Documentation

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

Definition at line 53 of file G4PhononReflection.cc.

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

Definition at line 57 of file G4PhononReflection.cc.

57 {;}

Member Function Documentation

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

Implements G4VDiscreteProcess.

Definition at line 66 of file G4PhononReflection.cc.

67  {
68  *condition = Forced;
69  return DBL_MAX;
70 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
G4VParticleChange * G4PhononReflection::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 80 of file G4PhononReflection.cc.

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

Here is the call graph for this function:


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