Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/dna/models/include/G4DNAPTBAugerModel.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 // Authors: S. Meylan and C. Villagrasa (IRSN, France)
 27 // Models come from
 28 // M. Bug et al, Rad. Phys and Chem. 130, 459-479 (2017)
 29 // 
 30 //
 31 // -------------------------------------------------------------------
 32 //
 33 // Geant4 Header G4DNAPTBAugerModel
 34 //  
 35 // -------------------------------------------------------------------
 36 //
 37 // Class description:
 38 // Implementation of atomic deexcitation 
 39 //
 40 // -------------------------------------------------------------------
 41 
 42 #ifndef G4DNAPTBAugerModel_h
 43 #define G4DNAPTBAugerModel_h 1
 44 
 45 #include "G4VAtomDeexcitation.hh"
 46 #include "G4AtomicShell.hh"
 47 #include "globals.hh"
 48 #include "G4DynamicParticle.hh"
 49 #include <vector>
 50 
 51 class G4AtomicTransitionManager;
 52 class G4VhShellCrossSection;
 53 class G4EmCorrections;
 54 class G4Material;
 55 
 56 /*!
 57  * \brief The G4DNAPTBAugerModel class
 58  * Implement the PTB Auger model
 59  */
 60 class G4DNAPTBAugerModel
 61 {  
 62 public: 
 63   
 64    /*!
 65    * \brief G4DNAPTBAugerModel
 66    * Constructor
 67    * \param modelName
 68    */
 69   G4DNAPTBAugerModel(const G4String &modelName);
 70 
 71   /*!
 72    * \brief ~G4DNAPTBAugerModel
 73    * Destructor
 74    */
 75   virtual ~G4DNAPTBAugerModel();
 76    
 77   G4DNAPTBAugerModel(G4DNAPTBAugerModel &) = delete;  // prevent copy-construction
 78   G4DNAPTBAugerModel & operator=(const G4DNAPTBAugerModel &right) = delete;  // prevent assignement
 79  
 80   /*!
 81    * \brief Initialise
 82    * Set the verbose value
 83    */
 84   virtual void Initialise();
 85 
 86   /*!
 87    * \brief SetCutForAugerElectrons
 88    * Set the cut for the auger electrons production
 89    * \param cut
 90    */
 91   void SetCutForAugerElectrons(G4double cut);
 92 
 93   /*!
 94    * \brief ComputeAugerEffect
 95    * Main method to be called by the ionisation model.
 96    * \param fvect
 97    * \param materialNameIni
 98    * \param bindingEnergy
 99    */
100   void ComputeAugerEffect(std::vector<G4DynamicParticle *> *fvect, const G4String& materialNameIni, G4double bindingEnergy);
101 
102 private:
103  
104   const G4String modelName; ///< name of the auger model
105 
106   G4int verboseLevel;
107   G4double minElectronEnergy;
108 
109   /*!
110    * \brief GenerateAugerWithRandomDirection
111    * Generates the auger particle
112    * \param fvect
113    * \param kineticEnergy
114    */
115   void GenerateAugerWithRandomDirection(std::vector<G4DynamicParticle*>* fvect, G4double kineticEnergy);
116 
117   /*!
118    * \brief CalculAugerEnergyFor
119    * \param atomId
120    * \return the auger particle energy
121    */
122   G4double CalculAugerEnergyFor(G4int atomId);
123 
124   /*!
125    * \brief DetermineIonisedAtom
126    * \param atomId
127    * \param materialName
128    * \param bindingEnergy
129    * \return the id of the chosen ionised atom
130    */
131   G4int DetermineIonisedAtom(G4int atomId, const G4String &materialName, G4double bindingEnergy);
132 
133 
134 };
135 
136 #endif
137 
138 
139 
140 
141