Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:14

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 //  G4eSingleCoulombScatteringModel.hh
0027 // -------------------------------------------------------------------
0028 //
0029 // GEANT4 Class header file
0030 //
0031 // File name:    G4eSingleCoulombScatteringModel
0032 //
0033 // Author:  Cristina Consolandi
0034 //
0035 // Creation date: 20.10.2011
0036 //
0037 // Class Description:
0038 //      Single Scattering model for electron-nuclei interaction.
0039 //      Suitable for high energy electrons and low scattering angles.
0040 //
0041 // Reference:
0042 //      M.J. Boschini et al.
0043 //      "Non Ionizing Energy Loss induced by Electrons in the Space Environment"
0044 //       Proc. of the 13th International Conference on Particle Physics and Advanced Technology
0045 //       (13th ICPPAT, Como 3-7/10/2011), World Scientific (Singapore).
0046 //      Available at: http://arxiv.org/abs/1111.4042v4
0047 //
0048 //
0049 // -------------------------------------------------------------------
0050 //
0051 
0052 #ifndef G4eSingleCoulombScatteringModel_h
0053 #define G4eSingleCoulombScatteringModel_h 1
0054 
0055 #include "G4VEmModel.hh"
0056 #include "globals.hh"
0057 #include "G4ScreeningMottCrossSection.hh"
0058 
0059 #include <vector>
0060 
0061 class G4ParticleChangeForGamma;
0062 class G4ParticleDefinition;
0063 class G4IonTable;
0064 class G4Nistmanager;
0065 
0066 class G4eSingleCoulombScatteringModel : public G4VEmModel
0067 {
0068 
0069 public:
0070 
0071   explicit G4eSingleCoulombScatteringModel(const G4String& nam = "eSingleCoulombScat");
0072 
0073   ~G4eSingleCoulombScatteringModel() override;
0074 
0075   void Initialise(const G4ParticleDefinition*, const G4DataVector&) final;
0076 
0077   void InitialiseLocal(const G4ParticleDefinition*,
0078                G4VEmModel* masterModel) final;
0079 
0080   G4double ComputeCrossSectionPerAtom(
0081                                 const G4ParticleDefinition*,
0082                 G4double kinEnergy,
0083                 G4double Z,
0084                 G4double A,
0085                 G4double cut,
0086                 G4double emax) final;
0087 
0088   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0089              const G4MaterialCutsCouple*,
0090              const G4DynamicParticle*,
0091              G4double tmin,
0092              G4double maxEnergy) final;
0093 
0094   inline void SetRecoilThreshold(G4double eth);
0095 
0096   void SetXSectionModel(const G4String& model);
0097 
0098   // hide assignment operator
0099   G4eSingleCoulombScatteringModel & operator=
0100   (const G4eSingleCoulombScatteringModel &right) = delete;
0101   G4eSingleCoulombScatteringModel(const G4eSingleCoulombScatteringModel&) = delete;
0102 
0103 private:
0104 
0105   inline void DefineMaterial(const G4MaterialCutsCouple*);
0106   inline void SetupParticle(const G4ParticleDefinition*);
0107 
0108   G4IonTable*               theIonTable;
0109   G4ParticleChangeForGamma* fParticleChange;
0110   G4NistManager*            fNistManager;
0111   G4ScreeningMottCrossSection* Mottcross;
0112   const G4ParticleDefinition* particle;
0113 
0114   const std::vector<G4double>* pCuts;
0115   const G4MaterialCutsCouple* currentCouple;
0116   const G4Material*           currentMaterial;
0117   const G4Element*            currentElement;
0118 
0119   G4int                     currentMaterialIndex;
0120   G4int                     FormFactor;
0121   G4int                     XSectionModel;
0122 
0123   G4double                  cosThetaMin;
0124   G4double                  recoilThreshold;
0125   G4double                  mass;
0126   G4double                  lowEnergyLimit;
0127 
0128 };
0129 
0130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0131 
0132 inline void
0133 G4eSingleCoulombScatteringModel::DefineMaterial(const G4MaterialCutsCouple* cup)
0134 {
0135   if(cup != currentCouple) {
0136     currentCouple = cup;
0137     currentMaterial = cup->GetMaterial();
0138     currentMaterialIndex = currentCouple->GetIndex();
0139   }
0140 }
0141 
0142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0143 
0144 inline void
0145 G4eSingleCoulombScatteringModel::SetupParticle(const G4ParticleDefinition* p)
0146 {
0147   if(p != particle) {
0148     particle = p;
0149     mass = particle->GetPDGMass();
0150     Mottcross->SetupParticle(p);
0151   }
0152 }
0153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0154 
0155 inline void G4eSingleCoulombScatteringModel::SetRecoilThreshold(G4double eth)
0156 {
0157   recoilThreshold = eth;
0158 }
0159 
0160 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
0161 
0162 #endif