Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VITRestDiscreteProcess.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// Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
28//
29// History:
30// -----------
31// 10 Oct 2011 M.Karamitros created
32//
33// -------------------------------------------------------------------
34
36
38 G4ProcessType aType)
39 : G4VITProcess(aName, aType)
40{
41 enableAlongStepDoIt = false;
42}
43
45
47 G4double previousStepSize,
49{
50 if ((previousStepSize < 0.0) || (fpState->theNumberOfInteractionLengthLeft
51 <= 0.0))
52 {
53 // beggining of tracking (or just after DoIt of this process)
55 }
56 else if (previousStepSize > 0.0)
57 {
58 // subtract NumberOfInteractionLengthLeft
60 }
61 else
62 {
63 // zero step
64 // DO NOTHING
65 }
66
67 // condition is set to "Not Forced"
69
70 // get mean free path
71 fpState->currentInteractionLength = GetMeanFreePath(track,
72 previousStepSize,
73 condition);
74
75 G4double value;
76 if (fpState->currentInteractionLength < DBL_MAX)
77 {
78 value = fpState->theNumberOfInteractionLengthLeft * (fpState->currentInteractionLength);
79 }
80 else
81 {
82 value = DBL_MAX;
83 }
84#ifdef G4VERBOSE
85 if (verboseLevel > 1)
86 {
87 G4cout << "G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength ";
88 G4cout << "[ " << GetProcessName() << "]" << G4endl;
90 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
91 G4cout << "InteractionLength= " << value / CLHEP::cm << "[cm] " << G4endl;
92 }
93#endif
94 return value;
95}
96
98 const G4Step&)
99{
100// reset NumberOfInteractionLengthLeft
102
103 return pParticleChange;
104}
105
108{
109 // beggining of tracking
111
112 // condition is set to "Not Forced"
114
115 // get mean life time
116 fpState->currentInteractionLength = GetMeanLifeTime(track, condition);
117
118#ifdef G4VERBOSE
119 if ((fpState->currentInteractionLength < 0.0) || (verboseLevel > 2))
120 {
121 G4cout << "G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength ";
122 G4cout << "[ " << GetProcessName() << "]" << G4endl;
123 track.GetDynamicParticle()->DumpInfo();
124 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
125 G4cout << "MeanLifeTime = " << fpState->currentInteractionLength / CLHEP::ns
126 << "[ns]" << G4endl;
127 }
128#endif
129
130 return fpState->theNumberOfInteractionLengthLeft
131 * (fpState->currentInteractionLength);
132}
133
135 const G4Step&)
136{
137// clear NumberOfInteractionLengthLeft
139
140 return pParticleChange;
141}
142
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
G4ProcessType
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void DumpInfo(G4int mode=0) const
const G4String & GetName() const
Definition: G4Material.hh:172
Definition: G4Step.hh:62
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
G4shared_ptr< G4ProcessState > fpState
virtual void ResetNumberOfInteractionLengthLeft()
virtual void ClearNumberOfInteractionLengthLeft()
G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &) override
G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *) override
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)=0
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
~G4VITRestDiscreteProcess() override
G4int verboseLevel
Definition: G4VProcess.hh:360
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:364
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
#define DBL_MAX
Definition: templates.hh:62