Geant4  10.02.p03
G4LowEnergyRayleigh Class Reference

#include <G4LowEnergyRayleigh.hh>

Inheritance diagram for G4LowEnergyRayleigh:
Collaboration diagram for G4LowEnergyRayleigh:

Public Member Functions

 G4LowEnergyRayleigh (const G4String &processName="LowEnRayleigh")
 
 ~G4LowEnergyRayleigh ()
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &photon)
 
G4VParticleChange * PostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
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 G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChange * AlongStepDoIt (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 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
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Private Member Functions

G4LowEnergyRayleighoperator= (const G4LowEnergyRayleigh &right)
 
 G4LowEnergyRayleigh (const G4LowEnergyRayleigh &)
 

Private Attributes

G4double lowEnergyLimit
 
G4double highEnergyLimit
 
G4RDVEMDataSetmeanFreePathTable
 
G4RDVEMDataSetformFactorData
 
G4RDVCrossSectionHandlercrossSectionHandler
 
const G4double intrinsicLowEnergyLimit
 
const G4double intrinsicHighEnergyLimit
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
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 60 of file G4LowEnergyRayleigh.hh.

Constructor & Destructor Documentation

◆ G4LowEnergyRayleigh() [1/2]

G4LowEnergyRayleigh::G4LowEnergyRayleigh ( const G4String processName = "LowEnRayleigh")

Definition at line 73 of file G4LowEnergyRayleigh.cc.

74  : G4VDiscreteProcess(processName),
75  lowEnergyLimit(250*eV),
76  highEnergyLimit(100*GeV),
79 {
82  {
83  G4Exception("G4LowEnergyRayleigh::G4LowEnergyRayleigh()",
84  "OutOfRange", FatalException,
85  "Energy limit outside intrinsic process validity range!");
86  }
87 
89 
90  G4RDVDataSetAlgorithm* ffInterpolation = new G4RDLogLogInterpolation;
91  G4String formFactorFile = "rayl/re-ff-";
92  formFactorData = new G4RDCompositeEMDataSet(ffInterpolation,1.,1.);
93  formFactorData->LoadData(formFactorFile);
94 
96 
97  if (verboseLevel > 0)
98  {
99  G4cout << GetProcessName() << " is created " << G4endl
100  << "Energy range: "
101  << lowEnergyLimit / keV << " keV - "
102  << highEnergyLimit / GeV << " GeV"
103  << G4endl;
104  }
105 }
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4double intrinsicLowEnergyLimit
G4RDVEMDataSet * formFactorData
G4RDVEMDataSet * meanFreePathTable
const G4double intrinsicHighEnergyLimit
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
static const double GeV
Definition: G4SIunits.hh:214
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double eV
Definition: G4SIunits.hh:212
G4RDVCrossSectionHandler * crossSectionHandler
#define G4endl
Definition: G4ios.hh:61
static const double keV
Definition: G4SIunits.hh:213
virtual G4bool LoadData(const G4String &fileName)=0
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ~G4LowEnergyRayleigh()

G4LowEnergyRayleigh::~G4LowEnergyRayleigh ( )

Definition at line 107 of file G4LowEnergyRayleigh.cc.

108 {
109  delete meanFreePathTable;
110  delete crossSectionHandler;
111  delete formFactorData;
112 }
G4RDVEMDataSet * formFactorData
G4RDVEMDataSet * meanFreePathTable
G4RDVCrossSectionHandler * crossSectionHandler

◆ G4LowEnergyRayleigh() [2/2]

G4LowEnergyRayleigh::G4LowEnergyRayleigh ( const G4LowEnergyRayleigh )
private

Member Function Documentation

◆ BuildPhysicsTable()

void G4LowEnergyRayleigh::BuildPhysicsTable ( const G4ParticleDefinition photon)
virtual

Reimplemented from G4VProcess.

Definition at line 114 of file G4LowEnergyRayleigh.cc.

115 {
116 
118  G4String crossSectionFile = "rayl/re-cs-";
119  crossSectionHandler->LoadData(crossSectionFile);
120 
121  delete meanFreePathTable;
123 }
G4RDVEMDataSet * meanFreePathTable
G4RDVEMDataSet * BuildMeanFreePathForMaterials(const G4DataVector *energyCuts=0)
G4RDVCrossSectionHandler * crossSectionHandler
void LoadData(const G4String &dataFile)
Here is the call graph for this function:

◆ DumpMeanFreePath()

G4double G4LowEnergyRayleigh::DumpMeanFreePath ( const G4Track &  aTrack,
G4double  previousStepSize,
G4ForceCondition *  condition 
)
inline

Definition at line 75 of file G4LowEnergyRayleigh.hh.

78  { return GetMeanFreePath(aTrack, previousStepSize, condition); }
G4double condition(const G4ErrorSymMatrix &m)
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
Here is the call graph for this function:

◆ GetMeanFreePath()

G4double G4LowEnergyRayleigh::GetMeanFreePath ( const G4Track &  aTrack,
G4double  previousStepSize,
G4ForceCondition *  condition 
)
protectedvirtual

Implements G4VDiscreteProcess.

Definition at line 202 of file G4LowEnergyRayleigh.cc.

205 {
206  const G4DynamicParticle* photon = track.GetDynamicParticle();
207  G4double energy = photon->GetKineticEnergy();
208  const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
209  size_t materialIndex = couple->GetIndex();
210 
211  G4double meanFreePath;
212  if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
213  else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
214  else meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
215  return meanFreePath;
216 }
G4RDVEMDataSet * meanFreePathTable
G4double GetKineticEnergy() const
double energy
Definition: plottest35.C:25
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsApplicable()

G4bool G4LowEnergyRayleigh::IsApplicable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4VProcess.

Definition at line 197 of file G4LowEnergyRayleigh.cc.

198 {
199  return ( &particle == G4Gamma::Gamma() );
200 }
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
Here is the call graph for this function:

◆ operator=()

G4LowEnergyRayleigh& G4LowEnergyRayleigh::operator= ( const G4LowEnergyRayleigh right)
private
Here is the caller graph for this function:

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 125 of file G4LowEnergyRayleigh.cc.

127 {
128 
129  aParticleChange.Initialize(aTrack);
130 
131  const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle();
132  G4double photonEnergy0 = incidentPhoton->GetKineticEnergy();
133 
134  if (photonEnergy0 <= lowEnergyLimit)
135  {
136  aParticleChange.ProposeTrackStatus(fStopAndKill);
137  aParticleChange.ProposeEnergy(0.);
138  aParticleChange.ProposeLocalEnergyDeposit(photonEnergy0);
139  return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
140  }
141 
142  // G4double e0m = photonEnergy0 / electron_mass_c2 ;
143  G4ParticleMomentum photonDirection0 = incidentPhoton->GetMomentumDirection();
144 
145  // Select randomly one element in the current material
146  const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
147  G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0);
148 
149  // Sample the angle of the scattered photon
150 
151  G4double wlPhoton = h_Planck*c_light/photonEnergy0;
152 
153  G4double gReject,x,dataFormFactor;
154  G4double randomFormFactor;
155  G4double cosTheta;
156  G4double sinTheta;
157  G4double fcostheta;
158 
159  do
160  {
161  do
162  {
163  cosTheta = 2. * G4UniformRand() - 1.;
164  fcostheta = ( 1. + cosTheta*cosTheta)/2.;
165  } while (fcostheta < G4UniformRand());
166 
167  G4double sinThetaHalf = std::sqrt((1. - cosTheta) / 2.);
168  x = sinThetaHalf / (wlPhoton/cm);
169  if (x > 1.e+005)
170  dataFormFactor = formFactorData->FindValue(x,Z-1);
171  else
172  dataFormFactor = formFactorData->FindValue(0.,Z-1);
173  randomFormFactor = G4UniformRand() * Z * Z;
174  sinTheta = std::sqrt(1. - cosTheta*cosTheta);
175  gReject = dataFormFactor * dataFormFactor;
176 
177  } while( gReject < randomFormFactor);
178 
179  // Scattered photon angles. ( Z - axis along the parent photon)
180  G4double phi = twopi * G4UniformRand() ;
181  G4double dirX = sinTheta*std::cos(phi);
182  G4double dirY = sinTheta*std::sin(phi);
183  G4double dirZ = cosTheta;
184 
185  // Update G4VParticleChange for the scattered photon
186  G4ThreeVector photonDirection1(dirX, dirY, dirZ);
187 
188  photonDirection1.rotateUz(photonDirection0);
189  aParticleChange.ProposeEnergy(photonEnergy0);
190  aParticleChange.ProposeMomentumDirection(photonDirection1);
191 
192  aParticleChange.SetNumberOfSecondaries(0);
193 
194  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
195 }
static const double cm
Definition: G4SIunits.hh:118
G4int SelectRandomAtom(const G4MaterialCutsCouple *couple, G4double e) const
G4RDVEMDataSet * formFactorData
float h_Planck
Definition: hepunit.py:263
int G4int
Definition: G4Types.hh:78
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:97
Float_t Z
static const double twopi
Definition: G4SIunits.hh:75
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
const G4ThreeVector & GetMomentumDirection() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4RDVCrossSectionHandler * crossSectionHandler
double G4double
Definition: G4Types.hh:76
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
float c_light
Definition: hepunit.py:257
Here is the call graph for this function:

Member Data Documentation

◆ crossSectionHandler

G4RDVCrossSectionHandler* G4LowEnergyRayleigh::crossSectionHandler
private

Definition at line 98 of file G4LowEnergyRayleigh.hh.

◆ formFactorData

G4RDVEMDataSet* G4LowEnergyRayleigh::formFactorData
private

Definition at line 96 of file G4LowEnergyRayleigh.hh.

◆ highEnergyLimit

G4double G4LowEnergyRayleigh::highEnergyLimit
private

Definition at line 93 of file G4LowEnergyRayleigh.hh.

◆ intrinsicHighEnergyLimit

const G4double G4LowEnergyRayleigh::intrinsicHighEnergyLimit
private

Definition at line 101 of file G4LowEnergyRayleigh.hh.

◆ intrinsicLowEnergyLimit

const G4double G4LowEnergyRayleigh::intrinsicLowEnergyLimit
private

Definition at line 100 of file G4LowEnergyRayleigh.hh.

◆ lowEnergyLimit

G4double G4LowEnergyRayleigh::lowEnergyLimit
private

Definition at line 92 of file G4LowEnergyRayleigh.hh.

◆ meanFreePathTable

G4RDVEMDataSet* G4LowEnergyRayleigh::meanFreePathTable
private

Definition at line 95 of file G4LowEnergyRayleigh.hh.


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