source: trunk/source/track/include/G4VParticleChange.icc@ 1199

Last change on this file since 1199 was 1058, checked in by garnier, 17 years ago

file release beta

File size: 9.5 KB
Line 
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: G4VParticleChange.icc,v 1.16 2007/03/25 22:54:52 kurasige Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30// remove obsolete methods of SetXXX 19 Sep, 04 H.Kurashige
31//----------------------------------------------------------------
32// Virtual methods for updating G4Step
33//
34
35inline G4Step* G4VParticleChange::UpdateStepInfo(G4Step* pStep)
36{
37 // Update the G4Step specific attributes
38 pStep->SetStepLength( theTrueStepLength );
39 pStep->AddTotalEnergyDeposit( theLocalEnergyDeposit );
40 pStep->AddNonIonizingEnergyDeposit( theNonIonizingEnergyDeposit );
41 pStep->SetControlFlag( theSteppingControlFlag );
42
43 if (theFirstStepInVolume) {pStep->SetFirstStepFlag();}
44 else {pStep->ClearFirstStepFlag();}
45 if (theLastStepInVolume) {pStep->SetLastStepFlag();}
46 else {pStep->ClearLastStepFlag();}
47
48 return pStep;
49}
50
51inline
52 G4Step* G4VParticleChange::UpdateStepForAtRest(G4Step* Step)
53{
54 if (!fSetParentWeightByProcess){
55 Step->GetPostStepPoint()->SetWeight( theParentWeight );
56 }
57 return UpdateStepInfo(Step);
58}
59
60inline
61 G4Step* G4VParticleChange::UpdateStepForAlongStep(G4Step* Step)
62{
63 if (!fSetParentWeightByProcess){
64 G4double newWeight= theParentWeight/(Step->GetPreStepPoint()->GetWeight())*(Step->GetPostStepPoint()->GetWeight());
65 Step->GetPostStepPoint()->SetWeight( newWeight );
66 }
67 return UpdateStepInfo(Step);
68}
69
70inline
71 G4Step* G4VParticleChange::UpdateStepForPostStep(G4Step* Step)
72{
73 if (!fSetParentWeightByProcess){
74 Step->GetPostStepPoint()->SetWeight( theParentWeight );
75 }
76 return UpdateStepInfo(Step);
77}
78
79//----------------------------------------------------------------
80// Set/Get inline functions
81//
82inline
83 G4Track* G4VParticleChange::GetSecondary(G4int anIndex) const
84{
85 return (*theListOfSecondaries)[anIndex];
86}
87
88inline
89 G4int G4VParticleChange::GetNumberOfSecondaries() const
90{
91 return theNumberOfSecondaries;
92}
93
94inline
95 void G4VParticleChange::ProposeTrackStatus(G4TrackStatus aStatus)
96{
97 theStatusChange = aStatus;
98}
99
100inline
101 G4TrackStatus G4VParticleChange::GetTrackStatus() const
102{
103 return theStatusChange;
104}
105
106inline
107G4SteppingControl G4VParticleChange::GetSteppingControl() const
108{
109 return theSteppingControlFlag;
110}
111
112inline
113void G4VParticleChange::ProposeSteppingControl(G4SteppingControl StepControlFlag)
114{
115 theSteppingControlFlag = StepControlFlag;
116}
117
118inline
119G4bool G4VParticleChange::GetFirstStepInVolume() const
120{
121 return theFirstStepInVolume;
122}
123
124inline
125G4bool G4VParticleChange::GetLastStepInVolume() const
126{
127 return theLastStepInVolume;
128}
129
130inline
131void G4VParticleChange::ProposeFirstStepInVolume(G4bool flag)
132{
133 theFirstStepInVolume = flag;
134}
135
136inline
137void G4VParticleChange::ProposeLastStepInVolume(G4bool flag)
138{
139 theLastStepInVolume = flag;
140}
141
142//----------------------------------------------------------------
143// Set/Get inline functions
144//
145
146inline
147 G4double G4VParticleChange::GetLocalEnergyDeposit() const
148{
149 return theLocalEnergyDeposit;
150}
151
152inline
153 void G4VParticleChange::ProposeLocalEnergyDeposit(G4double anEnergyPart)
154{
155 theLocalEnergyDeposit = anEnergyPart;
156}
157
158inline
159 G4double G4VParticleChange::GetNonIonizingEnergyDeposit() const
160{
161 return theNonIonizingEnergyDeposit;
162}
163
164inline
165 void G4VParticleChange::ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
166{
167 theNonIonizingEnergyDeposit = anEnergyPart;
168}
169
170inline
171 G4double G4VParticleChange::GetTrueStepLength() const
172{
173 return theTrueStepLength;
174}
175
176inline
177 void G4VParticleChange::ProposeTrueStepLength(G4double aLength)
178{
179 theTrueStepLength = aLength;
180}
181
182
183inline
184 void G4VParticleChange::SetVerboseLevel(G4int vLevel)
185{
186 verboseLevel = vLevel;
187}
188
189inline
190 G4int G4VParticleChange::GetVerboseLevel() const
191{
192 return verboseLevel;
193}
194
195inline
196 G4double G4VParticleChange::GetParentWeight() const
197{
198 return theParentWeight;
199}
200
201//----------------------------------------------------------------
202// inline functions for Initialization
203//
204
205inline
206 void G4VParticleChange::InitializeLocalEnergyDeposit(const G4Track&)
207{
208 // clear theLocalEnergyDeposited
209 theLocalEnergyDeposit = 0.0;
210 theNonIonizingEnergyDeposit = 0.0;
211}
212
213inline
214 void G4VParticleChange::InitializeSteppingControl(const G4Track& )
215{
216 // SteppingControlFlag
217 theSteppingControlFlag = NormalCondition;
218}
219
220inline
221 void G4VParticleChange::Clear()
222{
223 theNumberOfSecondaries = 0;
224 theFirstStepInVolume = false;
225 theLastStepInVolume = false;
226}
227
228//----------------------------------------------------------------
229// functions for Initialization
230//
231
232inline void G4VParticleChange::InitializeStatusChange(const G4Track& track)
233{
234 // set TrackStatus equal to the parent track's one
235 theStatusChange = track.GetTrackStatus();
236}
237
238inline void G4VParticleChange::InitializeParentWeight(const G4Track& track)
239{
240 // set the parent track's weight
241 theParentWeight = track.GetWeight();
242}
243
244inline void G4VParticleChange::InitializeTrueStepLength(const G4Track& track)
245{
246 // Reset theTrueStepLength
247 // !! Caution !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
248 theTrueStepLength = track.GetStep()->GetStepLength();
249 // !! TrueStepLength should be copied from G4Step not G4Track
250 // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
251}
252
253//----------------------------------------------------------------
254// methods for initialize
255inline
256void G4VParticleChange::InitializeStepInVolumeFlags(const G4Track& track)
257{
258 const G4Step* aStep = track.GetStep();
259 theFirstStepInVolume = aStep-> IsFirstStepInVolume();
260 theLastStepInVolume = aStep-> IsLastStepInVolume();
261}
262
263inline void G4VParticleChange::InitializeSecondaries(const G4Track&)
264{
265 // clear secondaries
266 if (theNumberOfSecondaries>0) {
267#ifdef G4VERBOSE
268 if (verboseLevel>0) {
269 G4cerr << "G4VParticleChange::Initialize() Warning ";
270 G4cerr << "theListOfSecondaries is not empty " << G4endl;
271 G4cerr << "All objects in theListOfSecondaries are destroyed!" << G4endl;
272 }
273#endif
274 for (G4int index= 0; index<theNumberOfSecondaries; index++){
275 if ( (*theListOfSecondaries)[index] ){
276 delete (*theListOfSecondaries)[index] ;
277 }
278 }
279 }
280 theNumberOfSecondaries = 0;
281}
282
283//----------------------------------------------------------------
284// methods for handling secondaries
285//
286
287inline void G4VParticleChange::SetNumberOfSecondaries(G4int totSecondaries)
288{
289 // check if tracks still exist in theListOfSecondaries
290 if (theNumberOfSecondaries>0) {
291#ifdef G4VERBOSE
292 if (verboseLevel>0) {
293 G4cerr << "G4VParticleChange::SetNumberOfSecondaries() Warning ";
294 G4cerr << "theListOfSecondaries is not empty ";
295 }
296#endif
297 for (G4int index= 0; index<theNumberOfSecondaries; index++){
298 if ( (*theListOfSecondaries)[index] ){
299 delete (*theListOfSecondaries)[index] ;
300 }
301 }
302 }
303 theNumberOfSecondaries = 0;
304 theSizeOftheListOfSecondaries = totSecondaries;
305
306 // Initialize ListOfSecondaries
307 theListOfSecondaries->Initialize(totSecondaries);
308}
309
310inline void G4VParticleChange::Initialize(const G4Track& track)
311{
312 InitializeStatusChange(track);
313 InitializeLocalEnergyDeposit(track);
314 InitializeSteppingControl(track);
315 InitializeTrueStepLength(track);
316 InitializeSecondaries(track);
317 InitializeParentWeight(track);
318 InitializeStepInVolumeFlags(track);
319}
320
321inline
322 void G4VParticleChange::ClearDebugFlag()
323{
324 debugFlag = false;
325}
326
327inline
328 void G4VParticleChange::SetDebugFlag()
329{
330 debugFlag = true;
331}
332
333inline
334 G4bool G4VParticleChange::GetDebugFlag() const
335{
336 return debugFlag;
337}
338
339inline
340 void G4VParticleChange::SetSecondaryWeightByProcess(G4bool flag)
341{
342 fSetSecondaryWeightByProcess = flag;
343}
344
345inline
346 G4bool G4VParticleChange::IsSecondaryWeightSetByProcess() const
347{
348 return fSetSecondaryWeightByProcess;
349}
350
351inline
352 void G4VParticleChange::ProposeParentWeight(G4double w)
353{
354 theParentWeight = w;
355}
356
357inline
358 void G4VParticleChange::SetParentWeightByProcess(G4bool flag)
359{
360 fSetParentWeightByProcess = flag;
361}
362
363inline
364 G4bool G4VParticleChange::IsParentWeightSetByProcess() const
365{
366 return fSetParentWeightByProcess;
367}
368
Note: See TracBrowser for help on using the repository browser.