Back to home page

EIC code displayed by LXR

 
 

    


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

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:     G4ICRU73QOModel
0033 //
0034 // Author:        Alexander Bagulya
0035 //
0036 // Creation date: 21.05.2010
0037 //
0038 // Modifications:
0039 //
0040 //
0041 // Class Description:
0042 //
0043 // Quantum Harmonic Oscillator Model for energy loss using atomic shell 
0044 // structure of atoms taking into account Q^2 (main for projectile charge Q), 
0045 // Q^3 and Q^4 terms for computation of energy loss due to binary collisions. 
0046 // Can be applied on heavy negatively charged particles for the energy interval 
0047 // 10 keV - 10 MeV scaled to the proton mass.
0048 //
0049 // Used data and formula of 
0050 // 1. G4QAOLowEnergyLoss class, S.Chauvie, P.Nieminen, M.G.Pia. IEEE Trans. 
0051 //    Nucl. Sci. 54 (2007) 578.
0052 // 2. ShellStrength and ShellEnergy from ICRU'73 Report 2005,
0053 // 3. Data for Ta (Z=73) from P.Sigmund, A.Shinner. Eur. Phys. J. D15 (2001) 
0054 //    165-172
0055 //
0056 // -------------------------------------------------------------------
0057 //
0058 
0059 #ifndef G4ICRU73QOModel_h
0060 #define G4ICRU73QOModel_h 1
0061 
0062 #include <CLHEP/Units/PhysicalConstants.h>
0063 
0064 #include "G4VEmModel.hh"
0065 #include "G4AtomicShells.hh"
0066 #include "G4DensityEffectData.hh"
0067 
0068 class G4ParticleChangeForLoss;
0069 
0070 class G4ICRU73QOModel : public G4VEmModel
0071 {
0072 
0073 public:
0074 
0075   explicit G4ICRU73QOModel(const G4ParticleDefinition* p = nullptr,
0076                            const G4String& nam = "ICRU73QO");
0077 
0078   ~G4ICRU73QOModel() = default;
0079 
0080   void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
0081 
0082   G4double ComputeCrossSectionPerElectron(
0083                                  const G4ParticleDefinition*,
0084                                  G4double kineticEnergy,
0085                                  G4double cutEnergy,
0086                                  G4double maxEnergy);
0087                                  
0088   G4double ComputeCrossSectionPerAtom(
0089                                  const G4ParticleDefinition*,
0090                                  G4double kineticEnergy,
0091                                  G4double Z, G4double A,
0092                                  G4double cutEnergy,
0093                                  G4double maxEnergy) override;
0094                                                                   
0095   G4double CrossSectionPerVolume(const G4Material*,
0096                                  const G4ParticleDefinition*,
0097                                  G4double kineticEnergy,
0098                                  G4double cutEnergy,
0099                                  G4double maxEnergy) override;
0100                                  
0101   G4double ComputeDEDXPerVolume(const G4Material*,
0102                 const G4ParticleDefinition*,
0103                 G4double kineticEnergy,
0104                 G4double) override;
0105 
0106   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
0107              const G4MaterialCutsCouple*,
0108              const G4DynamicParticle*,
0109              G4double tmin,
0110              G4double maxEnergy) override;
0111 
0112   // add correction to energy loss and compute non-ionizing energy loss
0113   void CorrectionsAlongStep(const G4MaterialCutsCouple*,
0114                 const G4DynamicParticle*,
0115                 const G4double& length,
0116                 G4double& eloss) override;
0117 
0118   // hide assignment operator
0119   G4ICRU73QOModel & operator=(const  G4ICRU73QOModel &right) = delete;
0120   G4ICRU73QOModel(const  G4ICRU73QOModel&) = delete;
0121  
0122 protected:
0123 
0124   G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
0125                   G4double kinEnergy) final;
0126 
0127 private:
0128 
0129   inline void SetParticle(const G4ParticleDefinition* p);
0130   inline void SetLowestKinEnergy(G4double val);
0131 
0132   G4double DEDX(const G4Material* material, G4double kineticEnergy);
0133 
0134   G4double DEDXPerElement(G4int Z, G4double kineticEnergy);
0135 
0136   // get number of shell, energy and oscillator strengths for material
0137   G4int GetNumberOfShells(G4int Z) const;
0138 
0139   G4double GetShellEnergy(G4int Z, G4int nbOfTheShell) const; 
0140   G4double GetOscillatorEnergy(G4int Z, G4int nbOfTheShell) const; 
0141   G4double GetShellStrength(G4int Z, G4int nbOfTheShell) const;
0142 
0143   // calculate stopping number for L's term
0144   G4double GetL0(G4double normEnergy) const;
0145   // terms in Z^2
0146   G4double GetL1(G4double normEnergy) const;
0147   // terms in Z^3
0148   G4double GetL2(G4double normEnergy) const;
0149   // terms in Z^4
0150   
0151   const G4ParticleDefinition* particle;
0152   G4ParticleDefinition*       theElectron;   
0153   G4ParticleChangeForLoss*    fParticleChange;
0154   G4DensityEffectData*        denEffData;
0155 
0156   G4double mass;
0157   G4double charge;
0158   G4double chargeSquare;
0159   G4double massRate;
0160   G4double ratio;
0161   G4double lowestKinEnergy;
0162 
0163   G4bool   isInitialised;
0164 
0165   // Z of element at now avaliable for the model
0166   static const G4int NQOELEM  = 26;
0167   static const G4int NQODATA  = 130;
0168   static const G4int ZElementAvailable[NQOELEM];
0169   
0170   // number, energy and oscillator strengths
0171   // for an harmonic oscillator model of material
0172   static const G4int startElemIndex[NQOELEM];
0173   static const G4int nbofShellsForElement[NQOELEM];
0174   static const G4double ShellEnergy[NQODATA];
0175   static const G4double SubShellOccupation[NQODATA];  // Z * ShellStrength
0176 
0177   G4int indexZ[100];
0178 
0179   //  variable for calculation of stopping number of L's term
0180   static const G4double L0[67][2];
0181   static const G4double L1[22][2];
0182   static const G4double L2[14][2];
0183   
0184   G4int sizeL0;
0185   G4int sizeL1;
0186   G4int sizeL2;
0187 
0188   static const G4double factorBethe[99];
0189   
0190 };
0191 
0192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0193 
0194 inline void G4ICRU73QOModel::SetParticle(const G4ParticleDefinition* p)
0195 {
0196   particle = p;
0197   mass = particle->GetPDGMass();
0198   charge = particle->GetPDGCharge()/CLHEP::eplus;
0199   chargeSquare = charge*charge;
0200   massRate     = mass/CLHEP::proton_mass_c2;
0201   ratio = CLHEP::electron_mass_c2/mass;
0202 }
0203 
0204 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0205 
0206 inline void G4ICRU73QOModel::SetLowestKinEnergy(G4double val)
0207 {
0208   lowestKinEnergy = val;
0209 }
0210 
0211 #endif