Geant4-11
G4VRestDiscreteProcess.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// G4VRestDiscreteProcess class implementation
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 "Randomize.hh"
37#include "G4Step.hh"
38#include "G4Track.hh"
39#include "G4MaterialTable.hh"
40#include "G4VParticleChange.hh"
41
42// --------------------------------------------------------------------
44 : G4VProcess("No Name Discrete Process")
45{
46 G4Exception("G4VRestDiscreteProcess::G4VRestDiscreteProcess", "ProcMan102",
47 JustWarning, "Default constructor is called");
48}
49
50// --------------------------------------------------------------------
53 : G4VProcess(aName, aType)
54{
55 enableAlongStepDoIt = false;
56}
57
58// --------------------------------------------------------------------
60{
61}
62
63// --------------------------------------------------------------------
66 : G4VProcess(right)
67{
68}
69
70// --------------------------------------------------------------------
73 G4double previousStepSize,
75{
76 if ( (previousStepSize < 0.0) || (theNumberOfInteractionLengthLeft<=0.0))
77 {
78 // beginning of tracking (or just after DoIt() of this process)
80 }
81 else if ( previousStepSize > 0.0)
82 {
83 // subtract NumberOfInteractionLengthLeft
85 }
86 else
87 {
88 // zero step
89 // DO NOTHING
90 }
91
92 // condition is set to "Not Forced"
94
95 // get mean free path
96 currentInteractionLength = GetMeanFreePath(track,previousStepSize,condition);
97
98 G4double value;
100 {
102 }
103 else
104 {
105 value = DBL_MAX;
106 }
107#ifdef G4VERBOSE
108 if (verboseLevel>1)
109 {
110 G4cout << "G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength() - ";
111 G4cout << "[ " << GetProcessName() << "]" << G4endl;
112 track.GetDynamicParticle()->DumpInfo();
113 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
114 G4cout << "InteractionLength= " << value/cm <<"[cm] " << G4endl;
115 }
116#endif
117 return value;
118}
119
120// --------------------------------------------------------------------
123 const G4Step& )
124{
126
127 return pParticleChange;
128}
129
130// --------------------------------------------------------------------
134{
135 // beginning of tracking
137
138 // condition is set to "Not Forced"
140
141 // get mean life time
143
144#ifdef G4VERBOSE
145 if ((currentInteractionLength <0.0) || (verboseLevel>2))
146 {
147 G4cout << "G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength() - ";
148 G4cout << "[ " << GetProcessName() << "]" << G4endl;
149 track.GetDynamicParticle()->DumpInfo();
150 G4cout << " in Material " << track.GetMaterial()->GetName() << G4endl;
151 G4cout << "MeanLifeTime = " << currentInteractionLength/ns << "[ns]"
152 << G4endl;
153 }
154#endif
155
157}
158
159// --------------------------------------------------------------------
162 const G4Step& )
163{
165
166 return pParticleChange;
167}
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
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
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:360
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:331
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
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 GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
#define DBL_MAX
Definition: templates.hh:62
#define ns
Definition: xmlparse.cc:614