Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/lowenergy/include/G4MicroElecElasticModel_new.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //
 27 // G4MicroElecElasticModel_new.hh, 2011/08/29 A.Valentin, M. Raine are with CEA [a]
 28 //                          2020/05/20 P. Caron, C. Inguimbert are with ONERA [b] 
 29 //                   Q. Gibaru is with CEA [a], ONERA [b] and CNES [c]
 30 //                   M. Raine and D. Lambert are with CEA [a]
 31 //
 32 // A part of this work has been funded by the French space agency(CNES[c])
 33 // [a] CEA, DAM, DIF - 91297 ARPAJON, France
 34 // [b] ONERA - DPHY, 2 avenue E.Belin, 31055 Toulouse, France
 35 // [c] CNES, 18 av.E.Belin, 31401 Toulouse CEDEX, France
 36 //
 37 // Based on the following publications
 38 //  - A.Valentin, M. Raine, 
 39 //    Inelastic cross-sections of low energy electrons in silicon
 40 //        for the simulation of heavy ion tracks with the Geant4-DNA toolkit,
 41 //        NSS Conf. Record 2010, pp. 80-85
 42 //             https://doi.org/10.1109/NSSMIC.2010.5873720
 43 //
 44 //      - A.Valentin, M. Raine, M.Gaillardin, P.Paillet
 45 //        Geant4 physics processes for microdosimetry simulation:
 46 //        very low energy electromagnetic models for electrons in Silicon,
 47 //             https://doi.org/10.1016/j.nimb.2012.06.007
 48 //        NIM B, vol. 288, pp. 66-73, 2012, part A
 49 //        heavy ions in Si, NIM B, vol. 287, pp. 124-129, 2012, part B
 50 //             https://doi.org/10.1016/j.nimb.2012.07.028
 51 //
 52 //  - M. Raine, M. Gaillardin, P. Paillet
 53 //        Geant4 physics processes for silicon microdosimetry simulation: 
 54 //        Improvements and extension of the energy-range validity up to 10 GeV/nucleon
 55 //        NIM B, vol. 325, pp. 97-100, 2014
 56 //             https://doi.org/10.1016/j.nimb.2014.01.014
 57 //
 58 //      - J. Pierron, C. Inguimbert, M. Belhaj, T. Gineste, J. Puech, M. Raine
 59 //        Electron emission yield for low energy electrons: 
 60 //        Monte Carlo simulation and experimental comparison for Al, Ag, and Si
 61 //        Journal of Applied Physics 121 (2017) 215107. 
 62 //               https://doi.org/10.1063/1.4984761
 63 //
 64 //      - P. Caron,
 65 //        Study of Electron-Induced Single-Event Upset in Integrated Memory Devices
 66 //        PHD, 16th October 2019
 67 //
 68 //  - Q.Gibaru, C.Inguimbert, P.Caron, M.Raine, D.Lambert, J.Puech, 
 69 //        Geant4 physics processes for microdosimetry and secondary electron emission simulation : 
 70 //        Extension of MicroElec to very low energies and new materials
 71 //        NIM B, 2020, in review.
 72 //
 73 //
 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 
 75 
 76 #ifndef G4MICROELECELASTICMODEL_NEW_HH
 77 #define G4MICROELECELASTICMODEL_NEW_HH 1
 78 
 79 #include <map>
 80 #include <CLHEP/Units/SystemOfUnits.h>
 81 #include "G4MicroElecMaterialStructure.hh"
 82 #include "G4MicroElecCrossSectionDataSet_new.hh"
 83 #include "G4VEmModel.hh"
 84 #include "G4Electron.hh"
 85 #include "G4ParticleChangeForGamma.hh"
 86 #include "G4LogLogInterpolation.hh"
 87 #include "G4ProductionCutsTable.hh"
 88 #include "G4NistManager.hh"
 89 
 90 class G4MicroElecElasticModel_new : public G4VEmModel
 91 {
 92 
 93 public:
 94   G4MicroElecElasticModel_new(const G4ParticleDefinition* p = 0, 
 95               const G4String& nam = "MicroElecElasticModel");
 96   ~G4MicroElecElasticModel_new() override;
 97   
 98   void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
 99 
100   G4double CrossSectionPerVolume(const G4Material* material,
101          const G4ParticleDefinition* p,
102          G4double ekin,
103          G4double emin,
104          G4double emax) override;
105 
106   G4double AcousticCrossSectionPerVolume(G4double ekin, G4double kbz, G4double rho,
107            G4double cs, G4double Aac, G4double Eac,
108            G4double prefactor);
109 
110   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
111        const G4MaterialCutsCouple*,
112        const G4DynamicParticle*,
113        G4double tmin,
114        G4double maxEnergy) override;
115 
116   void SetKillBelowThreshold (G4double threshold);     
117 
118   G4double GetKillBelowThreshold () { return killBelowEnergy; } 
119 
120   G4double DamageEnergy(G4double T,G4double A, G4double Z);
121 
122 protected:
123   G4ParticleChangeForGamma* fParticleChangeForGamma;
124   
125 private:
126 
127   G4MicroElecElasticModel_new & operator=(const  G4MicroElecElasticModel_new &right);
128   G4MicroElecElasticModel_new(const  G4MicroElecElasticModel_new&);
129 
130   // Final state
131   G4double Theta(G4ParticleDefinition * aParticleDefinition, G4double k, G4double integrDiff);
132   G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
133   G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
134   G4double LinLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
135   G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, 
136           G4double x11, G4double x12, G4double x21, G4double x22, 
137           G4double t1, G4double t2, G4double t, G4double e);
138 
139   G4double RandomizeCosTheta(G4double k);
140 
141   G4Material* nistSi = nullptr;
142   G4double killBelowEnergy;  
143   G4double lowEnergyLimit;  
144   G4double lowEnergyLimitOfModel;  
145   G4double highEnergyLimit; 
146   G4bool isInitialised;
147   G4int verboseLevel;
148   // Cross section
149   typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
150   MapFile tableFile;
151   typedef std::map<G4String,G4MicroElecCrossSectionDataSet_new*,std::less<G4String> > MapData;
152   //MapData tableData;
153   
154   typedef std::map<G4String, MapData*, std::less<G4String> > TCSMap;
155   TCSMap tableTCS;
156 
157   //Maps for multilayers
158   typedef std::map<G4double, std::map<G4double, G4double> > TriDimensionMap;
159 
160   typedef std::map<G4String, TriDimensionMap* > ThetaMap;
161   ThetaMap thetaDataStorage; //Storage of angles (cumulated)
162 
163   typedef std::map<G4String, std::vector<G4double>* > energyMap;
164   energyMap eIncidentEnergyStorage;
165 
166   typedef std::map<G4double, std::vector<G4double> > VecMap;
167 
168   typedef std::map<G4String, VecMap* > ProbaMap;
169   ProbaMap eProbaStorage; //Storage of probabilities for cumulated sections
170 
171   typedef std::map<G4String, G4MicroElecMaterialStructure*, std::less<G4String> > MapStructure;
172 
173   MapStructure tableMaterialsStructures; //Structures of all materials simulated
174 
175   G4MicroElecMaterialStructure* currentMaterialStructure = nullptr;
176   typedef std::map<G4String, G4double, std::less<G4String> > MapEnergy;
177   MapEnergy lowEnergyLimitTable;
178   MapEnergy highEnergyLimitTable;
179   MapEnergy workFunctionTable;
180       
181   G4bool killElectron, acousticModelEnabled;
182   G4String currentMaterialName;
183   G4bool isOkToBeInitialised;
184 
185 };
186 
187 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
188 
189 #endif
190