Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VITDiscreteProcess.cc
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// --------------------------------------------------------------
30// GEANT 4 class implementation file
31//
32// History: first implementation, based on object model of
33// 2nd December 1995, G.Cosmo
34// --------------------------------------------------------------
35// New Physics scheme 8 Jan. 1997 H.Kurahige
36// ------------------------------------------------------------
37
39#include "G4SystemOfUnits.hh"
40
41#include "G4Step.hh"
42#include "G4Track.hh"
43#include "G4MaterialTable.hh"
44#include "G4VParticleChange.hh"
45
47 G4VITProcess("No Name Discrete Process")
48{
49 G4Exception("G4VDiscreteProcess::G4VDiscreteProcess()",
50 "ProcMan102",
52 "Default constructor is called");
53}
54
55//------------------------------------------------------------------------------
56
58 G4ProcessType aType) :
59 G4VITProcess(aName, aType)
60{
61 enableAtRestDoIt = false;
62 enableAlongStepDoIt = false;
63}
64
65//------------------------------------------------------------------------------
66
68= default;
69
70//------------------------------------------------------------------------------
71
72G4VITDiscreteProcess::G4VITDiscreteProcess(G4VITDiscreteProcess& right)
73
74= default;
75
76//------------------------------------------------------------------------------
77
81 G4double previousStepSize,
83{
84 if((previousStepSize < 0.0)
85 || (fpState->theNumberOfInteractionLengthLeft <= 0.0))
86 {
87 // beginning of tracking (or just after DoIt of this process)
89 } else if(previousStepSize > 0.0)
90 {
91 // subtract NumberOfInteractionLengthLeft
93 } else
94 {
95 // zero step
96 // DO NOTHING
97 }
98
99 // condition is set to "Not Forced"
101
102 // get mean free path
103 fpState->currentInteractionLength = GetMeanFreePath(track,
104 previousStepSize,
105 condition);
106
107 G4double value;
108 if(fpState->currentInteractionLength < DBL_MAX)
109 {
110 value = fpState->theNumberOfInteractionLengthLeft
111 * fpState->currentInteractionLength;
112 } else
113 {
114 value = DBL_MAX;
115 }
116#ifdef G4VERBOSE
117 if(verboseLevel > 1)
118 {
119 G4cout << "G4VDiscreteProcess::PostStepGetPhysicalInteractionLength ";
120 G4cout << "[ " << GetProcessName() << "]" << G4endl;
121 track.GetDynamicParticle()->DumpInfo();
122 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
123 G4cout << "InteractionLength= " << value / cm << "[cm] " << G4endl;
124 }
125#endif
126 return value;
127}
128
129//------------------------------------------------------------------------------
130
132 const G4Step&)
133{
134// clear NumberOfInteractionLengthLeft
136
137 return pParticleChange;
138}
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4ForceCondition
@ NotForced
G4ProcessType
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void DumpInfo(G4int mode=0) const
const G4String & GetName() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
G4VITDiscreteProcess(const G4String &, G4ProcessType aType=fNotDefined)
~G4VITDiscreteProcess() override
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
G4shared_ptr< G4ProcessState > fpState
virtual void ClearNumberOfInteractionLengthLeft()
void ResetNumberOfInteractionLengthLeft() override
G4bool enableAtRestDoIt
G4int verboseLevel
G4bool enableAlongStepDoIt
G4VParticleChange * pParticleChange
const G4String & GetProcessName() const
#define DBL_MAX
Definition templates.hh:62