Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4Scheduler Class Reference

#include <G4Scheduler.hh>

Inheritance diagram for G4Scheduler:
Collaboration 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
 
G4ITTrackingManagerGetTrackingManager () const
 
- Public Member Functions inherited from G4VStateDependent
 G4VStateDependent (G4bool bottom=false)
 
virtual ~G4VStateDependent ()
 
G4int operator== (const G4VStateDependent &right) const
 
G4int operator!= (const G4VStateDependent &right) const
 

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 88 of file G4Scheduler.hh.

Constructor & Destructor Documentation

G4Scheduler::~G4Scheduler ( )
protectedvirtual

Definition at line 201 of file G4Scheduler.cc.

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

Here is the call graph for this function:

Member Function Documentation

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

Reimplemented from G4VScheduler.

Definition at line 310 of file G4Scheduler.hh.

311 {
312  if (fpUserTimeSteps == 0)
313  {
314  fpUserTimeSteps = new std::map<double, double>();
315  fUsePreDefinedTimeSteps = true;
316  }
317 
318  (*fpUserTimeSteps)[startingTime] = timeStep;
319 }
void G4Scheduler::AddWatchedTime ( double  time)
inline

Definition at line 177 of file G4Scheduler.hh.

178  {
179  fWatchedTimes.insert(time);
180  }
G4bool G4Scheduler::AreDefaultTimeStepsUsed ( )
inline

Definition at line 443 of file G4Scheduler.hh.

444 {
445  return (fUseDefaultTimeSteps == false && fUsePreDefinedTimeSteps == false);
446 }

Here is the caller graph for this function:

bool G4Scheduler::CanICarryOn ( )
protected

Definition at line 551 of file G4Scheduler.cc.

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

Here is the caller graph for this function:

void G4Scheduler::Clear ( )

Definition at line 217 of file G4Scheduler.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Scheduler::ClearList ( )

Definition at line 265 of file G4Scheduler.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Scheduler::DeleteInstance ( )
static

DeleteInstance should be used instead of the destructor

Definition at line 124 of file G4Scheduler.cc.

125 {
126  if(fgScheduler)
127  {
128  delete fgScheduler;
129  }
130 }
void G4Scheduler::DoProcess ( )
protected

Definition at line 603 of file G4Scheduler.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Scheduler::EndTracking ( )

Definition at line 1062 of file G4Scheduler.cc.

1063 {
1064  if(fRunning)
1065  {
1066  G4ExceptionDescription exceptionDescription;
1067  exceptionDescription
1068  << "End tracking is called while G4Scheduler is still running."
1069  << G4endl;
1070 
1071  G4Exception("G4Scheduler::EndTracking",
1072  "Scheduler017",
1074  exceptionDescription);
1075  }
1076 
1077  fTrackContainer.MergeSecondariesWithMainList();
1078 
1079  if (fTrackContainer.MainListsNOTEmpty())
1080  {
1081  G4TrackManyList* mainList = fTrackContainer.GetMainList();
1082  G4TrackManyList::iterator it = mainList->begin();
1083  G4TrackManyList::iterator end = mainList->end();
1084  for (; it != end; ++it)
1085  {
1086  fpTrackingManager->EndTrackingWOKill(*it);
1087  }
1088  }
1089 
1090  if (fTrackContainer.SecondaryListsNOTEmpty()) // should be empty
1091  {
1092  G4TrackManyList* secondaries = fTrackContainer.GetSecondariesList();
1093  G4TrackManyList::iterator it = secondaries->begin();
1094  G4TrackManyList::iterator end = secondaries->end();
1095 
1096  for (; it != end; ++it)
1097  {
1098  fpTrackingManager->EndTrackingWOKill(*it);
1099  }
1100  }
1101 }
void MergeSecondariesWithMainList()
G4TrackManyList * GetSecondariesList()
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void EndTrackingWOKill(G4Track *)
bool SecondaryListsNOTEmpty()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
G4TrackList * GetMainList(Key)

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Scheduler::FindUserPreDefinedTimeStep ( )
protected

Definition at line 1007 of file G4Scheduler.cc.

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

Here is the call graph for this function:

void G4Scheduler::ForceReinitialization ( )

Definition at line 1116 of file G4Scheduler.cc.

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

Here is the call graph for this function:

void G4Scheduler::GetCollisionType ( G4String interactionType)

Definition at line 1147 of file G4Scheduler.cc.

1148 {
1149  switch(fITStepStatus)
1150  {
1152  interactionType = "eInteractionWithMedium";
1153  break;
1155  interactionType = "eCollisionBetweenTracks";
1156  break;
1157  default:
1158  interactionType = "eCollisionBetweenTracks";
1159  break;
1160  }
1161 }

Here is the caller graph for this function:

G4double G4Scheduler::GetEndTime ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 341 of file G4Scheduler.hh.

342 {
343  return fEndTime;
344 }

Here is the caller graph for this function:

G4double G4Scheduler::GetGlobalTime ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 356 of file G4Scheduler.hh.

357 {
358  return fGlobalTime;
359 }

Here is the caller graph for this function:

G4ITGun * G4Scheduler::GetGun ( )
inline

Definition at line 428 of file G4Scheduler.hh.

429 {
430  return fpGun;
431 }
G4ITTrackingInteractivity * G4Scheduler::GetInteractivity ( )
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 418 of file G4Scheduler.hh.

419 {
420  return fpTrackingInteractivity;
421 }
double G4Scheduler::GetLimitingTimeStep ( ) const
virtual

Reimplemented from G4VScheduler.

Definition at line 941 of file G4Scheduler.cc.

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

Here is the caller graph for this function:

G4int G4Scheduler::GetMaxNbSteps ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 331 of file G4Scheduler.hh.

332 {
333  return fMaxSteps;
334 }

Here is the caller graph for this function:

double G4Scheduler::GetMaxTimeStep ( ) const
inline

Definition at line 189 of file G4Scheduler.hh.

190  {
191  return fMaxTimeStep;
192  }
int G4Scheduler::GetMaxZeroTimeAllowed ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 388 of file G4Scheduler.hh.

389 {
390  return fMaxNZeroTimeStepsAllowed;
391 }

Here is the caller graph for this function:

G4ITModelHandler * G4Scheduler::GetModelHandler ( )
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 293 of file G4Scheduler.hh.

294 {
295  return fpModelHandler;
296 }
G4int G4Scheduler::GetNbSteps ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 321 of file G4Scheduler.hh.

322 {
323  return fNbSteps;
324 }
double G4Scheduler::GetNextWatchedTime ( ) const

Definition at line 489 of file G4Scheduler.cc.

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

Here is the caller graph for this function:

size_t G4Scheduler::GetNTracks ( )
virtual

Definition at line 1141 of file G4Scheduler.cc.

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

Here is the call graph for this function:

G4double G4Scheduler::GetPreviousTimeStep ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 403 of file G4Scheduler.hh.

404 {
405  return fPreviousTimeStep;
406 }
G4double G4Scheduler::GetStartTime ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 336 of file G4Scheduler.hh.

337 {
338  return fStartTime;
339 }
G4ITStepStatus G4Scheduler::GetStatus ( ) const
inline

Definition at line 408 of file G4Scheduler.hh.

409 {
410  return fITStepStatus;
411 }
G4double G4Scheduler::GetTimeStep ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 346 of file G4Scheduler.hh.

347 {
348  return fTimeStep;
349 }
double G4Scheduler::GetTimeTolerance ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 398 of file G4Scheduler.hh.

399 {
400  return fTimeTolerance;
401 }

Here is the caller graph for this function:

G4ITTrackingManager* G4Scheduler::GetTrackingManager ( ) const
inline

Definition at line 194 of file G4Scheduler.hh.

195  {
196  return fpTrackingManager;
197  }
G4UserTimeStepAction * G4Scheduler::GetUserTimeStepAction ( ) const
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 367 of file G4Scheduler.hh.

368 {
369  return fpUserTimeStepAction;
370 }
int G4Scheduler::GetVerbose ( ) const
inline

Definition at line 377 of file G4Scheduler.hh.

378 {
379  return fVerbose;
380 }

Here is the caller graph for this function:

void G4Scheduler::Initialize ( )
virtual

Reimplemented from G4VScheduler.

Definition at line 282 of file G4Scheduler.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

G4Scheduler * G4Scheduler::Instance ( )
static

Definition at line 102 of file G4Scheduler.cc.

103 {
104  if(fgScheduler == 0) fgScheduler = new G4Scheduler();
105  return fgScheduler;
106 }

Here is the caller graph for this function:

bool G4Scheduler::IsInitialized ( )
inline

Definition at line 288 of file G4Scheduler.hh.

289 {
290  return fInitialized;
291 }

Here is the caller graph for this function:

bool G4Scheduler::IsRunning ( )
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 108 of file G4Scheduler.hh.

108 {return fRunning;}

Here is the caller graph for this function:

G4bool G4Scheduler::Notify ( G4ApplicationState  requestedState)
virtual

Implements G4VStateDependent.

Definition at line 109 of file G4Scheduler.cc.

110 {
111  if(requestedState == G4State_Quit)
112  {
113  if(fVerbose >= 4)
114  {
115  G4cout << "G4Scheduler received G4State_Quit" << G4endl;
116  }
117  Clear();
118  //DeleteInstance();
119  }
120  return true;
121 }
G4GLOB_DLL std::ostream G4cout
void Clear()
Definition: G4Scheduler.cc:217
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4Scheduler::PrintWhyDoYouStop ( )
protected

Definition at line 559 of file G4Scheduler.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Scheduler::Process ( )
virtual

Reimplemented from G4VScheduler.

Definition at line 377 of file G4Scheduler.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VScheduler.

Definition at line 275 of file G4Scheduler.cc.

276 {
277  fpModelHandler->RegisterModel(model, time);
278 }
void RegisterModel(G4VITStepModel *aModel, const G4double globalTime)

Here is the call graph for this function:

void G4Scheduler::Reset ( )
virtual

Reimplemented from G4VScheduler.

Definition at line 357 of file G4Scheduler.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Scheduler::SetDefaultTimeStep ( double  timeStep)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 351 of file G4Scheduler.hh.

352 {
353  fDefaultMinTimeStep = timeStep;
354 }
void G4Scheduler::SetEndTime ( const double  __endtime)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 298 of file G4Scheduler.hh.

299 {
300  fEndTime = __endtime;
301 }

Here is the caller graph for this function:

void G4Scheduler::SetGun ( G4ITGun gun)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 423 of file G4Scheduler.hh.

424 {
425  fpGun = gun;
426 }

Here is the caller graph for this function:

void G4Scheduler::SetInteractivity ( G4ITTrackingInteractivity interactivity)
virtual

Reimplemented from G4VScheduler.

Definition at line 1104 of file G4Scheduler.cc.

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

Here is the call graph for this function:

void G4Scheduler::SetMaxNbSteps ( G4int  maxSteps)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 326 of file G4Scheduler.hh.

327 {
328  fMaxSteps = maxSteps;
329 }

Here is the caller graph for this function:

void G4Scheduler::SetMaxTimeStep ( double  maxTimeStep)
inline

Definition at line 184 of file G4Scheduler.hh.

185  {
186  fMaxTimeStep = maxTimeStep;
187  }
void G4Scheduler::SetMaxZeroTimeAllowed ( int  maxTimeStepAllowed)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 383 of file G4Scheduler.hh.

384 {
385  fMaxNZeroTimeStepsAllowed = maxTimeStepAllowed;
386 }

Here is the caller graph for this function:

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

Reimplemented from G4VScheduler.

Definition at line 304 of file G4Scheduler.hh.

305 {
306  fUsePreDefinedTimeSteps = true;
307  fpUserTimeSteps = steps;
308 }
void G4Scheduler::SetTimeTolerance ( double  time)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 393 of file G4Scheduler.hh.

394 {
395  fTimeTolerance = time;
396 }

Here is the caller graph for this function:

void G4Scheduler::SetUserAction ( G4UserTimeStepAction userITAction)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 362 of file G4Scheduler.hh.

363 {
364  fpUserTimeStepAction = userITAction;
365 }
void G4Scheduler::SetVerbose ( int  verbose)
inlinevirtual

Reimplemented from G4VScheduler.

Definition at line 372 of file G4Scheduler.hh.

373 {
374  fVerbose = verbose;
375 }

Here is the caller graph for this function:

void G4Scheduler::Stepping ( )
protected

Definition at line 650 of file G4Scheduler.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Scheduler::Stop ( )
inline

Definition at line 413 of file G4Scheduler.hh.

414 {
415  fContinue = false;
416 }
void G4Scheduler::SynchronizeTracks ( )
protected

Definition at line 498 of file G4Scheduler.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4Scheduler::UseDefaultTimeSteps ( G4bool  flag)
inline

Definition at line 438 of file G4Scheduler.hh.

439 {
440  fUseDefaultTimeSteps = flag;
441 }

Here is the caller graph for this function:

void G4Scheduler::WhyDoYouStop ( )
inline

Definition at line 433 of file G4Scheduler.hh.

434 {
435  fWhyDoYouStop = true;
436 }

Here is the caller graph for this function:


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