Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4OpAbsorption.hh
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.hh 69576 2013-05-08 13:48:13Z gcosmo $
28 //
29 ////////////////////////////////////////////////////////////////////////
30 // Optical Photon Absorption Class Definition
31 ////////////////////////////////////////////////////////////////////////
32 //
33 // File: G4OpAbsorption.hh
34 // Description: Discrete Process -- Bulk 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 // 1999-10-29 add method and class descriptors
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 // mail: gum@triumf.ca
45 // magni@mi.infn.it
46 //
47 ////////////////////////////////////////////////////////////////////////
48 
49 #ifndef G4OpAbsorption_h
50 #define G4OpAbsorption_h 1
51 
52 /////////////
53 // Includes
54 /////////////
55 
56 #include "globals.hh"
57 #include "templates.hh"
58 #include "Randomize.hh"
59 #include "G4Step.hh"
60 #include "G4VDiscreteProcess.hh"
61 #include "G4DynamicParticle.hh"
62 #include "G4Material.hh"
63 #include "G4OpticalPhoton.hh"
64 
65 // Class Description:
66 // Discrete Process -- Bulk absorption of Optical Photons.
67 // Class inherits publicly from G4VDiscreteProcess
68 // Class Description - End:
69 
70 /////////////////////
71 // Class Definition
72 /////////////////////
73 
75 {
76 
77 public:
78 
79  ////////////////////////////////
80  // Constructors and Destructor
81  ////////////////////////////////
82 
83  G4OpAbsorption(const G4String& processName = "OpAbsorption",
84  G4ProcessType type = fOptical);
86 
87 private:
88 
90 
91  //////////////
92  // Operators
93  //////////////
94 
95  G4OpAbsorption& operator=(const G4OpAbsorption &right);
96 
97 public:
98 
99  ////////////
100  // Methods
101  ////////////
102 
103  G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
104  // Returns true -> 'is applicable' only for an optical photon.
105 
106  G4double GetMeanFreePath(const G4Track& aTrack,
107  G4double ,
108  G4ForceCondition* );
109  // Returns the absorption length for bulk absorption of optical
110  // photons in media with a specified attenuation length.
111 
112  G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
113  const G4Step& aStep);
114  // This is the method implementing bulk absorption of optical
115  // photons.
116 
117 };
118 
119 ////////////////////
120 // Inline methods
121 ////////////////////
122 
123 inline
125 {
126  return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
127 }
128 
129 #endif /* G4OpAbsorption_h */
G4OpAbsorption(const G4String &processName="OpAbsorption", G4ProcessType type=fOptical)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
bool G4bool
Definition: G4Types.hh:79
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4Step.hh:76
static G4OpticalPhoton * OpticalPhoton()
double G4double
Definition: G4Types.hh:76
G4ForceCondition
G4ProcessType
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)