Geant4-11
G4ITStepProcessor.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 "G4UImanager.hh"
37#include "G4ForceCondition.hh"
38#include "G4GPILSelection.hh"
40// #include "G4VSensitiveDetector.hh" // Include from 'hits/digi'
42#include "G4ParticleTable.hh"
45#include "G4IT.hh"
46#include "G4ITNavigator.hh" // Include from 'geometry'
47
48#include "G4VITProcess.hh"
49#include "G4VProcess.hh"
50#include "G4ITTransportation.hh"
51
52#include "G4ITTrackHolder.hh"
53
57
59
60#include <iomanip> // Include from 'system'
61#include <vector> // Include from 'system'
62
63using namespace std;
64
65static const size_t SizeOfSelectedDoItVector = 100;
66
67template<typename T>
68 inline bool IsInf(T value)
69 {
70 return std::numeric_limits<T>::has_infinity
71 && value == std::numeric_limits<T>::infinity();
72 }
73
74//______________________________________________________________________________
75
77{
78 fpVerbose = 0;
79 // fpUserSteppingAction = 0 ;
82 fpNavigator = 0;
83 kCarTolerance = -1.;
84 fInitialized = false;
88
91}
92
93//______________________________________________________________________________
94
95//G4ITStepProcessor::
98 fSelectedAtRestDoItVector(G4VITProcess::GetMaxProcessIndex(), 0),
99 fSelectedPostStepDoItVector(G4VITProcess::GetMaxProcessIndex(), 0)
100{
101 fPhysicalStep = -1.;
102 fPreviousStepSize = -1.;
103
104 fSafety = -1.;
105 fProposedSafety = -1.;
106 fEndpointSafety = -1;
107
109
111}
112
113//______________________________________________________________________________
114
115//G4ITStepProcessor::
118 fSelectedAtRestDoItVector(right.fSelectedAtRestDoItVector),
119 fSelectedPostStepDoItVector(right.fSelectedPostStepDoItVector)
120{
123
124 fSafety = right.fSafety;
127
128 fStepStatus = right.fStepStatus;
129
131}
132
133//______________________________________________________________________________
134
135//G4ITStepProcessor::
137//G4ITStepProcessor::
139{
140 if(this == &right) return *this;
141
146
149
150 fSafety = right.fSafety;
153
154 fStepStatus = right.fStepStatus;
155
157 return *this;
158}
159
160//______________________________________________________________________________
161
162//G4ITStepProcessor::
164{
165 ;
166}
167
168//______________________________________________________________________________
169
171{
172 std::map<const G4ParticleDefinition*, ProcessGeneralInfo*>::iterator it;
173
174 for(it = fProcessGeneralInfoMap.begin(); it != fProcessGeneralInfoMap.end();
175 it++)
176 {
177 if(it->second)
178 {
179 delete it->second;
180 it->second = 0;
181 }
182 }
183
185}
186
187//______________________________________________________________________________
188
190{
191 fInitialized = false;
193 Initialize();
194}
195
196//______________________________________________________________________________
197
199{
201 if(fInitialized) return;
202 // ActiveOnlyITProcess();
203
205 ->GetNavigatorForTracking());
206
208 kCarTolerance = 0.5
210
211 if(fpVerbose == 0)
212 {
214
215 if(interactivity)
216 {
217 fpVerbose = interactivity->GetSteppingVerbose();
219 }
220 }
221
223
224 fInitialized = true;
225}
226//______________________________________________________________________________
227
229{
230 if(fpStep)
231 {
233 delete fpStep;
234 }
235
236 if(fpSecondary) delete fpSecondary;
238 //G4ITTransportationManager::DeleteInstance();
239
240 // if(fpUserSteppingAction) delete fpUserSteppingAction;
241}
242//______________________________________________________________________________
243// should not be used
245{
246 fpVerbose = rhs.fpVerbose;
248
249 // fpUserSteppingAction = 0 ;
251 fpNavigator = 0;
252 fInitialized = false;
253
255 fInitialized = false;
257
262}
263
264//______________________________________________________________________________
265
267{
269}
270
271//______________________________________________________________________________
272
274{
276}
277
278//______________________________________________________________________________
279
281{
282 if(this == &rhs) return *this; // handle self assignment
283 //assignment operator
284 return *this;
285}
286
287//______________________________________________________________________________
288
290{
291 // Method not used for the time being
292#ifdef debug
293 G4cout<<"G4ITStepProcessor::CloneProcesses: is called"<<G4endl;
294#endif
295
298 ->GetIterator();
299
300 theParticleIterator->reset();
301 // TODO : Ne faire la boucle que sur les IT **** !!!
302 while((*theParticleIterator)())
303 {
304 G4ParticleDefinition* particle = theParticleIterator->value();
305 G4ProcessManager* pm = particle->GetProcessManager();
306
307 if(!pm)
308 {
309 G4cerr << "ERROR - G4ITStepProcessor::GetProcessNumber()" << G4endl<< " ProcessManager is NULL for particle = "
310 << particle->GetParticleName() << ", PDG_code = "
311 << particle->GetPDGEncoding() << G4endl;
312 G4Exception("G4ITStepProcessor::GetProcessNumber()", "ITStepProcessor0001",
313 FatalException, "Process Manager is not found.");
314 return;
315 }
316
318 }
319}
320
321//______________________________________________________________________________
322
324{
325 // Method not used for the time being
326 G4ProcessVector* processVector = processManager->GetProcessList();
327
328 G4VITProcess* itProcess = 0;
329 for(std::size_t i = 0; i < processVector->size(); ++i)
330 {
331 G4VProcess* base_process = (*processVector)[i];
332 itProcess = dynamic_cast<G4VITProcess*>(base_process);
333
334 if(!itProcess)
335 {
336 processManager->SetProcessActivation(base_process, false);
337 }
338 }
339}
340
341//______________________________________________________________________________
342
345{
346
347#ifdef debug
348 G4cout<<"G4ITStepProcessor::GetProcessNumber: is called track"<<G4endl;
349#endif
350 if(!pm)
351 {
352 G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl<< " ProcessManager is NULL for particle = "
353 << particle->GetParticleName() << ", PDG_code = "
354 << particle->GetPDGEncoding() << G4endl;
355 G4Exception("G4SteppingManager::GetProcessNumber()", "ITStepProcessor0002",
356 FatalException, "Process Manager is not found.");
357 return;
358 }
359
360 std::map<const G4ParticleDefinition*, ProcessGeneralInfo*>::iterator it =
361 fProcessGeneralInfoMap.find(particle);
362 if(it != fProcessGeneralInfoMap.end())
363 {
364 G4Exception("G4SteppingManager::SetupGeneralProcessInfo()",
365 "ITStepProcessor0003",
366 FatalException, "Process info already registered.");
367 return;
368 }
369
370 // here used as temporary
372
373 // AtRestDoits
378#ifdef debug
379 G4cout << "G4ITStepProcessor::GetProcessNumber: #ofAtRest="
381#endif
382
383 // AlongStepDoits
390#ifdef debug
391 G4cout << "G4ITStepProcessor::GetProcessNumber:#ofAlongStp="
393#endif
394
395 // PostStepDoits
401#ifdef debug
402 G4cout << "G4ITStepProcessor::GetProcessNumber: #ofPostStep="
404#endif
405
406 if (SizeOfSelectedDoItVector<fpProcessInfo->MAXofAtRestLoops ||
407 SizeOfSelectedDoItVector<fpProcessInfo->MAXofAlongStepLoops ||
408 SizeOfSelectedDoItVector<fpProcessInfo->MAXofPostStepLoops )
409 {
410 G4cerr << "ERROR - G4ITStepProcessor::GetProcessNumber()" << G4endl
411 << " SizeOfSelectedDoItVector= " << SizeOfSelectedDoItVector
412 << " ; is smaller then one of MAXofAtRestLoops= "
414 << " or MAXofAlongStepLoops= " << fpProcessInfo->MAXofAlongStepLoops
415 << " or MAXofPostStepLoops= " << fpProcessInfo->MAXofPostStepLoops << G4endl;
416 G4Exception("G4ITStepProcessor::GetProcessNumber()",
417 "ITStepProcessor0004", FatalException,
418 "The array size is smaller than the actual No of processes.");
419 }
420
424 {
425 G4ExceptionDescription exceptionDescription;
426 exceptionDescription << "No DoIt process found ";
427 G4Exception("G4ITStepProcessor::DoStepping","ITStepProcessor0005",
428 FatalErrorInArgument,exceptionDescription);
429 return;
430 }
431
434 {
438
440 {
441 G4ExceptionDescription exceptionDescription;
442 exceptionDescription << "No transportation process found ";
443 G4Exception("G4ITStepProcessor::SetupGeneralProcessInfo",
444 "ITStepProcessor0006",
445 FatalErrorInArgument,exceptionDescription);
446 }
447 }
449 // fpProcessInfo = 0;
450}
451
452//______________________________________________________________________________
453
455{
456 fpTrack = track;
457 if(fpTrack)
458 {
460 fpStep = const_cast<G4Step*>(fpTrack->GetStep());
461
462 if(fpITrack)
463 {
465 }
466 else
467 {
468 fpTrackingInfo = 0;
469 G4cerr << "Track ID : " << fpTrack->GetTrackID() << G4endl;
470
472 errMsg << "No IT pointer was attached to the track you try to process.";
473 G4Exception("G4ITStepProcessor::SetTrack",
474 "ITStepProcessor0007",
476 errMsg);
477 }
478 }
479 else
480 {
481 fpITrack = 0;
482 fpStep = 0;
483 }
484}
485//______________________________________________________________________________
486
488{
490 std::map<const G4ParticleDefinition*, ProcessGeneralInfo*>::iterator it =
491 fProcessGeneralInfoMap.find(particle);
492
493 if(it == fProcessGeneralInfoMap.end())
494 {
497 if(fpProcessInfo == 0)
498 {
499 G4ExceptionDescription exceptionDescription("...");
500 G4Exception("G4ITStepProcessor::GetProcessNumber",
501 "ITStepProcessor0008",
503 exceptionDescription);
504 return;
505 }
506 }
507 else
508 {
509 fpProcessInfo = it->second;
510 }
511}
512
513//______________________________________________________________________________
514
516{
520
523
526}
527
528//______________________________________________________________________________
529
531{
532 // Reset the secondary particles
536}
537
538//______________________________________________________________________________
539
541{
542 // Select the rest process which has the shortest time before
543 // it is invoked. In rest processes, GPIL()
544 // returns the time before a process occurs.
545 G4double lifeTime(DBL_MAX), shortestLifeTime (DBL_MAX);
546
548 shortestLifeTime = DBL_MAX;
549
550 unsigned int NofInactiveProc=0;
551
552 for( size_t ri=0; ri < fpProcessInfo->MAXofAtRestLoops; ri++ )
553 {
555 if (fpCurrentProcess== 0)
556 {
558 NofInactiveProc++;
559 continue;
560 } // NULL means the process is inactivated by a user on fly.
561
565
568
569 if(fCondition==Forced)
570 {
572 }
573 else
574 {
576 if(lifeTime < shortestLifeTime )
577 {
578 shortestLifeTime = lifeTime;
580 }
581 }
582 }
583
585
586// G4cout << " --> Selected at rest process : "
587// << (*fpProcessInfo->fpAtRestGetPhysIntVector)[fAtRestDoItProcTriggered]
588// ->GetProcessName()
589// << G4endl;
590
591 fTimeStep = shortestLifeTime;
592
593 // at least one process is necessary to destroy the particle
594 // exit with warning
595 if(NofInactiveProc==fpProcessInfo->MAXofAtRestLoops)
596 {
597 G4cerr << "ERROR - G4ITStepProcessor::InvokeAtRestDoItProcs()" << G4endl
598 << " No AtRestDoIt process is active!" << G4endl;
599 }
600}
601
602//_________________________________________________________________________
604{
606 G4TrackManyList::iterator it = mainList ->begin();
607 G4TrackManyList::iterator end = mainList ->end();
608
609 SetPreviousStepTime(previousTimeStep);
610
612
613 for (; it != end; )
614 {
615 G4Track * track = *it;
616
617#ifdef DEBUG
618 G4cout << "*CIL* " << GetIT(track)->GetName()
619 << " ID: " << track->GetTrackID()
620 << " at time : " << track->GetGlobalTime()
621 << G4endl;
622#endif
623
624 ++it;
626
628 }
629
630 return fILTimeStep;
631}
632
633//_________________________________________________________________________
634
636{
637 assert(fpTrack != 0);
638 if (fpTrack == 0)
639 {
641 return;
642 }
643
644 // assert(fpTrack->GetTrackStatus() != fStopAndKill);
645
647 {
648// trackContainer->GetMainList()->pop(fpTrack);
651 return;
652 }
653
654 if (IsInf(fTimeStep))
655 {
656 // G4cout << "!!!!!!!!!!!! IS INF " << track->GetTrackID() << G4endl;
658 return;
659 }
660 else if (fTimeStep < fILTimeStep - DBL_EPSILON)
661 {
662 // G4cout << "!!!!!!!!!!!! TEMPS DIFFERENTS "
663 // << track->GetTrackID() << G4endl;
664
666
668
669// G4cout << "Will set leading step to true for time :"
670// << SP -> GetInteractionTime() << " against fTimeStep : "
671// << G4BestUnit(fILTimeStep, "Time") << " the trackID is : " << track->GetTrackID()
672// << G4endl;
673
674// GetIT(fpTrack)->GetTrackingInfo()->SetLeadingStep(true);
676 }
677 else if(fabs(fILTimeStep - fTimeStep) < DBL_EPSILON )
678 {
679
680 // G4cout << "!!!!!!!!!!!! MEME TEMPS " << track->GetTrackID() << G4endl;
681 // G4cout << "Will set leading step to true for time :"
682 // << SP -> GetInteractionTime() << " against fTimeStep : "
683 // << fTimeStep << " the trackID is : " << track->GetTrackID()<< G4endl;
684// GetIT(fpTrack)->GetTrackingInfo()->SetLeadingStep(true);
686 }
687 // else
688 // {
689 // G4cout << "!!!! Bigger time : " << "currentTime : "<<fILTimeStep
690 // << " proposedTime : " << SP -> GetInteractionTime() << G4endl;
691 // }
692
694}
695
696//___________________________________________________________________________
697
699{
700 SetTrack(track);
702}
703
704//______________________________________________________________________________
705
707{
708 // DEBUG
709 // G4cout << "SetInitialStep for : " << fpITrack-> GetName() << G4endl;
710 //________________________________________________________
711 // Initialize geometry
712
714 {
715 //==========================================================================
716 // Create navigator state and Locate particle in geometry
717 //==========================================================================
718 /*
719 fpNavigator->NewNavigatorStateAndLocate(fpTrack->GetPosition(),
720 fpTrack->GetMomentumDirection());
721
722 fpITrack->GetTrackingInfo()->
723 SetNavigatorState(fpNavigator->GetNavigatorState());
724 */
725 fpNavigator->NewNavigatorState();
727 ->GetNavigatorState());
728
730 fpNavigator->LocateGlobalPointAndSetup(fpTrack->GetPosition(),
731 &direction,
732 false,
733 false); // was false, false
734
735 fpState->fTouchableHandle = fpNavigator->CreateTouchableHistory();
736
739 }
740 else
741 {
744
745 //==========================================================================
746 // Create OR set navigator state
747 //==========================================================================
748
750 {
751 fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()
754 ->GetNavigatorState());
755 }
756 else
757 {
758 fpNavigator->NewNavigatorState(*((G4TouchableHistory*) fpState
759 ->fTouchableHandle()));
761 ->GetNavigatorState());
762 }
763
764 G4VPhysicalVolume* oldTopVolume =
766
767 //==========================================================================
768 // Locate particle in geometry
769 //==========================================================================
770
771// G4VPhysicalVolume* newTopVolume =
772// fpNavigator->LocateGlobalPointAndSetup(
773// fpTrack->GetPosition(),
774// &fpTrack->GetMomentumDirection(),
775// true, false);
776
777 G4VPhysicalVolume* newTopVolume =
778 fpNavigator->ResetHierarchyAndLocate(fpTrack->GetPosition(),
781 ->GetTouchableHandle()()));
782
783 if(newTopVolume != oldTopVolume || oldTopVolume->GetRegularStructureId()
784 == 1)
785 {
786 fpState->fTouchableHandle = fpNavigator->CreateTouchableHistory();
789 }
790 }
791
793
794 //________________________________________________________
795 // If the primary track has 'Suspend' or 'PostponeToNextEvent' state,
796 // set the track state to 'Alive'.
799 {
801 }
802
803 //HoangTRAN: it's better to check the status here
804 if(fpTrack->GetTrackStatus() == fStopAndKill) return;
805
806 // If the primary track has 'zero' kinetic energy, set the track
807 // state to 'StopButAlive'.
808 if(fpTrack->GetKineticEnergy() <= 0.0)
809 {
811 }
812 //________________________________________________________
813 // Set vertex information of G4Track at here
814 if(fpTrack->GetCurrentStepNumber() == 0)
815 {
820 }
821 //________________________________________________________
822 // If track is already outside the world boundary, kill it
823 if(fpCurrentVolume == 0)
824 {
825 // If the track is a primary, stop processing
826 if(fpTrack->GetParentID() == 0)
827 {
828 G4cerr << "ERROR - G4ITStepProcessor::SetInitialStep()" << G4endl<< " Primary particle starting at - "
829 << fpTrack->GetPosition()
830 << " - is outside of the world volume." << G4endl;
831 G4Exception("G4ITStepProcessor::SetInitialStep()", "ITStepProcessor0011",
832 FatalException, "Primary vertex outside of the world!");
833 }
834
836 G4cout << "WARNING - G4ITStepProcessor::SetInitialStep()" << G4endl
837 << " Initial track position is outside world! - "
838 << fpTrack->GetPosition() << G4endl;
839 }
840 else
841 {
842 // Initial set up for attribues of 'Step'
844 }
845
847}
848//______________________________________________________________________________
849
851{
852
853 if(!fpStep)
854 {
855 // Create new Step and give it to the track
856 fpStep = new G4Step();
859
860 // Create new state and set it in the trackingInfo
863
864 SetupMembers();
866
868 }
869 else
870 {
871 SetupMembers();
872
874 /***
875 // Send G4Step information to Hit/Dig if the volume is sensitive
876 fpCurrentVolume = fpStep->GetPreStepPoint()->GetPhysicalVolume();
877 StepControlFlag = fpStep->GetControlFlag();
878 if( fpCurrentVolume != 0 && StepControlFlag != AvoidHitInvocation)
879 {
880 fpSensitive = fpStep->GetPreStepPoint()->
881 GetSensitiveDetector();
882
883 // if( fSensitive != 0 ) {
884 // fSensitive->Hit(fStep);
885 // }
886 }
887 ***/
888 // Store last PostStepPoint to PreStepPoint, and swap current and next
889 // volume information of G4Track. Reset total energy deposit in one Step.
892
893 //JA Set the volume before it is used (in DefineStepLength() for User Limit)
895 /*
896 G4cout << G4endl;
897 G4cout << "!!!!!!!!!!!!!!!!!!!!!!!!!!!" << G4endl;
898 G4cout << "PreStepPoint Volume : "
899 << fpCurrentVolume->GetName() << G4endl;
900 G4cout << "Track Touchable : "
901 << fpTrack->GetTouchableHandle()->GetVolume()->GetName() << G4endl;
902 G4cout << "Track NextTouchable : "
903 << fpTrack->GetNextTouchableHandle()->GetVolume()->GetName()
904 << G4endl;
905 */
906 // Reset the step's auxiliary points vector pointer
908
909 // Switch next touchable in track to current one
913
915 /*
916 G4VPhysicalVolume* oldTopVolume =
917 fpTrack->GetTouchableHandle()->GetVolume();
918 fpNavigator->SetNavigatorState(
919 fpITrack->GetTrackingInfo()->GetNavigatorState());
920
921 G4VPhysicalVolume* newTopVolume = fpNavigator->ResetHierarchyAndLocate(
922 fpTrack->GetPosition(), fpTrack->GetMomentumDirection(),
923 *((G4TouchableHistory*) fpTrack->GetTouchableHandle()()));
924
925 // G4VPhysicalVolume* newTopVolume=
926 // fpNavigator->LocateGlobalPointAndSetup(fpTrack->GetPosition(),
927 // &fpTrack->GetMomentumDirection(),
928 // true, false);
929
930 // G4cout << "New Top Volume : " << newTopVolume->GetName() << G4endl;
931
932 if (newTopVolume != oldTopVolume || oldTopVolume->GetRegularStructureId()
933 == 1)
934 {
935 fpState->fTouchableHandle = fpNavigator->CreateTouchableHistory();
936 fpTrack->SetTouchableHandle(fpState->fTouchableHandle);
937 fpTrack->SetNextTouchableHandle(fpState->fTouchableHandle);
938 }
939
940 */
942 //==========================================================================
943 // Only reset navigator state + reset volume hierarchy (internal)
944 // No need to relocate
945 //==========================================================================
946 fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()
948 }
949}
950
951//______________________________________________________________________________
952
953// ------------------------------------------------------------------------
954// Compute Interaction Length
955// ------------------------------------------------------------------------
957{
958
960
961#ifdef G4VERBOSE
962 // !!!!! Verbose
964#endif
965
966 G4TrackStatus trackStatus = fpTrack->GetTrackStatus();
967
968 if(trackStatus == fStopAndKill)
969 {
970 return;
971 }
972
973 if(trackStatus == fStopButAlive)
974 {
976 ->GetNavigatorState());
977 fpNavigator->ResetNavigatorState();
978 return GetAtRestIL();
979 }
980
981 // Find minimum Step length and corresponding time
982 // demanded by active disc./cont. processes
983
984 // ReSet the counter etc.
985 fpState->fPhysicalStep = DBL_MAX; // Initialize by a huge number
986 fPhysIntLength = DBL_MAX; // Initialize by a huge number
987
988 double proposedTimeStep = DBL_MAX;
989 G4VProcess* processWithPostStepGivenByTimeStep(0);
990
991 // GPIL for PostStep
994
995 // G4cout << "fpProcessInfo->MAXofPostStepLoops : "
996 // << fpProcessInfo->MAXofPostStepLoops
997 // << " mol : " << fpITrack -> GetName()
998 // << " id : " << fpTrack->GetTrackID()
999 // << G4endl;
1000
1001 for(size_t np = 0; np < fpProcessInfo->MAXofPostStepLoops; np++)
1002 {
1005 if(fpCurrentProcess == 0)
1006 {
1008 continue;
1009 } // NULL means the process is inactivated by a user on fly.
1010
1013 ->GetProcessID()));
1014
1015 // G4cout << "Is going to call : "
1016 // << fpCurrentProcess -> GetProcessName()
1017 // << G4endl;
1020 &fCondition);
1021
1022#ifdef G4VERBOSE
1023 // !!!!! Verbose
1025#endif
1026
1028 //fpCurrentProcess->SetProcessState(0);
1029
1030 switch(fCondition)
1031 {
1032 case ExclusivelyForced: // Will need special treatment
1036 break;
1037
1038 case Conditionally:
1039 // (fpState->fSelectedPostStepDoItVector)[np] = Conditionally;
1040 G4Exception("G4ITStepProcessor::DefinePhysicalStepLength()",
1041 "ITStepProcessor0008",
1043 "This feature is no more supported");
1044 break;
1045
1046 case Forced:
1048 break;
1049
1050 case StronglyForced:
1052 break;
1053
1054 default:
1056 break;
1057 }
1058
1060 {
1061 for(size_t nrest = np + 1; nrest < fpProcessInfo->MAXofPostStepLoops;
1062 nrest++)
1063 {
1065 }
1066 return; // Please note the 'return' at here !!!
1067 }
1068 else
1069 {
1070 if(fPhysIntLength < fpState->fPhysicalStep)
1071 {
1072 // To avoid checking whether the process is actually
1073 // proposing a time step, the returned time steps are
1074 // negative (just for tagging)
1076 {
1077 fPhysIntLength *= -1;
1078 if(fPhysIntLength < proposedTimeStep)
1079 {
1080 proposedTimeStep = fPhysIntLength;
1082 processWithPostStepGivenByTimeStep = fpCurrentProcess;
1083 }
1084 }
1085 else
1086 {
1091 }
1092 }
1093 }
1094 }
1095
1096 // GPIL for AlongStep
1098 G4double safetyProposedToAndByProcess = fpState->fProposedSafety;
1099
1100 for(size_t kp = 0; kp < fpProcessInfo->MAXofAlongStepLoops; kp++)
1101 {
1104 if(fpCurrentProcess == 0) continue;
1105 // NULL means the process is inactivated by a user on fly.
1106
1113 safetyProposedToAndByProcess,
1115
1116#ifdef G4VERBOSE
1117 // !!!!! Verbose
1119#endif
1120
1121 if(fPhysIntLength < fpState->fPhysicalStep)
1122 {
1124 // Should save PS and TS in IT
1125
1126 // Check if the process wants to be the GPIL winner. For example,
1127 // multi-scattering proposes Step limit, but won't be the winner.
1129 {
1132 }
1133
1134 // Transportation is assumed to be the last process in the vector
1135 if(kp == fpProcessInfo->MAXofAlongStepLoops - 1)
1136 {
1138
1139 if(!fpTransportation)
1140 {
1141 G4ExceptionDescription exceptionDescription;
1142 exceptionDescription << "No transportation process found ";
1143 G4Exception("G4ITStepProcessor::DoDefinePhysicalStepLength",
1144 "ITStepProcessor0009",
1146 exceptionDescription);
1147 }
1148
1150
1153 }
1154 }
1155 else
1156 {
1157 if(kp == fpProcessInfo->MAXofAlongStepLoops - 1)
1158 {
1160
1161 if(!fpTransportation)
1162 {
1163 G4ExceptionDescription exceptionDescription;
1164 exceptionDescription << "No transportation process found ";
1165 G4Exception("G4ITStepProcessor::DoDefinePhysicalStepLength",
1166 "ITStepProcessor0010",
1168 exceptionDescription);
1169 }
1170
1172 }
1173 }
1174
1175 // Handle PostStep processes sending back time steps rather than space length
1176 if(proposedTimeStep < fTimeStep)
1177 {
1178 if(fPostStepAtTimeDoItProcTriggered < fpProcessInfo->MAXofPostStepLoops)
1179 {
1181 {
1183 NotForced;
1184 // (fpState->fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] = InActivated;
1185
1187 fpStep->GetPostStepPoint()->SetProcessDefinedStep(processWithPostStepGivenByTimeStep);
1188
1189 fTimeStep = proposedTimeStep;
1190
1192 *fpStep,
1193 fTimeStep,
1195 }
1196 }
1197 }
1198 else
1199 {
1200 if(fPostStepDoItProcTriggered < fpProcessInfo->MAXofPostStepLoops)
1201 {
1203 {
1205 NotForced;
1206 }
1207 }
1208 }
1209
1210// fpCurrentProcess->SetProcessState(0);
1212
1213 // Make sure to check the safety, even if Step is not limited
1214 // by this process. J. Apostolakis, June 20, 1998
1215 //
1216 if(safetyProposedToAndByProcess < fpState->fProposedSafety)
1217 // proposedSafety keeps the smallest value:
1218 fpState->fProposedSafety = safetyProposedToAndByProcess;
1219 else
1220 // safetyProposedToAndByProcess always proposes a valid safety:
1221 safetyProposedToAndByProcess = fpState->fProposedSafety;
1222
1223 }
1224
1225 fpITrack->GetTrackingInfo()->SetNavigatorState(fpNavigator->GetNavigatorState());
1226 fpNavigator->ResetNavigatorState();
1227}
1228
1229//______________________________________________________________________________
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
@ InActivated
@ StronglyForced
@ NotForced
@ Conditionally
@ ExclusivelyForced
@ Forced
@ CandidateForSelection
bool IsInf(T value)
static const size_t SizeOfSelectedDoItVector
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:48
@ typeGPIL
@ typeDoIt
@ fGeomBoundary
Definition: G4StepStatus.hh:43
@ fWorldBoundary
Definition: G4StepStatus.hh:41
@ fUndefined
Definition: G4StepStatus.hh:55
@ fPostStepDoItProc
Definition: G4StepStatus.hh:49
@ fAlongStepDoItProc
Definition: G4StepStatus.hh:47
@ fExclusivelyForcedProc
Definition: G4StepStatus.hh:53
G4TrackStatus
@ fSuspend
@ fAlive
@ fStopAndKill
@ fStopButAlive
@ fPostponeToNextEvent
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define theParticleIterator
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
void Push(G4Track *)
G4TouchableHandle fTouchableHandle
G4SelectedPostStepDoItVector fSelectedPostStepDoItVector
G4SelectedAtRestDoItVector fSelectedAtRestDoItVector
G4ITStepProcessorState & operator=(const G4ITStepProcessorState &)
G4VITProcess * fpCurrentProcess
std::map< const G4ParticleDefinition *, ProcessGeneralInfo * > fProcessGeneralInfoMap
ProcessGeneralInfo * fpProcessInfo
G4ITStepProcessorState * fpState
G4StepPoint * fpPostStepPoint
G4ITNavigator * fpNavigator
void SetNavigator(G4ITNavigator *value)
G4ITTransportation * fpTransportation
G4GPILSelection fGPILSelection
G4ForceCondition fCondition
void SetTrack(G4Track *)
G4ITStepProcessor & operator=(const G4ITStepProcessor &other)
G4double ComputeInteractionLength(double previousTimeStep)
void SetPreviousStepTime(G4double)
G4StepPoint * fpPreStepPoint
G4TrackingInformation * fpTrackingInfo
void DefinePhysicalStepLength(G4Track *)
void SetupGeneralProcessInfo(G4ParticleDefinition *, G4ProcessManager *)
G4VITSteppingVerbose * fpVerbose
G4ITTrackingManager * fpTrackingManager
virtual void Initialize()
G4ITTrackHolder * fpTrackContainer
G4TrackVector * fpSecondary
G4ITLeadingTracks fLeadingTracks
size_t fPostStepAtTimeDoItProcTriggered
G4VPhysicalVolume * fpCurrentVolume
G4TrackList * GetMainList(Key)
static G4ITTrackHolder * Instance()
G4VITSteppingVerbose * GetSteppingVerbose()
G4ITTrackingInteractivity * GetInteractivity()
void StartTracking(G4Track *)
void EndTracking(G4Track *)
static G4ITTransportationManager * GetTransportationManager()
virtual void ComputeStep(const G4Track &, const G4Step &, const double timeStep, double &spaceStep)
G4TrackingInformation * GetTrackingInfo()
Definition: G4IT.hh:143
virtual const G4String & GetName() const =0
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
G4PTblDicIterator * GetIterator() const
static G4ParticleTable * GetParticleTable()
G4VProcess * SetProcessActivation(G4VProcess *aProcess, G4bool fActive)
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ProcessVector * GetProcessList() const
G4ProcessVector * GetPostStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
std::size_t entries() const
void SetProcessDefinedStep(const G4VProcess *aValue)
G4VPhysicalVolume * GetPhysicalVolume() const
Definition: G4Step.hh:62
void DeleteSecondaryVector()
void SetPointerToVectorOfAuxiliaryPoints(std::vector< G4ThreeVector > *vec)
void InitializeStep(G4Track *aValue)
void ResetTotalEnergyDeposit()
G4TrackVector * GetfSecondary()
void CopyPostToPreStepPoint()
G4StepPoint * GetPreStepPoint() const
G4TrackVector * NewSecondaryVector()
G4StepPoint * GetPostStepPoint() const
G4TrackStatus GetTrackStatus() const
void SetTrackStatus(const G4TrackStatus aTrackStatus)
G4int GetTrackID() const
void SetStep(const G4Step *aValue)
G4VPhysicalVolume * GetVolume() const
void SetVertexPosition(const G4ThreeVector &aValue)
const G4TouchableHandle & GetNextTouchableHandle() const
void SetVertexMomentumDirection(const G4ThreeVector &aValue)
G4VPhysicalVolume * GetNextVolume() const
void SetNextTouchableHandle(const G4TouchableHandle &apValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4int GetCurrentStepNumber() const
G4ParticleDefinition * GetDefinition() const
const G4TouchableHandle & GetTouchableHandle() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetStepLength() const
void SetVertexKineticEnergy(const G4double aValue)
G4int GetParentID() const
void SetLogicalVolumeAtVertex(const G4LogicalVolume *)
const G4Step * GetStep() const
G4shared_ptr< G4ProcessState_Lock > GetProcessState(size_t index)
void SetStepProcessorState(G4ITStepProcessorState_Lock *)
void SetNavigatorState(G4ITNavigatorState_Lock *)
G4ITStepProcessorState_Lock * GetStepProcessorState()
G4ITNavigatorState_Lock * GetNavigatorState() const
G4bool ProposesTimeStep() const
size_t GetProcessID() const
void ResetProcessState()
void SetProcessState(G4shared_ptr< G4ProcessState_Lock > aProcInfo)
G4double GetInteractionTimeLeft()
virtual void DPSLAlongStep()=0
virtual void DPSLStarted()=0
void SetStepProcessor(const G4ITStepProcessor *stepProcessor)
virtual void DPSLPostStep()=0
G4LogicalVolume * GetLogicalVolume() const
virtual G4int GetRegularStructureId() const =0
G4double PostStepGPIL(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
Definition: G4VProcess.hh:479
G4double AtRestGPIL(const G4Track &track, G4ForceCondition *condition)
Definition: G4VProcess.hh:472
G4double AlongStepGPIL(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
Definition: G4VProcess.hh:461
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34
G4ProcessVector * fpPostStepDoItVector
G4ProcessVector * fpPostStepGetPhysIntVector
G4ProcessVector * fpAlongStepGetPhysIntVector
G4ProcessVector * fpAlongStepDoItVector
G4ProcessVector * fpAtRestGetPhysIntVector
G4ITTransportation * fpTransportation
G4ProcessVector * fpAtRestDoItVector
#define DBL_EPSILON
Definition: templates.hh:66
#define DBL_MAX
Definition: templates.hh:62