Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/lowenergy/include/G4PenelopePhotoElectricModel.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 // Author: Luciano Pandola
 28 //
 29 // History:
 30 // -----------
 31 // 08 Jan 2010   L. Pandola   1st implementation. 
 32 // 25 May 2011   L. Pandola   Renamed (make v2008 as default Penelope)
 33 // 10 Jun 2011   L. Pandola   Migrated to the new AtomDeexcitation interface
 34 // 18 Sep 2013   L. Pandola   Migrated to the MT paradigm
 35 //
 36 // -------------------------------------------------------------------
 37 //
 38 // Class description:
 39 // Low Energy Electromagnetic Physics, Photo-electric effect
 40 // with Penelope Model, version 2008
 41 // -------------------------------------------------------------------
 42 
 43 #ifndef G4PENELOPEPHOTOELECTRICMODEL_HH
 44 #define G4PENELOPEPHOTOELECTRICMODEL_HH 1
 45 
 46 #include "globals.hh"
 47 #include "G4VEmModel.hh"
 48 #include "G4DataVector.hh"
 49 #include "G4ParticleChangeForGamma.hh"
 50 #include "G4AtomicTransitionManager.hh"
 51 #include "G4VAtomDeexcitation.hh"
 52 
 53 class G4ParticleDefinition;
 54 class G4DynamicParticle;
 55 class G4MaterialCutsCouple;
 56 class G4Material;
 57 
 58 class G4PenelopePhotoElectricModel : public G4VEmModel 
 59 {
 60 public:
 61   explicit G4PenelopePhotoElectricModel(const G4ParticleDefinition* p=nullptr,
 62           const G4String& processName ="PenPhotoElec");
 63   virtual ~G4PenelopePhotoElectricModel();
 64 
 65   void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
 66   void InitialiseLocal(const G4ParticleDefinition*,
 67                                G4VEmModel *masterModel) override;
 68 
 69   G4double ComputeCrossSectionPerAtom(
 70               const G4ParticleDefinition*,
 71               G4double kinEnergy,
 72               G4double Z,
 73               G4double A=0,
 74               G4double cut=0,
 75               G4double emax=DBL_MAX) override;
 76 
 77   void SampleSecondaries(std::vector<G4DynamicParticle*>*,
 78        const G4MaterialCutsCouple*,
 79        const G4DynamicParticle*,
 80        G4double tmin,
 81        G4double maxEnergy) override;
 82 
 83   void SetVerbosityLevel(G4int lev){fVerboseLevel = lev;};
 84   G4int GetVerbosityLevel(){return fVerboseLevel;};
 85 
 86   //testing purposes
 87   std::size_t GetNumberOfShellXS(G4int);
 88   G4double GetShellCrossSection(G4int Z,std::size_t shellID,G4double energy);
 89 
 90   G4PenelopePhotoElectricModel & operator=(const G4PenelopePhotoElectricModel &right) = delete;
 91   G4PenelopePhotoElectricModel(const G4PenelopePhotoElectricModel&) = delete;
 92 
 93 protected:
 94   G4ParticleChangeForGamma* fParticleChange;
 95   const G4ParticleDefinition* fParticle;
 96 
 97 private:
 98   void SetParticle(const G4ParticleDefinition*);
 99   G4double SampleElectronDirection(G4double energy);
100   void ReadDataFile(G4int Z);
101   std::size_t SelectRandomShell(G4int Z,G4double energy);
102   G4String WriteTargetShell(std::size_t shellID);
103 
104   //For each Z, the PhysicsTable contains nShell+1 physics vectors
105   //with log(E) vs. log(XS)
106   //Element [0] of the table is the total XS, element [iS] is the 
107   //partial cross section for shell iS-1
108   static const G4int fMaxZ =99;
109   static G4PhysicsTable* fLogAtomicShellXS[fMaxZ+1];
110 
111   G4VAtomDeexcitation*             fAtomDeexcitation;
112   const G4AtomicTransitionManager* fTransitionManager;
113 
114   //Intrinsic energy limits of the model: cannot be extended by the parent process
115   G4double fIntrinsicLowEnergyLimit;
116   G4double fIntrinsicHighEnergyLimit;
117   G4int fVerboseLevel;
118   G4bool fIsInitialised; 
119   //Used only for G4EmCalculator and Unit Tests
120   G4bool fLocalTable;
121 };
122 
123 #endif
124 
125