|
||||
File indexing completed on 2025-01-18 09:57:52
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 // Class: G4AdjointPhotoElectricModel 0028 // Author: L. Desorgher 0029 // Organisation: SpaceIT GmbH 0030 // 0031 // Model for the adjoint photo electric process. 0032 // Put a higher limit on the CS to avoid a high rate of Inverse Photo e-effect 0033 // at low energy. The very high adjoint CS of the reverse photo electric 0034 // reaction produce a high rate of reverse photo electric reaction in the inner 0035 // side of a shielding for eaxmple, the correction of this occurrence by weight 0036 // correction in the StepDoIt method is not statistically sufficient at small 0037 // energy. The problem is partially solved by setting a higher CS limit and 0038 // compensating it by an extra weight correction factor. However when coupling 0039 // it with other reverse processes the reverse photo-electric is still the 0040 // source of very occasional high weights that decrease the efficiency of the 0041 // computation. A way to solve this problemn is still needed but is difficult 0042 // to find as it happens in rare cases but does give a weight that is outside 0043 // the normal distribution. (Very Tricky!) 0044 // 0045 //////////////////////////////////////////////////////////////////////////////// 0046 0047 #ifndef G4AdjointPhotoElectricModel_h 0048 #define G4AdjointPhotoElectricModel_h 1 0049 0050 #include "globals.hh" 0051 #include "G4VEmAdjointModel.hh" 0052 0053 class G4AdjointPhotoElectricModel : public G4VEmAdjointModel 0054 { 0055 public: 0056 G4AdjointPhotoElectricModel(); 0057 ~G4AdjointPhotoElectricModel() override; 0058 0059 void SampleSecondaries(const G4Track& aTrack, G4bool isScatProjToProj, 0060 G4ParticleChange* fParticleChange) override; 0061 0062 G4double AdjointCrossSection(const G4MaterialCutsCouple* aCouple, 0063 G4double primEnergy, 0064 G4bool isScatProjToProj) override; 0065 0066 G4double AdjointCrossSectionPerAtom(const G4Element* anElement, 0067 G4double electronEnergy); 0068 0069 G4AdjointPhotoElectricModel(G4AdjointPhotoElectricModel&) = delete; 0070 G4AdjointPhotoElectricModel& operator=( 0071 const G4AdjointPhotoElectricModel& right) = delete; 0072 0073 protected: 0074 void CorrectPostStepWeight(G4ParticleChange* fParticleChange, 0075 G4double old_weight, G4double adjointPrimKinEnergy, 0076 G4double projectileKinEnergy, 0077 G4bool isScatProjToProj) override; 0078 0079 private: 0080 void DefineCurrentMaterialAndElectronEnergy( 0081 const G4MaterialCutsCouple* aCouple, G4double eEnergy); 0082 0083 G4double fShellProb[40][40]; 0084 G4double fXsec[40]; 0085 G4double fTotAdjointCS = 0.; 0086 G4double fFactorCSBiasing = 1.; 0087 G4double fPreStepAdjointCS = 0.; 0088 G4double fPostStepAdjointCS = 0.; 0089 G4double fCurrenteEnergy = 0.; 0090 0091 size_t fIndexElement = 0; 0092 }; 0093 0094 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |