CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CMSmplIonisationWithDeltaModel.h
Go to the documentation of this file.
1 //
2 // -------------------------------------------------------------------
3 //
4 // GEANT4 Class header file
5 //
6 //
7 // File name: CMSmplIonisationWithDeltaModel
8 //
9 // Author: Vladimir Ivanchenko copied from Geant4 10.5p01
10 //
11 // Creation date: 02.03.2019
12 //
13 // Class Description:
14 //
15 // Implementation of model of energy loss of the magnetic monopole
16 
17 // -------------------------------------------------------------------
18 //
19 
20 #ifndef CMSmplIonisationWithDeltaModel_h
21 #define CMSmplIonisationWithDeltaModel_h 1
22 
23 #include "G4VEmModel.hh"
24 #include "G4VEmFluctuationModel.hh"
25 #include <vector>
26 
27 class G4ParticleChangeForLoss;
28 
29 class CMSmplIonisationWithDeltaModel : public G4VEmModel, public G4VEmFluctuationModel
30 {
31 
32 public:
33 
34  explicit CMSmplIonisationWithDeltaModel(G4double mCharge,
35  const G4String& nam = "mplIonisationWithDelta");
36 
38 
39  void Initialise(const G4ParticleDefinition*,
40  const G4DataVector&) override;
41 
42  G4double ComputeDEDXPerVolume(const G4Material*,
43  const G4ParticleDefinition*,
44  G4double kineticEnergy,
45  G4double cutEnergy) override;
46 
47  virtual G4double ComputeCrossSectionPerElectron(
48  const G4ParticleDefinition*,
49  G4double kineticEnergy,
50  G4double cutEnergy,
51  G4double maxEnergy);
52 
54  const G4ParticleDefinition*,
55  G4double kineticEnergy,
56  G4double Z, G4double A,
57  G4double cutEnergy,
58  G4double maxEnergy) override;
59 
60  void SampleSecondaries(std::vector<G4DynamicParticle*>*,
61  const G4MaterialCutsCouple*,
62  const G4DynamicParticle*,
63  G4double tmin,
64  G4double maxEnergy) override;
65 
66 
67  G4double SampleFluctuations(const G4MaterialCutsCouple*,
68  const G4DynamicParticle*,
69  G4double tmax,
70  G4double length,
71  G4double meanLoss) override;
72 
73  G4double Dispersion(const G4Material*,
74  const G4DynamicParticle*,
75  G4double tmax,
76  G4double length) override;
77 
78  G4double MinEnergyCut(const G4ParticleDefinition*,
79  const G4MaterialCutsCouple* couple) override;
80 
81  void SetParticle(const G4ParticleDefinition* p);
82 
83 protected:
84 
85  G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
86  G4double kinEnergy) override;
87 
88 private:
89 
90  G4double ComputeDEDXAhlen(const G4Material* material, G4double bg2,
91  G4double cut);
92 
93  const G4ParticleDefinition* monopole;
94  G4ParticleDefinition* theElectron;
95  G4ParticleChangeForLoss* fParticleChange;
96 
97  G4double mass;
98  G4double magCharge;
99  G4double twoln10;
100  G4double betalow;
101  G4double betalim;
102  G4double beta2lim;
103  G4double bg2lim;
104  G4double chargeSquare;
105  G4double dedxlim;
106  G4int nmpl;
108 
109  static std::vector<G4double>* dedx0;
110 };
111 
112 #endif
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
virtual G4double ComputeCrossSectionPerElectron(const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
const double Z[kNumberCalorimeter]
G4double ComputeDEDXAhlen(const G4Material *material, G4double bg2, G4double cut)
CMSmplIonisationWithDeltaModel(G4double mCharge, const G4String &nam="mplIonisationWithDelta")
static std::vector< G4double > * dedx0
G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy) override
void SetParticle(const G4ParticleDefinition *p)
G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmax, G4double length, G4double meanLoss) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
static const double tmax[3]
G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kinEnergy) override
G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double tmax, G4double length) override
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy) override
G4double MinEnergyCut(const G4ParticleDefinition *, const G4MaterialCutsCouple *couple) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override