Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/optical/include/G4OpRayleigh.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //
 28 //
 29 ////////////////////////////////////////////////////////////////////////
 30 // Optical Photon Rayleigh Scattering Class Definition
 31 ////////////////////////////////////////////////////////////////////////
 32 //
 33 // File:        G4OpRayleigh.hh
 34 // Description: Discrete Process -- Rayleigh scattering of optical photons
 35 // Version:     1.0
 36 // Created:     1996-05-31
 37 // Author:      Juliet Armstrong
 38 // Updated:     2014-08-20 allow for more material types
 39 //              2005-07-28 add G4ProcessType to constructor
 40 //              1999-10-29 add method and class descriptors
 41 //              1997-04-09 by Peter Gumplinger
 42 //              > new physics/tracking scheme
 43 //
 44 ////////////////////////////////////////////////////////////////////////
 45 
 46 #ifndef G4OpRayleigh_h
 47 #define G4OpRayleigh_h 1
 48 
 49 #include "G4VDiscreteProcess.hh"
 50 #include "G4OpticalPhoton.hh"
 51 #include "G4PhysicsTable.hh"
 52 
 53 class G4OpRayleigh : public G4VDiscreteProcess
 54 {
 55  public:
 56   explicit G4OpRayleigh(const G4String& processName = "OpRayleigh",
 57                         G4ProcessType type          = fOptical);
 58   virtual ~G4OpRayleigh();
 59 
 60   virtual G4bool IsApplicable(
 61     const G4ParticleDefinition& aParticleType) override;
 62   // Returns true -> 'is applicable' only for an optical photon.
 63 
 64   virtual void BuildPhysicsTable(
 65     const G4ParticleDefinition& aParticleType) override;
 66   // Build thePhysicsTable at a right time
 67 
 68   virtual G4double GetMeanFreePath(const G4Track& aTrack, G4double,
 69                                    G4ForceCondition*) override;
 70   // Returns the mean free path for Rayleigh scattering
 71 
 72   virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
 73                                           const G4Step& aStep) override;
 74   // This is the method implementing Rayleigh scattering.
 75 
 76   virtual G4PhysicsTable* GetPhysicsTable() const;
 77   // Returns the address of the physics table.
 78 
 79   virtual void DumpPhysicsTable() const;
 80   // Prints the physics table.
 81 
 82   virtual void PreparePhysicsTable(const G4ParticleDefinition&) override;
 83   virtual void Initialise();
 84 
 85   void SetVerboseLevel(G4int);
 86 
 87  protected:
 88   G4PhysicsTable* thePhysicsTable;
 89 
 90  private:
 91   G4OpRayleigh(const G4OpRayleigh& right) = delete;
 92   G4OpRayleigh& operator=(const G4OpRayleigh& right) = delete;
 93 
 94   /// Calculates the mean free paths for a material as a function of
 95   /// photon energy
 96   G4PhysicsFreeVector* CalculateRayleighMeanFreePaths(
 97     const G4Material* material) const;
 98 
 99   size_t idx_rslength = 0;
100 };
101 
102 ////////////////////
103 // Inline methods
104 ////////////////////
105 
106 inline G4bool G4OpRayleigh::IsApplicable(
107   const G4ParticleDefinition& aParticleType)
108 {
109   return (&aParticleType == G4OpticalPhoton::OpticalPhoton());
110 }
111 
112 inline void G4OpRayleigh::DumpPhysicsTable() const
113 {
114   for(size_t i = 0; i < thePhysicsTable->entries(); ++i)
115   {
116     ((G4PhysicsFreeVector*) (*thePhysicsTable)[i])->DumpValues();
117   }
118 }
119 
120 inline G4PhysicsTable* G4OpRayleigh::GetPhysicsTable() const
121 {
122   return thePhysicsTable;
123 }
124 
125 #endif /* G4OpRayleigh_h */
126