Geant4-11
G4VParticleChange.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// G4VParticleChange inline methods implementation
27//
28// Author: Hisaya Kurashige, 23 March 1998
29// --------------------------------------------------------------------
30
31inline G4Track* G4VParticleChange::GetSecondary(G4int anIndex) const
32{
33 return (*theListOfSecondaries)[anIndex];
34}
35
36inline G4int G4VParticleChange::GetNumberOfSecondaries() const
37{
38 return theNumberOfSecondaries;
39}
40
41inline void G4VParticleChange::ProposeTrackStatus(G4TrackStatus aStatus)
42{
43 theStatusChange = aStatus;
44}
45
46inline G4TrackStatus G4VParticleChange::GetTrackStatus() const
47{
48 return theStatusChange;
49}
50
51inline G4SteppingControl G4VParticleChange::GetSteppingControl() const
52{
53 return theSteppingControlFlag;
54}
55
56inline void
57G4VParticleChange::ProposeSteppingControl(G4SteppingControl StepControlFlag)
58{
59 theSteppingControlFlag = StepControlFlag;
60}
61
62inline G4bool G4VParticleChange::GetFirstStepInVolume() const
63{
64 return theFirstStepInVolume;
65}
66
67inline G4bool G4VParticleChange::GetLastStepInVolume() const
68{
69 return theLastStepInVolume;
70}
71
72inline void G4VParticleChange::ProposeFirstStepInVolume(G4bool flag)
73{
74 theFirstStepInVolume = flag;
75}
76
77inline void G4VParticleChange::ProposeLastStepInVolume(G4bool flag)
78{
79 theLastStepInVolume = flag;
80}
81
82inline G4double G4VParticleChange::GetLocalEnergyDeposit() const
83{
84 return theLocalEnergyDeposit;
85}
86
87inline void G4VParticleChange::ProposeLocalEnergyDeposit(G4double anEnergyPart)
88{
89 theLocalEnergyDeposit = anEnergyPart;
90}
91
92inline G4double G4VParticleChange::GetNonIonizingEnergyDeposit() const
93{
94 return theNonIonizingEnergyDeposit;
95}
96
97inline void
98G4VParticleChange::ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
99{
100 theNonIonizingEnergyDeposit = anEnergyPart;
101}
102
103inline G4double G4VParticleChange::GetTrueStepLength() const
104{
105 return theTrueStepLength;
106}
107
108inline void G4VParticleChange::ProposeTrueStepLength(G4double aLength)
109{
110 theTrueStepLength = aLength;
111}
112
113inline void G4VParticleChange::SetVerboseLevel(G4int vLevel)
114{
115 verboseLevel = vLevel;
116}
117
118inline G4int G4VParticleChange::GetVerboseLevel() const
119{
120 return verboseLevel;
121}
122
123inline G4double G4VParticleChange::GetParentWeight() const
124{
125 return theParentWeight;
126}
127
128inline G4double G4VParticleChange::GetWeight() const
129{
130 return theParentWeight;
131}
132
133// ----------------------------------------------------------------
134// inline functions for Initialization
135//
136inline void G4VParticleChange::InitializeLocalEnergyDeposit(const G4Track&)
137{
138 // clear theLocalEnergyDeposited
139 theLocalEnergyDeposit = 0.0;
140 theNonIonizingEnergyDeposit = 0.0;
141}
142
143inline void G4VParticleChange::InitializeSteppingControl(const G4Track&)
144{
145 // SteppingControlFlag
146 theSteppingControlFlag = NormalCondition;
147}
148
149inline void G4VParticleChange::Clear()
150{
151 theNumberOfSecondaries = 0;
152 theFirstStepInVolume = false;
153 theLastStepInVolume = false;
154}
155
156inline void G4VParticleChange::InitializeStatusChange(const G4Track& track)
157{
158 // set TrackStatus equal to the parent track's one
159 theStatusChange = track.GetTrackStatus();
160}
161
162inline void G4VParticleChange::InitializeParentWeight(const G4Track& track)
163{
164 // set the parent track's weight
165 theParentWeight = track.GetWeight();
166 isParentWeightProposed = false;
167}
168
169inline void G4VParticleChange::InitializeParentGlobalTime(const G4Track& track)
170{
171 // set the parent track's global time at the pre-step point
172 theParentGlobalTime = track.GetStep()->GetPreStepPoint()->GetGlobalTime();
173}
174
175inline void G4VParticleChange::InitializeTrueStepLength(const G4Track& track)
176{
177 // Reset theTrueStepLength
178 // !! Caution !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
179 theTrueStepLength = track.GetStep()->GetStepLength();
180 // !! TrueStepLength should be copied from G4Step not G4Track
181 // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
182}
183
184inline void G4VParticleChange::InitializeStepInVolumeFlags(const G4Track& track)
185{
186 const G4Step* aStep = track.GetStep();
187 theFirstStepInVolume = aStep->IsFirstStepInVolume();
188 theLastStepInVolume = aStep->IsLastStepInVolume();
189}
190
191inline void G4VParticleChange::InitializeSecondaries(const G4Track&)
192{
193 // clear secondaries
194 if(theNumberOfSecondaries > 0)
195 {
196#ifdef G4VERBOSE
197 if(verboseLevel > 0)
198 {
199 G4cerr << "G4VParticleChange::Initialize() Warning ";
200 G4cerr << "theListOfSecondaries is not empty " << G4endl;
201 G4cerr << "All objects in theListOfSecondaries are destroyed!" << G4endl;
202 }
203#endif
204 for(G4int index = 0; index < theNumberOfSecondaries; index++)
205 {
206 if((*theListOfSecondaries)[index])
207 {
208 delete(*theListOfSecondaries)[index];
209 }
210 }
211 }
212 theNumberOfSecondaries = 0;
213}
214
215// ----------------------------------------------------------------
216// methods for handling secondaries
217//
218inline void G4VParticleChange::SetNumberOfSecondaries(G4int totSecondaries)
219{
220 // check if tracks still exist in theListOfSecondaries
221 if(theNumberOfSecondaries > 0)
222 {
223#ifdef G4VERBOSE
224 if(verboseLevel > 0)
225 {
226 G4cerr << "G4VParticleChange::SetNumberOfSecondaries() Warning ";
227 G4cerr << "theListOfSecondaries is not empty ";
228 }
229#endif
230 for(G4int index = 0; index < theNumberOfSecondaries; index++)
231 {
232 if((*theListOfSecondaries)[index])
233 {
234 delete(*theListOfSecondaries)[index];
235 }
236 }
237 }
238 theNumberOfSecondaries = 0;
239 theSizeOftheListOfSecondaries = totSecondaries;
240
241 // Initialize ListOfSecondaries
242 theListOfSecondaries->Initialize(totSecondaries);
243}
244
245inline void G4VParticleChange::Initialize(const G4Track& track)
246{
247 InitializeStatusChange(track);
248 InitializeLocalEnergyDeposit(track);
249 InitializeSteppingControl(track);
250 InitializeTrueStepLength(track);
251 InitializeSecondaries(track);
252 InitializeParentWeight(track);
253 InitializeParentGlobalTime(track);
254 InitializeStepInVolumeFlags(track);
255}
256
257inline void G4VParticleChange::ClearDebugFlag()
258{
259 debugFlag = false;
260}
261
262inline void G4VParticleChange::SetDebugFlag()
263{
264 debugFlag = true;
265}
266
267inline G4bool G4VParticleChange::GetDebugFlag() const
268{
269 return debugFlag;
270}
271
272inline void G4VParticleChange::SetSecondaryWeightByProcess(G4bool flag)
273{
274 fSetSecondaryWeightByProcess = flag;
275}
276
277inline G4bool G4VParticleChange::IsSecondaryWeightSetByProcess() const
278{
279 return fSetSecondaryWeightByProcess;
280}
281
282inline void G4VParticleChange::ProposeWeight(G4double w)
283{
284 theParentWeight = w;
285 isParentWeightProposed = true;
286}
287
288inline void G4VParticleChange::ProposeParentWeight(G4double w)
289{
290 ProposeWeight(w);
291}