00001 // 00002 // ******************************************************************** 00003 // * License and Disclaimer * 00004 // * * 00005 // * The Geant4 software is copyright of the Copyright Holders of * 00006 // * the Geant4 Collaboration. It is provided under the terms and * 00007 // * conditions of the Geant4 Software License, included in the file * 00008 // * LICENSE and available at http://cern.ch/geant4/license . These * 00009 // * include a list of copyright holders. * 00010 // * * 00011 // * Neither the authors of this software system, nor their employing * 00012 // * institutes,nor the agencies providing financial support for this * 00013 // * work make any representation or warranty, express or implied, * 00014 // * regarding this software system or assume any liability for its * 00015 // * use. Please see the license in the file LICENSE and URL above * 00016 // * for the full disclaimer and the limitation of liability. * 00017 // * * 00018 // * This code implementation is the result of the scientific and * 00019 // * technical work of the GEANT4 collaboration. * 00020 // * By using, copying, modifying or distributing the software (or * 00021 // * any work based on the software) you agree to acknowledge its * 00022 // * use in resulting scientific publications, and indicate your * 00023 // * acceptance of all terms of the Geant4 Software license. * 00024 // ******************************************************************** 00025 // 00026 // $Id: G4AdjointPhotoElectricModel.hh 69844 2013-05-16 09:19:33Z gcosmo $ 00027 // 00029 // Module: G4AdjointPhotoElectricModel 00030 // Author: L. Desorgher 00031 // Organisation: SpaceIT GmbH 00032 // Contract: ESA contract 21435/08/NL/AT 00033 // Customer: ESA/ESTEC 00035 // 00036 // CHANGE HISTORY 00037 // -------------- 00038 // ChangeHistory: 00039 // -1 September 2007 creation by L. Desorgher 00040 // 00041 // -January 2009. L. Desorgher 00042 // Put a higher limit on the CS to avoid a high rate of Inverse Photo e- effect at low energy. The very high adjoint CS of the reverse 00043 // photo electric reaction produce a high rate of reverse photo electric reaction in the inner side of a shielding for eaxmple, the correction of this occurence 00044 // by weight correction in the StepDoIt method is not statistically sufficient at small energy. The problem is partially solved by setting an higher CS limit 00045 // and compensating it by an extra weight correction factor. However when coupling it with other reverse processes the reverse photo-electric is still 00046 // the source of very occasional high weight that decrease the efficiency of the computation. A way to solve this problemn is still needed but is difficult 00047 // to find as it happens in rarea case but does give a weighrt that is outside the noemal distribution. (Very Tricky!) 00048 // 00049 // -October 2009 Correction of Element sampling. L. Desorgher 00050 // 00051 //------------------------------------------------------------- 00052 // Documentation: 00053 // Model for the adjoint photo electric process 00054 // 00055 #ifndef G4AdjointPhotoElectricModel_h 00056 #define G4AdjointPhotoElectricModel_h 1 00057 00058 00059 #include "globals.hh" 00060 #include "G4VEmAdjointModel.hh" 00061 #include "G4PEEffectModel.hh" 00062 class G4AdjointPhotoElectricModel: public G4VEmAdjointModel 00063 00064 { 00065 public: 00066 00067 G4AdjointPhotoElectricModel(); 00068 ~G4AdjointPhotoElectricModel(); 00069 00070 00071 00072 virtual void SampleSecondaries(const G4Track& aTrack, 00073 G4bool IsScatProjToProjCase, 00074 G4ParticleChange* fParticleChange); 00075 virtual G4double AdjointCrossSection(const G4MaterialCutsCouple* aCouple, 00076 G4double primEnergy, 00077 G4bool IsScatProjToProjCase); 00078 virtual G4double GetAdjointCrossSection(const G4MaterialCutsCouple* aCouple, 00079 G4double primEnergy, 00080 G4bool IsScatProjToProjCase); 00081 00082 G4double AdjointCrossSectionPerAtom(const G4Element* anElement,G4double electronEnergy); 00083 00084 00085 00086 inline void SetTheDirectPEEffectModel(G4PEEffectModel* aModel){theDirectPEEffectModel = aModel; 00087 DefineDirectEMModel(aModel);} 00088 00089 virtual void CorrectPostStepWeight(G4ParticleChange* fParticleChange, 00090 G4double old_weight, 00091 G4double adjointPrimKinEnergy, 00092 G4double projectileKinEnergy, 00093 G4bool IsScatProjToProjCase); 00094 00095 00096 private: 00097 G4double xsec[40]; 00098 G4double totAdjointCS; 00099 G4double totBiasedAdjointCS; 00100 G4double factorCSBiasing; 00101 G4double pre_step_AdjointCS; 00102 G4double post_step_AdjointCS; 00103 00104 00105 G4double shell_prob[40][40]; 00106 00107 00108 G4PEEffectModel* theDirectPEEffectModel; 00109 size_t index_element; 00110 G4double current_eEnergy; 00111 00112 00113 private: 00114 void DefineCurrentMaterialAndElectronEnergy(const G4MaterialCutsCouple* aCouple, 00115 G4double eEnergy); 00116 00117 }; 00118 00119 #endif