Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4VContinuousDiscreteProcess.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// $Id$
28//
29//
30// --------------------------------------------------------------
31// GEANT 4 class implementation file
32//
33// History: first implementation, based on object model of
34// 2nd December 1995, G.Cosmo
35// --------------------------------------------------------------
36// New Physics scheme 8 Jan. 1997 H.Kurahige
37// ------------------------------------------------------------
38
40#include "G4SystemOfUnits.hh"
41
42#include "G4Step.hh"
43#include "G4Track.hh"
44#include "G4MaterialTable.hh"
45#include "G4VParticleChange.hh"
46
48 :G4VProcess("No Name Discrete Process"),
49 valueGPILSelection(CandidateForSelection)
50{
51 G4Exception("G4VContinuousDiscreteProcess::G4VContinuousDiscreteProcess()","ProcMan102",
52 JustWarning,"Default constructor is called");
53}
54
56 : G4VProcess(aName, aType),
57 valueGPILSelection(CandidateForSelection)
58{
59 enableAtRestDoIt = false;
60}
61
63{
64}
65
67 : G4VProcess(right),
68 valueGPILSelection(right.valueGPILSelection)
69{
70}
71
73 const G4Track& track,
74 G4double previousStepSize,
76 )
77{
78 if ( (previousStepSize <=0.0) || (theNumberOfInteractionLengthLeft<=0.0)) {
79 // beggining of tracking (or just after DoIt of this process)
81 } else if ( previousStepSize > 0.0) {
82 // subtract NumberOfInteractionLengthLeft
84 } else {
85 // zero step
86 // DO NOTHING
87 }
88
89 // condition is set to "Not Forced"
91
92 // get mean free path
93 currentInteractionLength = GetMeanFreePath(track, previousStepSize, condition);
94
95 G4double value;
98 } else {
99 value = DBL_MAX;
100 }
101#ifdef G4VERBOSE
102 if (verboseLevel>1){
103 G4cout << "G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength ";
104 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
105 track.GetDynamicParticle()->DumpInfo();
106 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
107 G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl;
108 }
109#endif
110 return value;
111}
112
113
115 const G4Track& ,
116 const G4Step&
117 )
118{
119// clear NumberOfInteractionLengthLeft
121 return pParticleChange;
122}
123
125 const G4Track& ,
126 const G4Step&
127 )
128{
129// clear NumberOfInteractionLengthLeft
131 return pParticleChange;
132}
133
135 const G4Track& track,
136 G4double previousStepSize,
137 G4double currentMinimumStep,
138 G4double& currentSafety,
139 G4GPILSelection* selection
140 )
141{
142 // GPILSelection is set to defaule value of CandidateForSelection
143 valueGPILSelection = CandidateForSelection;
144
145 // get Step limit proposed by the process
146 G4double steplength = GetContinuousStepLimit(track,previousStepSize,currentMinimumStep, currentSafety);
147
148 // set return value for G4GPILSelection
149 *selection = valueGPILSelection;
150
151#ifdef G4VERBOSE
152 if (verboseLevel>1){
153 G4cout << "G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength ";
154 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
155 track.GetDynamicParticle()->DumpInfo();
156 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
157 G4cout << "IntractionLength= " << steplength/cm <<"[cm] " <<G4endl;
158 }
159#endif
160 return steplength ;
161}
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
G4ForceCondition
@ NotForced
G4GPILSelection
@ CandidateForSelection
G4ProcessType
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void DumpInfo(G4int mode=0) const
const G4String & GetName() const
Definition: G4Material.hh:177
Definition: G4Step.hh:78
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
virtual G4double GetContinuousStepLimit(const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4VContinuousDiscreteProcess(const G4String &, G4ProcessType aType=fNotDefined)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
G4double currentInteractionLength
Definition: G4VProcess.hh:297
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:92
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:418
void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
Definition: G4VProcess.cc:98
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define DBL_MAX
Definition: templates.hh:83