Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DynamicParticleMSC.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// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4DynamicParticleMSC
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 23.08.2024
37//
38// Modifications:
39//
40//
41// Class Description:
42//
43// This class manages the multiple scattering process for heavy charged
44// particles without any use of G4ParticleDefinition.
45//
46// -------------------------------------------------------------------
47//
48
49#ifndef G4DynamicParticleMSC_h
50#define G4DynamicParticleMSC_h 1
51
54
55class G4Material;
57
59{
60public:
61
63
64 ~G4DynamicParticleMSC() override;
65
66 // Step limit from AlongStep
68 const G4Track&,
69 G4double previousStepSize,
70 G4double currentMinimumStep,
71 G4double& currentSafety,
72 G4GPILSelection* selection) override;
73
74 // Step limit from cross section
76 const G4Track& track,
77 G4double previousStepSize,
78 G4ForceCondition* condition) override;
79
80 // AlongStep computations
81 G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&) override;
82
83 // This method is not used for tracking, it returns mean free path value
85
86 // This method is not used for tracking, it returns step limit
87 G4double GetContinuousStepLimit(const G4Track& track, G4double previousStepSize,
88 G4double currentMinimalStep, G4double& currentSafety) override;
89
90 // print description in html
91 void ProcessDescription(std::ostream&) const override;
92
93 // hide assignment operator
96
97private:
98
99 // all parameters are dynamic
100 void PreStepInitialisation(const G4Track&);
101
102 G4LossTableManager* lManager;
103 const G4Material* fMaterial{nullptr};
104
105 G4double fMass{0.0};
106 G4double fCharge{0.0};
107 G4double fEkinPreStep{0.0};
108 G4double fBeta{0.0};
109 G4double fZeff{0.0};
110
111 G4ThreeVector fNewDir{G4ThreeVector(0.0, 0.0, 0.0)};
112 G4ParticleChangeForMSC fParticleChange;
113};
114
115//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
116
117#endif
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4GPILSelection
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
void ProcessDescription(std::ostream &) const override
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
G4double GetMeanFreePath(const G4Track &, G4double, G4ForceCondition *condition) override
G4DynamicParticleMSC(const G4DynamicParticleMSC &)=delete
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4DynamicParticleMSC & operator=(const G4DynamicParticleMSC &right)=delete
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety) override
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
G4VContinuousDiscreteProcess(const G4String &, G4ProcessType aType=fNotDefined)