Geant4  10.00.p03
G4OpMieHG.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 //
27 //
28 // File G4OpMieHG.hh
29 // Description: Discrete Process -- Mie Scattering of Optical Photons
30 // Created: 2010-07-03
31 // Author: Xin Qian
32 // Based on work from Vlasios Vasileiou
33 //
34 // This subroutine will mimic the Mie scattering based on
35 // Henyey-Greenstein phase function
36 // Forward and backward angles are treated separately.
37 //
38 // mail: gum@triumf.ca
39 //
41 
42 #ifndef G4OpMieHG_h
43 #define G4OpMieHG_h 1
44 
45 #include "G4VDiscreteProcess.hh"
46 #include "G4OpticalPhoton.hh"
47 
49 {
50 
51 public:
52 
54  // Constructors and Destructor
56 
57  G4OpMieHG(const G4String& processName = "OpMieHG",
58  G4ProcessType type = fOptical);
59  ~G4OpMieHG();
60 
61 private:
62 
63  G4OpMieHG(const G4OpMieHG &right);
64 
66  // Operators
68 
69  G4OpMieHG& operator=(const G4OpMieHG &right);
70 
71 public:
72 
74  // Methods
76 
77  G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
78  // Returns true -> 'is applicable' only for an optical photon.
79 
80  G4double GetMeanFreePath(const G4Track& aTrack,
81  G4double,
83  // Return the mean free path of Mie scattering
84 
85  G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
86  const G4Step& aStep);
87  // This is the method implementing Mie scattering.
88 };
89 
90 
91 inline
93 {
94  return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
95 }
96 
97 #endif /* G4OpMieHG_h */
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4OpMieHG.cc:159
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4OpMieHG.cc:67
bool G4bool
Definition: G4Types.hh:79
Definition: G4Step.hh:76
static G4OpticalPhoton * OpticalPhoton()
G4OpMieHG(const G4String &processName="OpMieHG", G4ProcessType type=fOptical)
Definition: G4OpMieHG.cc:47
double G4double
Definition: G4Types.hh:76
~G4OpMieHG()
Definition: G4OpMieHG.cc:57
G4ForceCondition
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
Definition: G4OpMieHG.hh:92
G4OpMieHG & operator=(const G4OpMieHG &right)
G4ProcessType