|
||||
File indexing completed on 2025-01-18 09:58:55
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // 0027 //------------------ G4PhotoElectricEffect physics process ------------------ 0028 // by Michel Maire, 24 May 1996 0029 // 0030 // 12-06-96, Added SelectRandomAtom() method and new data member 0031 // for cumulative total cross section, by M.Maire 0032 // 21-06-96, SetCuts implementation, M.Maire 0033 // 17-09-96, Dynamic array PartialSumSigma 0034 // split ComputeBindingEnergy(), M.Maire 0035 // 08-01-97, crossection table + meanfreepath table, M.Maire 0036 // 13-03-97, adapted for the new physics scheme, M.Maire 0037 // 13-08-98, new methods SetBining() PrintInfo() 0038 // 17-11-98, use table of atomic shells in PostStepDoIt, mma 0039 // 06-01-99, Sandia crossSection below 50 keV, V.Grichine mma 0040 // 03-08-01, new methods Store/Retrieve PhysicsTable (mma) 0041 // 06-08-01, BuildThePhysicsTable() called from constructor (mma) 0042 // 19-09-01, come back to previous process name "phot" 0043 // 20-09-01, DoIt: fminimalEnergy = 1*eV (mma) 0044 // 01-10-01, come back to BuildPhysicsTable(const G4ParticleDefinition&) 0045 // 10-01-02, moved few function from icc to cc 0046 // 17-04-02, Keep only Sandia crossSections. Remove BuildPhysicsTables. 0047 // Simplify public interface (mma) 0048 // 29-04-02, Generate theta angle of the photoelectron from Sauter-Gavrila 0049 // distribution (mma) 0050 // 13-08-04, suppress icc file; make public ComputeCrossSectionPerAtom() (mma) 0051 // 21-04-05, Redesign - use G4VEmProcess interface (V.Ivanchenko) 0052 // 02-05-05, move ParticleChange actions in model (mma) 0053 // 04-05-05, Make class to be default (V.Ivanchenko) 0054 // 09-08-06, add SetModel(G4VEmModel*) (mma) 0055 // 12-09-06, move SetModel(G4VEmModel*) in G4VEmProcess (mma) 0056 // ----------------------------------------------------------------------------- 0057 0058 // class description 0059 // 0060 0061 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0062 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0063 0064 #ifndef G4PhotoElectricEffect_h 0065 #define G4PhotoElectricEffect_h 1 0066 0067 #include "globals.hh" 0068 #include "G4VEmProcess.hh" 0069 #include "G4Gamma.hh" 0070 0071 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0072 0073 class G4ParticleDefinition; 0074 class G4VEmModel; 0075 class G4MaterialCutsCouple; 0076 class G4DynamicParticle; 0077 0078 class G4PhotoElectricEffect : public G4VEmProcess 0079 0080 { 0081 public: // with description 0082 0083 explicit G4PhotoElectricEffect(const G4String& processName ="phot", 0084 G4ProcessType type = fElectromagnetic); 0085 0086 ~G4PhotoElectricEffect() override; 0087 0088 // true for Gamma only. 0089 G4bool IsApplicable(const G4ParticleDefinition&) final; 0090 0091 // print documentation in html format 0092 void ProcessDescription(std::ostream&) const override; 0093 0094 G4PhotoElectricEffect & operator=(const G4PhotoElectricEffect &right) = delete; 0095 G4PhotoElectricEffect(const G4PhotoElectricEffect&) = delete; 0096 0097 protected: 0098 0099 virtual void InitialiseProcess(const G4ParticleDefinition*) override; 0100 0101 private: 0102 0103 G4bool isInitialised = false; 0104 }; 0105 0106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0107 0108 #endif 0109
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |