|
||||
File indexing completed on 2025-01-18 09:57:58
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 // ------------------------------------------------------------------- 0028 // 0029 // GEANT4 Class header file 0030 // 0031 // 0032 // File name: G4BetheHeitler5DModel 0033 // 0034 // Authors: 0035 // Igor Semeniouk and Denis Bernard, 0036 // LLR, Ecole polytechnique & CNRS/IN2P3, 91128 Palaiseau, France 0037 // 0038 // Modifications: 0039 // 27-10-17 New class (IgS) 0040 // 19-01-18 version that calculates the pdf in the same way as in the fortran 0041 // version (Denis Bernard) 0042 // 04-06-18 Performance optimization of the final state sampling (M. Novak) 0043 // 10-04-19 CLHEP for boost and rotation, remove local functions (IgS) 0044 // 02-09-19 Set base class to be G4PairProductionRelModel that can provide 0045 // accurate sections now from threshold to very high energies 0046 // including the LPM effect. (M. Novak) 0047 // 14-10-19 Muon's pair genaration in SampleSecondaries 0048 // 0049 // Class Description: 0050 // 0051 // Implementation of gamma conversion to e+e- in the field of a nucleus 0052 // 0053 0054 // ------------------------------------------------------------------- 0055 // 0056 0057 #ifndef G4BetheHeitler5DModel_h 0058 #define G4BetheHeitler5DModel_h 1 0059 0060 #include "G4PairProductionRelModel.hh" 0061 0062 class G4IonTable; 0063 0064 class G4BetheHeitler5DModel : public G4PairProductionRelModel 0065 { 0066 0067 public: 0068 0069 explicit G4BetheHeitler5DModel(const G4ParticleDefinition* p = nullptr, 0070 const G4String& nam = "BetheHeitler5D"); 0071 0072 ~G4BetheHeitler5DModel() override; 0073 0074 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override; 0075 0076 void SampleSecondaries(std::vector<G4DynamicParticle*>* fvect, 0077 const G4MaterialCutsCouple* couple, 0078 const G4DynamicParticle* aDynamicGamma, 0079 G4double, G4double) override; 0080 0081 inline void SetVerbose(G4int val) { fVerbose = val; } 0082 0083 // Only e+, e+ or mu+, mu- pairs supported 0084 void SetLeptonPair(const G4ParticleDefinition* p1, 0085 const G4ParticleDefinition* p2); 0086 0087 // hide assignment operator 0088 G4BetheHeitler5DModel& operator=(const G4BetheHeitler5DModel& right) = delete; 0089 G4BetheHeitler5DModel(const G4BetheHeitler5DModel&) = delete; 0090 0091 private: 0092 0093 G4double MaxDiffCrossSection(const G4double* par, G4double eZ, 0094 G4double e, G4double loge) const; 0095 0096 inline void SetConversionMode(G4int to) { fConvMode = to; } 0097 0098 G4IonTable* theIonTable; 0099 0100 const G4ParticleDefinition* fLepton1; 0101 const G4ParticleDefinition* fLepton2; 0102 0103 const G4ParticleDefinition* fTheMuPlus; 0104 const G4ParticleDefinition* fTheMuMinus; 0105 0106 G4int fVerbose; 0107 G4int fConversionType; 0108 G4int fConvMode; 0109 G4bool iraw; 0110 }; 0111 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |