Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LivermoreIonisationModel.hh
Go to the documentation of this file.
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// on base of G4LowEnergyIonisation developed by A.Forti and V.Ivanchenko
29//
30// History:
31// -----------
32// 12 Jan 2009 L. Pandola 1st implementation. Migration from EM process
33// to EM model. Physics is unchanged.
34// 23 Oct 2009 L. Pandola remove un-necessary methods to manage atomic
35// deexcitation (done by G4VEmModel)
36// 01 Jun 2011 V Ivanchenko general cleanup - all old deexcitation code removed
37// 04 Jul 2011 L Pandola removed unused private member
38//
39// -------------------------------------------------------------------
40//
41// Class description:
42// Low Energy Electromagnetic Physics, e- ionisation
43// with Livermore Model
44// -------------------------------------------------------------------
45
46#ifndef G4LIVERMOREIONISATIONMODEL_HH
47#define G4LIVERMOREIONISATIONMODEL_HH 1
48
49#include "G4VEmModel.hh"
50#include "globals.hh"
51
56
58{
59public:
61 const G4String& processName = "LowEnergyIoni");
63
64 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
65
66
68 G4double kinEnergy,
69 G4double Z,
70 G4double A=0,
71 G4double cut=0,
72 G4double emax=DBL_MAX) override;
73
74 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
76 const G4DynamicParticle*,
77 G4double tmin,
78 G4double maxEnergy) override;
79
82 G4double kineticEnergy,
83 G4double cutEnergy) override;
84
85 void SetVerboseLevel(G4int vl) {verboseLevel = vl;};
86 G4int GetVerboseLevel(){return verboseLevel;};
87
90
91protected:
93
94private:
95 G4eIonisationCrossSectionHandler* crossSectionHandler;
96 G4VEnergySpectrum* energySpectrum;
97 G4AtomicTransitionManager* transitionManager;
98
99 //Intrinsic energy limits of the model: cannot be extended by the parent process
100 G4double fIntrinsicLowEnergyLimit;
101 G4double fIntrinsicHighEnergyLimit;
102 G4int verboseLevel;
103 G4bool isInitialised;
104};
105
106#endif
107
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4ParticleChangeForLoss * fParticleChange
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4LivermoreIonisationModel(const G4LivermoreIonisationModel &)=delete
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
G4LivermoreIonisationModel & operator=(const G4LivermoreIonisationModel &right)=delete
#define DBL_MAX
Definition: templates.hh:62