Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4mplIonisationWithDeltaModel.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// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class header file
31//
32//
33// File name: G4mplIonisationWithDeltaModel
34//
35// Author: Vladimir Ivanchenko
36//
37// Creation date: 06.09.2005
38//
39// Modifications:
40// 12.08.2007 ComputeDEDXAhlen function added (M. Vladymyrov)
41//
42// Class Description:
43//
44// Implementation of model of energy loss of the magnetic monopole
45
46// -------------------------------------------------------------------
47//
48
49#ifndef G4mplIonisationWithDeltaModel_h
50#define G4mplIonisationWithDeltaModel_h 1
51
52#include "G4VEmModel.hh"
54
56
58{
59
60public:
61
63 const G4String& nam = "mplIonisationWithDelta");
64
66
67 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
68
71 G4double kineticEnergy,
72 G4double cutEnergy);
73
76 G4double kineticEnergy,
77 G4double cutEnergy,
78 G4double maxEnergy);
79
82 G4double kineticEnergy,
83 G4double Z, G4double A,
84 G4double cutEnergy,
85 G4double maxEnergy);
86
87 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
89 const G4DynamicParticle*,
90 G4double tmin,
91 G4double maxEnergy);
92
93
95 const G4DynamicParticle*,
96 G4double& tmax,
97 G4double& length,
98 G4double& meanLoss);
99
100 virtual G4double Dispersion(const G4Material*,
101 const G4DynamicParticle*,
102 G4double& tmax,
103 G4double& length);
104
105 void SetParticle(const G4ParticleDefinition* p);
106
107protected:
108
110 G4double kinEnergy);
111
112private:
113
114 G4double ComputeDEDXAhlen(const G4Material* material, G4double bg2, G4double cut);
115
116 // hide assignment operator
119
120 const G4ParticleDefinition* monopole;
121 G4ParticleDefinition* theElectron;
122 G4ParticleChangeForLoss* fParticleChange;
123
124 G4double mass;
125 G4double magCharge;
126 G4double twoln10;
127 G4double betalow;
128 G4double betalim;
129 G4double beta2lim;
130 G4double bg2lim;
131 G4double chargeSquare;
132 G4double dedxlim;
133 G4int nmpl;
134 G4double pi_hbarc2_over_mc2;
135};
136
137#endif
138
139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy)
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double &tmax, G4double &length)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
void SetParticle(const G4ParticleDefinition *p)
virtual G4double SampleFluctuations(const G4Material *, const G4DynamicParticle *, G4double &tmax, G4double &length, G4double &meanLoss)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)