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

#include <G4VLowEnergyDiscretePhotonProcess.hh>

+ Inheritance diagram for G4VLowEnergyDiscretePhotonProcess:

Public Member Functions

 G4VLowEnergyDiscretePhotonProcess (const G4String &processName, const G4String &aCrossSectionFileName, const G4String &aScatterFileName, G4VDataSetAlgorithm *aScatterInterpolation, G4double aLowEnergyLimit, G4double aHighEnergyLimit)
 
virtual ~G4VLowEnergyDiscretePhotonProcess (void)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &particleDefinition)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &photon)
 
G4double GetLowEnergyLimit (void) const
 
G4double GetHighEnergyLimit (void) const
 
- Public Member Functions inherited from G4VLowEnergyTestableDiscreteProcess
 G4VLowEnergyTestableDiscreteProcess (const G4String &name)
 
virtual ~G4VLowEnergyTestableDiscreteProcess ()
 
G4double DumpMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
- 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 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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Member Functions

virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
const G4VCrossSectionHandlerGetCrossSectionHandler (void) const
 
const G4VEMDataSetGetMeanFreePathTable (void) const
 
const G4VEMDataSetGetScatterFunctionData (void) const
 
G4ThreeVector GetPhotonPolarization (const G4DynamicParticle &photon)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- 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 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 62 of file G4VLowEnergyDiscretePhotonProcess.hh.

Constructor & Destructor Documentation

◆ G4VLowEnergyDiscretePhotonProcess()

G4VLowEnergyDiscretePhotonProcess::G4VLowEnergyDiscretePhotonProcess ( const G4String processName,
const G4String aCrossSectionFileName,
const G4String aScatterFileName,
G4VDataSetAlgorithm aScatterInterpolation,
G4double  aLowEnergyLimit,
G4double  aHighEnergyLimit 
)

Definition at line 57 of file G4VLowEnergyDiscretePhotonProcess.cc.

63 :
65 lowEnergyLimit(aLowEnergyLimit),
66 highEnergyLimit(aHighEnergyLimit),
67 crossSectionFileName(aCrossSectionFileName),
68 meanFreePathTable(0)
69{
70 crossSectionHandler = new G4CrossSectionHandler();
71 scatterFunctionData = new G4CompositeEMDataSet(aScatterInterpolation, 1., 1.);
72 scatterFunctionData->LoadData(aScatterFileName);
73
74 if (verboseLevel > 0)
75 {
76 G4cout << GetProcessName() << " is created " << G4endl
77 << "Energy range: "
78 << lowEnergyLimit / keV << " keV - "
79 << highEnergyLimit / GeV << " GeV"
80 << G4endl;
81 }
82}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
virtual G4bool LoadData(const G4String &fileName)=0
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ ~G4VLowEnergyDiscretePhotonProcess()

G4VLowEnergyDiscretePhotonProcess::~G4VLowEnergyDiscretePhotonProcess ( void  )
virtual

Definition at line 86 of file G4VLowEnergyDiscretePhotonProcess.cc.

87{
88 if (meanFreePathTable)
89 delete meanFreePathTable;
90
91 delete crossSectionHandler;
92 delete scatterFunctionData;
93}

Member Function Documentation

◆ BuildPhysicsTable()

void G4VLowEnergyDiscretePhotonProcess::BuildPhysicsTable ( const G4ParticleDefinition photon)
virtual

Reimplemented from G4VProcess.

Definition at line 106 of file G4VLowEnergyDiscretePhotonProcess.cc.

107{
108 crossSectionHandler->Clear();
109 crossSectionHandler->LoadData(crossSectionFileName);
110
111 if (meanFreePathTable)
112 delete meanFreePathTable;
113 meanFreePathTable=crossSectionHandler->BuildMeanFreePathForMaterials();
114}
G4VEMDataSet * BuildMeanFreePathForMaterials(const G4DataVector *energyCuts=0)
void LoadData(const G4String &dataFile)

◆ GetCrossSectionHandler()

const G4VCrossSectionHandler * G4VLowEnergyDiscretePhotonProcess::GetCrossSectionHandler ( void  ) const
inlineprotected

Definition at line 129 of file G4VLowEnergyDiscretePhotonProcess.hh.

129{return crossSectionHandler;}

◆ GetHighEnergyLimit()

G4double G4VLowEnergyDiscretePhotonProcess::GetHighEnergyLimit ( void  ) const
inline

Definition at line 110 of file G4VLowEnergyDiscretePhotonProcess.hh.

110{return highEnergyLimit;}

◆ GetLowEnergyLimit()

G4double G4VLowEnergyDiscretePhotonProcess::GetLowEnergyLimit ( void  ) const
inline

Definition at line 106 of file G4VLowEnergyDiscretePhotonProcess.hh.

106{return lowEnergyLimit;}

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 120 of file G4VLowEnergyDiscretePhotonProcess.cc.

121{
122 G4double photonEnergy;
123 photonEnergy = aTrack.GetDynamicParticle()->GetKineticEnergy();
124
125 if (photonEnergy < lowEnergyLimit)
126 return DBL_MAX;
127
128 if (photonEnergy > highEnergyLimit)
129 photonEnergy=highEnergyLimit;
130
131 size_t materialIndex;
132 materialIndex = aTrack.GetMaterialCutsCouple()->GetIndex();
133
134 return meanFreePathTable->FindValue(photonEnergy, materialIndex);
135}
double G4double
Definition: G4Types.hh:64
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
#define DBL_MAX
Definition: templates.hh:83

◆ GetMeanFreePathTable()

const G4VEMDataSet * G4VLowEnergyDiscretePhotonProcess::GetMeanFreePathTable ( void  ) const
inlineprotected

Definition at line 133 of file G4VLowEnergyDiscretePhotonProcess.hh.

133{return meanFreePathTable;}

◆ GetPhotonPolarization()

G4ThreeVector G4VLowEnergyDiscretePhotonProcess::GetPhotonPolarization ( const G4DynamicParticle photon)
protected

Definition at line 141 of file G4VLowEnergyDiscretePhotonProcess.cc.

142{
143 G4ThreeVector photonMomentumDirection;
144 G4ThreeVector photonPolarization;
145
146 photonPolarization = photon.GetPolarization();
147 photonMomentumDirection = photon.GetMomentumDirection();
148
149 if ((!photonPolarization.isOrthogonal(photonMomentumDirection, 1e-6)) || photonPolarization.mag()==0.)
150 {
151 // if |photonPolarization|==0. or |photonPolarization * photonDirection0| > 1e-6 * |photonPolarization ^ photonDirection0|
152 // then polarization is choosen randomly.
153
154 G4ThreeVector e1(photonMomentumDirection.orthogonal().unit());
155 G4ThreeVector e2(photonMomentumDirection.cross(e1).unit());
156
157 G4double angle(G4UniformRand() * twopi);
158
159 e1*=std::cos(angle);
160 e2*=std::sin(angle);
161
162 photonPolarization=e1+e2;
163 }
164 else if (photonPolarization.howOrthogonal(photonMomentumDirection) != 0.)
165 {
166 // if |photonPolarization * photonDirection0| != 0.
167 // then polarization is made orthonormal;
168
169 photonPolarization=photonPolarization.perpPart(photonMomentumDirection);
170 }
171
172 return photonPolarization.unit();
173}
@ photon
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
bool isOrthogonal(const Hep3Vector &v, double epsilon=tolerance) const
Definition: SpaceVector.cc:237
Hep3Vector perpPart() const
double mag() const
double howOrthogonal(const Hep3Vector &v) const
Definition: SpaceVector.cc:219

◆ GetScatterFunctionData()

const G4VEMDataSet * G4VLowEnergyDiscretePhotonProcess::GetScatterFunctionData ( void  ) const
inlineprotected

Definition at line 137 of file G4VLowEnergyDiscretePhotonProcess.hh.

137{return scatterFunctionData;}

◆ IsApplicable()

G4bool G4VLowEnergyDiscretePhotonProcess::IsApplicable ( const G4ParticleDefinition particleDefinition)
virtual

Reimplemented from G4VProcess.

Definition at line 99 of file G4VLowEnergyDiscretePhotonProcess.cc.

100{
101 return (&particleDefinition)==G4Gamma::Gamma();
102}
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86

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