ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4QAOLowEnergyLoss.hh
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4QAOLowEnergyLoss.hh
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 // ------------------------------------------------------------
28 // GEANT 4 class header file
29 //
30 // History: New Implementation
31 //
32 // ---------- G4QAOLowEnergyLoss physics process -------
33 // by Stephane Chauvie, 21 May 2000
34 //
35 // Modified:
36 // 16/09/2000 S. Chauvie Oscillator for all materials
37 // 23/05/2000 MGP Made compliant to design
38 // 01/06/2001 V.Ivanchenko replace names by Z
39 //
40 // Class description:
41 // Quantal Harmonic Oscillator Model for energy loss of low energy antiprotons
42 // Further documentation available from http://www.ge.infn.it/geant4/lowE
43 
44 // ------------------------------------------------------------
45 
46 
47 #ifndef G4QAOLowEnergyLoss_hh
48 #define G4QAOLowEnergyLoss_hh 1
49 
50 #include "G4VLowEnergyModel.hh"
51 #include "globals.hh"
52 
54 {
55 public:
56 
58 
60 
62  const G4Material* material) const;
63  // returns the higher limit for model validity
64 
66  const G4Material* material) const;
67  // returns the lower limit for model validity
68 
69  G4double HighEnergyLimit(const G4ParticleDefinition* aParticle) const;
70  // returns the higher limit for model validity
71 
72  G4double LowEnergyLimit(const G4ParticleDefinition* aParticle) const;
73  // returns the lower limit for model validity
74 
76  const G4Material* material) const;
77  // returns true if the model is applicable at that energy for
78  // that particle for that material
79 
80  G4bool IsInCharge(const G4ParticleDefinition* aParticle,
81  const G4Material* material) const;
82  // returns true if the model is applicable at that energy for
83  // that particle for that material
84 
85  G4double TheValue(const G4DynamicParticle* particle,
86  const G4Material* material);
87  // returns the energy loss via the quantal harmonic oscillator model
88 
89  G4double TheValue(const G4ParticleDefinition* aParticle,
90  const G4Material* material,
91  G4double kineticEnergy);
92  // returns the energy loss via the quantal harmonic oscillator model
93 
94 private:
95 
96  G4double EnergyLoss(const G4Material* material,
97  G4double kineticEnergy,
98  G4double zParticle) const;
99  // returns the energy loss via the quantal harmonic oscillator model
100 
101  // get number of shell, energy and oscillator strenghts for material
102  G4int GetNumberOfShell(const G4Material* material) const;
103 
104  G4double GetShellEnergy(const G4Material* material,G4int nbOfTheShell) const;
105  G4double GetOscillatorEnergy(const G4Material* material,G4int nbOfTheShell) const;
106  G4double GetShellStrength(const G4Material* material,G4int nbOfTheShell) const;
107  G4double GetOccupationNumber(G4int Z, G4int ShellNb) const;
108 
109  // calculate stopping number for L's term
110  G4double GetL0(G4double normEnergy) const;
111  // terms in Z^2
112  G4double GetL1(G4double normEnergy) const;
113  // terms in Z^3
114  G4double GetL2(G4double normEnergy) const;
115  // terms in Z^4
116 
117  // Z of element at now avaliable for the model
118  static const G4int materialAvailable[6];
119 
120  // number, energy and oscillator strenghts
121  // for an harmonic oscillator model of material
122  static const G4int nbofShellForMaterial[6];
123  static const G4double alShellEnergy[3];
124  static const G4double alShellStrength[3];
125  static const G4double siShellEnergy[3];
126  static const G4double siShellStrength[3];
127  static const G4double cuShellEnergy[4];
128  static const G4double cuShellStrength[4];
129  static const G4double taShellEnergy[6];
130  static const G4double taShellStrength[6];
131  static const G4double auShellEnergy[6];
132  static const G4double auShellStrength[6];
133  static const G4double ptShellEnergy[6];
134  static const G4double ptShellStrength[6];
135 
137 
138  // variable for calculation of stopping number of L's term
139  static const G4double L0[67][2];
140  static const G4double L1[22][2];
141  static const G4double L2[14][2];
142  static const G4int nbOfElectronPerSubShell[1540];
143  static const G4int fNumberOfShells[101];
144 
148 
149 };
150 
151 #endif