G4OpAbsorption Class Reference

#include <G4OpAbsorption.hh>

Inheritance diagram for G4OpAbsorption:

G4VDiscreteProcess G4VProcess

Public Member Functions

 G4OpAbsorption (const G4String &processName="OpAbsorption", G4ProcessType type=fOptical)
 ~G4OpAbsorption ()
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)

Detailed Description

Definition at line 74 of file G4OpAbsorption.hh.


Constructor & Destructor Documentation

G4OpAbsorption::G4OpAbsorption ( const G4String processName = "OpAbsorption",
G4ProcessType  type = fOptical 
)

Definition at line 73 of file G4OpAbsorption.cc.

References fOpAbsorption, G4cout, G4endl, G4VProcess::GetProcessName(), G4VProcess::SetProcessSubType(), and G4VProcess::verboseLevel.

00074               : G4VDiscreteProcess(processName, type)
00075 {
00076         if (verboseLevel>0) {
00077            G4cout << GetProcessName() << " is created " << G4endl;
00078         }
00079 
00080         SetProcessSubType(fOpAbsorption);
00081 }

G4OpAbsorption::~G4OpAbsorption (  ) 

Definition at line 91 of file G4OpAbsorption.cc.

00091 {}


Member Function Documentation

G4double G4OpAbsorption::GetMeanFreePath ( const G4Track aTrack,
G4double  ,
G4ForceCondition  
) [virtual]

Implements G4VDiscreteProcess.

Definition at line 117 of file G4OpAbsorption.cc.

References DBL_MAX, G4Track::GetDynamicParticle(), G4Track::GetMaterial(), G4Material::GetMaterialPropertiesTable(), and G4DynamicParticle::GetTotalMomentum().

00120 {
00121         const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
00122         const G4Material* aMaterial = aTrack.GetMaterial();
00123 
00124         G4double thePhotonMomentum = aParticle->GetTotalMomentum();
00125 
00126         G4MaterialPropertiesTable* aMaterialPropertyTable;
00127         G4MaterialPropertyVector* AttenuationLengthVector;
00128         
00129         G4double AttenuationLength = DBL_MAX;
00130 
00131         aMaterialPropertyTable = aMaterial->GetMaterialPropertiesTable();
00132 
00133         if ( aMaterialPropertyTable ) {
00134            AttenuationLengthVector = aMaterialPropertyTable->
00135                                                 GetProperty("ABSLENGTH");
00136            if ( AttenuationLengthVector ){
00137              AttenuationLength = AttenuationLengthVector->
00138                                          Value(thePhotonMomentum);
00139            }
00140            else {
00141 //             G4cout << "No Absorption length specified" << G4endl;
00142            }
00143         } 
00144         else {
00145 //           G4cout << "No Absorption length specified" << G4endl;
00146         }
00147 
00148         return AttenuationLength;
00149 }

G4bool G4OpAbsorption::IsApplicable ( const G4ParticleDefinition aParticleType  )  [inline, virtual]

Reimplemented from G4VProcess.

Definition at line 124 of file G4OpAbsorption.hh.

References G4OpticalPhoton::OpticalPhoton().

00125 {
00126    return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
00127 }

G4VParticleChange * G4OpAbsorption::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
) [virtual]

Reimplemented from G4VDiscreteProcess.

Definition at line 101 of file G4OpAbsorption.cc.

References G4VProcess::aParticleChange, fStopAndKill, G4cout, G4endl, G4ParticleChange::Initialize(), G4VDiscreteProcess::PostStepDoIt(), G4VParticleChange::ProposeTrackStatus(), and G4VProcess::verboseLevel.

00102 {
00103         aParticleChange.Initialize(aTrack);
00104 
00105         aParticleChange.ProposeTrackStatus(fStopAndKill);
00106 
00107         if (verboseLevel>0) {
00108            G4cout << "\n** Photon absorbed! **" << G4endl;
00109         }
00110         return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
00111 }


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