G4VLowEnergyDiscretePhotonProcess Class Reference

#include <G4VLowEnergyDiscretePhotonProcess.hh>

Inheritance diagram for G4VLowEnergyDiscretePhotonProcess:

G4VLowEnergyTestableDiscreteProcess G4VDiscreteProcess G4VProcess

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

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)

Detailed Description

Definition at line 62 of file G4VLowEnergyDiscretePhotonProcess.hh.


Constructor & Destructor Documentation

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.

References G4cout, G4endl, G4VProcess::GetProcessName(), G4VEMDataSet::LoadData(), and G4VProcess::verboseLevel.

00063   :
00064   G4VLowEnergyTestableDiscreteProcess(processName),
00065   lowEnergyLimit(aLowEnergyLimit),
00066   highEnergyLimit(aHighEnergyLimit),
00067   crossSectionFileName(aCrossSectionFileName),
00068   meanFreePathTable(0)
00069 {
00070   crossSectionHandler = new G4CrossSectionHandler();
00071   scatterFunctionData = new G4CompositeEMDataSet(aScatterInterpolation, 1., 1.);
00072   scatterFunctionData->LoadData(aScatterFileName);
00073  
00074   if (verboseLevel > 0)
00075     {
00076       G4cout << GetProcessName() << " is created " << G4endl
00077              << "Energy range: "
00078              << lowEnergyLimit / keV << " keV - "
00079              << highEnergyLimit / GeV << " GeV"
00080              << G4endl;
00081     }
00082 }

G4VLowEnergyDiscretePhotonProcess::~G4VLowEnergyDiscretePhotonProcess ( void   )  [virtual]

Definition at line 86 of file G4VLowEnergyDiscretePhotonProcess.cc.

00087 {
00088   if (meanFreePathTable)
00089     delete meanFreePathTable;
00090  
00091   delete crossSectionHandler;
00092   delete scatterFunctionData;
00093 }


Member Function Documentation

void G4VLowEnergyDiscretePhotonProcess::BuildPhysicsTable ( const G4ParticleDefinition photon  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 106 of file G4VLowEnergyDiscretePhotonProcess.cc.

References G4VCrossSectionHandler::BuildMeanFreePathForMaterials(), G4VCrossSectionHandler::Clear(), and G4VCrossSectionHandler::LoadData().

00107 {
00108   crossSectionHandler->Clear();
00109   crossSectionHandler->LoadData(crossSectionFileName);
00110 
00111   if (meanFreePathTable)
00112     delete meanFreePathTable; 
00113   meanFreePathTable=crossSectionHandler->BuildMeanFreePathForMaterials();
00114 }

const G4VCrossSectionHandler* G4VLowEnergyDiscretePhotonProcess::GetCrossSectionHandler ( void   )  const [inline, protected]

Definition at line 129 of file G4VLowEnergyDiscretePhotonProcess.hh.

00129 {return crossSectionHandler;}

G4double G4VLowEnergyDiscretePhotonProcess::GetHighEnergyLimit ( void   )  const [inline]

Definition at line 110 of file G4VLowEnergyDiscretePhotonProcess.hh.

00110 {return highEnergyLimit;}

G4double G4VLowEnergyDiscretePhotonProcess::GetLowEnergyLimit ( void   )  const [inline]

Definition at line 106 of file G4VLowEnergyDiscretePhotonProcess.hh.

00106 {return lowEnergyLimit;}

G4double G4VLowEnergyDiscretePhotonProcess::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
) [protected, virtual]

Implements G4VDiscreteProcess.

Definition at line 120 of file G4VLowEnergyDiscretePhotonProcess.cc.

References DBL_MAX, G4VEMDataSet::FindValue(), G4Track::GetDynamicParticle(), G4MaterialCutsCouple::GetIndex(), G4DynamicParticle::GetKineticEnergy(), and G4Track::GetMaterialCutsCouple().

00121 {
00122   G4double photonEnergy;
00123   photonEnergy = aTrack.GetDynamicParticle()->GetKineticEnergy();
00124 
00125   if (photonEnergy < lowEnergyLimit)
00126     return DBL_MAX;
00127  
00128   if (photonEnergy > highEnergyLimit)
00129     photonEnergy=highEnergyLimit;
00130  
00131   size_t materialIndex;
00132   materialIndex = aTrack.GetMaterialCutsCouple()->GetIndex(); 
00133 
00134   return meanFreePathTable->FindValue(photonEnergy, materialIndex);
00135 }

const G4VEMDataSet* G4VLowEnergyDiscretePhotonProcess::GetMeanFreePathTable ( void   )  const [inline, protected]

Definition at line 133 of file G4VLowEnergyDiscretePhotonProcess.hh.

00133 {return meanFreePathTable;}

G4ThreeVector G4VLowEnergyDiscretePhotonProcess::GetPhotonPolarization ( const G4DynamicParticle photon  )  [protected]

Definition at line 141 of file G4VLowEnergyDiscretePhotonProcess.cc.

References G4UniformRand, and photon.

00142 {
00143   G4ThreeVector photonMomentumDirection;
00144   G4ThreeVector photonPolarization;
00145 
00146   photonPolarization = photon.GetPolarization(); 
00147   photonMomentumDirection = photon.GetMomentumDirection();
00148 
00149   if ((!photonPolarization.isOrthogonal(photonMomentumDirection, 1e-6)) || photonPolarization.mag()==0.)
00150     {
00151       // if |photonPolarization|==0. or |photonPolarization * photonDirection0| > 1e-6 * |photonPolarization ^ photonDirection0|
00152       // then polarization is choosen randomly.
00153   
00154       G4ThreeVector e1(photonMomentumDirection.orthogonal().unit());
00155       G4ThreeVector e2(photonMomentumDirection.cross(e1).unit());
00156   
00157       G4double angle(G4UniformRand() * twopi);
00158   
00159       e1*=std::cos(angle);
00160       e2*=std::sin(angle);
00161   
00162       photonPolarization=e1+e2;
00163     }
00164   else if (photonPolarization.howOrthogonal(photonMomentumDirection) != 0.)
00165     {
00166       // if |photonPolarization * photonDirection0| != 0.
00167       // then polarization is made orthonormal;
00168   
00169       photonPolarization=photonPolarization.perpPart(photonMomentumDirection);
00170     }
00171  
00172   return photonPolarization.unit();
00173 }

const G4VEMDataSet* G4VLowEnergyDiscretePhotonProcess::GetScatterFunctionData ( void   )  const [inline, protected]

Definition at line 137 of file G4VLowEnergyDiscretePhotonProcess.hh.

00137 {return scatterFunctionData;}

G4bool G4VLowEnergyDiscretePhotonProcess::IsApplicable ( const G4ParticleDefinition particleDefinition  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 99 of file G4VLowEnergyDiscretePhotonProcess.cc.

References G4Gamma::Gamma().

00100 {
00101   return (&particleDefinition)==G4Gamma::Gamma(); 
00102 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:50 2013 for Geant4 by  doxygen 1.4.7