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
G4ITStepProcessor2.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//
27// Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
28//
29// History:
30// -----------
31// 10 Oct 2011 M.Karamitros created
32//
33// -------------------------------------------------------------------
34
35#include "G4ITStepProcessor.hh"
36#include "G4LossTableManager.hh"
37#include "G4EnergyLossTables.hh"
38#include "G4ProductionCuts.hh"
40#include "G4VITProcess.hh"
42#include "G4IT.hh"
44#include "G4ITTransportation.hh"
45
46#include "G4ITNavigator.hh" // Include from 'geometry'
47
50
51#include "G4ITTrackHolder.hh"
52#include "G4ITReaction.hh"
53
55
56//#define DEBUG_MEM 1
57
58#ifdef DEBUG_MEM
59#include "G4MemStat.hh"
60using namespace G4MemStat;
62#endif
63
65{
66 // Now Store the secondaries from ParticleChange to SecondaryList
67 G4Track* tempSecondaryTrack;
68
69 for(G4int DSecLoop = 0; DSecLoop < fpParticleChange->GetNumberOfSecondaries();
70 DSecLoop++)
71 {
72 tempSecondaryTrack = fpParticleChange->GetSecondary(DSecLoop);
73
74 if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
75 {
76 ApplyProductionCut(tempSecondaryTrack);
77 }
78
79 // Set parentID
80 tempSecondaryTrack->SetParentID(fpTrack->GetTrackID());
81
82 // Set the process pointer which created this track
83 tempSecondaryTrack->SetCreatorProcess(fpCurrentProcess);
84
85 // If this 2ndry particle has 'zero' kinetic energy, make sure
86 // it invokes a rest process at the beginning of the tracking
87 if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN)
88 {
89 G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
90 if (pm->GetAtRestProcessVector()->entries()>0)
91 {
92 tempSecondaryTrack->SetTrackStatus( fStopButAlive );
93 fpSecondary->push_back( tempSecondaryTrack );
94 fN2ndariesAtRestDoIt++;
95 }
96 else
97 {
98 delete tempSecondaryTrack;
99 }
100 }
101 else
102 {
103 fpSecondary->push_back( tempSecondaryTrack );
104 counter++;
105 }
106 } //end of loop on secondary
107}
108
109//_________________________________________________________________________
110
111void G4ITStepProcessor::DoIt(double timeStep)
112
113// Call the process having the min step length or just propagate the track on the given time step
114
115// If the track is "leading the step" (ie one of its process has been selected
116// as the one having the minimum time step over all tracks and processes),
117// it will undergo its selected processes. Otherwise, it will just propagate the track
118// on the given time step.
119
120{
121 if(fpVerbose) fpVerbose->DoItStarted();
122
123 G4TrackManyList* mainList = fpTrackContainer->GetMainList();
124 G4TrackManyList::iterator it = mainList->end();
125 it--;
126 std::size_t initialSize = mainList->size();
127
128// G4cout << "initialSize = " << initialSize << G4endl;
129
130 for(std::size_t i = 0 ; i < initialSize ; ++i)
131 {
132
133// G4cout << "i = " << i << G4endl;
134
135 G4Track* track = *it;
136 if (!track)
137 {
138 G4ExceptionDescription exceptionDescription;
139 exceptionDescription << "No track was pop back the main track list.";
140 G4Exception("G4ITStepProcessor::DoIt", "NO_TRACK",
141 FatalException, exceptionDescription);
142 }
143 // G4TrackManyList::iterator next_it (it);
144 // next_it--;
145 // it = next_it;
146
147 it--;
148 // Must be called before EndTracking(track)
149 // Otherwise the iterator will point to the list of killed tracks
150
151 if(track->GetTrackStatus() == fStopAndKill)
152 {
153 fpTrackingManager->EndTracking(track);
154// G4cout << GetIT(track)->GetName() << G4endl;
155// G4cout << " ************************ CONTINUE ********************" << G4endl;
156 continue;
157 }
158
159#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DOIT)
160 MemStat mem_first, mem_second, mem_diff;
161 mem_first = MemoryUsage();
162#endif
163
164 Stepping(track, timeStep);
165
166#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DOIT)
167 MemStat mem_intermediaire = MemoryUsage();
168 mem_diff = mem_intermediaire-mem_first;
169 G4cout << "\t\t >> || MEM || In DoIT with track "
170 << track->GetTrackID() << ", diff is : " << mem_diff << G4endl;
171#endif
172
174
175#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DOIT)
176 mem_second = MemoryUsage();
177 mem_diff = mem_second-mem_first;
178 G4cout << "\t >> || MEM || In DoIT with track "
179 << track->GetTrackID()
180 << ", diff is : " << mem_diff << G4endl;
181#endif
182 }
183
184
185 fpTrackContainer->MergeSecondariesWithMainList();
186 fpTrackContainer->KillTracks(); // (18-06-15 : MK) Remove it ?
187 fLeadingTracks.Reset();
188}
189
190//_________________________________________________________________________
191
193{
194 if (!fpTrack)
195 {
197 return;
198 }
199
200 G4TrackStatus status = fpTrack->GetTrackStatus();
201
202 switch (status)
203 {
204 case fAlive:
205 case fStopButAlive:
206 case fSuspend:
208 default:
210 break;
211
212 case fStopAndKill:
215// G4TrackList::Pop(fpTrack);
216 fpTrackingManager->EndTracking(fpTrack);
217// fTrackContainer->PushToKill(fpTrack);
218 break;
219
222 if (fpSecondary)
223 {
224 for (std::size_t i = 0; i < fpSecondary->size(); ++i)
225 {
226 delete (*fpSecondary)[i];
227 }
228 fpSecondary->clear();
229 }
230// G4TrackList::Pop(fpTrack);
231 fpTrackingManager->EndTracking(fpTrack);
232// fTrackContainer->PushToKill(fpTrack);
233 break;
234 }
235
237}
238
239//_________________________________________________________________________
240
242{
243 if (!fpSecondary || fpSecondary->empty())
244 {
245 // DEBUG
246 // G4cout << "NO SECONDARIES !!! " << G4endl;
247 return;
248 }
249
250 // DEBUG
251 // G4cout << "There are secondaries : "<< secondaries -> size() << G4endl ;
252
253 G4TrackVector::iterator secondaries_i = fpSecondary->begin();
254
255 for (; secondaries_i != fpSecondary->end(); ++secondaries_i)
256 {
257 G4Track* secondary = *secondaries_i;
258 fpTrackContainer->_PushTrack(secondary);
259 }
260}
261
262//______________________________________________________________________________
263
264void G4ITStepProcessor::Stepping(G4Track* track, const double & timeStep)
265{
266
267#ifdef DEBUG_MEM
268 MemStat mem_first, mem_second, mem_diff;
269#endif
270
271#ifdef DEBUG_MEM
272 mem_first = MemoryUsage();
273#endif
274
276
277#ifdef DEBUG_MEM
278 MemStat mem_intermediaire = MemoryUsage();
279 mem_diff = mem_intermediaire-mem_first;
280 G4cout << "\t\t\t >> || MEM || After CleanProcessor " << track->GetTrackID() << ", diff is : " << mem_diff << G4endl;
281#endif
282
283 if(track == 0) return; // maybe put an exception here
284 fTimeStep = timeStep;
285 SetTrack(track);
286 DoStepping();
287}
288//______________________________________________________________________________
289
290// ************************************************************************
291// Stepping
292// ************************************************************************
294{
295 SetupMembers();
296
297#ifdef DEBUG_MEM
298 MemStat mem_first, mem_second, mem_diff;
299#endif
300
301#ifdef DEBUG_MEM
302 mem_first = MemoryUsage();
303#endif
304
305#ifdef G4VERBOSE
306 if(fpVerbose) fpVerbose->PreStepVerbose(fpTrack);
307#endif
308
309 if(!fpProcessInfo)
310 {
311 G4ExceptionDescription exceptionDescription;
312 exceptionDescription << "No process info found for particle :"
313 << fpTrack->GetDefinition()->GetParticleName();
314 G4Exception("G4ITStepProcessor::DoStepping",
315 "ITStepProcessor0012",
317 exceptionDescription);
318 return;
319 }
320// else if(fpTrack->GetTrackStatus() == fStopAndKill)
321// {
322// fpState->fStepStatus = fUndefined;
323// return;
324// }
325
326 if(fpProcessInfo->MAXofPostStepLoops == 0 &&
327 fpProcessInfo->MAXofAlongStepLoops == 0
328 && fpProcessInfo->MAXofAtRestLoops == 0)
329 {/*
330 G4ExceptionDescription exceptionDescription;
331 exceptionDescription << "No process was found for particle :"
332 << fpTrack->GetDefinition()->GetParticleName();
333 G4Exception("G4ITStepProcessor::DoStepping",
334 "ITStepProcessorNoProcess",
335 JustWarning,
336 exceptionDescription);
337
338 fpTrack->SetTrackStatus(fStopAndKill);
339 fpState->fStepStatus = fUndefined;*/
340 return;
341 }
342
343 //--------
344 // Prelude
345 //--------
346#ifdef G4VERBOSE
347 // !!!!! Verbose
348 if(fpVerbose) fpVerbose->NewStep();
349#endif
350
351 //---------------------------------
352 // AtRestStep, AlongStep and PostStep Processes
353 //---------------------------------
354
355 fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()->GetNavigatorState());
356// fpNavigator->ResetHierarchyAndLocate( fpTrack->GetPosition(),
357// fpTrack->GetMomentumDirection(),
358// *((G4TouchableHistory*)fpTrack->GetTouchableHandle()()) );
359// fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()->GetNavigatorState());
360 // We reset the navigator state before checking for AtRest
361 // in case a AtRest processe would use a navigator info
362
363#ifdef DEBUG_MEM
364 MemStat mem_intermediaire = MemoryUsage();
365 mem_diff = mem_intermediaire-mem_first;
366 G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After dealing with navigator with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
367#endif
368
369 if(fpTrack->GetTrackStatus() == fStopButAlive)
370 {
371 if(fpProcessInfo->MAXofAtRestLoops > 0 && fpProcessInfo->fpAtRestDoItVector
372 != 0) // second condition to make coverity happy
373 {
374 //-----------------
375 // AtRestStepDoIt
376 //-----------------
378 fpState->fStepStatus = fAtRestDoItProc;
379 fpStep->GetPostStepPoint()->SetStepStatus(fpState->fStepStatus);
380
381#ifdef G4VERBOSE
382 // !!!!! Verbose
383 if(fpVerbose) fpVerbose->AtRestDoItInvoked();
384#endif
385
386 }
387 // Make sure the track is killed
388 // fpTrack->SetTrackStatus(fStopAndKill);
389 }
390 else // if(fTimeStep > 0.) // Bye, because PostStepIL can return 0 => time =0
391 {
392 if(fpITrack == 0)
393 {
394 G4ExceptionDescription exceptionDescription;
395 exceptionDescription << " !!! TrackID : " << fpTrack->GetTrackID()
396 << G4endl<< " !!! Track status : "<< fpTrack->GetTrackStatus() << G4endl
397 << " !!! Particle Name : "<< fpTrack -> GetDefinition() -> GetParticleName() << G4endl
398 << "No G4ITStepProcessor::fpITrack found" << G4endl;
399
400 G4Exception("G4ITStepProcessor::DoStepping",
401 "ITStepProcessor0013",
403 exceptionDescription);
404 return; // to make coverity happy
405 }
406
407 if(fpITrack->GetTrackingInfo()->IsLeadingStep() == false)
408 {
409 // In case the track has NOT the minimum step length
410 // Given the final step time, the transportation
411 // will compute the final position of the particle
413 fpStep->GetPostStepPoint()->SetProcessDefinedStep(fpTransportation);
415 }
416
417#ifdef DEBUG_MEM
418 mem_intermediaire = MemoryUsage();
419 mem_diff = mem_intermediaire-mem_first;
420 G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After FindTransportationStep() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
421#endif
422
423 // Store the Step length (geometrical length) to G4Step and G4Track
424 fpTrack->SetStepLength(fpState->fPhysicalStep);
425 fpStep->SetStepLength(fpState->fPhysicalStep);
426
427 G4double GeomStepLength = fpState->fPhysicalStep;
428
429 // Store StepStatus to PostStepPoint
430 fpStep->GetPostStepPoint()->SetStepStatus(fpState->fStepStatus);
431
432 // Invoke AlongStepDoIt
434
435#ifdef DEBUG_MEM
436 mem_intermediaire = MemoryUsage();
437 mem_diff = mem_intermediaire-mem_first;
438 G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After InvokeAlongStepDoItProcs() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
439#endif
440
441#ifdef G4VERBOSE
442 // !!!!! Verbose
443 if(fpVerbose) fpVerbose->AlongStepDoItAllDone();
444#endif
445
446 // Update track by taking into account all changes by AlongStepDoIt
447 // fpStep->UpdateTrack(); // done in InvokeAlongStepDoItProcs
448
449 // Update safety after invocation of all AlongStepDoIts
450 fpState->fEndpointSafOrigin = fpPostStepPoint->GetPosition();
451
452 fpState->fEndpointSafety =
453 std::max(fpState->fProposedSafety - GeomStepLength, kCarTolerance);
454
455 fpStep->GetPostStepPoint()->SetSafety(fpState->fEndpointSafety);
456
457 if(GetIT(fpTrack)->GetTrackingInfo()->IsLeadingStep())
458 {
459 // Invoke PostStepDoIt including G4ITTransportation::PSDI
461
462#ifdef DEBUG_MEM
463 mem_intermediaire = MemoryUsage();
464 mem_diff = mem_intermediaire-mem_first;
465 G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After InvokePostStepDoItProcs() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
466#endif
467#ifdef G4VERBOSE
468 // !!!!! Verbose
469 if(fpVerbose) fpVerbose->StepInfoForLeadingTrack();
470#endif
471 }
472 else
473 {
474 // Only invoke transportation and all other forced processes
476 fpStep->GetPostStepPoint()->SetProcessDefinedStep(fpTransportation);
477
478#ifdef DEBUG_MEM
479 mem_intermediaire = MemoryUsage();
480 mem_diff = mem_intermediaire-mem_first;
481 G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After InvokeTransportationProc() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
482#endif
483 }
484
485#ifdef G4VERBOSE
486 // !!!!! Verbose
487 if(fpVerbose) fpVerbose->PostStepDoItAllDone();
488#endif
489 }
490
491 fpNavigator->ResetNavigatorState();
492
493#ifdef DEBUG_MEM
494 mem_intermediaire = MemoryUsage();
495 mem_diff = mem_intermediaire-mem_first;
496 G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After fpNavigator->SetNavigatorState with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
497#endif
498
499 //-------
500 // Finale
501 //-------
502
503 // Update 'TrackLength' and remeber the Step length of the current Step
504 fpTrack->AddTrackLength(fpStep->GetStepLength());
506
507//#ifdef G4VERBOSE
508// // !!!!! Verbose
509// if(fpVerbose) fpVerbose->StepInfo();
510//#endif
511
512#ifdef G4VERBOSE
513 if(fpVerbose) fpVerbose->PostStepVerbose(fpTrack);
514#endif
515
516// G4cout << " G4ITStepProcessor::DoStepping -- " <<fpTrack->GetTrackID() << " tps = " << fpTrack->GetGlobalTime() << G4endl;
517
518 // Send G4Step information to Hit/Dig if the volume is sensitive
519 /***
520 fpCurrentVolume = fpStep->GetPreStepPoint()->GetPhysicalVolume();
521 StepControlFlag = fpStep->GetControlFlag();
522
523 if( fpCurrentVolume != 0 && StepControlFlag != AvoidHitInvocation)
524 {
525 fpSensitive = fpStep->GetPreStepPoint()->
526 GetSensitiveDetector();
527 if( fpSensitive != 0 )
528 {
529 fpSensitive->Hit(fpStep);
530 }
531 }
532
533 User intervention process.
534 if( fpUserSteppingAction != 0 )
535 {
536 fpUserSteppingAction->UserSteppingAction(fpStep);
537 }
538 G4UserSteppingAction* regionalAction
539 = fpStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetRegion()
540 ->GetRegionalSteppingAction();
541 if( regionalAction ) regionalAction->UserSteppingAction(fpStep);
542 ***/
543 fpTrackingManager->AppendStep(fpTrack, fpStep);
544 // Stepping process finish. Return the value of the StepStatus.
545
546#ifdef DEBUG_MEM
547 MemStat mem_intermediaire = MemoryUsage();
548 mem_diff = mem_intermediaire-mem_first;
549 G4cout << "\t\t\t >> || MEM || End of DoStepping() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
550#endif
551
552 // return fpState->fStepStatus;
553}
554
555//______________________________________________________________________________
556
557// ************************************************************************
558// AtRestDoIt
559// ************************************************************************
560
562{
563 fpStep->SetStepLength(0.); //the particle has stopped
564 fpTrack->SetStepLength(0.);
565
566 G4SelectedAtRestDoItVector& selectedAtRestDoItVector =
568
569 // invoke selected process
570 for(std::size_t np = 0; np < fpProcessInfo->MAXofAtRestLoops; ++np)
571 {
572 //
573 // Note: DoItVector has inverse order against GetPhysIntVector
574 // and SelectedAtRestDoItVector.
575 //
576 if(selectedAtRestDoItVector[fpProcessInfo->MAXofAtRestLoops - np - 1] != InActivated)
577 {
578 fpCurrentProcess =
579 (G4VITProcess*) (*fpProcessInfo->fpAtRestDoItVector)[(G4int)np];
580
581// G4cout << " Invoke : "
582// << fpCurrentProcess->GetProcessName()
583// << G4endl;
584
585 // if(fpVerbose)
586 // {
587 // fpVerbose->AtRestDoItOneByOne();
588 // }
589
590 fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess
591 ->GetProcessID()));
592 fpParticleChange = fpCurrentProcess->AtRestDoIt(*fpTrack, *fpStep);
593 fpCurrentProcess->ResetProcessState();
594
595 // Set the current process as a process which defined this Step length
596 fpStep->GetPostStepPoint()->SetProcessDefinedStep(fpCurrentProcess);
597
598 // Update Step
599 fpParticleChange->UpdateStepForAtRest(fpStep);
600
601 // Now Store the secondaries from ParticleChange to SecondaryList
602 DealWithSecondaries(fN2ndariesAtRestDoIt);
603
604 // Set the track status according to what the process defined
605 // if kinetic energy >0, otherwise set fStopButAlive
606 fpTrack->SetTrackStatus(fpParticleChange->GetTrackStatus());
607
608 // clear ParticleChange
609 fpParticleChange->Clear();
610
611 } //if(fSelectedAtRestDoItVector[np] != InActivated){
612 } //for(std::size_t np=0; np < MAXofAtRestLoops; ++np){
613 fpStep->UpdateTrack();
614
615 // Modification par rapport au transport standard :
616 // fStopAndKill doit etre propose par le modele
617 // sinon d autres processus AtRest seront appeles
618 // au pas suivant
619 // fpTrack->SetTrackStatus(fStopAndKill);
620}
621
622//______________________________________________________________________________
623
624// ************************************************************************
625// AlongStepDoIt
626// ************************************************************************
627
629{
630
631#ifdef DEBUG_MEM
632 MemStat mem_first, mem_second, mem_diff;
633#endif
634
635#ifdef DEBUG_MEM
636 mem_first = MemoryUsage();
637#endif
638
639 // If the current Step is defined by a 'ExclusivelyForced'
640 // PostStepDoIt, then don't invoke any AlongStepDoIt
641 if(fpState->fStepStatus == fExclusivelyForcedProc)
642 {
643 return; // Take note 'return' at here !!!
644 }
645
646 // Invoke the all active continuous processes
647 for(std::size_t ci = 0; ci < fpProcessInfo->MAXofAlongStepLoops; ++ci)
648 {
649 fpCurrentProcess =
650 (G4VITProcess*) (*fpProcessInfo->fpAlongStepDoItVector)[(G4int)ci];
651 if(fpCurrentProcess == 0) continue;
652 // NULL means the process is inactivated by a user on fly.
653
654 fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess
655 ->GetProcessID()));
656 fpParticleChange = fpCurrentProcess->AlongStepDoIt(*fpTrack, *fpStep);
657
658#ifdef DEBUG_MEM
659 MemStat mem_intermediaire = MemoryUsage();
660 mem_diff = mem_intermediaire-mem_first;
661 G4cout << "\t\t\t >> || MEM || After calling AlongStepDoIt for " << fpCurrentProcess->GetProcessName() << " and track "<< fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
662#endif
663
664// fpCurrentProcess->SetProcessState(0);
665 fpCurrentProcess->ResetProcessState();
666 // Update the PostStepPoint of Step according to ParticleChange
667
668 fpParticleChange->UpdateStepForAlongStep(fpStep);
669
670#ifdef G4VERBOSE
671 // !!!!! Verbose
672 if(fpVerbose) fpVerbose->AlongStepDoItOneByOne();
673#endif
674
675 // Now Store the secondaries from ParticleChange to SecondaryList
676 DealWithSecondaries(fN2ndariesAlongStepDoIt);
677
678 // Set the track status according to what the process defined
679 // if kinetic energy >0, otherwise set fStopButAlive
680 fpTrack->SetTrackStatus(fpParticleChange->GetTrackStatus());
681
682 // clear ParticleChange
683 fpParticleChange->Clear();
684 }
685
686#ifdef DEBUG_MEM
687 MemStat mem_intermediaire = MemoryUsage();
688 mem_diff = mem_intermediaire-mem_first;
689 G4cout << "\t\t\t >> || MEM || After looping on processes with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
690#endif
691
692 fpStep->UpdateTrack();
693
694 G4TrackStatus fNewStatus = fpTrack->GetTrackStatus();
695
696 if(fNewStatus == fAlive && fpTrack->GetKineticEnergy() <= DBL_MIN)
697 {
698 // G4cout << "G4ITStepProcessor::InvokeAlongStepDoItProcs : Track will be killed" << G4endl;
699 if(fpProcessInfo->MAXofAtRestLoops>0) fNewStatus = fStopButAlive;
700 else fNewStatus = fStopAndKill;
701 fpTrack->SetTrackStatus( fNewStatus );
702 }
703
704}
705
706//______________________________________________________________________________
707
708// ************************************************************************
709// PostStepDoIt
710// ************************************************************************
711
713{
714 std::size_t _MAXofPostStepLoops = fpProcessInfo->MAXofPostStepLoops;
715 G4SelectedPostStepDoItVector& selectedPostStepDoItVector = fpState
717 G4StepStatus& stepStatus = fpState->fStepStatus;
718
719 // Invoke the specified discrete processes
720 for(std::size_t np = 0; np < _MAXofPostStepLoops; ++np)
721 {
722 //
723 // Note: DoItVector has inverse order against GetPhysIntVector
724 // and SelectedPostStepDoItVector.
725 //
726 G4int Cond = selectedPostStepDoItVector[_MAXofPostStepLoops
727 - np - 1];
728 if(Cond != InActivated)
729 {
730 if(((Cond == NotForced) && (stepStatus == fPostStepDoItProc)) ||
731 ((Cond == Forced) && (stepStatus != fExclusivelyForcedProc))
732 ||
733 // ((Cond == Conditionally) && (stepStatus == fAlongStepDoItProc)) ||
734 ((Cond == ExclusivelyForced) && (stepStatus == fExclusivelyForcedProc))
735 || ((Cond == StronglyForced)))
736 {
737
738 InvokePSDIP(np);
739 }
740 } //if(*fSelectedPostStepDoItVector(np)........
741
742 // Exit from PostStepLoop if the track has been killed,
743 // but extra treatment for processes with Strongly Forced flag
744 if(fpTrack->GetTrackStatus() == fStopAndKill)
745 {
746 for(std::size_t np1 = np + 1; np1 < _MAXofPostStepLoops; ++np1)
747 {
748 G4int Cond2 = selectedPostStepDoItVector[_MAXofPostStepLoops
749 - np1 - 1];
750 if(Cond2 == StronglyForced)
751 {
752 InvokePSDIP(np1);
753 }
754 }
755 break;
756 }
757 } //for(std::size_t np=0; np < MAXofPostStepLoops; ++np){
758}
759
760//______________________________________________________________________________
761
763{
764 fpCurrentProcess = (G4VITProcess*) (*fpProcessInfo->fpPostStepDoItVector)[(G4int)np];
765
766 fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess
767 ->GetProcessID()));
768 fpParticleChange = fpCurrentProcess->PostStepDoIt(*fpTrack, *fpStep);
769// fpCurrentProcess->SetProcessState(0);
770 fpCurrentProcess->ResetProcessState();
771
772 // Update PostStepPoint of Step according to ParticleChange
773 fpParticleChange->UpdateStepForPostStep(fpStep);
774
775#ifdef G4VERBOSE
776 // !!!!! Verbose
777 if(fpVerbose) fpVerbose->PostStepDoItOneByOne();
778#endif
779
780 // Update G4Track according to ParticleChange after each PostStepDoIt
781 fpStep->UpdateTrack();
782
783 // Update safety after each invocation of PostStepDoIts
785
786 // Now Store the secondaries from ParticleChange to SecondaryList
787 DealWithSecondaries(fN2ndariesPostStepDoIt);
788
789 // Set the track status according to what the process defined
790 fpTrack->SetTrackStatus(fpParticleChange->GetTrackStatus());
791
792 // clear ParticleChange
793 fpParticleChange->Clear();
794}
795
796//______________________________________________________________________________
797
798// ************************************************************************
799// Transport on a given time
800// ************************************************************************
801
803{
804 double physicalStep(0.);
805
806 fpTransportation = fpProcessInfo->fpTransportation;
807 // dynamic_cast<G4ITTransportation*>((fpProcessInfo->fpAlongStepGetPhysIntVector)[MAXofAlongStepLoops-1]);
808
809 if(!fpTrack)
810 {
811 G4ExceptionDescription exceptionDescription;
812 exceptionDescription << "No G4ITStepProcessor::fpTrack found";
813 G4Exception("G4ITStepProcessor::FindTransportationStep",
814 "ITStepProcessor0013",
816 exceptionDescription);
817 return;
818
819 }
820 if(!fpITrack)
821 {
822 G4ExceptionDescription exceptionDescription;
823 exceptionDescription << "No G4ITStepProcessor::fITrack";
824 G4Exception("G4ITStepProcessor::FindTransportationStep",
825 "ITStepProcessor0014",
827 exceptionDescription);
828 return;
829 }
830 if(!(fpITrack->GetTrack()))
831 {
832 G4ExceptionDescription exceptionDescription;
833 exceptionDescription << "No G4ITStepProcessor::fITrack->GetTrack()";
834 G4Exception("G4ITStepProcessor::FindTransportationStep",
835 "ITStepProcessor0015",
837 exceptionDescription);
838 return;
839 }
840
841 if(fpTransportation)
842 {
843 fpTransportation->SetProcessState(fpTrackingInfo->GetProcessState(fpTransportation
844 ->GetProcessID()));
845 fpTransportation->ComputeStep(*fpTrack, *fpStep, fTimeStep, physicalStep);
846
847// fpTransportation->SetProcessState(0);
848 fpTransportation->ResetProcessState();
849 }
850
851 if(physicalStep >= DBL_MAX)
852 {
853// G4cout << "---- 2) Setting stop and kill for " << GetIT(fpTrack)->GetName() << G4endl;
854 fpTrack -> SetTrackStatus(fStopAndKill);
855 return;
856 }
857
858 fpState->fPhysicalStep = physicalStep;
859}
860
861//______________________________________________________________________________
862
864{
865 std::size_t _MAXofPostStepLoops = fpProcessInfo->MAXofPostStepLoops;
866 G4SelectedPostStepDoItVector& selectedPostStepDoItVector = fpState
868 G4StepStatus& stepStatus = fpState->fStepStatus;
869
870 // Invoke the specified discrete processes
871 for(std::size_t np = 0; np < _MAXofPostStepLoops; ++np)
872 {
873 //
874 // Note: DoItVector has inverse order against GetPhysIntVector
875 // and SelectedPostStepDoItVector.
876 //
877 G4int Cond = selectedPostStepDoItVector[_MAXofPostStepLoops - np - 1];
878 if(Cond != InActivated)
879 {
880 if(((Cond == Forced) && (stepStatus != fExclusivelyForcedProc)) ||
881 // ((Cond == Conditionally) && (stepStatus == fAlongStepDoItProc)) ||
882 ((Cond == ExclusivelyForced) && (stepStatus == fExclusivelyForcedProc))
883 || ((Cond == StronglyForced)))
884 {
885
886 InvokePSDIP(np);
887 }
888 } //if(Cond != InActivated)
889
890 // Exit from PostStepLoop if the track has been killed,
891 // but extra treatment for processes with Strongly Forced flag
892 if(fpTrack->GetTrackStatus() == fStopAndKill)
893 {
894 for(std::size_t np1 = np + 1; np1 < _MAXofPostStepLoops; ++np1)
895 {
896 G4int Cond2 = selectedPostStepDoItVector[_MAXofPostStepLoops - np1 - 1];
897 if(Cond2 == StronglyForced)
898 {
899 InvokePSDIP(np1);
900 }
901 }
902 break;
903 }
904 }
905}
906
907//______________________________________________________________________________
908
909// ************************************************************************
910// Apply cuts
911// ************************************************************************
912
914{
915 G4bool tBelowCutEnergyAndSafety = false;
916 G4int tPtclIdx = G4ProductionCuts::GetIndex(aSecondary->GetDefinition());
917 if(tPtclIdx < 0)
918 {
919 return;
920 }
921 G4ProductionCutsTable* tCutsTbl =
923 G4int tCoupleIdx = tCutsTbl->GetCoupleIndex(fpPreStepPoint
924 ->GetMaterialCutsCouple());
925 G4double tProdThreshold =
926 (*(tCutsTbl->GetEnergyCutsVector(tPtclIdx)))[tCoupleIdx];
927 if(aSecondary->GetKineticEnergy() < tProdThreshold)
928 {
929 tBelowCutEnergyAndSafety = true;
930 if(std::abs(aSecondary->GetDynamicParticle()->GetCharge()) > DBL_MIN)
931 {
932 G4double currentRange
934 aSecondary->GetKineticEnergy(),
935 fpPreStepPoint->GetMaterialCutsCouple());
936 tBelowCutEnergyAndSafety = (currentRange < CalculateSafety() );
937 }
938 }
939
940 if(tBelowCutEnergyAndSafety)
941 {
942 if(!(aSecondary->IsGoodForTracking()))
943 {
944 // Add kinetic energy to the total energy deposit
945 fpStep->AddTotalEnergyDeposit(aSecondary->GetKineticEnergy());
946 aSecondary->SetKineticEnergy(0.0);
947 }
948 }
949}
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ InActivated
@ StronglyForced
@ NotForced
@ ExclusivelyForced
@ Forced
class std::vector< int, std::allocator< int > > G4SelectedPostStepDoItVector
class std::vector< int, std::allocator< int > > G4SelectedAtRestDoItVector
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:48
G4StepStatus
Definition: G4StepStatus.hh:40
@ fPostStepDoItProc
Definition: G4StepStatus.hh:49
@ fAtRestDoItProc
Definition: G4StepStatus.hh:45
@ fExclusivelyForcedProc
Definition: G4StepStatus.hh:53
G4TrackStatus
@ fKillTrackAndSecondaries
@ fSuspend
@ fAlive
@ fStopAndKill
@ fStopButAlive
@ fPostponeToNextEvent
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetCharge() const
void RemoveReactionSet(G4Track *track)
static G4ITReactionSet * Instance()
G4ThreeVector fEndpointSafOrigin
G4SelectedPostStepDoItVector fSelectedPostStepDoItVector
G4SelectedAtRestDoItVector fSelectedAtRestDoItVector
void DealWithSecondaries(G4int &)
void SetTrack(G4Track *)
void Stepping(G4Track *, const double &)
void DoIt(double timeStep)
void ApplyProductionCut(G4Track *)
void InvokePSDIP(std::size_t)
void _PushTrack(G4Track *track)
G4TrackList * GetMainList(Key)
void MergeSecondariesWithMainList()
void AppendStep(G4Track *track, G4Step *step)
void EndTracking(G4Track *)
virtual void ComputeStep(const G4Track &, const G4Step &, const double timeStep, double &spaceStep)
G4TrackingInformation * GetTrackingInfo()
Definition: G4IT.hh:143
G4Track * GetTrack()
Definition: G4IT.hh:218
static G4LossTableManager * Instance()
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
size_t size() const
G4ProcessManager * GetProcessManager() const
G4bool GetApplyCutsFlag() const
const G4String & GetParticleName() const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t entries() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4int GetCoupleIndex(const G4MaterialCutsCouple *aCouple) const
static G4int GetIndex(const G4String &name)
void SetSafety(const G4double aValue)
void SetStepStatus(const G4StepStatus aValue)
void SetProcessDefinedStep(const G4VProcess *aValue)
const G4ThreeVector & GetPosition() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void UpdateTrack()
void SetStepLength(G4double value)
void AddTotalEnergyDeposit(G4double value)
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
G4TrackStatus GetTrackStatus() const
void SetTrackStatus(const G4TrackStatus aTrackStatus)
void SetStepLength(G4double value)
G4int GetTrackID() const
G4ParticleDefinition * GetDefinition() const
void AddTrackLength(const G4double aValue)
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
void IncrementCurrentStepNumber()
void SetKineticEnergy(const G4double aValue)
G4bool IsGoodForTracking() const
void SetParentID(const G4int aValue)
void SetCreatorProcess(const G4VProcess *aValue)
G4shared_ptr< G4ProcessState_Lock > GetProcessState(size_t index)
G4ITNavigatorState_Lock * GetNavigatorState() const
size_t GetProcessID() const
void ResetProcessState()
void SetProcessState(G4shared_ptr< G4ProcessState_Lock > aProcInfo)
virtual void StepInfoForLeadingTrack()=0
virtual void PostStepDoItAllDone()=0
virtual void DoItStarted()=0
virtual void AlongStepDoItOneByOne()=0
virtual void PreStepVerbose(G4Track *)=0
virtual void PostStepDoItOneByOne()=0
virtual void PostStepVerbose(G4Track *)=0
virtual void AtRestDoItInvoked()=0
virtual void NewStep()=0
virtual void AlongStepDoItAllDone()=0
virtual G4Step * UpdateStepForAlongStep(G4Step *Step)
virtual G4Step * UpdateStepForAtRest(G4Step *Step)
virtual G4Step * UpdateStepForPostStep(G4Step *Step)
G4int GetNumberOfSecondaries() const
G4Track * GetSecondary(G4int anIndex) const
G4TrackStatus GetTrackStatus() const
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)=0
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)=0
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
G4ProcessVector * fpPostStepDoItVector
std::size_t MAXofAlongStepLoops
std::size_t MAXofAtRestLoops
G4ProcessVector * fpAlongStepDoItVector
G4ITTransportation * fpTransportation
std::size_t MAXofPostStepLoops
G4ProcessVector * fpAtRestDoItVector
#define DBL_MIN
Definition: templates.hh:54
#define DBL_MAX
Definition: templates.hh:62