Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/lowenergy/include/G4QAOLowEnergyLoss.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 // ------------------------------------------------------------
 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 
 53 class G4QAOLowEnergyLoss : public G4VLowEnergyModel
 54 {
 55 public:   
 56   explicit G4QAOLowEnergyLoss(const G4String& name);   
 57   ~G4QAOLowEnergyLoss();
 58     
 59   G4double HighEnergyLimit(const G4ParticleDefinition* aParticle,
 60                            const G4Material* material) const override;
 61   // returns the higher limit for model validity
 62 
 63   G4double LowEnergyLimit(const G4ParticleDefinition* aParticle,
 64                           const G4Material* material) const override;
 65   // returns the lower limit for model validity
 66  
 67   G4double HighEnergyLimit(const G4ParticleDefinition* aParticle) const override;
 68   // returns the higher limit for model validity
 69 
 70   G4double LowEnergyLimit(const G4ParticleDefinition* aParticle) const override;
 71   // returns the lower limit for model validity
 72  
 73   G4bool IsInCharge(const G4DynamicParticle* particle,
 74         const G4Material* material) const override;
 75   // returns true if the model is applicable at that energy for
 76   // that particle for that material
 77 
 78   G4bool IsInCharge(const G4ParticleDefinition* aParticle,
 79         const G4Material* material) const override;
 80   // returns true if the model is applicable at that energy for
 81   // that particle for that material
 82   
 83   G4double TheValue(const G4DynamicParticle* particle,
 84                const G4Material* material) override;
 85   // returns the energy loss via the quantal harmonic oscillator model 
 86   
 87   G4double TheValue(const G4ParticleDefinition* aParticle,
 88                     const G4Material* material,
 89                                 G4double kineticEnergy) override;
 90   // returns the energy loss via the quantal harmonic oscillator model 
 91 
 92 private:  
 93   G4double EnergyLoss(const G4Material* material,
 94                             G4double kineticEnergy,
 95                             G4double zParticle) const;
 96   // returns the energy loss via the quantal harmonic oscillator model 
 97    
 98   // get number of shell, energy and oscillator strengths for material
 99   G4int GetNumberOfShell(const G4Material* material) const;
100 
101   G4double GetShellEnergy(const G4Material* material,G4int nbOfTheShell) const; 
102   G4double GetOscillatorEnergy(const G4Material* material,G4int nbOfTheShell) const; 
103   G4double GetShellStrength(const G4Material* material,G4int nbOfTheShell) const;
104   G4double GetOccupationNumber(G4int Z, G4int ShellNb) const;
105 
106   // calculate stopping number for L's term
107   G4double GetL0(G4double normEnergy) const;
108   // terms in Z^2
109   G4double GetL1(G4double normEnergy) const;
110   // terms in Z^3
111   G4double GetL2(G4double normEnergy) const;
112   // terms in Z^4
113     
114   // number, energy and oscillator strengths
115   // for an harmonic oscillator model of material
116   static const G4int nbofShellForMaterial[6];
117   static const G4double alShellEnergy[3];
118   static const G4double alShellStrength[3];
119   static const G4double siShellEnergy[3];
120   static const G4double siShellStrength[3];
121   static const G4double cuShellEnergy[4];
122   static const G4double cuShellStrength[4];
123   static const G4double taShellEnergy[6];
124   static const G4double taShellStrength[6];
125   static const G4double auShellEnergy[6];
126   static const G4double auShellStrength[6];
127   static const G4double ptShellEnergy[6];
128   static const G4double ptShellStrength[6];
129   //  variable for calculation of stopping number of L's term
130   static const G4double L0[67][2];
131   static const G4double L1[22][2];
132   static const G4double L2[14][2];
133   static const G4int nbOfElectronPerSubShell[1540];
134   static const G4int fNumberOfShells[101];
135 
136   // Z of element at now avaliable for the model
137   static const G4int materialAvailable[6];
138   
139   G4int numberOfMaterials;
140   G4int sizeL0;
141   G4int sizeL1;
142   G4int sizeL2;
143   
144 }; 
145 
146 #endif
147