Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Scheduler Class Reference

#include <G4Scheduler.hh>

+ Inheritance diagram for G4Scheduler:

Public Member Functions

virtual G4bool Notify (G4ApplicationState requestedState)
 
virtual void RegisterModel (G4VITStepModel *, double)
 
void Initialize ()
 
void ForceReinitialization ()
 
bool IsInitialized ()
 
bool IsRunning ()
 
void Reset ()
 
void Process ()
 
void ClearList ()
 
void SetGun (G4ITGun *)
 
G4ITGunGetGun ()
 
void Stop ()
 
void Clear ()
 
void EndTracking ()
 
void SetEndTime (const double)
 
void SetTimeTolerance (double)
 
double GetTimeTolerance () const
 
void SetMaxZeroTimeAllowed (int)
 
int GetMaxZeroTimeAllowed () const
 
G4ITModelHandlerGetModelHandler ()
 
void SetTimeSteps (std::map< double, double > *)
 
void AddTimeStep (double, double)
 
void SetDefaultTimeStep (double)
 
double GetLimitingTimeStep () const
 
G4int GetNbSteps () const
 
void SetMaxNbSteps (G4int)
 
G4int GetMaxNbSteps () const
 
G4double GetStartTime () const
 
G4double GetEndTime () const
 
virtual G4double GetTimeStep () const
 
G4double GetPreviousTimeStep () const
 
G4double GetGlobalTime () const
 
void SetUserAction (G4UserTimeStepAction *)
 
G4UserTimeStepActionGetUserTimeStepAction () const
 
void UseDefaultTimeSteps (G4bool)
 
G4bool AreDefaultTimeStepsUsed ()
 
G4ITStepStatus GetStatus () const
 
void SetVerbose (int)
 
int GetVerbose () const
 
void WhyDoYouStop ()
 
void SetInteractivity (G4ITTrackingInteractivity *)
 
G4ITTrackingInteractivityGetInteractivity ()
 
virtual size_t GetNTracks ()
 
void GetCollisionType (G4String &interactionType)
 
void AddWatchedTime (double time)
 
double GetNextWatchedTime () const
 
void SetMaxTimeStep (double maxTimeStep)
 
double GetMaxTimeStep () const
 
virtual void Initialize ()
 
virtual void Reset ()
 
virtual void SetVerbose (int)
 
virtual void SetGun (G4ITGun *)
 
virtual void Process ()
 
virtual G4bool IsRunning ()
 
virtual G4ITModelHandlerGetModelHandler ()
 
virtual void RegisterModel (G4VITStepModel *, double)
 
virtual void SetEndTime (const double)
 
virtual void SetTimeTolerance (double)
 
virtual double GetTimeTolerance () const
 
virtual void SetMaxZeroTimeAllowed (int)
 
virtual int GetMaxZeroTimeAllowed () const
 
virtual void SetTimeSteps (std::map< double, double > *)
 
virtual void AddTimeStep (double, double)
 
virtual void SetDefaultTimeStep (double)
 
virtual double GetLimitingTimeStep () const
 
virtual G4int GetNbSteps () const
 
virtual void SetMaxNbSteps (G4int)
 
virtual G4int GetMaxNbSteps () const
 
virtual G4double GetStartTime () const
 
virtual G4double GetEndTime () const
 
virtual G4double GetTimeStep () const
 
virtual G4double GetPreviousTimeStep () const
 
virtual G4double GetGlobalTime () const
 
virtual void SetUserAction (G4UserTimeStepAction *)
 
virtual G4UserTimeStepActionGetUserTimeStepAction () const
 
virtual void SetInteractivity (G4ITTrackingInteractivity *)
 
virtual G4ITTrackingInteractivityGetInteractivity ()
 
- Public Member Functions inherited from G4VStateDependent
 G4VStateDependent (G4bool bottom=false)
 
virtual ~G4VStateDependent ()
 
G4bool operator== (const G4VStateDependent &right) const
 
G4bool operator!= (const G4VStateDependent &right) const
 
virtual G4bool Notify (G4ApplicationState requestedState)=0
 

Static Public Member Functions

static G4SchedulerInstance ()
 
static void DeleteInstance ()
 
- Static Public Member Functions inherited from G4VScheduler
static G4VSchedulerInstance ()
 

Protected Member Functions

virtual ~G4Scheduler ()
 
void DoProcess ()
 
void SynchronizeTracks ()
 
void Stepping ()
 
void FindUserPreDefinedTimeStep ()
 
bool CanICarryOn ()
 
void PrintWhyDoYouStop ()
 
- Protected Member Functions inherited from G4VScheduler
 G4VScheduler ()
 
virtual ~G4VScheduler ()
 

Detailed Description

G4Scheduler synchronizes (in time) track stepping

Definition at line 87 of file G4Scheduler.hh.

Constructor & Destructor Documentation

◆ ~G4Scheduler()

G4Scheduler::~G4Scheduler ( )
protectedvirtual

Definition at line 200 of file G4Scheduler.cc.

201{
202
203 if(fpMessenger) // is used as a flag to know whether the manager was cleared
204 {
205 Clear();
206 }
207
208 fgScheduler = 0;
209
210// if (fVerbose >= 1)
211// {
212// G4cout << "G4Scheduler is being deleted. Bye :) !" << G4endl;
213// }
214}
void Clear()
Definition: G4Scheduler.cc:216

Member Function Documentation

◆ AddTimeStep()

void G4Scheduler::AddTimeStep ( double  startingTime,
double  timeStep 
)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 304 of file G4Scheduler.hh.

305{
306 if (fpUserTimeSteps == 0)
307 {
308 fpUserTimeSteps = new std::map<double, double>();
309 fUsePreDefinedTimeSteps = true;
310 }
311
312 (*fpUserTimeSteps)[startingTime] = timeStep;
313}

◆ AddWatchedTime()

void G4Scheduler::AddWatchedTime ( double  time)
inline

Definition at line 176 of file G4Scheduler.hh.

177 {
178 fWatchedTimes.insert(time);
179 }

◆ AreDefaultTimeStepsUsed()

G4bool G4Scheduler::AreDefaultTimeStepsUsed ( )
inline

Definition at line 437 of file G4Scheduler.hh.

438{
439 return (fUseDefaultTimeSteps == false && fUsePreDefinedTimeSteps == false);
440}

Referenced by G4SchedulerMessenger::GetCurrentValue().

◆ CanICarryOn()

bool G4Scheduler::CanICarryOn ( )
protected

Definition at line 550 of file G4Scheduler.cc.

551{
552 return fGlobalTime < fEndTime && (fMaxSteps == -1 ? true : fNbSteps < fMaxSteps)
553 && fContinue == true;
554}

Referenced by SynchronizeTracks().

◆ Clear()

void G4Scheduler::Clear ( )

Definition at line 216 of file G4Scheduler.cc.

217{
218// if (fVerbose) G4cout << "*** G4Scheduler is being cleared ***" << G4endl;
219
220 if(fpMessenger)
221 {
222 delete fpMessenger;
223 fpMessenger = 0;
224 }
225 if(fpStepProcessor)
226 {
227 delete fpStepProcessor;
228 fpStepProcessor = 0;
229 }
230 if(fpModelProcessor)
231 {
232 delete fpModelProcessor;
233 fpModelProcessor = 0;
234 }
235
237 ClearList();
238 if(fpTrackingManager)
239 {
240 delete fpTrackingManager;
241 fpTrackingManager = 0;
242 }
243
244 if(fReactionSet)
245 {
246 delete fReactionSet;
247 fReactionSet = 0;
248 }
249
250 if(fpModelHandler)
251 {
252 delete fpModelHandler;
253 fpModelHandler = 0;
254 }
255
256 //* DEBUG
257 //* assert(G4StateManager::GetStateManager()->
258 //* DeregisterDependent(this) == true);
259
260}
static G4ITTypeManager * Instance()
Definition: G4ITType.cc:57
void ReleaseRessource()
Definition: G4ITType.cc:82
void ClearList()
Definition: G4Scheduler.cc:264

Referenced by Notify(), and ~G4Scheduler().

◆ ClearList()

void G4Scheduler::ClearList ( )

Definition at line 264 of file G4Scheduler.cc.

265{
266// if (fNbTracks == 0) return;
267
268 fTrackContainer.Clear();
269
271}
static void DeleteInstance()

Referenced by Clear(), and Process().

◆ DeleteInstance()

void G4Scheduler::DeleteInstance ( )
static

DeleteInstance should be used instead of the destructor

Definition at line 123 of file G4Scheduler.cc.

124{
125 if(fgScheduler)
126 {
127 delete fgScheduler;
128 }
129}

◆ DoProcess()

void G4Scheduler::DoProcess ( )
protected

Definition at line 602 of file G4Scheduler.cc.

604{
605 if(fpUserTimeStepAction) fpUserTimeStepAction->NewStage();
606
607#if defined (DEBUG_MEM) && defined (DEBUG_MEM_STEPPING)
608 MemStat mem_first, mem_second, mem_diff;
609#endif
610
611#if defined (DEBUG_MEM) && defined (DEBUG_MEM_STEPPING)
612 mem_first = MemoryUsage();
613#endif
614
615 while (fGlobalTime < fStopTime
616 && fTrackContainer.MainListsNOTEmpty()
617 && (fMaxSteps == -1 ? true : fNbSteps < fMaxSteps)
618 && fContinue == true)
619 {
620// G4cout << "Mainlist size : " << fTrackContainer.GetMainList()->size()
621// << G4endl;
622
623 Stepping();
624
625#if defined (DEBUG_MEM) && defined (DEBUG_MEM_STEPPING)
626 mem_second = MemoryUsage();
627 mem_diff = mem_second-mem_first;
628 G4cout << "\t || MEM || After step " << fNbSteps << ", diff is : "
629 << mem_diff << G4endl;
630#endif
631 }
632
634
635#if defined (DEBUG_MEM) && defined (DEBUG_MEM_STEPPING)
636 mem_second = MemoryUsage();
637 mem_diff = mem_second-mem_first;
638 G4cout << "\t || MEM || After stepping, diff is : " << mem_diff << G4endl;
639#endif
640
641#ifdef G4VERBOSE
642 if(fVerbose > 2)
643 G4cout << "*** G4Scheduler has finished processing a track list at time : "
644 << G4BestUnit(fGlobalTime, "Time") << G4endl;
645#endif
646}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void Stepping()
Definition: G4Scheduler.cc:649
void PrintWhyDoYouStop()
Definition: G4Scheduler.cc:558
MemStat MemoryUsage()
Definition: G4MemStat.cc:55

Referenced by SynchronizeTracks().

◆ EndTracking()

void G4Scheduler::EndTracking ( )

Definition at line 1061 of file G4Scheduler.cc.

1062{
1063 if(fRunning)
1064 {
1065 G4ExceptionDescription exceptionDescription;
1066 exceptionDescription
1067 << "End tracking is called while G4Scheduler is still running."
1068 << G4endl;
1069
1070 G4Exception("G4Scheduler::EndTracking",
1071 "Scheduler017",
1073 exceptionDescription);
1074 }
1075
1076 fTrackContainer.MergeSecondariesWithMainList();
1077
1078 if (fTrackContainer.MainListsNOTEmpty())
1079 {
1080 G4TrackManyList* mainList = fTrackContainer.GetMainList();
1081 G4TrackManyList::iterator it = mainList->begin();
1082 G4TrackManyList::iterator end = mainList->end();
1083 for (; it != end; ++it)
1084 {
1085 fpTrackingManager->EndTrackingWOKill(*it);
1086 }
1087 }
1088
1089 if (fTrackContainer.SecondaryListsNOTEmpty()) // should be empty
1090 {
1091 G4TrackManyList* secondaries = fTrackContainer.GetSecondariesList();
1092 G4TrackManyList::iterator it = secondaries->begin();
1093 G4TrackManyList::iterator end = secondaries->end();
1094
1095 for (; it != end; ++it)
1096 {
1097 fpTrackingManager->EndTrackingWOKill(*it);
1098 }
1099 }
1100}
@ 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
G4TrackList * GetMainList(Key)
void MergeSecondariesWithMainList()
G4TrackManyList * GetSecondariesList()
bool SecondaryListsNOTEmpty()
void EndTrackingWOKill(G4Track *)

Referenced by Process().

◆ FindUserPreDefinedTimeStep()

void G4Scheduler::FindUserPreDefinedTimeStep ( )
protected

Definition at line 1006 of file G4Scheduler.cc.

1007{
1008
1009 if(fpUserTimeSteps == 0)
1010 {
1011 G4ExceptionDescription exceptionDescription;
1012 exceptionDescription
1013 << "You are asking to use user defined steps but you did not give any.";
1014 G4Exception("G4Scheduler::FindUserPreDefinedTimeStep",
1015 "Scheduler004",
1017 exceptionDescription);
1018 return; // makes coverity happy
1019 }
1020 map<double, double>::iterator fpUserTimeSteps_i =
1021 fpUserTimeSteps->upper_bound(fGlobalTime);
1022 map<double, double>::iterator fpUserTimeSteps_low = fpUserTimeSteps
1023 ->lower_bound(fGlobalTime);
1024
1025 // DEBUG
1026 // G4cout << "fGlobalTime : " << G4BestUnit(fGlobalTime,"Time") << G4endl;
1027 // G4cout << "fpUserTimeSteps_i : "
1028 // <<"<"<<G4BestUnit(fpUserTimeSteps_i->first,"Time")<<", "
1029 // << G4BestUnit(fpUserTimeSteps_i->second,"Time")<<">"
1030 // << "\t fpUserTimeSteps_low : "
1031 // <<"<"<<G4BestUnit(fpUserTimeSteps_low->first,"Time")<<", "
1032 // << G4BestUnit(fpUserTimeSteps_low->second,"Time")<<">"
1033 // << G4endl;
1034
1035 if(fpUserTimeSteps_i == fpUserTimeSteps->end())
1036 {
1037 fpUserTimeSteps_i--;
1038 }
1039 else if(fabs(fGlobalTime - fpUserTimeSteps_low->first) < fTimeTolerance)
1040 {
1041 // Case : fGlobalTime = X picosecond
1042 // and fpUserTimeSteps_low->first = X picosecond
1043 // but the precision is not good enough
1044 fpUserTimeSteps_i = fpUserTimeSteps_low;
1045 }
1046 else if(fpUserTimeSteps_i == fpUserTimeSteps_low)
1047 {
1048 // "Normal" cases
1049 fpUserTimeSteps_i--;
1050 }
1051 else
1052 {
1053 fpUserTimeSteps_i = fpUserTimeSteps_low;
1054 }
1055
1056 fDefinedMinTimeStep = fpUserTimeSteps_i->second;
1057}

◆ ForceReinitialization()

void G4Scheduler::ForceReinitialization ( )

Definition at line 1115 of file G4Scheduler.cc.

1116{
1117 fInitialized = false;
1118 Initialize();
1119}
void Initialize()
Definition: G4Scheduler.cc:281

◆ GetCollisionType()

void G4Scheduler::GetCollisionType ( G4String interactionType)

Definition at line 1146 of file G4Scheduler.cc.

1147{
1148 switch(fITStepStatus)
1149 {
1151 interactionType = "eInteractionWithMedium";
1152 break;
1154 interactionType = "eCollisionBetweenTracks";
1155 break;
1156 default:
1157 interactionType = "eCollisionBetweenTracks";
1158 break;
1159 }
1160}
@ eInteractionWithMedium
@ eCollisionBetweenTracks

Referenced by Stepping().

◆ GetEndTime()

G4double G4Scheduler::GetEndTime ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 335 of file G4Scheduler.hh.

336{
337 return fEndTime;
338}

Referenced by G4DNAIRT::G4DNAIRT(), G4SchedulerMessenger::GetCurrentValue(), and G4DNAIRT::Initialize().

◆ GetGlobalTime()

G4double G4Scheduler::GetGlobalTime ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 350 of file G4Scheduler.hh.

351{
352 return fGlobalTime;
353}

Referenced by G4ITTrackHolder::_PushTrack(), G4DNAIRT::MakeReaction(), and G4DNAIRT::Sampling().

◆ GetGun()

G4ITGun * G4Scheduler::GetGun ( )
inline

Definition at line 422 of file G4Scheduler.hh.

423{
424 return fpGun;
425}

◆ GetInteractivity()

G4ITTrackingInteractivity * G4Scheduler::GetInteractivity ( )
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 412 of file G4Scheduler.hh.

413{
414 return fpTrackingInteractivity;
415}

◆ GetLimitingTimeStep()

double G4Scheduler::GetLimitingTimeStep ( ) const
virtual

Reimplemented from G4VScheduler.

Definition at line 940 of file G4Scheduler.cc.

941{
942 if (fpUserTimeSteps == 0) return fDefaultMinTimeStep;
943 if (fabs(fGlobalTime - fUserUpperTimeLimit) < fTimeTolerance)
944 return fDefinedMinTimeStep;
945
946 map<double, double>::const_iterator it_fpUserTimeSteps_i = fpUserTimeSteps
947 ->upper_bound(fGlobalTime);
948 map<double, double>::const_iterator it_fpUserTimeSteps_low = fpUserTimeSteps
949 ->lower_bound(fGlobalTime);
950
951 // DEBUG
952 // G4cout << "fGlobalTime : " << G4BestUnit(fGlobalTime,"Time")
953 // << G4endl;
954 // G4cout << "fpUserTimeSteps_i : "
955 // <<"<"<<G4BestUnit(it_fpUserTimeSteps->first,"Time")
956 // <<", "<< G4BestUnit(it_fpUserTimeSteps->second,"Time")<<">"
957 // << "\t fpUserTimeSteps_low : "
958 // <<"<"<<G4BestUnit(fpUserTimeSteps_low->first,"Time")<<", "*
959 // << G4BestUnit(fpUserTimeSteps_low->second,"Time")<<">"
960 // << G4endl;
961
962 if (it_fpUserTimeSteps_i == fpUserTimeSteps->end())
963 {
964 it_fpUserTimeSteps_i--;
965 fUserUpperTimeLimit = fStopTime;
966 }
967 else if (fabs(fGlobalTime - it_fpUserTimeSteps_low->first) < fTimeTolerance)
968 {
969 // Case : fGlobalTime = X picosecond
970 // and fpUserTimeSteps_low->first = X picosecond
971 // but the precision is not good enough
972 it_fpUserTimeSteps_i = it_fpUserTimeSteps_low;
973 map<double, double>::const_iterator tmp_it = it_fpUserTimeSteps_low;
974 ++tmp_it;
975 if (tmp_it == fpUserTimeSteps->end())
976 {
977 fUserUpperTimeLimit = fStopTime;
978 }
979 else
980 {
981 fUserUpperTimeLimit = tmp_it->first;
982 }
983 }
984 else if (it_fpUserTimeSteps_i == it_fpUserTimeSteps_low)
985 {
986 // "Normal" cases
987 fUserUpperTimeLimit = it_fpUserTimeSteps_i->first;
988// it_fpUserTimeSteps_i++;
989// G4cout << "Global time = " << fGlobalTime << G4endl;
990// G4cout << "Is begin = "
991// << (it_fpUserTimeSteps_i == fpUserTimeSteps->begin())<< G4endl;
992
993 if(it_fpUserTimeSteps_i != fpUserTimeSteps->begin()) it_fpUserTimeSteps_i--;
994 }
995 else
996 {
997 fUserUpperTimeLimit = it_fpUserTimeSteps_i->first;
998 it_fpUserTimeSteps_i = it_fpUserTimeSteps_low;
999 }
1000
1001 return it_fpUserTimeSteps_i->second;
1002}

Referenced by Stepping().

◆ GetMaxNbSteps()

G4int G4Scheduler::GetMaxNbSteps ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 325 of file G4Scheduler.hh.

326{
327 return fMaxSteps;
328}

Referenced by G4SchedulerMessenger::GetCurrentValue().

◆ GetMaxTimeStep()

double G4Scheduler::GetMaxTimeStep ( ) const
inline

Definition at line 188 of file G4Scheduler.hh.

189 {
190 return fMaxTimeStep;
191 }

◆ GetMaxZeroTimeAllowed()

int G4Scheduler::GetMaxZeroTimeAllowed ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 382 of file G4Scheduler.hh.

383{
384 return fMaxNZeroTimeStepsAllowed;
385}

Referenced by G4SchedulerMessenger::GetCurrentValue().

◆ GetModelHandler()

G4ITModelHandler * G4Scheduler::GetModelHandler ( )
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 287 of file G4Scheduler.hh.

288{
289 return fpModelHandler;
290}

◆ GetNbSteps()

G4int G4Scheduler::GetNbSteps ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 315 of file G4Scheduler.hh.

316{
317 return fNbSteps;
318}

◆ GetNextWatchedTime()

double G4Scheduler::GetNextWatchedTime ( ) const

Definition at line 488 of file G4Scheduler.cc.

489{
490 std::set<double>::const_iterator up = fWatchedTimes.upper_bound(fGlobalTime);
491 if(up == fWatchedTimes.end()) return DBL_MAX;
492 return *up;
493}
#define DBL_MAX
Definition: templates.hh:62

Referenced by SynchronizeTracks().

◆ GetNTracks()

size_t G4Scheduler::GetNTracks ( )
virtual

Definition at line 1140 of file G4Scheduler.cc.

1141{
1142 return fTrackContainer.GetNTracks();
1143}
virtual size_t GetNTracks()

◆ GetPreviousTimeStep()

G4double G4Scheduler::GetPreviousTimeStep ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 397 of file G4Scheduler.hh.

398{
399 return fPreviousTimeStep;
400}

◆ GetStartTime()

G4double G4Scheduler::GetStartTime ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 330 of file G4Scheduler.hh.

331{
332 return fStartTime;
333}

Referenced by G4DNAIRT::G4DNAIRT(), and G4DNAIRT::Initialize().

◆ GetStatus()

G4ITStepStatus G4Scheduler::GetStatus ( ) const
inline

Definition at line 402 of file G4Scheduler.hh.

403{
404 return fITStepStatus;
405}

◆ GetTimeStep()

G4double G4Scheduler::GetTimeStep ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 340 of file G4Scheduler.hh.

341{
342 return fTimeStep;
343}

◆ GetTimeTolerance()

double G4Scheduler::GetTimeTolerance ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 392 of file G4Scheduler.hh.

393{
394 return fTimeTolerance;
395}

Referenced by G4ITTrackHolder::_PushTrack(), and G4SchedulerMessenger::GetCurrentValue().

◆ GetUserTimeStepAction()

G4UserTimeStepAction * G4Scheduler::GetUserTimeStepAction ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 361 of file G4Scheduler.hh.

362{
363 return fpUserTimeStepAction;
364}

◆ GetVerbose()

int G4Scheduler::GetVerbose ( ) const
inline

Definition at line 371 of file G4Scheduler.hh.

372{
373 return fVerbose;
374}

Referenced by G4SchedulerMessenger::GetCurrentValue().

◆ Initialize()

void G4Scheduler::Initialize ( )
virtual

Reimplemented from G4VScheduler.

Definition at line 281 of file G4Scheduler.cc.

282{
283 if(fpStepProcessor)
284 {
285 delete fpStepProcessor;
286 }
287 if(fpModelProcessor)
288 {
289 delete fpModelProcessor;
290 }
291 // if(fpMasterModelProcessor)
292 // {
293 // delete fpMasterModelProcessor;
294 // }
295
296 //______________________________________________________________
297
298 fpModelProcessor = new G4ITModelProcessor();
299 fpModelProcessor->SetModelHandler(fpModelHandler);
300 fpModelProcessor->SetTrackingManager(fpTrackingManager);
301
302 // fpMasterModelProcessor = new G4ITModelProcessor();
303 // fpMasterModelProcessor->SetModelHandler(fpModelHandler);
304 // fpModelProcessor = fpMasterModelProcessor;
305
306 //______________________________________________________________
307
308 fpStepProcessor = new G4ITStepProcessor();
309 fpStepProcessor->SetTrackingManager(fpTrackingManager);
310
311 fpTrackingManager->SetInteractivity(fpTrackingInteractivity);
312
313 // fpMasterStepProcessor = new G4ITStepProcessor();
314 // fpMasterStepProcessor->SetTrackingManager(fpTrackingManager);
315 // fpStepProcessor = fpMasterStepProcessor ;
316 //______________________________________________________________
317
318 if(fUsePreDefinedTimeSteps)
319 {
320 if(fpUserTimeSteps == 0) // Extra checking, is it necessary ?
321 {
322 G4ExceptionDescription exceptionDescription;
323 exceptionDescription
324 << "You are asking to use user defined steps but you did not give any.";
325 G4Exception("G4Scheduler::FindUserPreDefinedTimeStep",
326 "Scheduler004",
328 exceptionDescription);
329 return; // makes coverity happy
330 }
331 }
332
333// if (fComputeTimeStep)
334// {
335// if (fpModelProcessor == 0)
336// {
337// G4ExceptionDescription exceptionDescription;
338// exceptionDescription
339// << "There is no G4ITModelProcessor to handle IT reaction. ";
340// exceptionDescription
341// << "You probably did not initialize the G4Scheduler. ";
342// exceptionDescription
343// << "Just do G4Scheduler::Instance()->Initialize(); ";
344// exceptionDescription << " but only after initializing the run manager.";
345// G4Exception("G4Scheduler::CalculateMinStep", "ITScheduler005",
346// FatalErrorInArgument, exceptionDescription);
347// return; // makes coverity happy
348// }
349// }
350
351 fInitialized = true;
352}
void SetModelHandler(G4ITModelHandler *)
void SetTrackingManager(G4ITTrackingManager *trackingManager)
void SetTrackingManager(G4ITTrackingManager *trackMan)
void SetInteractivity(G4ITTrackingInteractivity *)

Referenced by ForceReinitialization(), G4DNAChemistryManager::InitializeThread(), Process(), and G4SchedulerMessenger::SetNewValue().

◆ Instance()

◆ IsInitialized()

bool G4Scheduler::IsInitialized ( )
inline

Definition at line 282 of file G4Scheduler.hh.

283{
284 return fInitialized;
285}

Referenced by G4SchedulerMessenger::GetCurrentValue().

◆ IsRunning()

bool G4Scheduler::IsRunning ( )
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 107 of file G4Scheduler.hh.

107{return fRunning;}

◆ Notify()

G4bool G4Scheduler::Notify ( G4ApplicationState  requestedState)
virtual

Implements G4VStateDependent.

Definition at line 108 of file G4Scheduler.cc.

109{
110 if(requestedState == G4State_Quit)
111 {
112 if(fVerbose >= 4)
113 {
114 G4cout << "G4Scheduler received G4State_Quit" << G4endl;
115 }
116 Clear();
117 //DeleteInstance();
118 }
119 return true;
120}
@ G4State_Quit

◆ PrintWhyDoYouStop()

void G4Scheduler::PrintWhyDoYouStop ( )
protected

Definition at line 558 of file G4Scheduler.cc.

559{
560#ifdef G4VERBOSE
561 if(fWhyDoYouStop)
562 {
563 G4cout << "G4Scheduler has reached a stage: it might be"
564 " a transition or the end"
565 << G4endl;
566
567 G4bool normalStop = false;
568
569 if(fGlobalTime >= fStopTime)
570 {
571 G4cout << "== G4Scheduler: I stop because I reached the stop time : "
572 << G4BestUnit(fStopTime,"Time") << " =="<< G4endl;
573 normalStop = true;
574 }
575 if(fTrackContainer.MainListsNOTEmpty() == false) // is empty
576 {
577 G4cout << "G4Scheduler: I stop because the current main list of tracks "
578 "is empty"
579 << G4endl;
580 normalStop = true;
581 }
582 if(fMaxSteps == -1 ? false : fNbSteps >= fMaxSteps)
583 {
584 G4cout << "G4Scheduler: I stop because I reached the maximum allowed "
585 "number of steps=" << fMaxSteps
586 << G4endl;
587 normalStop = true;
588 }
589 if(fContinue && normalStop == false)
590 {
591 G4cout << "G4Scheduler: It might be that I stop because "
592 "I have been told so. You may check "
593 "member fContinue and usage of the method G4Scheduler::Stop()."
594 << G4endl;
595 }
596 }
597#endif
598}
bool G4bool
Definition: G4Types.hh:86

Referenced by DoProcess().

◆ Process()

void G4Scheduler::Process ( )
virtual

Reimplemented from G4VScheduler.

Definition at line 376 of file G4Scheduler.cc.

377{
378
379#ifdef G4VERBOSE
380 if(fVerbose)
381 {
382 G4cout << "*** G4Scheduler starts processing " << G4endl;
383 if(fVerbose > 2)
384 G4cout << "___________________________________________"
385 "___________________________" << G4endl;
386 }
387#endif
388
389 if (fInitialized == false) Initialize();
390
391 // fpTrackingManager->Initialize();
392 fpModelProcessor->Initialize();
393 fpStepProcessor->Initialize();
394
395 // TODO
396 // fpMasterModelProcessor->Initialize();
397 // fpMasterStepProcessor->Initialize();
398
399 if (fpGun) fpGun->DefineTracks();
400
401 if (fpTrackingInteractivity) fpTrackingInteractivity->Initialize();
402
403 // ___________________
404 fRunning = true;
405 Reset();
406
407 if (fpUserTimeStepAction) fpUserTimeStepAction->StartProcessing();
408
409#ifdef G4VERBOSE
410 G4bool trackFound = false;
411 G4IosFlagsSaver iosfs(G4cout);
412 G4cout.precision(5);
413#endif
414
415 //===========================================================================
416 // By default, before the G4Scheduler is launched, the tracks are pushed to
417 // the delayed lists
418 //===========================================================================
419
420 if(fTrackContainer.DelayListsNOTEmpty())
421 {
422 fStartTime = fTrackContainer.GetNextTime();
423#ifdef G4VERBOSE
424 trackFound = true;
425 G4Timer localtimer;
426 if(fVerbose>1)
427 {
428 localtimer.Start();
429 }
430#endif
432#ifdef G4VERBOSE
433 if(fVerbose>1)
434 {
435 localtimer.Stop();
436 G4cout << "G4Scheduler: process time= "<< localtimer << G4endl;
437 }
438#endif
439 }
440
441// //---------------------------------
442// // TODO: This could be removed ?
443// if(fTrackContainer.MainListsNOTEmpty()
444// && (fMaxSteps == -1 ? true : fNbSteps < fMaxSteps)
445// && fGlobalTime < fEndTime
446// && fContinue == true)
447//{
448//#ifdef G4VERBOSE
449// trackFound = true;
450//#endif
451// DoProcess();
452//}
453// //---------------------------------
454
455#ifdef G4VERBOSE
456 if(fVerbose)
457 {
458 if(trackFound)
459 {
460 G4cout << "*** G4Scheduler ends at time : "
461 << G4BestUnit(fGlobalTime,"Time") << G4endl;
462 G4cout << "___________________________________" << G4endl;
463 }
464 else
465 {
466 G4cout << "*** G4Scheduler did not start because no "
467 "track was found to be processed"<< G4endl;
468 G4cout << "___________________________________" << G4endl;
469 }
470 }
471#endif
472
473 fRunning = false;
474
475 if (fpUserTimeStepAction) fpUserTimeStepAction->EndProcessing();
476
477 // ___________________
478 EndTracking();
479 ClearList();
480
481 Reset();
482
483 if (fpTrackingInteractivity) fpTrackingInteractivity->Finalize();
484}
virtual void DefineTracks()
Definition: G4ITGun.hh:56
virtual void Initialize()
void SynchronizeTracks()
Definition: G4Scheduler.cc:497
void EndTracking()
void Reset()
Definition: G4Scheduler.cc:356
void Stop()
void Start()
virtual void StartProcessing()

Referenced by G4DNAChemistryManager::Run(), and G4SchedulerMessenger::SetNewValue().

◆ RegisterModel()

void G4Scheduler::RegisterModel ( G4VITStepModel model,
double  time 
)
virtual

Reimplemented from G4VScheduler.

Definition at line 274 of file G4Scheduler.cc.

275{
276 fpModelHandler->RegisterModel(model, time);
277}
void RegisterModel(G4VITStepModel *pModel, G4double globalTime)

◆ Reset()

void G4Scheduler::Reset ( )
virtual

Reimplemented from G4VScheduler.

Definition at line 356 of file G4Scheduler.cc.

357{
358 fStartTime = 0;
359 fUserUpperTimeLimit = -1;
360 fTimeStep = DBL_MAX;
361 fTSTimeStep = DBL_MAX;
362 fILTimeStep = DBL_MAX;
363 fPreviousTimeStep = DBL_MAX;
364 fGlobalTime = -1;
365 fInteractionStep = true;
366 fITStepStatus = eUndefined;
367 fZeroTimeCount = 0;
368
369 fNbSteps = 0;
370 fContinue = true;
371 // fReactingTracks.clear();
372 fReactionSet->CleanAllReaction();
373}
@ eUndefined
void CleanAllReaction()

Referenced by Process().

◆ SetDefaultTimeStep()

void G4Scheduler::SetDefaultTimeStep ( double  timeStep)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 345 of file G4Scheduler.hh.

346{
347 fDefaultMinTimeStep = timeStep;
348}

◆ SetEndTime()

void G4Scheduler::SetEndTime ( const double  __endtime)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 292 of file G4Scheduler.hh.

293{
294 fEndTime = __endtime;
295}

Referenced by G4SchedulerMessenger::SetNewValue().

◆ SetGun()

void G4Scheduler::SetGun ( G4ITGun gun)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 417 of file G4Scheduler.hh.

418{
419 fpGun = gun;
420}

Referenced by G4DNAChemistryManager::SetGun().

◆ SetInteractivity()

void G4Scheduler::SetInteractivity ( G4ITTrackingInteractivity interactivity)
virtual

Reimplemented from G4VScheduler.

Definition at line 1103 of file G4Scheduler.cc.

1104{
1105 fpTrackingInteractivity = interactivity;
1106 if(fpTrackingManager)
1107 {
1108 fpTrackingManager->SetInteractivity(fpTrackingInteractivity);
1109 }
1110
1111 //G4MIWorkspace::GetWorldWorkspace()->SetTrackingInteractivity(interactivity);
1112}

◆ SetMaxNbSteps()

void G4Scheduler::SetMaxNbSteps ( G4int  maxSteps)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 320 of file G4Scheduler.hh.

321{
322 fMaxSteps = maxSteps;
323}

Referenced by G4SchedulerMessenger::SetNewValue().

◆ SetMaxTimeStep()

void G4Scheduler::SetMaxTimeStep ( double  maxTimeStep)
inline

Definition at line 183 of file G4Scheduler.hh.

184 {
185 fMaxTimeStep = maxTimeStep;
186 }

◆ SetMaxZeroTimeAllowed()

void G4Scheduler::SetMaxZeroTimeAllowed ( int  maxTimeStepAllowed)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 377 of file G4Scheduler.hh.

378{
379 fMaxNZeroTimeStepsAllowed = maxTimeStepAllowed;
380}

Referenced by G4SchedulerMessenger::SetNewValue().

◆ SetTimeSteps()

void G4Scheduler::SetTimeSteps ( std::map< double, double > *  steps)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 298 of file G4Scheduler.hh.

299{
300 fUsePreDefinedTimeSteps = true;
301 fpUserTimeSteps = steps;
302}

◆ SetTimeTolerance()

void G4Scheduler::SetTimeTolerance ( double  time)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 387 of file G4Scheduler.hh.

388{
389 fTimeTolerance = time;
390}

Referenced by G4SchedulerMessenger::SetNewValue().

◆ SetUserAction()

void G4Scheduler::SetUserAction ( G4UserTimeStepAction userITAction)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 356 of file G4Scheduler.hh.

357{
358 fpUserTimeStepAction = userITAction;
359}

◆ SetVerbose()

void G4Scheduler::SetVerbose ( int  verbose)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 366 of file G4Scheduler.hh.

367{
368 fVerbose = verbose;
369}

Referenced by G4SchedulerMessenger::SetNewValue().

◆ Stepping()

void G4Scheduler::Stepping ( )
protected

Definition at line 649 of file G4Scheduler.cc.

650{
651 fTimeStep = fMaxTimeStep;
652
653 fTSTimeStep = DBL_MAX;
654 fILTimeStep = DBL_MAX;
655
656 fInteractionStep = false;
657 fReachedUserTimeLimit = false;
658
659 fITStepStatus = eUndefined;
660
661 // Start of step
662#ifdef G4VERBOSE
663 if (fVerbose > 2)
664 {
665#ifdef USE_COLOR
666 G4cout << LIGHT_RED;
667#endif
668 G4cout << "*** Start Of Step N°" << fNbSteps + 1 << " ***" << G4endl;
669 G4cout << "Current Global time : " << G4BestUnit(fGlobalTime, "Time")
670 <<G4endl;
671#ifdef USE_COLOR
673#endif
674 }
675#endif
676
677#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
678 MemStat mem_first, mem_second, mem_diff;
679#endif
680
681#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
682 mem_first = MemoryUsage();
683#endif
684
685 fDefinedMinTimeStep = GetLimitingTimeStep();
686
687 if (fUsePreDefinedTimeSteps)
688 {
689#ifdef G4VERBOSE
690 if (fVerbose > 2)
691 {
692#ifdef USE_COLOR
693 G4cout << LIGHT_RED;
694#endif
695 G4cout << "*** At time : " << G4BestUnit(fGlobalTime, "Time")
696 << " the chosen user time step is : "
697 << G4BestUnit(fDefinedMinTimeStep, "Time") << " ***" << G4endl;
698#ifdef USE_COLOR
700#endif
701 }
702#endif
703 }
704
705 if (fpModelProcessor->GetComputeTimeStep()) // fComputeTimeStep)
706 {
707 fTSTimeStep = fpModelProcessor->CalculateMinTimeStep(fGlobalTime,
708 fDefinedMinTimeStep);
709 // => at least N (N = nb of tracks) loops
710 }
711 else if(fUseDefaultTimeSteps)
712 {
713 fTSTimeStep = fDefinedMinTimeStep;
714 }
715
716#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
717 mem_second = MemoryUsage();
718 mem_diff = mem_second-mem_first;
719 G4cout << "|| MEM || After computing TS, diff is : " << mem_diff << G4endl;
720#endif
721
722#ifdef G4VERBOSE
723 if (fVerbose > 2)
724 {
725#ifdef USE_COLOR
726 G4cout << LIGHT_RED;
727#endif
728 G4cout << "*** Time stepper returned : " << G4BestUnit(fTSTimeStep, "Time")
729 << " ***" << G4endl;
730#ifdef USE_COLOR
732#endif
733 }
734#endif
735
736#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
737 mem_first = MemoryUsage();
738#endif
739
740 // Call IL even if fTSTimeStep == 0
741 // if fILTimeStep == 0 give the priority to DoIt processes
742
743 fILTimeStep =
744 fpStepProcessor->ComputeInteractionLength(fPreviousTimeStep);
745 // => at least N loops
746 // All process returns the physical step of interaction
747 // The transportation process calculates the corresponding
748 // time step
749
750#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
751 mem_second = MemoryUsage();
752 mem_diff = mem_second-mem_first;
753 G4cout << "|| MEM || After IL, diff is : " << mem_diff << G4endl;
754#endif
755
756#ifdef G4VERBOSE
757 if (fVerbose > 2)
758 {
759#ifdef USE_COLOR
760 G4cout << LIGHT_RED;
761#endif
762 G4cout << "*** The minimum time returned by the processes is : "
763 << G4BestUnit(fILTimeStep, "Time") << " ***" << G4endl;
764#ifdef USE_COLOR
766#endif
767 }
768#endif
769
770#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
771 mem_first = MemoryUsage();
772#endif
773
774 if (fILTimeStep <= fTSTimeStep)
775 // Give the priority to the IL
776 {
777 fInteractionStep = true;
778 fReactionSet->CleanAllReaction();
779 fTimeStep = fILTimeStep;
780 fITStepStatus = eInteractionWithMedium;
781 fpStepProcessor->PrepareLeadingTracks();
782 }
783 else
784 {
785 fInteractionStep = false;
786 fpStepProcessor->ResetLeadingTracks();
787 fTimeStep = fTSTimeStep;
788 fITStepStatus = eCollisionBetweenTracks;
789 }
790
791 if (fGlobalTime + fTimeStep > fStopTime)
792 // This check is done at every time step
793 {
794 fTimeStep = fStopTime - fGlobalTime;
795 fITStepStatus = eInteractionWithMedium; // ie: transportation
796 fInteractionStep = true;
797 fReactionSet->CleanAllReaction();
798 fpStepProcessor->ResetLeadingTracks();
799 }
800
801 if (fTimeStep == 0) // < fTimeTolerance)
802 {
803 ++fZeroTimeCount;
804 if (fZeroTimeCount >= fMaxNZeroTimeStepsAllowed)
805 {
806 G4ExceptionDescription exceptionDescription;
807
808 exceptionDescription << "Too many zero time steps were detected. ";
809 exceptionDescription << "The simulation is probably stuck. ";
810 exceptionDescription
811 << "The maximum number of zero time steps is currently : "
812 << fMaxNZeroTimeStepsAllowed;
813 exceptionDescription << ".";
814
815 G4Exception("G4Scheduler::Stepping",
816 "SchedulerNullTimeSteps",
818 exceptionDescription);
819 }
820 }
821 else
822 {
823 fZeroTimeCount = 0;
824 }
825
826 fReachedUserTimeLimit =
827 ((fTimeStep <= fDefinedMinTimeStep) || ((fTimeStep > fDefinedMinTimeStep)
828 && fabs(fTimeStep - fDefinedMinTimeStep) < fTimeTolerance)) ?
829 true : false;
830
831 if (fpUserTimeStepAction) fpUserTimeStepAction->UserPreTimeStepAction();
832 // TODO: pre/post
833
834#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
835 mem_second = MemoryUsage();
836 mem_diff = mem_second-mem_first;
837 G4cout << "|| MEM || After LeadingTracks and UserPreTimeStepAction: "
838 << mem_diff << G4endl;
839#endif
840
841#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
842 mem_first = MemoryUsage();
843#endif
844
845
846 fGlobalTime += fTimeStep;
847
848 // if fTSTimeStep > 0 => still need to call the transportation process
849 // if fILTimeStep < fTSTimeStep => call only DoIt processes, no reactions
850 // if fILTimeStep == fTSTimeStep => give the priority to the DoIt processes
851 if (fTSTimeStep > 0 || fILTimeStep <= fTSTimeStep)
852 {
853 // G4cout << "Will call DoIT" << G4endl;
854 fpStepProcessor->DoIt(fTimeStep);
855 // fTrackContainer.MergeSecondariesWithMainList();
856 // fTrackContainer.KillTracks(); // remove ?
857 }
858 // else
859 // {
860 // G4cout << "fTSTimeStep : " << fTSTimeStep << G4endl;
861 // G4cout << "fILTimeStep : " << fILTimeStep << G4endl;
862 // }
863
864#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
865 mem_second = MemoryUsage();
866 mem_diff = mem_second-mem_first;
867 G4cout << "|| MEM || After DoIT, diff is : " << mem_diff << G4endl;
868#endif
869
870#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
871 mem_first = MemoryUsage();
872#endif
873
874 fpModelProcessor->ComputeTrackReaction(fITStepStatus,
875 fGlobalTime,
876 fTimeStep,
877 fPreviousTimeStep,
878 fReachedUserTimeLimit,
879 fTimeTolerance,
880 fpUserTimeStepAction,
881 fVerbose);
882
883 ++fNbSteps;
884
885 if (fpUserTimeStepAction)
886 {
887 fpUserTimeStepAction->UserPostTimeStepAction();
888 }
889
890 fPreviousTimeStep = fTimeStep;
891
892#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
893 mem_second = MemoryUsage();
894 mem_diff = mem_second-mem_first;
895 G4cout << "|| MEM || After computing reactions + UserPostTimeStepAction, "
896 "diff is : " << mem_diff << G4endl;
897#endif
898
899 // End of step
900#ifdef G4VERBOSE
901 if (fVerbose >= 2)
902 {
903#ifdef USE_COLOR
904 G4cout << LIGHT_RED;
905#endif
906
907 G4String interactionType;
908 GetCollisionType(interactionType);
909
910 std::stringstream finalOutput;
911
912 finalOutput << "*** End of step N°" << fNbSteps
913 << "\t T_i= " << G4BestUnit(fGlobalTime-fTimeStep, "Time")
914 << "\t dt= " << G4BestUnit(fTimeStep, "Time")
915 << "\t T_f= " << G4BestUnit(fGlobalTime, "Time")
916 << "\t " << interactionType
917 << G4endl;
918
919 if(fVerbose>2)
920 {
921 if(fReachedUserTimeLimit)
922 {
923 finalOutput << "It has also reached the user time limit" << G4endl;
924 }
925 finalOutput << "_______________________________________________________________"
926 "_______"<< G4endl;
927 }
928
929 G4cout << finalOutput.str();
930
931#ifdef USE_COLOR
933#endif
934 }
935#endif
936
937}
#define LIGHT_RED
Definition: G4Scheduler.cc:83
#define RESET_COLOR
Definition: G4Scheduler.cc:86
G4double CalculateMinTimeStep(G4double currentGlobalTime, G4double definedMinTimeStep)
void ComputeTrackReaction(G4ITStepStatus fITStepStatus, G4double fGlobalTime, G4double currentTimeStep, G4double previousTimeStep, G4bool reachedUserTimeLimit, G4double fTimeTolerance, G4UserTimeStepAction *fpUserTimeStepAction, G4int fVerbose)
bool GetComputeTimeStep() const
G4double ComputeInteractionLength(double previousTimeStep)
void DoIt(double timeStep)
double GetLimitingTimeStep() const
Definition: G4Scheduler.cc:940
void GetCollisionType(G4String &interactionType)
virtual void UserPostTimeStepAction()
virtual void UserPreTimeStepAction()

Referenced by DoProcess().

◆ Stop()

void G4Scheduler::Stop ( )
inline

Definition at line 407 of file G4Scheduler.hh.

408{
409 fContinue = false;
410}

Referenced by G4DNAIRTMoleculeEncounterStepper::CalculateMinTimeStep().

◆ SynchronizeTracks()

void G4Scheduler::SynchronizeTracks ( )
protected

Definition at line 497 of file G4Scheduler.cc.

498{
499// if(fTrackContainer.WaitingListsNOTEmpty())
500// {
501// G4ExceptionDescription exceptionDescription;
502// exceptionDescription
503// << "There is a waiting track list (fpWaitingList != 0).";
504// exceptionDescription
505// << " When G4Scheduler::SynchronizeTracks() is called, ";
506// exceptionDescription
507// << "no more tracks should remain in the fpWaitingList.";
508// G4Exception("G4Scheduler::SynchronizeTracks", "ITScheduler002",
509// FatalErrorInArgument, exceptionDescription);
510// }
511
512 // Backup main list and time feature
513 // Reminder : the main list here, should
514 // have the biggest global time
515 // fTrackContainer.MoveMainToWaitingList();
516 // TODO: not yet supported
517
518 fTmpGlobalTime = fGlobalTime;
519 //fTmpEndTime = fEndTime;
520
521 fGlobalTime = fTrackContainer.GetNextTime();
522 G4double tmpGlobalTime = fGlobalTime;
523
524 double nextWatchedTime = -1;
525 bool carryOn = true;
526
527 while(fTrackContainer.MergeNextTimeToMainList(tmpGlobalTime) && carryOn)
528 {
529// assert(tmpGlobalTime == fGlobalTime);
530 fStopTime = min(fTrackContainer.GetNextTime(), fEndTime);
531 while((nextWatchedTime = GetNextWatchedTime()) < fTrackContainer.GetNextTime()
532 && (carryOn = CanICarryOn()))
533 {
534 fStopTime = min(nextWatchedTime, fEndTime);
535 DoProcess();
536 }
537
538 carryOn = CanICarryOn();
539
540 if(nextWatchedTime > fEndTime && carryOn)
541 {
542 fStopTime = min(fTrackContainer.GetNextTime(), fEndTime);
543 DoProcess();
544 }
545 }
546}
double G4double
Definition: G4Types.hh:83
bool MergeNextTimeToMainList(double &time)
double GetNextWatchedTime() const
Definition: G4Scheduler.cc:488
bool CanICarryOn()
Definition: G4Scheduler.cc:550
void DoProcess()
Definition: G4Scheduler.cc:602
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

Referenced by Process().

◆ UseDefaultTimeSteps()

void G4Scheduler::UseDefaultTimeSteps ( G4bool  flag)
inline

Definition at line 432 of file G4Scheduler.hh.

433{
434 fUseDefaultTimeSteps = flag;
435}

Referenced by G4SchedulerMessenger::SetNewValue().

◆ WhyDoYouStop()

void G4Scheduler::WhyDoYouStop ( )
inline

Definition at line 427 of file G4Scheduler.hh.

428{
429 fWhyDoYouStop = true;
430}

Referenced by G4SchedulerMessenger::SetNewValue().


The documentation for this class was generated from the following files: