Geant4 Cross Reference |
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 // G4IonCoulombScatteringModel.hh 27 // ------------------------------------------------------------------- 28 // 29 // GEANT4 Class header file 30 // 31 // File name: G4IonCoulombScatteringModel 32 // 33 // Author: Cristina Consolandi 34 // 35 // Creation date: 05.10.2010 from G4eCoulombScatteringModel 36 // & G4CoulombScatteringModel 37 // 38 // 39 // Class Description: 40 // Single Scattering Model for 41 // for protons, alpha and heavy Ions 42 // 43 // Reference: 44 // M.J. Boschini et al. "Nuclear and Non-Ionizing Energy-Loss 45 // for Coulomb ScatteredParticles from Low Energy up to Relativistic 46 // Regime in Space Radiation Environment" 47 // Accepted for publication in the Proceedings of the ICATPP Conference 48 // on Cosmic Rays for Particle and Astroparticle Physics, Villa Olmo, 7-8 49 // October, 2010, to be published by World Scientific (Singapore). 50 // 51 // Available for downloading at: 52 // http://arxiv.org/abs/1011.4822 53 // 54 // ------------------------------------------------------------------- 55 // 56 57 #ifndef G4IonCoulombScatteringModel_h 58 #define G4IonCoulombScatteringModel_h 1 59 60 #include "G4VEmModel.hh" 61 #include "globals.hh" 62 #include "G4NistManager.hh" 63 #include "G4IonCoulombCrossSection.hh" 64 65 #include <vector> 66 67 class G4ParticleChangeForGamma; 68 class G4ParticleDefinition; 69 class G4IonTable; 70 71 class G4IonCoulombScatteringModel : public G4VEmModel 72 { 73 public: 74 75 explicit G4IonCoulombScatteringModel(const G4String& nam = 76 "IonCoulombScattering"); 77 78 ~G4IonCoulombScatteringModel() override; 79 80 void Initialise(const G4ParticleDefinition*, const G4DataVector&) final; 81 82 G4double ComputeCrossSectionPerAtom( 83 const G4ParticleDefinition*, 84 G4double kinEnergy, 85 G4double Z, 86 G4double A, 87 G4double cut, 88 G4double emax) final; 89 90 void SampleSecondaries(std::vector<G4DynamicParticle*>*, 91 const G4MaterialCutsCouple*, 92 const G4DynamicParticle*, 93 G4double tmin, 94 G4double maxEnergy) final; 95 96 97 98 inline void SetRecoilThreshold(G4double eth); 99 100 inline void SetHeavyIonCorr(G4int b); 101 102 inline G4int GetHeavyIonCorr(); 103 104 // hide assignment operator 105 G4IonCoulombScatteringModel & operator= 106 (const G4IonCoulombScatteringModel &right) = delete; 107 G4IonCoulombScatteringModel(const G4IonCoulombScatteringModel&) = delete; 108 109 private: 110 111 inline void DefineMaterial(const G4MaterialCutsCouple*); 112 113 inline void SetupParticle(const G4ParticleDefinition*); 114 115 G4IonTable* theIonTable; 116 G4ParticleChangeForGamma* fParticleChange; 117 G4NistManager* fNistManager; 118 G4IonCoulombCrossSection* ioncross; 119 120 const std::vector<G4double>* pCuts; 121 const G4MaterialCutsCouple* currentCouple; 122 const G4Material* currentMaterial; 123 const G4Element* currentElement; 124 G4int currentMaterialIndex; 125 126 G4int heavycorr; 127 128 G4double cosThetaMin; 129 G4double recoilThreshold; 130 131 // projectile 132 const G4ParticleDefinition* particle; 133 const G4ParticleDefinition* theProton; 134 G4double mass; 135 136 }; 137 138 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 139 140 inline void 141 G4IonCoulombScatteringModel::DefineMaterial(const G4MaterialCutsCouple* cup) 142 { 143 if(cup != currentCouple) { 144 currentCouple = cup; 145 currentMaterial = cup->GetMaterial(); 146 currentMaterialIndex = currentCouple->GetIndex(); 147 } 148 } 149 150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 151 152 inline 153 void G4IonCoulombScatteringModel::SetupParticle(const G4ParticleDefinition* p) 154 { 155 if(p != particle) { 156 particle = p; 157 mass = particle->GetPDGMass(); 158 ioncross->SetupParticle(p); 159 } 160 } 161 162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 163 164 inline void G4IonCoulombScatteringModel::SetRecoilThreshold(G4double eth) 165 { 166 recoilThreshold = eth; 167 } 168 169 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 170 171 inline void G4IonCoulombScatteringModel::SetHeavyIonCorr(G4int b) 172 { 173 heavycorr = b; 174 } 175 176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 177 178 inline G4int G4IonCoulombScatteringModel::GetHeavyIonCorr() 179 { 180 return heavycorr; 181 } 182 183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 184 185 #endif 186