Geant4 11.3.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)
 
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 168 of file G4Scheduler.cc.

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

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

302{
303 if (fpUserTimeSteps == nullptr) {
304 fpUserTimeSteps = new std::map<G4double, G4double>();
305 fUsePreDefinedTimeSteps = true;
306 }
307
308 (*fpUserTimeSteps)[startingTime] = timeStep;
309}

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

432{
433 return (!fUseDefaultTimeSteps && !fUsePreDefinedTimeSteps);
434}

◆ CanICarryOn()

G4bool G4Scheduler::CanICarryOn ( )
protected

Definition at line 408 of file G4Scheduler.cc.

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

Referenced by SynchronizeTracks().

◆ Clear()

void G4Scheduler::Clear ( )

Definition at line 177 of file G4Scheduler.cc.

178{
179 if (fpMessenger != nullptr) {
180 delete fpMessenger;
181 fpMessenger = nullptr;
182 }
183 if (fpStepProcessor != nullptr) {
184 delete fpStepProcessor;
185 fpStepProcessor = nullptr;
186 }
187 if (fpModelProcessor != nullptr) {
188 delete fpModelProcessor;
189 fpModelProcessor = nullptr;
190 }
191
193 ClearList();
194 if (fpTrackingManager != nullptr) {
195 delete fpTrackingManager;
196 fpTrackingManager = nullptr;
197 }
198
199 if (fReactionSet != nullptr) {
200 delete fReactionSet;
201 fReactionSet = nullptr;
202 }
203
204 if (fpModelHandler != nullptr) {
205 delete fpModelHandler;
206 fpModelHandler = nullptr;
207 }
208}
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 212 of file G4Scheduler.cc.

213{
214 fTrackContainer.Clear();
216}
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 455 of file G4Scheduler.cc.

456{
457 if (fpUserTimeStepAction != nullptr) fpUserTimeStepAction->NewStage();
458
459#if defined(DEBUG_MEM) && defined(DEBUG_MEM_STEPPING)
460 MemStat mem_first, mem_second, mem_diff;
461#endif
462
463#if defined(DEBUG_MEM) && defined(DEBUG_MEM_STEPPING)
464 mem_first = MemoryUsage();
465#endif
466
467 while (fGlobalTime < fStopTime && fTrackContainer.MainListsNOTEmpty()
468 && (fMaxSteps == -1 ? true : fNbSteps < fMaxSteps) && fContinue)
469 {
470 Stepping();
471
472#if defined(DEBUG_MEM) && defined(DEBUG_MEM_STEPPING)
473 mem_second = MemoryUsage();
474 mem_diff = mem_second - mem_first;
475 G4cout << "\t || MEM || After step " << fNbSteps << ", diff is : " << mem_diff << G4endl;
476#endif
477 }
478
480
481#if defined(DEBUG_MEM) && defined(DEBUG_MEM_STEPPING)
482 mem_second = MemoryUsage();
483 mem_diff = mem_second - mem_first;
484 G4cout << "\t || MEM || After stepping, diff is : " << mem_diff << G4endl;
485#endif
486
487#ifdef G4VERBOSE
488 if (fVerbose > 2)
489 G4cout << "*** G4Scheduler has finished processing a track list at time : "
490 << G4BestUnit(fGlobalTime, "Time") << G4endl;
491#endif
492}
#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 819 of file G4Scheduler.cc.

820{
821 if (fRunning) {
822 G4ExceptionDescription exceptionDescription;
823 exceptionDescription << "End tracking is called while G4Scheduler is still running." << G4endl;
824
825 G4Exception("G4Scheduler::EndTracking", "Scheduler017", FatalErrorInArgument,
826 exceptionDescription);
827 }
828
829 while (fTrackContainer.DelayListsNOTEmpty()) {
830 auto nextTime = fTrackContainer.GetNextTime();
831 fTrackContainer.MergeNextTimeToMainList(nextTime);
832 }
833
834 fTrackContainer.MergeSecondariesWithMainList();
835
836 if (fTrackContainer.MainListsNOTEmpty()) {
837 G4TrackManyList* mainList = fTrackContainer.GetMainList();
838 G4TrackManyList::iterator it = mainList->begin();
839 G4TrackManyList::iterator end = mainList->end();
840 for (; it != end; ++it) {
841 fpTrackingManager->EndTrackingWOKill(*it);
842 }
843 }
844
845 if (fTrackContainer.SecondaryListsNOTEmpty()) // should be empty
846 {
847 G4TrackManyList* secondaries = fTrackContainer.GetSecondariesList();
848 G4TrackManyList::iterator it = secondaries->begin();
849 G4TrackManyList::iterator end = secondaries->end();
850
851 for (; it != end; ++it) {
852 fpTrackingManager->EndTrackingWOKill(*it);
853 }
854 }
855}
@ 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 785 of file G4Scheduler.cc.

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

◆ ForceReinitialization()

void G4Scheduler::ForceReinitialization ( )

Definition at line 867 of file G4Scheduler.cc.

868{
869 fInitialized = false;
870 Initialize();
871}
void Initialize() override

◆ GetCollisionType()

void G4Scheduler::GetCollisionType ( G4String & interactionType)

Definition at line 879 of file G4Scheduler.cc.

880{
881 switch (fITStepStatus) {
883 interactionType = "eInteractionWithMedium";
884 break;
886 interactionType = "eCollisionBetweenTracks";
887 break;
888 default:
889 interactionType = "eCollisionBetweenTracks";
890 break;
891 }
892}
@ eInteractionWithMedium
@ eCollisionBetweenTracks

Referenced by Stepping().

◆ GetEndTime()

G4double G4Scheduler::GetEndTime ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 331 of file G4Scheduler.hh.

332{
333 return fEndTime;
334}

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

◆ GetGlobalTime()

G4double G4Scheduler::GetGlobalTime ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 346 of file G4Scheduler.hh.

347{
348 return fGlobalTime;
349}

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

◆ GetGun()

G4ITGun * G4Scheduler::GetGun ( )
inline

Definition at line 416 of file G4Scheduler.hh.

417{
418 return fpGun;
419}

◆ GetInteractivity()

G4ITTrackingInteractivity * G4Scheduler::GetInteractivity ( )
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 406 of file G4Scheduler.hh.

407{
408 return fpTrackingInteractivity;
409}

◆ GetLimitingTimeStep()

G4double G4Scheduler::GetLimitingTimeStep ( ) const
overridevirtual

Reimplemented from G4VScheduler.

Definition at line 749 of file G4Scheduler.cc.

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

Referenced by Stepping().

◆ GetMaxNbSteps()

G4int G4Scheduler::GetMaxNbSteps ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 321 of file G4Scheduler.hh.

322{
323 return fMaxSteps;
324}

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

377{
378 return fMaxNZeroTimeStepsAllowed;
379}

◆ GetModelHandler()

G4ITModelHandler * G4Scheduler::GetModelHandler ( )
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 285 of file G4Scheduler.hh.

286{
287 return fpModelHandler;
288}

◆ GetNbSteps()

G4int G4Scheduler::GetNbSteps ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 311 of file G4Scheduler.hh.

312{
313 return fNbSteps;
314}

◆ GetNextWatchedTime()

G4double G4Scheduler::GetNextWatchedTime ( ) const

Definition at line 368 of file G4Scheduler.cc.

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

Referenced by SynchronizeTracks().

◆ GetNTracks()

size_t G4Scheduler::GetNTracks ( )
virtual

Definition at line 873 of file G4Scheduler.cc.

874{
875 return fTrackContainer.GetNTracks();
876}

Referenced by Stepping().

◆ GetPreviousTimeStep()

G4double G4Scheduler::GetPreviousTimeStep ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 391 of file G4Scheduler.hh.

392{
393 return fPreviousTimeStep;
394}

◆ GetScavengerMaterial()

G4VScavengerMaterial * G4Scheduler::GetScavengerMaterial ( ) const
inline

◆ GetStartTime()

G4double G4Scheduler::GetStartTime ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 326 of file G4Scheduler.hh.

327{
328 return fStartTime;
329}

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

◆ GetStatus()

G4ITStepStatus G4Scheduler::GetStatus ( ) const
inline

Definition at line 396 of file G4Scheduler.hh.

397{
398 return fITStepStatus;
399}

◆ GetTimeStep()

G4double G4Scheduler::GetTimeStep ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 336 of file G4Scheduler.hh.

337{
338 return fTimeStep;
339}

◆ GetTimeTolerance()

G4double G4Scheduler::GetTimeTolerance ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 386 of file G4Scheduler.hh.

387{
388 return fTimeTolerance;
389}

Referenced by G4ITTrackHolder::_PushTrack().

◆ GetUserTimeStepAction()

G4UserTimeStepAction * G4Scheduler::GetUserTimeStepAction ( ) const
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 356 of file G4Scheduler.hh.

357{
358 return fpUserTimeStepAction;
359}

◆ GetVerbose()

G4int G4Scheduler::GetVerbose ( ) const
inline

Definition at line 366 of file G4Scheduler.hh.

367{
368 return fVerbose;
369}

◆ Initialize()

void G4Scheduler::Initialize ( )
overridevirtual

Reimplemented from G4VScheduler.

Definition at line 226 of file G4Scheduler.cc.

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

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

◆ Instance()

◆ IsInitialized()

G4bool G4Scheduler::IsInitialized ( )
inline

Definition at line 280 of file G4Scheduler.hh.

281{
282 return fInitialized;
283}

◆ 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 415 of file G4Scheduler.cc.

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

Referenced by DoProcess().

◆ Process()

void G4Scheduler::Process ( )
overridevirtual

Reimplemented from G4VScheduler.

Definition at line 272 of file G4Scheduler.cc.

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

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

◆ Reset()

void G4Scheduler::Reset ( )
overridevirtual

Reimplemented from G4VScheduler.

Definition at line 253 of file G4Scheduler.cc.

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

Referenced by Process().

◆ ResetScavenger()

void G4Scheduler::ResetScavenger ( bool value)
inline

Definition at line 436 of file G4Scheduler.hh.

437{
438 fResetScavenger = value;
439}

◆ SetDefaultTimeStep()

void G4Scheduler::SetDefaultTimeStep ( G4double timeStep)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 341 of file G4Scheduler.hh.

342{
343 fDefaultMinTimeStep = timeStep;
344}

◆ SetEndTime()

void G4Scheduler::SetEndTime ( const G4double __endtime)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 290 of file G4Scheduler.hh.

291{
292 fEndTime = __endtime;
293}

◆ SetGun()

void G4Scheduler::SetGun ( G4ITGun * gun)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 411 of file G4Scheduler.hh.

412{
413 fpGun = gun;
414}

Referenced by G4DNAChemistryManager::SetGun().

◆ SetInteractivity()

void G4Scheduler::SetInteractivity ( G4ITTrackingInteractivity * interactivity)
overridevirtual

Reimplemented from G4VScheduler.

Definition at line 858 of file G4Scheduler.cc.

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

◆ SetMaxNbSteps()

void G4Scheduler::SetMaxNbSteps ( G4int maxSteps)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 316 of file G4Scheduler.hh.

317{
318 fMaxSteps = maxSteps;
319}

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

372{
373 fMaxNZeroTimeStepsAllowed = maxTimeStepAllowed;
374}

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

296{
297 fUsePreDefinedTimeSteps = true;
298 fpUserTimeSteps = steps;
299}

◆ SetTimeTolerance()

void G4Scheduler::SetTimeTolerance ( G4double time)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 381 of file G4Scheduler.hh.

382{
383 fTimeTolerance = time;
384}

◆ SetUserAction()

void G4Scheduler::SetUserAction ( G4UserTimeStepAction * userITAction)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 351 of file G4Scheduler.hh.

352{
353 fpUserTimeStepAction = userITAction;
354}

◆ SetVerbose()

void G4Scheduler::SetVerbose ( G4int verbose)
inlineoverridevirtual

Reimplemented from G4VScheduler.

Definition at line 361 of file G4Scheduler.hh.

362{
363 fVerbose = verbose;
364}

◆ Stepping()

void G4Scheduler::Stepping ( )
protected

Definition at line 495 of file G4Scheduler.cc.

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

402{
403 fContinue = false;
404}

Referenced by G4DNAIRTMoleculeEncounterStepper::CalculateMinTimeStep().

◆ SynchronizeTracks()

void G4Scheduler::SynchronizeTracks ( )
protected

Definition at line 379 of file G4Scheduler.cc.

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

427{
428 fUseDefaultTimeSteps = flag;
429}

◆ WhyDoYouStop()

void G4Scheduler::WhyDoYouStop ( )
inline

Definition at line 421 of file G4Scheduler.hh.

422{
423 fWhyDoYouStop = true;
424}

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