Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4Step.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// G4Step class implementation
27//
28// Authors:
29// Katsuya Amako (e-mail: Katsuya.Amako@kek.jp)
30// Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp)
31// Revisions:
32// Hisaya Kurashige, 1998-2007
33// --------------------------------------------------------------------
34
35#include "G4Step.hh"
36
37// --------------------------------------------------------------------
39{
40 fpPreStepPoint = new G4StepPoint();
41 fpPostStepPoint = new G4StepPoint();
42
43 secondaryInCurrentStep = new std::vector<const G4Track*>;
44}
45
46// --------------------------------------------------------------------
48{
49 delete fpPreStepPoint;
50 delete fpPostStepPoint;
51
52 secondaryInCurrentStep->clear();
53 delete secondaryInCurrentStep;
54
55 if(fSecondary != nullptr)
56 {
57 fSecondary->clear();
58 delete fSecondary;
59 }
60}
61
62// --------------------------------------------------------------------
64 : fTotalEnergyDeposit(right.fTotalEnergyDeposit)
65 , fNonIonizingEnergyDeposit(right.fNonIonizingEnergyDeposit)
66 , fStepLength(right.fStepLength)
67 , fpTrack(right.fpTrack)
68 , fpSteppingControlFlag(right.fpSteppingControlFlag)
69 , fFirstStepInVolume(right.fFirstStepInVolume)
70 , fLastStepInVolume(right.fLastStepInVolume)
71 , nSecondaryByLastStep(right.nSecondaryByLastStep)
72 , secondaryInCurrentStep(right.secondaryInCurrentStep)
73 , fpVectorOfAuxiliaryPointsPointer(right.fpVectorOfAuxiliaryPointsPointer)
74{
75 if(right.fpPreStepPoint != nullptr)
76 {
77 fpPreStepPoint = new G4StepPoint(*(right.fpPreStepPoint));
78 }
79 else
80 {
81 fpPreStepPoint = new G4StepPoint();
82 }
83 if(right.fpPostStepPoint != nullptr)
84 {
85 fpPostStepPoint = new G4StepPoint(*(right.fpPostStepPoint));
86 }
87 else
88 {
89 fpPostStepPoint = new G4StepPoint();
90 }
91
92 if(right.fSecondary != nullptr)
93 {
94 fSecondary = new G4TrackVector(*(right.fSecondary));
95 }
96 else
97 {
98 fSecondary = new G4TrackVector();
99 }
100
101 // secondaryInCurrentStep is cleared
102 secondaryInCurrentStep = new std::vector<const G4Track*>;
103}
104
105// --------------------------------------------------------------------
107{
108 if(this != &right)
109 {
112 fStepLength = right.fStepLength;
113 fpTrack = right.fpTrack;
114 fpSteppingControlFlag = right.fpSteppingControlFlag;
115 fFirstStepInVolume = right.fFirstStepInVolume;
116 fLastStepInVolume = right.fLastStepInVolume;
117 nSecondaryByLastStep = right.nSecondaryByLastStep;
118 secondaryInCurrentStep = right.secondaryInCurrentStep;
119 fpVectorOfAuxiliaryPointsPointer = right.fpVectorOfAuxiliaryPointsPointer;
120
121 delete fpPreStepPoint;
122
123 if(right.fpPreStepPoint != nullptr)
124 {
125 fpPreStepPoint = new G4StepPoint(*(right.fpPreStepPoint));
126 }
127 else
128 {
129 fpPreStepPoint = new G4StepPoint();
130 }
131
132 delete fpPostStepPoint;
133
134 if(right.fpPostStepPoint != nullptr)
135 {
136 fpPostStepPoint = new G4StepPoint(*(right.fpPostStepPoint));
137 }
138 else
139 {
140 fpPostStepPoint = new G4StepPoint();
141 }
142
143 if(fSecondary != nullptr)
144 {
145 fSecondary->clear();
146 delete fSecondary;
147 }
148 if(right.fSecondary != nullptr)
149 {
150 fSecondary = new G4TrackVector(*(right.fSecondary));
151 }
152 else
153 {
154 fSecondary = new G4TrackVector();
155 }
156
157 // secondaryInCurrentStep is not copied
158 if(secondaryInCurrentStep != nullptr)
159 {
160 secondaryInCurrentStep->clear();
161 delete secondaryInCurrentStep;
162 }
163 secondaryInCurrentStep = new std::vector<const G4Track*>;
164 }
165 return *this;
166}
167
168// --------------------------------------------------------------------
170{
171 static G4ThreadLocal G4bool isFirstTime = true;
172 if(isFirstTime)
173 {
174 isFirstTime = false;
175#ifdef G4VERBOSE
176 G4Exception("G4Step::GetDeltaMomentum()", "Warning", JustWarning,
177 "This method is obsolete and will be removed soon");
178#endif
179 }
180
181 return fpPostStepPoint->GetMomentum() - fpPreStepPoint->GetMomentum();
182}
183
184// --------------------------------------------------------------------
186{
187 static G4ThreadLocal G4bool isFirstTime = true;
188 if(isFirstTime)
189 {
190 isFirstTime = false;
191#ifdef G4VERBOSE
192 G4Exception("G4Step::GetDeltaEnergy()", "Warning", JustWarning,
193 "This method is obsolete and will be removed soon");
194#endif
195 }
196
197 return fpPostStepPoint->GetKineticEnergy() -
198 fpPreStepPoint->GetKineticEnergy();
199}
200
201// --------------------------------------------------------------------
202const std::vector<const G4Track*>* G4Step::GetSecondaryInCurrentStep() const
203{
204 secondaryInCurrentStep->clear();
205 std::size_t nSecondary = fSecondary->size();
206 for(std::size_t i = nSecondaryByLastStep; i < nSecondary; ++i)
207 {
208 secondaryInCurrentStep->push_back((*fSecondary)[i]);
209 }
210 return secondaryInCurrentStep;
211}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::vector< G4Track * > G4TrackVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
G4ThreeVector GetMomentum() const
G4double GetKineticEnergy() const
Definition: G4Step.hh:62
~G4Step()
Definition: G4Step.cc:47
G4double GetDeltaEnergy() const
Definition: G4Step.cc:185
G4ThreeVector GetDeltaMomentum() const
Definition: G4Step.cc:169
G4Step()
Definition: G4Step.cc:38
G4Step & operator=(const G4Step &)
Definition: G4Step.cc:106
const std::vector< const G4Track * > * GetSecondaryInCurrentStep() const
Definition: G4Step.cc:202
G4double fNonIonizingEnergyDeposit
Definition: G4Step.hh:187
G4double fTotalEnergyDeposit
Definition: G4Step.hh:184
#define G4ThreadLocal
Definition: tls.hh:77