Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4PenelopeIonisationModel.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//
29// History:
30// -----------
31// 30 Mar 2010 L. Pandola 1st implementation.
32// 25 May 2011 L. Pandola Renamed (make v2008 as default Penelope)
33// 09 Mar 2012 L. Pandola Moved the management and calculation of
34// cross sections to a separate class
35// 07 Oct 2013 L. Pandola Migration to MT
36// 23 Jun 2015 L. Pandola Added private member to store the PIXE flag
37//
38// -------------------------------------------------------------------
39//
40// Class description:
41// Low Energy Electromagnetic Physics, e+ and e- ionisation
42// with Penelope Model, version 2008
43// -------------------------------------------------------------------
44
45#ifndef G4PENELOPEIONISATIONMODEL_HH
46#define G4PENELOPEIONISATIONMODEL_HH 1
47
48#include "globals.hh"
49#include "G4VEmModel.hh"
50#include "G4DataVector.hh"
53
59class G4Material;
64
66{
67public:
68 explicit G4PenelopeIonisationModel(const G4ParticleDefinition* p=nullptr,
69 const G4String& processName ="PenIoni");
71
72 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
74 G4VEmModel*) override;
75
76 //*This is a dummy method. Never inkoved by the tracking, it just issues
77 //*a warning if one tries to get Cross Sections per Atom via the
78 //*G4EmCalculator.
84 G4double) override;
85
88 theParticle,
89 G4double kineticEnergy,
90 G4double cutEnergy,
91 G4double maxEnergy = DBL_MAX) override;
92
93 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
95 const G4DynamicParticle*,
96 G4double tmin,
97 G4double maxEnergy) override;
98
101 G4double kineticEnergy,
102 G4double cutEnergy) override;
103
104 // Min cut in kinetic energy allowed by the model
106 const G4MaterialCutsCouple*) override;
107
108 void SetVerbosityLevel(G4int lev){fVerboseLevel = lev;};
109 G4int GetVerbosityLevel(){return fVerboseLevel;};
110
113
114protected:
117
118private:
119 void SetParticle(const G4ParticleDefinition*);
120 void SampleFinalStateElectron(const G4Material*,
121 G4double cutEnergy,
122 G4double kineticEnergy);
123 void SampleFinalStatePositron(const G4Material*,
124 G4double cutEnergy,
125 G4double kineticEnergy);
126
127 G4PenelopeOscillatorManager* fOscManager;
128 G4PenelopeIonisationXSHandler* fCrossSectionHandler;
129 G4VAtomDeexcitation* fAtomDeexcitation;
130
131 G4double fKineticEnergy1;
132 G4double fCosThetaPrimary;
133 G4double fEnergySecondary;
134 G4double fCosThetaSecondary;
135
136 //Intrinsic energy limits of the model: cannot be extended by the parent process
137 G4double fIntrinsicLowEnergyLimit;
138 G4double fIntrinsicHighEnergyLimit;
139
140 G4int fVerboseLevel;
141 G4int fTargetOscillator;
142 size_t fNBins;
143 G4bool fIsInitialised;
144 G4bool fPIXEflag;
145 //Used only for G4EmCalculator and Unit Tests
146 G4bool fLocalTable;
147};
148
149#endif
150
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4PenelopeIonisationModel & operator=(const G4PenelopeIonisationModel &right)=delete
G4ParticleChangeForLoss * fParticleChange
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double, G4double, G4double, G4double, G4double) override
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *theParticle, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy=DBL_MAX) override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *) override
const G4ParticleDefinition * fParticle
G4PenelopeIonisationModel(const G4PenelopeIonisationModel &)=delete
#define DBL_MAX
Definition: templates.hh:62