Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VRestContinuousDiscreteProcess.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// G4VRestContinuousDiscreteProcess
27//
28// Authors:
29// - 2 December 1995, G.Cosmo - First implementation, based on object model
30// - 8 January 1997, H.Kurashige - New Physics scheme
31// --------------------------------------------------------------------
32
34#include "G4SystemOfUnits.hh"
35
36#include "G4Step.hh"
37#include "G4Track.hh"
38#include "G4MaterialTable.hh"
39#include "G4VParticleChange.hh"
40
41// --------------------------------------------------------------------
43 : G4VProcess("No Name Discrete Process")
44{
45 G4Exception("G4VRestContinuousDiscreteProcess::G4VRestContinuousDiscreteProcess()",
46 "ProcMan102", JustWarning, "Default constructor is called");
47}
48
49// --------------------------------------------------------------------
50G4VRestContinuousDiscreteProcess::
51G4VRestContinuousDiscreteProcess(const G4String& aName, G4ProcessType aType)
52 : G4VProcess(aName, aType)
53{
54}
55
56// --------------------------------------------------------------------
60
61// --------------------------------------------------------------------
62G4VRestContinuousDiscreteProcess::
63G4VRestContinuousDiscreteProcess(G4VRestContinuousDiscreteProcess& right)
64 : G4VProcess(right),
65 valueGPILSelection(right.valueGPILSelection)
66{
67}
68
69// --------------------------------------------------------------------
73{
74 // beginning of tracking
76
77 // condition is set to "Not Forced"
79
80 // get mean life time
82
85
86#ifdef G4VERBOSE
87 if ((currentInteractionLength <0.0) || (verboseLevel>2))
88 {
91 G4cout << "G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength() - ";
92 G4cout << "[ " << GetProcessName() << "]" << G4endl;
94 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
95 G4cout << "MeanLifeTime = " << t << " [ns]"
96 << G4endl;
97 }
98#endif
99
100 return time;
101}
102
103// --------------------------------------------------------------------
112
113// --------------------------------------------------------------------
116 G4double previousStepSize,
117 G4double currentMinimumStep,
118 G4double& currentSafety,
119 G4GPILSelection* selection )
120{
121 // GPILSelection is set to defaule value of CandidateForSelection
122 valueGPILSelection = CandidateForSelection;
123
124 // get Step limit proposed by the process
125 G4double steplength = GetContinuousStepLimit(track, previousStepSize,
126 currentMinimumStep, currentSafety);
127
128 // set return value for G4GPILSelection
129 *selection = valueGPILSelection;
130
131#ifdef G4VERBOSE
132 if (verboseLevel>1)
133 {
134 G4cout << "G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength() - ";
135 G4cout << "[ " << GetProcessName() << "]" << G4endl;
136 track.GetDynamicParticle()->DumpInfo();
137 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
138 G4cout << "IntractionLength= " << steplength/cm <<"[cm] " << G4endl;
139 }
140#endif
141 return steplength;
142}
143
144// --------------------------------------------------------------------
151
152// --------------------------------------------------------------------
155 G4double previousStepSize,
157{
158 if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0))
159 {
160 // beginning of tracking (or just after DoIt() of this process)
162 }
163 else if ( previousStepSize > 0.0)
164 {
165 // subtract NumberOfInteractionLengthLeft
167 }
168 else
169 {
170 // zero step
171 // DO NOTHING
172 }
173
174 // condition is set to "Not Forced"
176
177 // get mean free path
178 currentInteractionLength = GetMeanFreePath(track,previousStepSize,condition);
179
180
181 G4double value;
183 {
185 }
186 else
187 {
188 value = DBL_MAX;
189 }
190#ifdef G4VERBOSE
191 if (verboseLevel>1)
192 {
193 G4cout << "G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength() - ";
194 G4cout << "[ " << GetProcessName() << "]" << G4endl;
195 track.GetDynamicParticle()->DumpInfo();
196 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
197 G4cout << "InteractionLength= " << value/cm <<"[cm] " << G4endl;
198 }
199#endif
200 return value;
201}
202
203// --------------------------------------------------------------------
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4ForceCondition
@ NotForced
G4GPILSelection
@ CandidateForSelection
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
G4double currentInteractionLength
void SubtractNumberOfInteractionLengthLeft(G4double prevStepSize)
virtual void ResetNumberOfInteractionLengthLeft()
Definition G4VProcess.cc:80
void ClearNumberOfInteractionLengthLeft()
G4int verboseLevel
G4double theNumberOfInteractionLengthLeft
G4VParticleChange * pParticleChange
const G4String & GetProcessName() const
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)=0
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual G4double GetContinuousStepLimit(const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
G4VRestContinuousDiscreteProcess(const G4String &, G4ProcessType aType=fNotDefined)
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
#define DBL_MAX
Definition templates.hh:62