Geant4 11.4.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

 G4Scheduler (const G4Scheduler &)=delete
G4Scheduleroperator= (const G4Scheduler &)=delete
G4bool Notify (G4ApplicationState requestedState) override
void RegisterModel (G4VITStepModel *, G4double) override
void Initialize () override
void ForceReinitialization ()
G4bool IsInitialized ()
G4bool IsRunning () override
void Reset () override
void Process () override
void ClearList ()
void SetGun (G4ITGun *) override
G4ITGunGetGun ()
void Stop ()
void Clear ()
void EndTracking ()
void SetEndTime (const G4double) override
void SetTimeTolerance (G4double) override
G4double GetTimeTolerance () const override
void SetMaxZeroTimeAllowed (G4int) override
G4int GetMaxZeroTimeAllowed () const override
G4ITModelHandlerGetModelHandler () override
void SetTimeSteps (std::map< G4double, G4double > *) override
void AddTimeStep (G4double, G4double) override
void SetDefaultTimeStep (G4double) override
G4double GetLimitingTimeStep () const override
G4int GetNbSteps () const override
void SetMaxNbSteps (G4int) override
G4int GetMaxNbSteps () const override
G4double GetStartTime () const override
G4double GetEndTime () const override
G4double GetTimeStep () const override
G4double GetPreviousTimeStep () const override
G4double GetGlobalTime () const override
void SetUserAction (G4UserTimeStepAction *) override
G4UserTimeStepActionGetUserTimeStepAction () const override
void UseDefaultTimeSteps (G4bool)
G4bool AreDefaultTimeStepsUsed ()
G4ITStepStatus GetStatus () const
void SetVerbose (G4int) override
G4int GetVerbose () const
void WhyDoYouStop ()
void SetInteractivity (G4ITTrackingInteractivity *) override
G4ITTrackingInteractivityGetInteractivity () override
virtual size_t GetNTracks ()
void GetCollisionType (G4String &interactionType)
void AddWatchedTime (G4double time)
G4double GetNextWatchedTime () const
void SetMaxTimeStep (G4double maxTimeStep)
G4double GetMaxTimeStep () const
G4VScavengerMaterialGetScavengerMaterial () const
void SetScavengerMaterial (std::unique_ptr< G4VScavengerMaterial > scavengerMaterial)
G4bool IsInteractionStep ()
void SetInteractionStep (G4bool InteractionStep)
void ResetScavenger (bool)
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 void NotifyDeletion (const G4Event *)
virtual void NotifyDeletion (const G4Run *)

Static Public Member Functions

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

Protected Member Functions

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

Detailed Description

G4Scheduler synchronizes (in time) track stepping

Definition at line 86 of file G4Scheduler.hh.

Constructor & Destructor Documentation

◆ ~G4Scheduler()

G4Scheduler::~G4Scheduler ( )
overrideprotected

Definition at line 169 of file G4Scheduler.cc.

170{
171 if (fpMessenger != nullptr) // is used as a flag to know whether the manager was cleared
172 {
173 Clear();
174 }
175 fgScheduler = nullptr;
176}

◆ G4Scheduler()

G4Scheduler::G4Scheduler ( const G4Scheduler & )
delete

Referenced by G4Scheduler(), Instance(), and operator=().

Member Function Documentation

◆ AddTimeStep()

void G4Scheduler::AddTimeStep ( G4double startingTime,
G4double timeStep )
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 308 of file G4Scheduler.hh.

309{
310 if (fpUserTimeSteps == nullptr) {
311 fpUserTimeSteps = new std::map<G4double, G4double>();
312 fUsePreDefinedTimeSteps = true;
313 }
314
315 (*fpUserTimeSteps)[startingTime] = timeStep;
316}

◆ AddWatchedTime()

void G4Scheduler::AddWatchedTime ( G4double time)
inline

Definition at line 176 of file G4Scheduler.hh.

176{ fWatchedTimes.insert(time); }

◆ AreDefaultTimeStepsUsed()

G4bool G4Scheduler::AreDefaultTimeStepsUsed ( )
inline

Definition at line 438 of file G4Scheduler.hh.

439{
440 return (!fUseDefaultTimeSteps && !fUsePreDefinedTimeSteps);
441}

◆ CanICarryOn()

G4bool G4Scheduler::CanICarryOn ( )
protected

Definition at line 410 of file G4Scheduler.cc.

411{
412 return fGlobalTime < fEndTime && (fMaxSteps == -1 ? true : fNbSteps < fMaxSteps) && fContinue;
413}

Referenced by SynchronizeTracks().

◆ Clear()

void G4Scheduler::Clear ( )

Definition at line 178 of file G4Scheduler.cc.

179{
180 if (fpMessenger != nullptr) {
181 delete fpMessenger;
182 fpMessenger = nullptr;
183 }
184 if (fpStepProcessor != nullptr) {
185 delete fpStepProcessor;
186 fpStepProcessor = nullptr;
187 }
188 if (fpModelProcessor != nullptr) {
189 delete fpModelProcessor;
190 fpModelProcessor = nullptr;
191 }
192
194 ClearList();
195 if (fpTrackingManager != nullptr) {
196 delete fpTrackingManager;
197 fpTrackingManager = nullptr;
198 }
199
200 if (fReactionSet != nullptr) {
201 delete fReactionSet;
202 fReactionSet = nullptr;
203 }
204
205 if (fpModelHandler != nullptr) {
206 delete fpModelHandler;
207 fpModelHandler = nullptr;
208 }
209}
static G4ITTypeManager * Instance()
Definition G4ITType.cc:57
void ReleaseRessource()
Definition G4ITType.cc:82
void ClearList()

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

◆ ClearList()

void G4Scheduler::ClearList ( )

Definition at line 213 of file G4Scheduler.cc.

214{
215 fTrackContainer.Clear();
217}
static void DeleteInstance()

Referenced by Clear(), and Process().

◆ DeleteInstance()

void G4Scheduler::DeleteInstance ( )
static

DeleteInstance should be used instead of the destructor

Definition at line 110 of file G4Scheduler.cc.

111{
112 delete fgScheduler;
113}

◆ DoProcess()

void G4Scheduler::DoProcess ( )
protected

Definition at line 457 of file G4Scheduler.cc.

458{
459 if (fpUserTimeStepAction != nullptr) fpUserTimeStepAction->NewStage();
460
461#if defined(DEBUG_MEM) && defined(DEBUG_MEM_STEPPING)
462 MemStat mem_first, mem_second, mem_diff;
463#endif
464
465#if defined(DEBUG_MEM) && defined(DEBUG_MEM_STEPPING)
466 mem_first = MemoryUsage();
467#endif
468
469 while (fGlobalTime < fStopTime && fTrackContainer.MainListsNOTEmpty()
470 && (fMaxSteps == -1 ? true : fNbSteps < fMaxSteps) && fContinue)
471 {
472 Stepping();
473
474#if defined(DEBUG_MEM) && defined(DEBUG_MEM_STEPPING)
475 mem_second = MemoryUsage();
476 mem_diff = mem_second - mem_first;
477 G4cout << "\t || MEM || After step " << fNbSteps << ", diff is : " << mem_diff << G4endl;
478#endif
479 }
480
482
483#if defined(DEBUG_MEM) && defined(DEBUG_MEM_STEPPING)
484 mem_second = MemoryUsage();
485 mem_diff = mem_second - mem_first;
486 G4cout << "\t || MEM || After stepping, diff is : " << mem_diff << G4endl;
487#endif
488
489#ifdef G4VERBOSE
490 if (fVerbose > 2)
491 G4cout << "*** G4Scheduler has finished processing a track list at time : "
492 << G4BestUnit(fGlobalTime, "Time") << G4endl;
493#endif
494}
#define G4BestUnit(a, b)
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void Stepping()
void PrintWhyDoYouStop()
MemStat MemoryUsage()
Definition G4MemStat.cc:55

Referenced by SynchronizeTracks().

◆ EndTracking()

void G4Scheduler::EndTracking ( )

Definition at line 821 of file G4Scheduler.cc.

822{
823 if (fRunning) {
824 G4ExceptionDescription exceptionDescription;
825 exceptionDescription << "End tracking is called while G4Scheduler is still running." << G4endl;
826
827 G4Exception("G4Scheduler::EndTracking", "Scheduler017", FatalErrorInArgument,
828 exceptionDescription);
829 }
830
831 while (fTrackContainer.DelayListsNOTEmpty()) {
832 auto nextTime = fTrackContainer.GetNextTime();
833 fTrackContainer.MergeNextTimeToMainList(nextTime);
834 }
835
836 fTrackContainer.MergeSecondariesWithMainList();
837
838 if (fTrackContainer.MainListsNOTEmpty()) {
839 G4TrackManyList* mainList = fTrackContainer.GetMainList();
840 G4TrackManyList::iterator it = mainList->begin();
841 G4TrackManyList::iterator end = mainList->end();
842 for (; it != end; ++it) {
843 fpTrackingManager->EndTrackingWOKill(*it);
844 }
845 }
846
847 if (fTrackContainer.SecondaryListsNOTEmpty()) // should be empty
848 {
849 G4TrackManyList* secondaries = fTrackContainer.GetSecondariesList();
850 G4TrackManyList::iterator it = secondaries->begin();
851 G4TrackManyList::iterator end = secondaries->end();
852
853 for (; it != end; ++it) {
854 fpTrackingManager->EndTrackingWOKill(*it);
855 }
856 }
857}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ManyFastLists< G4Track > G4TrackManyList
G4ManyFastLists_iterator< G4Track > iterator

Referenced by Process().

◆ FindUserPreDefinedTimeStep()

void G4Scheduler::FindUserPreDefinedTimeStep ( )
protected

Definition at line 787 of file G4Scheduler.cc.

788{
789 if (fpUserTimeSteps == nullptr) {
790 G4ExceptionDescription exceptionDescription;
791 exceptionDescription << "You are asking to use user defined steps but you did not give any.";
792 G4Exception("G4Scheduler::FindUserPreDefinedTimeStep", "Scheduler004", FatalErrorInArgument,
793 exceptionDescription);
794 return; // makes coverity happy
795 }
796 auto fpUserTimeSteps_i = fpUserTimeSteps->upper_bound(fGlobalTime);
797 auto fpUserTimeSteps_low = fpUserTimeSteps->lower_bound(fGlobalTime);
798
799 if (fpUserTimeSteps_i == fpUserTimeSteps->end()) {
800 fpUserTimeSteps_i--;
801 }
802 else if (fabs(fGlobalTime - fpUserTimeSteps_low->first) < fTimeTolerance) {
803 // Case : fGlobalTime = X picosecond
804 // and fpUserTimeSteps_low->first = X picosecond
805 // but the precision is not good enough
806 fpUserTimeSteps_i = fpUserTimeSteps_low;
807 }
808 else if (fpUserTimeSteps_i == fpUserTimeSteps_low) {
809 // "Normal" cases
810 fpUserTimeSteps_i--;
811 }
812 else {
813 fpUserTimeSteps_i = fpUserTimeSteps_low;
814 }
815
816 fDefinedMinTimeStep = fpUserTimeSteps_i->second;
817}

◆ ForceReinitialization()

void G4Scheduler::ForceReinitialization ( )

Definition at line 869 of file G4Scheduler.cc.

870{
871 fInitialized = false;
872 Initialize();
873}
void Initialize() override

◆ GetCollisionType()

void G4Scheduler::GetCollisionType ( G4String & interactionType)

Definition at line 881 of file G4Scheduler.cc.

882{
883 switch (fITStepStatus) {
885 interactionType = "eInteractionWithMedium";
886 break;
888 interactionType = "eCollisionBetweenTracks";
889 break;
890 default:
891 interactionType = "eCollisionBetweenTracks";
892 break;
893 }
894}
@ eInteractionWithMedium
@ eCollisionBetweenTracks

Referenced by Stepping().

◆ GetEndTime()

G4double G4Scheduler::GetEndTime ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 338 of file G4Scheduler.hh.

339{
340 return fEndTime;
341}

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

◆ GetGlobalTime()

G4double G4Scheduler::GetGlobalTime ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 353 of file G4Scheduler.hh.

354{
355 return fGlobalTime;
356}

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

◆ GetGun()

G4ITGun * G4Scheduler::GetGun ( )
inline

Definition at line 423 of file G4Scheduler.hh.

424{
425 return fpGun;
426}

◆ GetInteractivity()

G4ITTrackingInteractivity * G4Scheduler::GetInteractivity ( )
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 413 of file G4Scheduler.hh.

414{
415 return fpTrackingInteractivity;
416}

◆ GetLimitingTimeStep()

G4double G4Scheduler::GetLimitingTimeStep ( ) const
overridevirtual

Reimplemented from G4VScheduler.

Definition at line 751 of file G4Scheduler.cc.

752{
753 if (fpUserTimeSteps == nullptr) return fDefaultMinTimeStep;
754 if (fabs(fGlobalTime - fUserUpperTimeLimit) < fTimeTolerance) return fDefinedMinTimeStep;
755
756 auto it_fpUserTimeSteps_i = fpUserTimeSteps->upper_bound(fGlobalTime);
757 auto it_fpUserTimeSteps_low = fpUserTimeSteps->lower_bound(fGlobalTime);
758
759 if (it_fpUserTimeSteps_i == fpUserTimeSteps->end()) {
760 it_fpUserTimeSteps_i--;
761 fUserUpperTimeLimit = fStopTime;
762 }
763 else if (fabs(fGlobalTime - it_fpUserTimeSteps_low->first) < fTimeTolerance) {
764 it_fpUserTimeSteps_i = it_fpUserTimeSteps_low;
765 auto tmp_it = it_fpUserTimeSteps_low;
766 ++tmp_it;
767 if (tmp_it == fpUserTimeSteps->end()) {
768 fUserUpperTimeLimit = fStopTime;
769 }
770 else {
771 fUserUpperTimeLimit = tmp_it->first;
772 }
773 }
774 else if (it_fpUserTimeSteps_i == it_fpUserTimeSteps_low) {
775 fUserUpperTimeLimit = it_fpUserTimeSteps_i->first;
776 if (it_fpUserTimeSteps_i != fpUserTimeSteps->begin()) it_fpUserTimeSteps_i--;
777 }
778 else {
779 fUserUpperTimeLimit = it_fpUserTimeSteps_i->first;
780 it_fpUserTimeSteps_i = it_fpUserTimeSteps_low;
781 }
782 return it_fpUserTimeSteps_i->second;
783}

Referenced by Stepping().

◆ GetMaxNbSteps()

G4int G4Scheduler::GetMaxNbSteps ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 328 of file G4Scheduler.hh.

329{
330 return fMaxSteps;
331}

◆ GetMaxTimeStep()

G4double G4Scheduler::GetMaxTimeStep ( ) const
inline

Definition at line 182 of file G4Scheduler.hh.

182{ return fMaxTimeStep; }

◆ GetMaxZeroTimeAllowed()

G4int G4Scheduler::GetMaxZeroTimeAllowed ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 383 of file G4Scheduler.hh.

384{
385 return fMaxNZeroTimeStepsAllowed;
386}

◆ GetModelHandler()

G4ITModelHandler * G4Scheduler::GetModelHandler ( )
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 292 of file G4Scheduler.hh.

293{
294 return fpModelHandler;
295}

◆ GetNbSteps()

G4int G4Scheduler::GetNbSteps ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 318 of file G4Scheduler.hh.

319{
320 return fNbSteps;
321}

◆ GetNextWatchedTime()

G4double G4Scheduler::GetNextWatchedTime ( ) const

Definition at line 369 of file G4Scheduler.cc.

370{
371 auto up = fWatchedTimes.upper_bound(fGlobalTime);
372 if (up == fWatchedTimes.end()) {
373 return DBL_MAX;
374 }
375 return *up;
376}
#define DBL_MAX
Definition templates.hh:62

Referenced by SynchronizeTracks().

◆ GetNTracks()

size_t G4Scheduler::GetNTracks ( )
virtual

Definition at line 875 of file G4Scheduler.cc.

876{
877 return fTrackContainer.GetNTracks();
878}

Referenced by Stepping().

◆ GetPreviousTimeStep()

G4double G4Scheduler::GetPreviousTimeStep ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 398 of file G4Scheduler.hh.

399{
400 return fPreviousTimeStep;
401}

◆ GetScavengerMaterial()

◆ GetStartTime()

G4double G4Scheduler::GetStartTime ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 333 of file G4Scheduler.hh.

334{
335 return fStartTime;
336}

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

◆ GetStatus()

G4ITStepStatus G4Scheduler::GetStatus ( ) const
inline

Definition at line 403 of file G4Scheduler.hh.

404{
405 return fITStepStatus;
406}

◆ GetTimeStep()

G4double G4Scheduler::GetTimeStep ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 343 of file G4Scheduler.hh.

344{
345 return fTimeStep;
346}

◆ GetTimeTolerance()

G4double G4Scheduler::GetTimeTolerance ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 393 of file G4Scheduler.hh.

394{
395 return fTimeTolerance;
396}

Referenced by G4ITTrackHolder::_PushTrack().

◆ GetUserTimeStepAction()

G4UserTimeStepAction * G4Scheduler::GetUserTimeStepAction ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 363 of file G4Scheduler.hh.

364{
365 return fpUserTimeStepAction;
366}

◆ GetVerbose()

G4int G4Scheduler::GetVerbose ( ) const
inline

Definition at line 373 of file G4Scheduler.hh.

374{
375 return fVerbose;
376}

◆ Initialize()

void G4Scheduler::Initialize ( )
overridevirtual

Reimplemented from G4VScheduler.

Definition at line 227 of file G4Scheduler.cc.

228{
229 delete fpStepProcessor;
230 delete fpModelProcessor;
231
232 fpModelProcessor = new G4ITModelProcessor();
233 fpModelProcessor->SetModelHandler(fpModelHandler);
234 fpModelProcessor->SetTrackingManager(fpTrackingManager);
235 fpStepProcessor = new G4ITStepProcessor();
236 fpStepProcessor->SetTrackingManager(fpTrackingManager);
237 fpTrackingManager->SetInteractivity(fpTrackingInteractivity);
238 if (fUsePreDefinedTimeSteps) {
239 if (fpUserTimeSteps == nullptr) // Extra checking, is it necessary ?
240 {
241 G4ExceptionDescription exceptionDescription;
242 exceptionDescription << "You are asking to use user defined steps but you did not give any.";
243 G4Exception("G4Scheduler::FindUserPreDefinedTimeStep", "Scheduler004", FatalErrorInArgument,
244 exceptionDescription);
245 return; // makes coverity happy
246 }
247 }
248
249 fInitialized = true;
250}

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

◆ Instance()

◆ IsInitialized()

G4bool G4Scheduler::IsInitialized ( )
inline

Definition at line 287 of file G4Scheduler.hh.

288{
289 return fInitialized;
290}

◆ IsInteractionStep()

G4bool G4Scheduler::IsInteractionStep ( )
inline

Definition at line 189 of file G4Scheduler.hh.

189 {
190 return fInteractionStep;
191 }

◆ IsRunning()

G4bool G4Scheduler::IsRunning ( )
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 107 of file G4Scheduler.hh.

107{ return fRunning; }

◆ Notify()

G4bool G4Scheduler::Notify ( G4ApplicationState requestedState)
overridevirtual

Implements G4VStateDependent.

Definition at line 98 of file G4Scheduler.cc.

99{
100 if (requestedState == G4State_Quit) {
101 if (fVerbose >= 4) {
102 G4cout << "G4Scheduler received G4State_Quit" << G4endl;
103 }
104 Clear();
105 }
106 return true;
107}
@ G4State_Quit

◆ operator=()

G4Scheduler & G4Scheduler::operator= ( const G4Scheduler & )
delete

◆ PrintWhyDoYouStop()

void G4Scheduler::PrintWhyDoYouStop ( )
protected

Definition at line 417 of file G4Scheduler.cc.

418{
419#ifdef G4VERBOSE
420 if (fWhyDoYouStop) {
421 G4cout << "G4Scheduler has reached a stage: it might be"
422 " a transition or the end"
423 << G4endl;
424
425 G4bool normalStop = false;
426
427 if (fGlobalTime >= fStopTime) {
428 G4cout << "== G4Scheduler: I stop because I reached the stop time : "
429 << G4BestUnit(fStopTime, "Time") << " ==" << G4endl;
430 normalStop = true;
431 }
432 if (!fTrackContainer.MainListsNOTEmpty()) // is empty
433 {
434 G4cout << "G4Scheduler: I stop because the current main list of tracks "
435 "is empty"
436 << G4endl;
437 normalStop = true;
438 }
439 if (fMaxSteps == -1 ? false : fNbSteps >= fMaxSteps) {
440 G4cout << "G4Scheduler: I stop because I reached the maximum allowed "
441 "number of steps="
442 << fMaxSteps << G4endl;
443 normalStop = true;
444 }
445 if (fContinue && !normalStop) {
446 G4cout << "G4Scheduler: It might be that I stop because "
447 "I have been told so. You may check "
448 "member fContinue and usage of the method G4Scheduler::Stop()."
449 << G4endl;
450 }
451 }
452#endif
453}
bool G4bool
Definition G4Types.hh:86

Referenced by DoProcess().

◆ Process()

void G4Scheduler::Process ( )
overridevirtual

Reimplemented from G4VScheduler.

Definition at line 273 of file G4Scheduler.cc.

274{
275#ifdef G4VERBOSE
276 if (fVerbose != 0) {
277 G4cout << "*** G4Scheduler starts processing " << G4endl;
278 if (fVerbose > 2)
279 G4cout << "___________________________________________"
280 "___________________________"
281 << G4endl;
282 }
283#endif
284
285 if (!fInitialized) {
286 Initialize();
287 }
288 fpModelProcessor->Initialize();
289 fpStepProcessor->Initialize();
290
291 if (fpGun != nullptr) fpGun->DefineTracks();
292
293 if (fpTrackingInteractivity != nullptr) fpTrackingInteractivity->Initialize();
294
295 // ___________________
296 fRunning = true;
297 Reset();
298
299 if (fResetScavenger) {
300 if (fpUserScavenger != nullptr) {
301 fpUserScavenger->Reset();
302 }
303 }
304
305 if (fpUserTimeStepAction != nullptr) {
306 fpUserTimeStepAction->StartProcessing();
307 }
308
309#ifdef G4VERBOSE
310 G4bool trackFound = false;
311 G4IosFlagsSaver iosfs(G4cout);
312 G4cout.precision(5);
313#endif
314
315 //===========================================================================
316 // By default, before the G4Scheduler is launched, the tracks are pushed to
317 // the delayed lists
318 //===========================================================================
319
320 if (fTrackContainer.DelayListsNOTEmpty()) {
321 fStartTime = fTrackContainer.GetNextTime();
322#ifdef G4VERBOSE
323 trackFound = true;
324 G4Timer localtimer;
325 if (fVerbose > 1) {
326 localtimer.Start();
327 }
328#endif
330#ifdef G4VERBOSE
331 if (fVerbose > 1) {
332 localtimer.Stop();
333 G4cout << "G4Scheduler: process time= " << localtimer << G4endl;
334 }
335#endif
336 }
337#ifdef G4VERBOSE
338 if (fVerbose != 0) {
339 if (trackFound) {
340 G4cout << "*** G4Scheduler ends at time : " << G4BestUnit(fGlobalTime, "Time") << G4endl;
341 G4cout << "___________________________________" << G4endl;
342 }
343 else {
344 G4cout << "*** G4Scheduler did not start because no "
345 "track was found to be processed"
346 << G4endl;
347 G4cout << "___________________________________" << G4endl;
348 }
349 }
350#endif
351
352 fRunning = false;
353
354 if (fpUserTimeStepAction != nullptr) {
355 fpUserTimeStepAction->EndProcessing();
356 }
357 // ___________________
358 EndTracking();
359 ClearList();
360 Reset();
361
362 if (fpTrackingInteractivity != nullptr) {
363 fpTrackingInteractivity->Finalize();
364 }
365}
void SynchronizeTracks()
void EndTracking()
void Reset() override
void Stop()
void Start()

Referenced by G4DNAChemistryManager::Run().

◆ RegisterModel()

void G4Scheduler::RegisterModel ( G4VITStepModel * model,
G4double time )
overridevirtual

Reimplemented from G4VScheduler.

Definition at line 220 of file G4Scheduler.cc.

221{
222 fpModelHandler->RegisterModel(model, time);
223}

◆ Reset()

void G4Scheduler::Reset ( )
overridevirtual

Reimplemented from G4VScheduler.

Definition at line 254 of file G4Scheduler.cc.

255{
256 fStartTime = 0;
257 fUserUpperTimeLimit = -1;
258 fTimeStep = DBL_MAX;
259 fTSTimeStep = DBL_MAX;
260 fILTimeStep = DBL_MAX;
261 fPreviousTimeStep = DBL_MAX;
262 fGlobalTime = -1;
263 fInteractionStep = true;
264 fITStepStatus = eUndefined;
265 fZeroTimeCount = 0;
266
267 fNbSteps = 0;
268 fContinue = true;
269 fReactionSet->CleanAllReaction();
270}
@ eUndefined

Referenced by Process().

◆ ResetScavenger()

void G4Scheduler::ResetScavenger ( bool value)
inline

Definition at line 443 of file G4Scheduler.hh.

444{
445 fResetScavenger = value;
446}

◆ SetDefaultTimeStep()

void G4Scheduler::SetDefaultTimeStep ( G4double timeStep)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 348 of file G4Scheduler.hh.

349{
350 fDefaultMinTimeStep = timeStep;
351}

◆ SetEndTime()

void G4Scheduler::SetEndTime ( const G4double __endtime)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 297 of file G4Scheduler.hh.

298{
299 fEndTime = __endtime;
300}

◆ SetGun()

void G4Scheduler::SetGun ( G4ITGun * gun)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 418 of file G4Scheduler.hh.

419{
420 fpGun = gun;
421}

Referenced by G4DNAChemistryManager::SetGun().

◆ SetInteractionStep()

void G4Scheduler::SetInteractionStep ( G4bool InteractionStep)
inline

Definition at line 192 of file G4Scheduler.hh.

192 {
193 fInteractionStep = InteractionStep;
194 }

Referenced by G4DNAScavengerProcess::PostStepDoIt().

◆ SetInteractivity()

void G4Scheduler::SetInteractivity ( G4ITTrackingInteractivity * interactivity)
overridevirtual

Reimplemented from G4VScheduler.

Definition at line 860 of file G4Scheduler.cc.

861{
862 fpTrackingInteractivity = interactivity;
863 if (fpTrackingManager != nullptr) {
864 fpTrackingManager->SetInteractivity(fpTrackingInteractivity);
865 }
866}

◆ SetMaxNbSteps()

void G4Scheduler::SetMaxNbSteps ( G4int maxSteps)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 323 of file G4Scheduler.hh.

324{
325 fMaxSteps = maxSteps;
326}

◆ SetMaxTimeStep()

void G4Scheduler::SetMaxTimeStep ( G4double maxTimeStep)
inline

Definition at line 180 of file G4Scheduler.hh.

180{ fMaxTimeStep = maxTimeStep; }

◆ SetMaxZeroTimeAllowed()

void G4Scheduler::SetMaxZeroTimeAllowed ( G4int maxTimeStepAllowed)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 378 of file G4Scheduler.hh.

379{
380 fMaxNZeroTimeStepsAllowed = maxTimeStepAllowed;
381}

◆ SetScavengerMaterial()

void G4Scheduler::SetScavengerMaterial ( std::unique_ptr< G4VScavengerMaterial > scavengerMaterial)
inline

Definition at line 185 of file G4Scheduler.hh.

186 {
187 fpUserScavenger = std::move(scavengerMaterial);
188 }

◆ SetTimeSteps()

void G4Scheduler::SetTimeSteps ( std::map< G4double, G4double > * steps)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 302 of file G4Scheduler.hh.

303{
304 fUsePreDefinedTimeSteps = true;
305 fpUserTimeSteps = steps;
306}

◆ SetTimeTolerance()

void G4Scheduler::SetTimeTolerance ( G4double time)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 388 of file G4Scheduler.hh.

389{
390 fTimeTolerance = time;
391}

◆ SetUserAction()

void G4Scheduler::SetUserAction ( G4UserTimeStepAction * userITAction)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 358 of file G4Scheduler.hh.

359{
360 fpUserTimeStepAction = userITAction;
361}

◆ SetVerbose()

void G4Scheduler::SetVerbose ( G4int verbose)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 368 of file G4Scheduler.hh.

369{
370 fVerbose = verbose;
371}

◆ Stepping()

void G4Scheduler::Stepping ( )
protected

Definition at line 497 of file G4Scheduler.cc.

498{
499 fTimeStep = fMaxTimeStep;
500
501 fTSTimeStep = DBL_MAX;
502 fILTimeStep = DBL_MAX;
503
504 //fInteractionStep = false;
505 fReachedUserTimeLimit = false;
506
507 fITStepStatus = eUndefined;
508
509 // Start of step
510#ifdef G4VERBOSE
511 if (fVerbose > 2) {
512# ifdef USE_COLOR
513 G4cout << LIGHT_RED;
514# endif
515 G4cout << "*** Start Of Step N°" << fNbSteps + 1 << " species number : " << GetNTracks()
516 << " ***" << G4endl;
517 G4cout << "Current Global time : " << G4BestUnit(fGlobalTime, "Time") << G4endl;
518# ifdef USE_COLOR
520# endif
521 }
522#endif
523
524#if defined(DEBUG_MEM) && defined(DEBUG_MEM_DETAILED_STEPPING)
525 MemStat mem_first, mem_second, mem_diff;
526#endif
527
528#if defined(DEBUG_MEM) && defined(DEBUG_MEM_DETAILED_STEPPING)
529 mem_first = MemoryUsage();
530#endif
531
532 fDefinedMinTimeStep = GetLimitingTimeStep();
533
534 if (fUsePreDefinedTimeSteps) {
535#ifdef G4VERBOSE
536 if (fVerbose > 2) {
537# ifdef USE_COLOR
538 G4cout << LIGHT_RED;
539# endif
540 G4cout << "*** At time : " << G4BestUnit(fGlobalTime, "Time")
541 << " the chosen user time step is : " << G4BestUnit(fDefinedMinTimeStep, "Time")
542 << " ***" << G4endl;
543# ifdef USE_COLOR
545# endif
546 }
547#endif
548 }
549
550 if (fpModelProcessor->GetComputeTimeStep()) // fComputeTimeStep)
551 {
552 fTSTimeStep = fpModelProcessor->CalculateMinTimeStep(fGlobalTime, fDefinedMinTimeStep);
553 // => at least N (N = nb of tracks) loops
554 }
555 else if (fUseDefaultTimeSteps) {
556 fTSTimeStep = fDefinedMinTimeStep;
557 }
558
559#if defined(DEBUG_MEM) && defined(DEBUG_MEM_DETAILED_STEPPING)
560 mem_second = MemoryUsage();
561 mem_diff = mem_second - mem_first;
562 G4cout << "|| MEM || After computing TS, diff is : " << mem_diff << G4endl;
563#endif
564
565#ifdef G4VERBOSE
566 if (fVerbose > 2) {
567# ifdef USE_COLOR
568 G4cout << LIGHT_RED;
569# endif
570 G4cout << "*** Time stepper returned : " << G4BestUnit(fTSTimeStep, "Time") << " ***" << G4endl;
571# ifdef USE_COLOR
573# endif
574 }
575#endif
576
577#if defined(DEBUG_MEM) && defined(DEBUG_MEM_DETAILED_STEPPING)
578 mem_first = MemoryUsage();
579#endif
580
581 // Call IL even if fTSTimeStep == 0
582 // if fILTimeStep == 0 give the priority to DoIt processes
583
584 fILTimeStep = fpStepProcessor->ComputeInteractionLength(fPreviousTimeStep);
585 // => at least N loops
586 // All process returns the physical step of interaction
587 // The transportation process calculates the corresponding
588 // time step
589
590#if defined(DEBUG_MEM) && defined(DEBUG_MEM_DETAILED_STEPPING)
591 mem_second = MemoryUsage();
592 mem_diff = mem_second - mem_first;
593 G4cout << "|| MEM || After IL, diff is : " << mem_diff << G4endl;
594#endif
595
596#ifdef G4VERBOSE
597 if (fVerbose > 2) {
598# ifdef USE_COLOR
599 G4cout << LIGHT_RED;
600# endif
601 G4cout << "*** The minimum time returned by the processes is : "
602 << G4BestUnit(fILTimeStep, "Time") << " ***" << G4endl;
603# ifdef USE_COLOR
605# endif
606 }
607#endif
608
609#if defined(DEBUG_MEM) && defined(DEBUG_MEM_DETAILED_STEPPING)
610 mem_first = MemoryUsage();
611#endif
612
613 if (fILTimeStep <= fTSTimeStep)
614 // Give the priority to the IL
615 {
616 fInteractionStep = true;
617 //fReactionSet->CleanAllReaction();
618 fTimeStep = fILTimeStep;
619 fITStepStatus = eInteractionWithMedium;
620 fpStepProcessor->PrepareLeadingTracks();
621 }
622 else {
623 fInteractionStep = false;
624 fpStepProcessor->ResetLeadingTracks();
625 fTimeStep = fTSTimeStep;
626 fITStepStatus = eCollisionBetweenTracks;
627 }
628
629 if (fGlobalTime + fTimeStep > fStopTime)
630 // This check is done at every time step
631 {
632 fTimeStep = fStopTime - fGlobalTime;
633 fITStepStatus = eInteractionWithMedium; // ie: transportation
634 fInteractionStep = true;
635 fReactionSet->CleanAllReaction();
636 fpStepProcessor->ResetLeadingTracks();
637 }
638
639 if (fTimeStep == 0) // < fTimeTolerance)
640 {
641 ++fZeroTimeCount;
642 if (fZeroTimeCount >= fMaxNZeroTimeStepsAllowed) {
643 G4ExceptionDescription exceptionDescription;
644
645 exceptionDescription << "Too many zero time steps were detected. ";
646 exceptionDescription << "The simulation is probably stuck. ";
647 exceptionDescription << "The maximum number of zero time steps is currently : "
648 << fMaxNZeroTimeStepsAllowed;
649 exceptionDescription << ".";
650
651 G4Exception("G4Scheduler::Stepping", "SchedulerNullTimeSteps", FatalErrorInArgument,
652 exceptionDescription);
653 }
654 }
655 else {
656 fZeroTimeCount = 0;
657 }
658
659 fReachedUserTimeLimit = (fTimeStep <= fDefinedMinTimeStep)
660 || ((fTimeStep > fDefinedMinTimeStep)
661 && fabs(fTimeStep - fDefinedMinTimeStep) < fTimeTolerance);
662
663 if (fpUserTimeStepAction != nullptr) fpUserTimeStepAction->UserPreTimeStepAction();
664 // TODO: pre/post
665
666#if defined(DEBUG_MEM) && defined(DEBUG_MEM_DETAILED_STEPPING)
667 mem_second = MemoryUsage();
668 mem_diff = mem_second - mem_first;
669 G4cout << "|| MEM || After LeadingTracks and UserPreTimeStepAction: " << mem_diff << G4endl;
670#endif
671
672#if defined(DEBUG_MEM) && defined(DEBUG_MEM_DETAILED_STEPPING)
673 mem_first = MemoryUsage();
674#endif
675
676 fGlobalTime += fTimeStep;
677
678 // if fTSTimeStep > 0 => still need to call the transportation process
679 // if fILTimeStep < fTSTimeStep => call only DoIt processes, no reactions
680 // if fILTimeStep == fTSTimeStep => give the priority to the DoIt processes
681 if (fTSTimeStep > 0 || fILTimeStep <= fTSTimeStep) {
682 fpStepProcessor->DoIt(fTimeStep);
683 }
684#if defined(DEBUG_MEM) && defined(DEBUG_MEM_DETAILED_STEPPING)
685 mem_second = MemoryUsage();
686 mem_diff = mem_second - mem_first;
687 G4cout << "|| MEM || After DoIT, diff is : " << mem_diff << G4endl;
688#endif
689
690#if defined(DEBUG_MEM) && defined(DEBUG_MEM_DETAILED_STEPPING)
691 mem_first = MemoryUsage();
692#endif
693
694 fpModelProcessor->ComputeTrackReaction(fITStepStatus, fGlobalTime, fTimeStep, fPreviousTimeStep,
695 fReachedUserTimeLimit, fTimeTolerance,
696 fpUserTimeStepAction, fVerbose);
697
698 ++fNbSteps;
699
700 if (fpUserTimeStepAction != nullptr) {
701 fpUserTimeStepAction->UserPostTimeStepAction();
702 }
703
704 fPreviousTimeStep = fTimeStep;
705
706#if defined(DEBUG_MEM) && defined(DEBUG_MEM_DETAILED_STEPPING)
707 mem_second = MemoryUsage();
708 mem_diff = mem_second - mem_first;
709 G4cout << "|| MEM || After computing reactions + UserPostTimeStepAction, "
710 "diff is : "
711 << mem_diff << G4endl;
712#endif
713
714 // End of step
715#ifdef G4VERBOSE
716 if (fVerbose >= 2) {
717# ifdef USE_COLOR
718 G4cout << LIGHT_RED;
719# endif
720
721 G4String interactionType;
722 GetCollisionType(interactionType);
723
724 std::stringstream finalOutput;
725
726 finalOutput << "*** End of step N°" << fNbSteps
727 << "\t T_i= " << G4BestUnit(fGlobalTime - fTimeStep, "Time")
728 << "\t dt= " << G4BestUnit(fTimeStep, "Time")
729 << "\t T_f= " << G4BestUnit(fGlobalTime, "Time") << "\t " << interactionType
730 << G4endl;
731
732 if (fVerbose > 2) {
733 if (fReachedUserTimeLimit) {
734 finalOutput << "It has also reached the user time limit" << G4endl;
735 }
736 finalOutput << "_______________________________________________________________"
737 "_______"
738 << G4endl;
739 }
740
741 G4cout << finalOutput.str();
742
743# ifdef USE_COLOR
745# endif
746 }
747#endif
748}
#define LIGHT_RED
#define RESET_COLOR
virtual size_t GetNTracks()
void GetCollisionType(G4String &interactionType)
G4double GetLimitingTimeStep() const override

Referenced by DoProcess().

◆ Stop()

void G4Scheduler::Stop ( )
inline

Definition at line 408 of file G4Scheduler.hh.

409{
410 fContinue = false;
411}

Referenced by G4DNAIRTMoleculeEncounterStepper::CalculateMinTimeStep().

◆ SynchronizeTracks()

void G4Scheduler::SynchronizeTracks ( )
protected

Definition at line 380 of file G4Scheduler.cc.

381{
382 fTmpGlobalTime = fGlobalTime;
383 fGlobalTime = fTrackContainer.GetNextTime();
384 G4double tmpGlobalTime = fGlobalTime;
385 G4double nextWatchedTime = -1;
386 G4bool carryOn = true;
387 while (fTrackContainer.MergeNextTimeToMainList(tmpGlobalTime) && carryOn) {
388 if (tmpGlobalTime != fGlobalTime) {
389 fGlobalTime = tmpGlobalTime;
390 };
391 fStopTime = min(fTrackContainer.GetNextTime(), fEndTime);
392 while ((nextWatchedTime = GetNextWatchedTime()) < fTrackContainer.GetNextTime()
393 && (carryOn = CanICarryOn()))
394 {
395 fStopTime = min(nextWatchedTime, fEndTime);
396 DoProcess();
397 }
398
399 carryOn = CanICarryOn();
400
401 if (nextWatchedTime > fEndTime && carryOn) {
402 fStopTime = min(fTrackContainer.GetNextTime(), fEndTime);
403 DoProcess();
404 }
405 }
406}
double G4double
Definition G4Types.hh:83
G4bool CanICarryOn()
G4double GetNextWatchedTime() const
void DoProcess()
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 433 of file G4Scheduler.hh.

434{
435 fUseDefaultTimeSteps = flag;
436}

◆ WhyDoYouStop()

void G4Scheduler::WhyDoYouStop ( )
inline

Definition at line 428 of file G4Scheduler.hh.

429{
430 fWhyDoYouStop = true;
431}

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