Geant4  10.01.p02
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 // $Id: G4OpAbsorption.cc 69576 2013-05-08 13:48:13Z gcosmo $
28 //
30 // Optical Photon Absorption Class Implementation
32 //
33 // File: G4OpAbsorption.cc
34 // Description: Discrete Process -- Absorption of Optical Photons
35 // Version: 1.0
36 // Created: 1996-05-21
37 // Author: Juliet Armstrong
38 // Updated: 2005-07-28 - add G4ProcessType to constructor
39 // 2000-09-18 by Peter Gumplinger
40 // > comment out warning - "No Absorption length specified"
41 // 1997-04-09 by Peter Gumplinger
42 // > new physics/tracking scheme
43 // 1998-08-25 by Stefano Magni
44 // > Change process to use G4MaterialPropertiesTables
45 // 1998-09-03 by Peter Gumplinger
46 // > Protect G4MaterialPropertyVector* AttenuationLengthVector
47 // mail: gum@triumf.ca
48 // magni@mi.infn.it
49 //
51 
52 #include "G4ios.hh"
53 #include "G4OpProcessSubType.hh"
54 
55 #include "G4OpAbsorption.hh"
56 
58 // Class Implementation
60 
62  // Operators
64 
65 // G4OpAbsorption::operator=(const G4OpAbsorption &right)
66 // {
67 // }
68 
70  // Constructors
72 
74  : G4VDiscreteProcess(processName, type)
75 {
76  if (verboseLevel>0) {
77  G4cout << GetProcessName() << " is created " << G4endl;
78  }
79 
81 }
82 
83 // G4OpAbsorption::G4OpAbsorption(const G4OpAbsorpton &right)
84 // {
85 // }
86 
88  // Destructors
90 
92 
94  // Methods
96 
97 // PostStepDoIt
98 // -------------
99 //
101 G4OpAbsorption::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
102 {
103  aParticleChange.Initialize(aTrack);
104 
106 
107  if (verboseLevel>0) {
108  G4cout << "\n** Photon absorbed! **" << G4endl;
109  }
110  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
111 }
112 
113 
114 // GetMeanFreePath
115 // ---------------
116 //
118  G4double ,
120 {
121  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
122  const G4Material* aMaterial = aTrack.GetMaterial();
123 
124  G4double thePhotonMomentum = aParticle->GetTotalMomentum();
125 
126  G4MaterialPropertiesTable* aMaterialPropertyTable;
127  G4MaterialPropertyVector* AttenuationLengthVector;
128 
129  G4double AttenuationLength = DBL_MAX;
130 
131  aMaterialPropertyTable = aMaterial->GetMaterialPropertiesTable();
132 
133  if ( aMaterialPropertyTable ) {
134  AttenuationLengthVector = aMaterialPropertyTable->
135  GetProperty("ABSLENGTH");
136  if ( AttenuationLengthVector ){
137  AttenuationLength = AttenuationLengthVector->
138  Value(thePhotonMomentum);
139  }
140  else {
141 // G4cout << "No Absorption length specified" << G4endl;
142  }
143  }
144  else {
145 // G4cout << "No Absorption length specified" << G4endl;
146  }
147 
148  return AttenuationLength;
149 }
G4OpAbsorption(const G4String &processName="OpAbsorption", G4ProcessType type=fOptical)
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4DynamicParticle * GetDynamicParticle() const
G4double GetTotalMomentum() const
G4GLOB_DLL std::ostream G4cout
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Step.hh:76
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4Material * GetMaterial() const
virtual void Initialize(const G4Track &)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
#define DBL_MAX
Definition: templates.hh:83
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4ProcessType
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)