Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4Step.icc
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 inline methods 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
35inline G4StepPoint* G4Step::GetPreStepPoint() const
36{
37 return fpPreStepPoint;
38}
39
40inline void G4Step::SetPreStepPoint(G4StepPoint* value)
41{
42 delete fpPreStepPoint;
43 fpPreStepPoint = value;
44}
45
46inline G4StepPoint* G4Step::GetPostStepPoint() const
47{
48 return fpPostStepPoint;
49}
50
51inline void G4Step::SetPostStepPoint(G4StepPoint* value)
52{
53 delete fpPostStepPoint;
54 fpPostStepPoint = value;
55}
56
57inline G4double G4Step::GetStepLength() const
58{
59 return fStepLength;
60}
61
62inline void G4Step::SetStepLength(G4double value)
63{
64 fStepLength = value;
65}
66
67inline G4ThreeVector G4Step::GetDeltaPosition() const
68{
69 return fpPostStepPoint->GetPosition() - fpPreStepPoint->GetPosition();
70}
71
72inline G4double G4Step::GetDeltaTime() const
73{
74 return fpPostStepPoint->GetLocalTime() - fpPreStepPoint->GetLocalTime();
75}
76
77inline G4double G4Step::GetTotalEnergyDeposit() const
78{
79 return fTotalEnergyDeposit;
80}
81
82inline void G4Step::SetTotalEnergyDeposit(G4double value)
83{
84 fTotalEnergyDeposit = value;
85}
86
87inline G4double G4Step::GetNonIonizingEnergyDeposit() const
88{
89 return fNonIonizingEnergyDeposit;
90}
91
92inline void G4Step::SetNonIonizingEnergyDeposit(G4double value)
93{
94 fNonIonizingEnergyDeposit = value;
95}
96
97inline void G4Step::AddTotalEnergyDeposit(G4double value)
98{
99 fTotalEnergyDeposit += value;
100}
101
102inline void G4Step::ResetTotalEnergyDeposit()
103{
104 fTotalEnergyDeposit = 0.;
105 fNonIonizingEnergyDeposit = 0.;
106}
107
108inline void G4Step::AddNonIonizingEnergyDeposit(G4double value)
109{
110 fNonIonizingEnergyDeposit += value;
111}
112
113inline void G4Step::ResetNonIonizingEnergyDeposit()
114{
115 fNonIonizingEnergyDeposit = 0.;
116}
117
118inline void G4Step::SetControlFlag(G4SteppingControl value)
119{
120 fpSteppingControlFlag = value;
121}
122
123inline G4SteppingControl G4Step::GetControlFlag() const
124{
125 return fpSteppingControlFlag;
126}
127
128inline void G4Step::CopyPostToPreStepPoint()
129{
130 // This method is called at the beginning of each step
131 *(fpPreStepPoint) = *(fpPostStepPoint);
132 fpPostStepPoint->SetStepStatus(fUndefined);
133
134 // store number of secondaries
135 nSecondaryByLastStep = fSecondary->size();
136}
137
138//-------------------------------------------------------------
139// To implement bi-directional association between G4Step and
140// and G4Track, a combined usage of 'forward declaration' and
141// 'include' is necessary.
142//-------------------------------------------------------------
143#include "G4Track.hh"
144
145inline G4Track* G4Step::GetTrack() const
146{
147 return fpTrack;
148}
149
150inline void G4Step::SetTrack(G4Track* value)
151{
152 fpTrack = value;
153}
154
155inline void G4Step::InitializeStep(G4Track* aValue)
156{
157 fpTrack = aValue;
158 fpTrack->SetStepLength(0.);
159
160 // Initialize G4Step attributes
161 fStepLength = 0.;
162 fTotalEnergyDeposit = 0.;
163 fNonIonizingEnergyDeposit = 0.;
164
165 nSecondaryByLastStep = 0;
166
167 // Initialize G4StepPoint attributes.
168 // To avoid the circular dependency between G4Track, G4Step
169 // and G4StepPoint, G4Step has to manage the copy actions.
170 //
171 fpPreStepPoint->SetSafety(0.0);
172 fpPreStepPoint->SetProcessDefinedStep(0);
173 fpPreStepPoint->SetStepStatus(fUndefined);
174
175 // G4Track properties
176
177 fpPreStepPoint->SetWeight(fpTrack->GetWeight());
178 fpPreStepPoint->SetPosition(fpTrack->GetPosition());
179 fpPreStepPoint->SetLocalTime(fpTrack->GetLocalTime());
180 fpPreStepPoint->SetGlobalTime(fpTrack->GetGlobalTime());
181
182 // G4Track properties (via G4DynamicParticle)
183
184 auto pParticle = fpTrack->GetDynamicParticle();
185
186 fpPreStepPoint->SetMass(pParticle->GetMass());
187 fpPreStepPoint->SetCharge(pParticle->GetCharge());
188 fpPreStepPoint->SetProperTime(pParticle->GetProperTime());
189 fpPreStepPoint->SetKineticEnergy(pParticle->GetKineticEnergy());
190 fpPreStepPoint->SetMomentumDirection(pParticle->GetMomentumDirection());
191 fpPreStepPoint->SetPolarization(pParticle->GetPolarization());
192
193 // G4Track properties (via G4LogicalVolume)
194
195 auto lvol = fpTrack->GetTouchable()->GetVolume()->GetLogicalVolume();
196
197 fpPreStepPoint->SetTouchableHandle(fpTrack->GetTouchableHandle());
198
199 fpPreStepPoint->SetMaterial(lvol->GetMaterial());
200 fpPreStepPoint->SetMaterialCutsCouple(lvol->GetMaterialCutsCouple());
201 fpPreStepPoint->SetSensitiveDetector(lvol->GetSensitiveDetector());
202
203 // SetVelocity() must be called after SetMaterial() for fpPreStepPoint.
204 fpPreStepPoint->SetVelocity(fpTrack->CalculateVelocity());
205
206 (*fpPostStepPoint) = (*fpPreStepPoint);
207}
208
209inline void G4Step::UpdateTrack()
210{
211 // To avoid the circular dependency between G4Track, G4Step
212 // and G4StepPoint, G4Step has to manage the update actions.
213 // position, time
214 fpTrack->SetPosition(fpPostStepPoint->GetPosition());
215 fpTrack->SetGlobalTime(fpPostStepPoint->GetGlobalTime());
216 fpTrack->SetLocalTime(fpPostStepPoint->GetLocalTime());
217 // mass, charge, proper time
218 G4DynamicParticle* pParticle =
219 (G4DynamicParticle*) (fpTrack->GetDynamicParticle());
220 pParticle->SetMass(fpPostStepPoint->GetMass());
221 pParticle->SetCharge(fpPostStepPoint->GetCharge());
222 pParticle->SetProperTime(fpPostStepPoint->GetProperTime());
223 // energy, momentum, polarization
224 pParticle->SetMomentumDirection(fpPostStepPoint->GetMomentumDirection());
225 pParticle->SetKineticEnergy(fpPostStepPoint->GetKineticEnergy());
226 pParticle->SetPolarization(fpPostStepPoint->GetPolarization());
227 // step length
228 fpTrack->SetStepLength(fStepLength);
229 // NextTouchable is updated
230 // (G4Track::Touchable points touchable of Pre-StepPoint)
231 fpTrack->SetNextTouchableHandle(fpPostStepPoint->GetTouchableHandle());
232 fpTrack->SetWeight(fpPostStepPoint->GetWeight());
233
234 // set velocity
235 fpTrack->SetVelocity(fpPostStepPoint->GetVelocity());
236}
237
238inline std::size_t G4Step::GetNumberOfSecondariesInCurrentStep() const
239{
240 return fSecondary->size() - nSecondaryByLastStep;
241}
242
243inline const G4TrackVector* G4Step::GetSecondary() const
244{
245 return fSecondary;
246}
247
248inline G4TrackVector* G4Step::GetfSecondary()
249{
250 return fSecondary;
251}
252
253inline void G4Step::SetSecondary(G4TrackVector* value)
254{
255 fSecondary = value;
256}
257
258inline G4TrackVector* G4Step::NewSecondaryVector()
259{
260 fSecondary = new G4TrackVector();
261 return fSecondary;
262}
263
264inline void G4Step::DeleteSecondaryVector()
265{
266 if(fSecondary != nullptr)
267 {
268 fSecondary->clear();
269 delete fSecondary;
270 fSecondary = nullptr;
271 }
272}
273
274inline G4bool G4Step::IsFirstStepInVolume() const
275{
276 return fFirstStepInVolume;
277}
278
279inline G4bool G4Step::IsLastStepInVolume() const
280{
281 return fLastStepInVolume;
282}
283
284inline void G4Step::SetFirstStepFlag()
285{
286 fFirstStepInVolume = true;
287}
288
289inline void G4Step::ClearFirstStepFlag()
290{
291 fFirstStepInVolume = false;
292}
293
294inline void G4Step::SetLastStepFlag()
295{
296 fLastStepInVolume = true;
297}
298
299inline void G4Step::ClearLastStepFlag()
300{
301 fLastStepInVolume = false;
302}
303
304inline void
305G4Step::SetPointerToVectorOfAuxiliaryPoints(std::vector<G4ThreeVector>* vec)
306{
307 fpVectorOfAuxiliaryPointsPointer = vec;
308}
309
310inline std::vector<G4ThreeVector>*
311G4Step::GetPointerToVectorOfAuxiliaryPoints() const
312{
313 return fpVectorOfAuxiliaryPointsPointer;
314}