Geant4-11
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// --------------------------------------------------------------------
52 : G4VProcess(aName, aType)
53{
54}
55
56// --------------------------------------------------------------------
58{
59}
60
61// --------------------------------------------------------------------
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
83#ifdef G4VERBOSE
84 if ((currentInteractionLength <0.0) || (verboseLevel>2))
85 {
86 G4cout << "G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength() - ";
87 G4cout << "[ " << GetProcessName() << "]" << G4endl;
89 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
90 G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]"
91 << G4endl;
92 }
93#endif
94
96}
97
98// --------------------------------------------------------------------
101 const G4Step& )
102{
104
105 return pParticleChange;
106}
107
108// --------------------------------------------------------------------
111 G4double previousStepSize,
112 G4double currentMinimumStep,
113 G4double& currentSafety,
114 G4GPILSelection* selection )
115{
116 // GPILSelection is set to defaule value of CandidateForSelection
118
119 // get Step limit proposed by the process
120 G4double steplength = GetContinuousStepLimit(track, previousStepSize,
121 currentMinimumStep, currentSafety);
122
123 // set return value for G4GPILSelection
124 *selection = valueGPILSelection;
125
126#ifdef G4VERBOSE
127 if (verboseLevel>1)
128 {
129 G4cout << "G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength() - ";
130 G4cout << "[ " << GetProcessName() << "]" << G4endl;
131 track.GetDynamicParticle()->DumpInfo();
132 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
133 G4cout << "IntractionLength= " << steplength/cm <<"[cm] " << G4endl;
134 }
135#endif
136 return steplength;
137}
138
139// --------------------------------------------------------------------
142 const G4Step& )
143{
144 return pParticleChange;
145}
146
147// --------------------------------------------------------------------
150 G4double previousStepSize,
152{
153 if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0))
154 {
155 // beginning of tracking (or just after DoIt() of this process)
157 }
158 else if ( previousStepSize > 0.0)
159 {
160 // subtract NumberOfInteractionLengthLeft
162 }
163 else
164 {
165 // zero step
166 // DO NOTHING
167 }
168
169 // condition is set to "Not Forced"
171
172 // get mean free path
173 currentInteractionLength = GetMeanFreePath(track,previousStepSize,condition);
174
175
176 G4double value;
178 {
180 }
181 else
182 {
183 value = DBL_MAX;
184 }
185#ifdef G4VERBOSE
186 if (verboseLevel>1)
187 {
188 G4cout << "G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength() - ";
189 G4cout << "[ " << GetProcessName() << "]" << G4endl;
190 track.GetDynamicParticle()->DumpInfo();
191 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
192 G4cout << "InteractionLength= " << value/cm <<"[cm] " << G4endl;
193 }
194#endif
195 return value;
196}
197
198// --------------------------------------------------------------------
201 const G4Step& )
202{
204
205 return pParticleChange;
206}
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4ForceCondition
@ NotForced
G4GPILSelection
@ CandidateForSelection
G4ProcessType
static constexpr double cm
Definition: G4SIunits.hh:99
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:173
Definition: G4Step.hh:62
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
G4double currentInteractionLength
Definition: G4VProcess.hh:335
void SubtractNumberOfInteractionLengthLeft(G4double prevStepSize)
Definition: G4VProcess.hh:524
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:80
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:424
G4int verboseLevel
Definition: G4VProcess.hh:356
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:331
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
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 &)
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
#define DBL_MAX
Definition: templates.hh:62
#define ns
Definition: xmlparse.cc:614