Geant4-11
G4OpAbsorption.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
29// Optical Photon Absorption Class Implementation
31//
32// File: G4OpAbsorption.cc
33// Description: Discrete Process -- Absorption of Optical Photons
34// Version: 1.0
35// Created: 1996-05-21
36// Author: Juliet Armstrong
37// Updated: 2005-07-28 - add G4ProcessType to constructor
38// 2000-09-18 by Peter Gumplinger
39// > comment out warning - "No Absorption length specified"
40// 1997-04-09 by Peter Gumplinger
41// > new physics/tracking scheme
42// 1998-08-25 by Stefano Magni
43// > Change process to use G4MaterialPropertiesTables
44// 1998-09-03 by Peter Gumplinger
45// > Protect G4MaterialPropertyVector* AttenuationLengthVector
46//
48
49#include "G4ios.hh"
50#include "G4OpProcessSubType.hh"
52
53#include "G4OpAbsorption.hh"
54
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 : G4VDiscreteProcess(processName, type)
58{
59 Initialise();
60 if(verboseLevel > 0)
61 {
62 G4cout << GetProcessName() << " is created " << G4endl;
63 }
65}
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72{
73 Initialise();
74}
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
78{
79 SetVerboseLevel(G4OpticalParameters::Instance()->GetAbsorptionVerboseLevel());
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84 const G4Step& aStep)
85{
87
88 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
89 G4double thePhotonMomentum = aParticle->GetTotalMomentum();
90
93
94 if(verboseLevel > 1)
95 {
96 G4cout << "\n** OpAbsorption: Photon absorbed! **" << G4endl;
97 }
98 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
99}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
104{
105 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
108 G4double attLength = DBL_MAX;
109
110 if(MPT)
111 {
113 if(attVector)
114 {
115 attLength =
116 attVector->Value(aParticle->GetTotalMomentum(), idx_absorption);
117 }
118 }
119
120 return attLength;
121}
122
123//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125{
126 verboseLevel = verbose;
128}
G4ForceCondition
@ fOpAbsorption
G4ProcessType
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetTotalMomentum() const
G4MaterialPropertyVector * GetProperty(const char *key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
void SetVerboseLevel(G4int)
G4OpAbsorption(const G4String &processName="OpAbsorption", G4ProcessType type=fOptical)
virtual ~G4OpAbsorption()
virtual void Initialise()
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
static G4OpticalParameters * Instance()
void SetAbsorptionVerboseLevel(G4int)
virtual void Initialize(const G4Track &)
G4double Value(const G4double energy, std::size_t &lastidx) const
Definition: G4Step.hh:62
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
#define DBL_MAX
Definition: templates.hh:62