Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Scheduler.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr)
28//
29// History:
30// -----------
31// 10 Oct 2011 M.Karamitros created
32//
33// -------------------------------------------------------------------
34
35#include <G4AllITFinder.hh>
36#include <G4Scheduler.hh>
38
39#include "G4SystemOfUnits.hh"
40#include "G4ITModelProcessor.hh"
41#include "G4ITStepProcessor.hh"
42#include "G4IT.hh"
43#include "G4ITReactionChange.hh"
44#include "G4ITModelHandler.hh"
45#include "G4VITStepModel.hh"
50#include "G4UnitsTable.hh"
51#include "G4ITStepStatus.hh"
52#include "G4ITGun.hh"
53#include "G4StateManager.hh"
54#include "G4Timer.hh"
55#include "G4IosFlagsSaver.hh"
56#include <sstream>
57
58//#include "G4Phenomenon.hh"
59//#include "G4MIWorkspace.hh"
60
61//#define DEBUG_MEM 1
62#define DEBUG_MEM_STEPPING
63#define DEBUG_MEM_DETAILED_STEPPING
64//#define DEBUG_MEM_DOIT
65
66#ifdef DEBUG_MEM
67#include "G4MemStat.hh"
68using namespace G4MemStat;
70#endif
71
72//COLOR FOR DEBUGING
73//#define USE_COLOR 1
74
75#ifdef USE_COLOR
76#define RED "\033[0;31m"
77#define LIGHT_RED "\33[1;31m"
78#define GREEN "\033[32;40m"
79#define GREEN_ON_BLUE "\033[1;32;44m"
80#define RESET_COLOR "\033[0m"
81#else
82#define RED ""
83#define LIGHT_RED ""
84#define GREEN ""
85#define GREEN_ON_BLUE ""
86#define RESET_COLOR ""
87#endif
88
89using namespace std;
90
91G4ThreadLocal G4Scheduler* G4Scheduler::fgScheduler(nullptr);
92
93template<typename T>
94 inline G4bool IsInf(T value)
95 {
96 return std::numeric_limits<T>::has_infinity
97 && value == std::numeric_limits<T>::infinity();
98 }
99//_________________________________________________________________________
100
102{
103 if(fgScheduler == nullptr) fgScheduler = new G4Scheduler();
104 return fgScheduler;
105}
106//_________________________________________________________________________
107
109{
110 if(requestedState == G4State_Quit)
111 {
112 if(fVerbose >= 4)
113 {
114 G4cout << "G4Scheduler received G4State_Quit" << G4endl;
115 }
116 Clear();
117 //DeleteInstance();
118 }
119 return true;
120}
121//_________________________________________________________________________
122
124{
125
126
127 delete fgScheduler;
128
129}
130//_________________________________________________________________________
131
132G4Scheduler::G4Scheduler() :
133
134 fTrackContainer((G4ITTrackHolder&) *G4ITTrackHolder::Instance())
135{
136 Create();
137}
138
139void G4Scheduler::Create()
140{
141 fUseDefaultTimeSteps = true;
142 fUserUpperTimeLimit = -1;
143 fpGun = nullptr;
144 fContinue = true;
145// fpMainList = 0;
146// fpWaitingList = 0;
147 fpTrackingInteractivity = nullptr;
148
149 fITStepStatus = eUndefined;
150
151 fpUserTimeSteps = nullptr;
152
153 fTimeStep = DBL_MAX;
154 fTSTimeStep = DBL_MAX;
155 fILTimeStep = DBL_MAX;
156 fPreviousTimeStep = DBL_MAX;
157
158 fZeroTimeCount = 0;
159 fMaxNZeroTimeStepsAllowed = 10000;
160
161 fStartTime = 0;
162 fTimeTolerance = 1 * picosecond;
163 fEndTime = 1 * microsecond;
164 fGlobalTime = -1;
165 fInteractionStep = true;
166 fUsePreDefinedTimeSteps = false;
167
168 fDefaultMinTimeStep = 1 * picosecond;
169
170 fpStepProcessor = nullptr;
171 fpModelProcessor = nullptr;
172
173 fNbSteps = 0;
174 fMaxSteps = -1;
175
176 fRunning = false;
177 fInitialized = false;
178
179 fpUserTimeStepAction = nullptr;
180 fpModelHandler = new G4ITModelHandler();
181 fpTrackingManager = new G4ITTrackingManager();
182
183 fVerbose = 0;
184 fWhyDoYouStop = false;
185 fDefinedMinTimeStep = -1.;
186 fReachedUserTimeLimit = false;
187 fStopTime = -1.;
188 fTmpGlobalTime = -1.;
189
190 fpMessenger = new G4SchedulerMessenger(this);
191
192 fReactionSet = G4ITReactionSet::Instance();
193 fMaxTimeStep = DBL_MAX;
194
195 //hoang add
196 fResetScavenger = true;//true by default
197
199}
200
201//_________________________________________________________________________
202
204{
205
206 if(fpMessenger != nullptr) // is used as a flag to know whether the manager was cleared
207 {
208 Clear();
209 }
210
211 fgScheduler = nullptr;
212
213// if (fVerbose >= 1)
214// {
215// G4cout << "G4Scheduler is being deleted. Bye :) !" << G4endl;
216// }
217}
218
220{
221// if (fVerbose) G4cout << "*** G4Scheduler is being cleared ***" << G4endl;
222
223 if(fpMessenger != nullptr)
224 {
225 delete fpMessenger;
226 fpMessenger = nullptr;
227 }
228 if(fpStepProcessor != nullptr)
229 {
230 delete fpStepProcessor;
231 fpStepProcessor = nullptr;
232 }
233 if(fpModelProcessor != nullptr)
234 {
235 delete fpModelProcessor;
236 fpModelProcessor = nullptr;
237 }
238
240 ClearList();
241 if(fpTrackingManager != nullptr)
242 {
243 delete fpTrackingManager;
244 fpTrackingManager = nullptr;
245 }
246
247 if(fReactionSet != nullptr)
248 {
249 delete fReactionSet;
250 fReactionSet = nullptr;
251 }
252
253 if(fpModelHandler != nullptr)
254 {
255 delete fpModelHandler;
256 fpModelHandler = nullptr;
257 }
258
259 //* DEBUG
260 //* assert(G4StateManager::GetStateManager()->
261 //* DeregisterDependent(this) == true);
262
263}
264
265//_________________________________________________________________________
266
268{
269// if (fNbTracks == 0) return;
270
271 fTrackContainer.Clear();
272
274}
275
276//_________________________________________________________________________
278{
279 fpModelHandler->RegisterModel(model, time);
280}
281
282//_________________________________________________________________________
283
285{
286
287
288 delete fpStepProcessor;
289
290
291
292 delete fpModelProcessor;
293
294 // if(fpMasterModelProcessor)
295 // {
296 // delete fpMasterModelProcessor;
297 // }
298
299 //______________________________________________________________
300
301 fpModelProcessor = new G4ITModelProcessor();
302 fpModelProcessor->SetModelHandler(fpModelHandler);
303 fpModelProcessor->SetTrackingManager(fpTrackingManager);
304
305 // fpMasterModelProcessor = new G4ITModelProcessor();
306 // fpMasterModelProcessor->SetModelHandler(fpModelHandler);
307 // fpModelProcessor = fpMasterModelProcessor;
308
309 //______________________________________________________________
310
311 fpStepProcessor = new G4ITStepProcessor();
312 fpStepProcessor->SetTrackingManager(fpTrackingManager);
313
314 fpTrackingManager->SetInteractivity(fpTrackingInteractivity);
315
316 // fpMasterStepProcessor = new G4ITStepProcessor();
317 // fpMasterStepProcessor->SetTrackingManager(fpTrackingManager);
318 // fpStepProcessor = fpMasterStepProcessor ;
319 //______________________________________________________________
320
321 if(fUsePreDefinedTimeSteps)
322 {
323 if(fpUserTimeSteps == nullptr) // Extra checking, is it necessary ?
324 {
325 G4ExceptionDescription exceptionDescription;
326 exceptionDescription
327 << "You are asking to use user defined steps but you did not give any.";
328 G4Exception("G4Scheduler::FindUserPreDefinedTimeStep",
329 "Scheduler004",
331 exceptionDescription);
332 return; // makes coverity happy
333 }
334 }
335
336// if (fComputeTimeStep)
337// {
338// if (fpModelProcessor == 0)
339// {
340// G4ExceptionDescription exceptionDescription;
341// exceptionDescription
342// << "There is no G4ITModelProcessor to handle IT reaction. ";
343// exceptionDescription
344// << "You probably did not initialize the G4Scheduler. ";
345// exceptionDescription
346// << "Just do G4Scheduler::Instance()->Initialize(); ";
347// exceptionDescription << " but only after initializing the run manager.";
348// G4Exception("G4Scheduler::CalculateMinStep", "ITScheduler005",
349// FatalErrorInArgument, exceptionDescription);
350// return; // makes coverity happy
351// }
352// }
353
354 fInitialized = true;
355}
356
357//_________________________________________________________________________
358
360{
361 fStartTime = 0;
362 fUserUpperTimeLimit = -1;
363 fTimeStep = DBL_MAX;
364 fTSTimeStep = DBL_MAX;
365 fILTimeStep = DBL_MAX;
366 fPreviousTimeStep = DBL_MAX;
367 fGlobalTime = -1;
368 fInteractionStep = true;
369 fITStepStatus = eUndefined;
370 fZeroTimeCount = 0;
371
372 fNbSteps = 0;
373 fContinue = true;
374 // fReactingTracks.clear();
375 fReactionSet->CleanAllReaction();
376}
377//_________________________________________________________________________
378
380{
381
382#ifdef G4VERBOSE
383 if(fVerbose != 0)
384 {
385 G4cout << "*** G4Scheduler starts processing " << G4endl;
386 if(fVerbose > 2)
387 G4cout << "___________________________________________"
388 "___________________________" << G4endl;
389 }
390#endif
391
392 if (!fInitialized) Initialize();
393
394 // fpTrackingManager->Initialize();
395 fpModelProcessor->Initialize();
396 fpStepProcessor->Initialize();
397
398 // TODO
399 // fpMasterModelProcessor->Initialize();
400 // fpMasterStepProcessor->Initialize();
401
402 if (fpGun != nullptr) fpGun->DefineTracks();
403
404 if (fpTrackingInteractivity != nullptr) fpTrackingInteractivity->Initialize();
405
406 // ___________________
407 fRunning = true;
408 Reset();
409
410 //hoang 7/12/2020
411
412 if(fResetScavenger) {
413 if (fpUserScavenger != nullptr) {
414 fpUserScavenger->Reset();
415 }
416 }
417
418 if (fpUserTimeStepAction != nullptr)
419 {
420 fpUserTimeStepAction->StartProcessing();
421 }
422
423#ifdef G4VERBOSE
424 G4bool trackFound = false;
425 G4IosFlagsSaver iosfs(G4cout);
426 G4cout.precision(5);
427#endif
428
429 //===========================================================================
430 // By default, before the G4Scheduler is launched, the tracks are pushed to
431 // the delayed lists
432 //===========================================================================
433
434 if(fTrackContainer.DelayListsNOTEmpty())
435 {
436 fStartTime = fTrackContainer.GetNextTime();
437#ifdef G4VERBOSE
438 trackFound = true;
439 G4Timer localtimer;
440 if(fVerbose>1)
441 {
442 localtimer.Start();
443 }
444#endif
446#ifdef G4VERBOSE
447 if(fVerbose>1)
448 {
449 localtimer.Stop();
450 G4cout << "G4Scheduler: process time= "<< localtimer << G4endl;
451 }
452#endif
453 }
454
455// //---------------------------------
456// // TODO: This could be removed ?
457// if(fTrackContainer.MainListsNOTEmpty()
458// && (fMaxSteps == -1 ? true : fNbSteps < fMaxSteps)
459// && fGlobalTime < fEndTime
460// && fContinue == true)
461//{
462//#ifdef G4VERBOSE
463// trackFound = true;
464//#endif
465// DoProcess();
466//}
467// //---------------------------------
468
469#ifdef G4VERBOSE
470 if(fVerbose != 0)
471 {
472 if(trackFound)
473 {
474 G4cout << "*** G4Scheduler ends at time : "
475 << G4BestUnit(fGlobalTime,"Time") << G4endl;
476 G4cout << "___________________________________" << G4endl;
477 }
478 else
479 {
480 G4cout << "*** G4Scheduler did not start because no "
481 "track was found to be processed"<< G4endl;
482 G4cout << "___________________________________" << G4endl;
483 }
484 }
485#endif
486
487 fRunning = false;
488
489 if (fpUserTimeStepAction != nullptr) fpUserTimeStepAction->EndProcessing();
490
491 // ___________________
492 EndTracking();
493 ClearList();
494
495 Reset();
496
497 if (fpTrackingInteractivity != nullptr) fpTrackingInteractivity->Finalize();
498}
499
500//_________________________________________________________________________
501
503{
504 auto up = fWatchedTimes.upper_bound(fGlobalTime);
505 if(up == fWatchedTimes.end()) return DBL_MAX;
506 return *up;
507}
508
509//_________________________________________________________________________
510
512{
513// if(fTrackContainer.WaitingListsNOTEmpty())
514// {
515// G4ExceptionDescription exceptionDescription;
516// exceptionDescription
517// << "There is a waiting track list (fpWaitingList != 0).";
518// exceptionDescription
519// << " When G4Scheduler::SynchronizeTracks() is called, ";
520// exceptionDescription
521// << "no more tracks should remain in the fpWaitingList.";
522// G4Exception("G4Scheduler::SynchronizeTracks", "ITScheduler002",
523// FatalErrorInArgument, exceptionDescription);
524// }
525
526 // Backup main list and time feature
527 // Reminder : the main list here, should
528 // have the biggest global time
529 // fTrackContainer.MoveMainToWaitingList();
530 // TODO: not yet supported
531
532 fTmpGlobalTime = fGlobalTime;
533 //fTmpEndTime = fEndTime;
534
535 fGlobalTime = fTrackContainer.GetNextTime();
536 G4double tmpGlobalTime = fGlobalTime;
537
538 G4double nextWatchedTime = -1;
539 G4bool carryOn = true;
540
541 while(fTrackContainer.MergeNextTimeToMainList(tmpGlobalTime) && carryOn)
542 {
543 if(tmpGlobalTime != fGlobalTime)
544 {
545 fGlobalTime = tmpGlobalTime;
546 };
547 fStopTime = min(fTrackContainer.GetNextTime(), fEndTime);
548 while((nextWatchedTime = GetNextWatchedTime()) < fTrackContainer.GetNextTime()
549 && (carryOn = CanICarryOn()))
550 {
551 fStopTime = min(nextWatchedTime, fEndTime);
552 DoProcess();
553 }
554
555 carryOn = CanICarryOn();
556
557 if(nextWatchedTime > fEndTime && carryOn)
558 {
559 fStopTime = min(fTrackContainer.GetNextTime(), fEndTime);
560 DoProcess();
561 }
562 }
563}
564
565//_________________________________________________________________________
566
568{
569 return fGlobalTime < fEndTime && (fMaxSteps == -1 ? true : fNbSteps < fMaxSteps)
570 && fContinue;
571}
572
573//_________________________________________________________________________
574
576{
577#ifdef G4VERBOSE
578 if(fWhyDoYouStop)
579 {
580 G4cout << "G4Scheduler has reached a stage: it might be"
581 " a transition or the end"
582 << G4endl;
583
584 G4bool normalStop = false;
585
586 if(fGlobalTime >= fStopTime)
587 {
588 G4cout << "== G4Scheduler: I stop because I reached the stop time : "
589 << G4BestUnit(fStopTime,"Time") << " =="<< G4endl;
590 normalStop = true;
591 }
592 if(!fTrackContainer.MainListsNOTEmpty()) // is empty
593 {
594 G4cout << "G4Scheduler: I stop because the current main list of tracks "
595 "is empty"
596 << G4endl;
597 normalStop = true;
598 }
599 if(fMaxSteps == -1 ? false : fNbSteps >= fMaxSteps)
600 {
601 G4cout << "G4Scheduler: I stop because I reached the maximum allowed "
602 "number of steps=" << fMaxSteps
603 << G4endl;
604 normalStop = true;
605 }
606 if(fContinue && !normalStop)
607 {
608 G4cout << "G4Scheduler: It might be that I stop because "
609 "I have been told so. You may check "
610 "member fContinue and usage of the method G4Scheduler::Stop()."
611 << G4endl;
612 }
613 }
614#endif
615}
616
617//_________________________________________________________________________
618
620// We split it from the Process() method to avoid repeating code in SynchronizeTracks
621{
622 if(fpUserTimeStepAction != nullptr) fpUserTimeStepAction->NewStage();
623
624#if defined (DEBUG_MEM) && defined (DEBUG_MEM_STEPPING)
625 MemStat mem_first, mem_second, mem_diff;
626#endif
627
628#if defined (DEBUG_MEM) && defined (DEBUG_MEM_STEPPING)
629 mem_first = MemoryUsage();
630#endif
631
632 while (fGlobalTime < fStopTime
633 && fTrackContainer.MainListsNOTEmpty()
634 && (fMaxSteps == -1 ? true : fNbSteps < fMaxSteps)
635 && fContinue)
636 {
637// G4cout << "Mainlist size : " << fTrackContainer.GetMainList()->size()
638// << G4endl;
639
640 Stepping();
641
642#if defined (DEBUG_MEM) && defined (DEBUG_MEM_STEPPING)
643 mem_second = MemoryUsage();
644 mem_diff = mem_second-mem_first;
645 G4cout << "\t || MEM || After step " << fNbSteps << ", diff is : "
646 << mem_diff << G4endl;
647#endif
648 }
649
651
652#if defined (DEBUG_MEM) && defined (DEBUG_MEM_STEPPING)
653 mem_second = MemoryUsage();
654 mem_diff = mem_second-mem_first;
655 G4cout << "\t || MEM || After stepping, diff is : " << mem_diff << G4endl;
656#endif
657
658#ifdef G4VERBOSE
659 if(fVerbose > 2)
660 G4cout << "*** G4Scheduler has finished processing a track list at time : "
661 << G4BestUnit(fGlobalTime, "Time") << G4endl;
662#endif
663}
664//_________________________________________________________________________
665
667{
668 fTimeStep = fMaxTimeStep;
669
670 fTSTimeStep = DBL_MAX;
671 fILTimeStep = DBL_MAX;
672
673 fInteractionStep = false;
674 fReachedUserTimeLimit = false;
675
676 fITStepStatus = eUndefined;
677
678 // Start of step
679#ifdef G4VERBOSE
680 if (fVerbose > 2)
681 {
682#ifdef USE_COLOR
683 G4cout << LIGHT_RED;
684#endif
685 G4cout << "*** Start Of Step N°" << fNbSteps + 1 << " ***" << G4endl;
686 G4cout << "Current Global time : " << G4BestUnit(fGlobalTime, "Time")
687 <<G4endl;
688#ifdef USE_COLOR
690#endif
691 }
692#endif
693
694#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
695 MemStat mem_first, mem_second, mem_diff;
696#endif
697
698#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
699 mem_first = MemoryUsage();
700#endif
701
702 fDefinedMinTimeStep = GetLimitingTimeStep();
703
704 if (fUsePreDefinedTimeSteps)
705 {
706#ifdef G4VERBOSE
707 if (fVerbose > 2)
708 {
709#ifdef USE_COLOR
710 G4cout << LIGHT_RED;
711#endif
712 G4cout << "*** At time : " << G4BestUnit(fGlobalTime, "Time")
713 << " the chosen user time step is : "
714 << G4BestUnit(fDefinedMinTimeStep, "Time") << " ***" << G4endl;
715#ifdef USE_COLOR
717#endif
718 }
719#endif
720 }
721
722 if (fpModelProcessor->GetComputeTimeStep()) // fComputeTimeStep)
723 {
724 fTSTimeStep = fpModelProcessor->CalculateMinTimeStep(fGlobalTime,
725 fDefinedMinTimeStep);
726 // => at least N (N = nb of tracks) loops
727 }
728 else if(fUseDefaultTimeSteps)
729 {
730 fTSTimeStep = fDefinedMinTimeStep;
731 }
732
733#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
734 mem_second = MemoryUsage();
735 mem_diff = mem_second-mem_first;
736 G4cout << "|| MEM || After computing TS, diff is : " << mem_diff << G4endl;
737#endif
738
739#ifdef G4VERBOSE
740 if (fVerbose > 2)
741 {
742#ifdef USE_COLOR
743 G4cout << LIGHT_RED;
744#endif
745 G4cout << "*** Time stepper returned : " << G4BestUnit(fTSTimeStep, "Time")
746 << " ***" << G4endl;
747#ifdef USE_COLOR
749#endif
750 }
751#endif
752
753#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
754 mem_first = MemoryUsage();
755#endif
756
757 // Call IL even if fTSTimeStep == 0
758 // if fILTimeStep == 0 give the priority to DoIt processes
759
760 fILTimeStep = fpStepProcessor->ComputeInteractionLength(fPreviousTimeStep);
761 // => at least N loops
762 // All process returns the physical step of interaction
763 // The transportation process calculates the corresponding
764 // time step
765
766#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
767 mem_second = MemoryUsage();
768 mem_diff = mem_second-mem_first;
769 G4cout << "|| MEM || After IL, diff is : " << mem_diff << G4endl;
770#endif
771
772#ifdef G4VERBOSE
773 if (fVerbose > 2)
774 {
775#ifdef USE_COLOR
776 G4cout << LIGHT_RED;
777#endif
778 G4cout << "*** The minimum time returned by the processes is : "
779 << G4BestUnit(fILTimeStep, "Time") << " ***" << G4endl;
780#ifdef USE_COLOR
782#endif
783 }
784#endif
785
786#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
787 mem_first = MemoryUsage();
788#endif
789
790 if (fILTimeStep <= fTSTimeStep)
791 // Give the priority to the IL
792 {
793 fInteractionStep = true;
794 fReactionSet->CleanAllReaction();
795 fTimeStep = fILTimeStep;
796 fITStepStatus = eInteractionWithMedium;
797 fpStepProcessor->PrepareLeadingTracks();
798 }
799 else
800 {
801 fInteractionStep = false;
802 fpStepProcessor->ResetLeadingTracks();
803 fTimeStep = fTSTimeStep;
804 fITStepStatus = eCollisionBetweenTracks;
805 }
806
807 if (fGlobalTime + fTimeStep > fStopTime)
808 // This check is done at every time step
809 {
810 fTimeStep = fStopTime - fGlobalTime;
811 fITStepStatus = eInteractionWithMedium; // ie: transportation
812 fInteractionStep = true;
813 fReactionSet->CleanAllReaction();
814 fpStepProcessor->ResetLeadingTracks();
815 }
816
817 if (fTimeStep == 0) // < fTimeTolerance)
818 {
819 ++fZeroTimeCount;
820 if (fZeroTimeCount >= fMaxNZeroTimeStepsAllowed)
821 {
822 G4ExceptionDescription exceptionDescription;
823
824 exceptionDescription << "Too many zero time steps were detected. ";
825 exceptionDescription << "The simulation is probably stuck. ";
826 exceptionDescription
827 << "The maximum number of zero time steps is currently : "
828 << fMaxNZeroTimeStepsAllowed;
829 exceptionDescription << ".";
830
831 G4Exception("G4Scheduler::Stepping",
832 "SchedulerNullTimeSteps",
834 exceptionDescription);
835 }
836 }
837 else
838 {
839 fZeroTimeCount = 0;
840 }
841
842 fReachedUserTimeLimit =
843 (fTimeStep <= fDefinedMinTimeStep) || ((fTimeStep > fDefinedMinTimeStep)
844 && fabs(fTimeStep - fDefinedMinTimeStep) < fTimeTolerance);
845
846 if (fpUserTimeStepAction != nullptr) fpUserTimeStepAction->UserPreTimeStepAction();
847 // TODO: pre/post
848
849#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
850 mem_second = MemoryUsage();
851 mem_diff = mem_second-mem_first;
852 G4cout << "|| MEM || After LeadingTracks and UserPreTimeStepAction: "
853 << mem_diff << G4endl;
854#endif
855
856#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
857 mem_first = MemoryUsage();
858#endif
859
860
861 fGlobalTime += fTimeStep;
862
863 // if fTSTimeStep > 0 => still need to call the transportation process
864 // if fILTimeStep < fTSTimeStep => call only DoIt processes, no reactions
865 // if fILTimeStep == fTSTimeStep => give the priority to the DoIt processes
866 if (fTSTimeStep > 0 || fILTimeStep <= fTSTimeStep)
867 {
868 // G4cout << "Will call DoIT" << G4endl;
869 fpStepProcessor->DoIt(fTimeStep);
870 // fTrackContainer.MergeSecondariesWithMainList();
871 // fTrackContainer.KillTracks(); // remove ?
872 }
873 // else
874 // {
875 // G4cout << "fTSTimeStep : " << fTSTimeStep << G4endl;
876 // G4cout << "fILTimeStep : " << fILTimeStep << G4endl;
877 // }
878
879#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
880 mem_second = MemoryUsage();
881 mem_diff = mem_second-mem_first;
882 G4cout << "|| MEM || After DoIT, diff is : " << mem_diff << G4endl;
883#endif
884
885#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
886 mem_first = MemoryUsage();
887#endif
888
889 fpModelProcessor->ComputeTrackReaction(fITStepStatus,
890 fGlobalTime,
891 fTimeStep,
892 fPreviousTimeStep,
893 fReachedUserTimeLimit,
894 fTimeTolerance,
895 fpUserTimeStepAction,
896 fVerbose);
897
898 ++fNbSteps;
899
900 if (fpUserTimeStepAction != nullptr)
901 {
902 fpUserTimeStepAction->UserPostTimeStepAction();
903 }
904
905 fPreviousTimeStep = fTimeStep;
906
907#if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
908 mem_second = MemoryUsage();
909 mem_diff = mem_second-mem_first;
910 G4cout << "|| MEM || After computing reactions + UserPostTimeStepAction, "
911 "diff is : " << mem_diff << G4endl;
912#endif
913
914 // End of step
915#ifdef G4VERBOSE
916 if (fVerbose >= 2)
917 {
918#ifdef USE_COLOR
919 G4cout << LIGHT_RED;
920#endif
921
922 G4String interactionType;
923 GetCollisionType(interactionType);
924
925 std::stringstream finalOutput;
926
927 finalOutput << "*** End of step N°" << fNbSteps
928 << "\t T_i= " << G4BestUnit(fGlobalTime-fTimeStep, "Time")
929 << "\t dt= " << G4BestUnit(fTimeStep, "Time")
930 << "\t T_f= " << G4BestUnit(fGlobalTime, "Time")
931 << "\t " << interactionType
932 << G4endl;
933
934 if(fVerbose>2)
935 {
936 if(fReachedUserTimeLimit)
937 {
938 finalOutput << "It has also reached the user time limit" << G4endl;
939 }
940 finalOutput << "_______________________________________________________________"
941 "_______"<< G4endl;
942 }
943
944 G4cout << finalOutput.str();
945
946#ifdef USE_COLOR
948#endif
949 }
950#endif
951
952}
953//_________________________________________________________________________
954
956{
957 if (fpUserTimeSteps == nullptr) return fDefaultMinTimeStep;
958 if (fabs(fGlobalTime - fUserUpperTimeLimit) < fTimeTolerance)
959 return fDefinedMinTimeStep;
960
961 auto it_fpUserTimeSteps_i = fpUserTimeSteps
962 ->upper_bound(fGlobalTime);
963 auto it_fpUserTimeSteps_low = fpUserTimeSteps
964 ->lower_bound(fGlobalTime);
965
966 // DEBUG
967 // G4cout << "fGlobalTime : " << G4BestUnit(fGlobalTime,"Time")
968 // << G4endl;
969 // G4cout << "fpUserTimeSteps_i : "
970 // <<"<"<<G4BestUnit(it_fpUserTimeSteps->first,"Time")
971 // <<", "<< G4BestUnit(it_fpUserTimeSteps->second,"Time")<<">"
972 // << "\t fpUserTimeSteps_low : "
973 // <<"<"<<G4BestUnit(fpUserTimeSteps_low->first,"Time")<<", "*
974 // << G4BestUnit(fpUserTimeSteps_low->second,"Time")<<">"
975 // << G4endl;
976
977 if (it_fpUserTimeSteps_i == fpUserTimeSteps->end())
978 {
979 it_fpUserTimeSteps_i--;
980 fUserUpperTimeLimit = fStopTime;
981 }
982 else if (fabs(fGlobalTime - it_fpUserTimeSteps_low->first) < fTimeTolerance)
983 {
984 // Case : fGlobalTime = X picosecond
985 // and fpUserTimeSteps_low->first = X picosecond
986 // but the precision is not good enough
987 it_fpUserTimeSteps_i = it_fpUserTimeSteps_low;
988 auto tmp_it = it_fpUserTimeSteps_low;
989 ++tmp_it;
990 if (tmp_it == fpUserTimeSteps->end())
991 {
992 fUserUpperTimeLimit = fStopTime;
993 }
994 else
995 {
996 fUserUpperTimeLimit = tmp_it->first;
997 }
998 }
999 else if (it_fpUserTimeSteps_i == it_fpUserTimeSteps_low)
1000 {
1001 // "Normal" cases
1002 fUserUpperTimeLimit = it_fpUserTimeSteps_i->first;
1003// it_fpUserTimeSteps_i++;
1004// G4cout << "Global time = " << fGlobalTime << G4endl;
1005// G4cout << "Is begin = "
1006// << (it_fpUserTimeSteps_i == fpUserTimeSteps->begin())<< G4endl;
1007
1008 if(it_fpUserTimeSteps_i != fpUserTimeSteps->begin()) it_fpUserTimeSteps_i--;
1009 }
1010 else
1011 {
1012 fUserUpperTimeLimit = it_fpUserTimeSteps_i->first;
1013 it_fpUserTimeSteps_i = it_fpUserTimeSteps_low;
1014 }
1015
1016 return it_fpUserTimeSteps_i->second;
1017}
1018
1019//_________________________________________________________________________
1020
1022{
1023
1024 if(fpUserTimeSteps == nullptr)
1025 {
1026 G4ExceptionDescription exceptionDescription;
1027 exceptionDescription
1028 << "You are asking to use user defined steps but you did not give any.";
1029 G4Exception("G4Scheduler::FindUserPreDefinedTimeStep",
1030 "Scheduler004",
1032 exceptionDescription);
1033 return; // makes coverity happy
1034 }
1035 auto fpUserTimeSteps_i =
1036 fpUserTimeSteps->upper_bound(fGlobalTime);
1037 auto fpUserTimeSteps_low = fpUserTimeSteps
1038 ->lower_bound(fGlobalTime);
1039
1040 // DEBUG
1041 // G4cout << "fGlobalTime : " << G4BestUnit(fGlobalTime,"Time") << G4endl;
1042 // G4cout << "fpUserTimeSteps_i : "
1043 // <<"<"<<G4BestUnit(fpUserTimeSteps_i->first,"Time")<<", "
1044 // << G4BestUnit(fpUserTimeSteps_i->second,"Time")<<">"
1045 // << "\t fpUserTimeSteps_low : "
1046 // <<"<"<<G4BestUnit(fpUserTimeSteps_low->first,"Time")<<", "
1047 // << G4BestUnit(fpUserTimeSteps_low->second,"Time")<<">"
1048 // << G4endl;
1049
1050 if(fpUserTimeSteps_i == fpUserTimeSteps->end())
1051 {
1052 fpUserTimeSteps_i--;
1053 }
1054 else if(fabs(fGlobalTime - fpUserTimeSteps_low->first) < fTimeTolerance)
1055 {
1056 // Case : fGlobalTime = X picosecond
1057 // and fpUserTimeSteps_low->first = X picosecond
1058 // but the precision is not good enough
1059 fpUserTimeSteps_i = fpUserTimeSteps_low;
1060 }
1061 else if(fpUserTimeSteps_i == fpUserTimeSteps_low)
1062 {
1063 // "Normal" cases
1064 fpUserTimeSteps_i--;
1065 }
1066 else
1067 {
1068 fpUserTimeSteps_i = fpUserTimeSteps_low;
1069 }
1070
1071 fDefinedMinTimeStep = fpUserTimeSteps_i->second;
1072}
1073
1074//_________________________________________________________________________
1075
1077{
1078 if(fRunning)
1079 {
1080 G4ExceptionDescription exceptionDescription;
1081 exceptionDescription
1082 << "End tracking is called while G4Scheduler is still running."
1083 << G4endl;
1084
1085 G4Exception("G4Scheduler::EndTracking",
1086 "Scheduler017",
1088 exceptionDescription);
1089 }
1090
1091 fTrackContainer.MergeSecondariesWithMainList();
1092
1093 if (fTrackContainer.MainListsNOTEmpty())
1094 {
1095 G4TrackManyList* mainList = fTrackContainer.GetMainList();
1096 G4TrackManyList::iterator it = mainList->begin();
1097 G4TrackManyList::iterator end = mainList->end();
1098 for (; it != end; ++it)
1099 {
1100 fpTrackingManager->EndTrackingWOKill(*it);
1101 }
1102 }
1103
1104 if (fTrackContainer.SecondaryListsNOTEmpty()) // should be empty
1105 {
1106 G4TrackManyList* secondaries = fTrackContainer.GetSecondariesList();
1107 G4TrackManyList::iterator it = secondaries->begin();
1108 G4TrackManyList::iterator end = secondaries->end();
1109
1110 for (; it != end; ++it)
1111 {
1112 fpTrackingManager->EndTrackingWOKill(*it);
1113 }
1114 }
1115}
1116
1117//_________________________________________________________________________
1119{
1120 fpTrackingInteractivity = interactivity;
1121 if(fpTrackingManager != nullptr)
1122 {
1123 fpTrackingManager->SetInteractivity(fpTrackingInteractivity);
1124 }
1125
1126 //G4MIWorkspace::GetWorldWorkspace()->SetTrackingInteractivity(interactivity);
1127}
1128
1129//_________________________________________________________________________
1131{
1132 fInitialized = false;
1133 Initialize();
1134}
1135
1137{
1138 return fTrackContainer.GetNTracks();
1139}
1140//_________________________________________________________________________
1141
1143{
1144 switch(fITStepStatus)
1145 {
1147 interactionType = "eInteractionWithMedium";
1148 break;
1150 interactionType = "eCollisionBetweenTracks";
1151 break;
1152 default:
1153 interactionType = "eCollisionBetweenTracks";
1154 break;
1155 }
1156}
G4ApplicationState
@ G4State_Quit
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
@ eInteractionWithMedium
@ eCollisionBetweenTracks
@ eUndefined
G4bool IsInf(T value)
#define LIGHT_RED
#define RESET_COLOR
#define G4BestUnit(a, b)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static void DeleteInstance()
virtual void DefineTracks()
Definition G4ITGun.hh:56
void RegisterModel(G4VITStepModel *pModel, G4double globalTime)
void SetModelHandler(G4ITModelHandler *)
void SetTrackingManager(G4ITTrackingManager *trackingManager)
G4double CalculateMinTimeStep(G4double currentGlobalTime, G4double definedMinTimeStep)
void ComputeTrackReaction(G4ITStepStatus fITStepStatus, G4double fGlobalTime, G4double currentTimeStep, G4double previousTimeStep, G4bool reachedUserTimeLimit, G4double fTimeTolerance, G4UserTimeStepAction *fpUserTimeStepAction, G4int fVerbose)
bool GetComputeTimeStep() const
static G4ITReactionSet * Instance()
void SetTrackingManager(G4ITTrackingManager *trackMan)
G4double ComputeInteractionLength(double previousTimeStep)
void DoIt(double timeStep)
virtual void Initialize()
size_t GetNTracks() override
G4TrackList * GetMainList(Key)
void MergeSecondariesWithMainList()
G4TrackManyList * GetSecondariesList()
bool MergeNextTimeToMainList(double &time)
void EndTrackingWOKill(G4Track *)
void SetInteractivity(G4ITTrackingInteractivity *)
void ReserveRessource()
Definition G4ITType.cc:76
static G4ITTypeManager * Instance()
Definition G4ITType.cc:57
void ReleaseRessource()
Definition G4ITType.cc:82
void SynchronizeTracks()
G4bool CanICarryOn()
void EndTracking()
void FindUserPreDefinedTimeStep()
void ForceReinitialization()
G4bool Notify(G4ApplicationState requestedState) override
G4double GetNextWatchedTime() const
void Reset() override
void Stepping()
static G4Scheduler * Instance()
void Process() override
void ClearList()
void Initialize() override
virtual size_t GetNTracks()
void GetCollisionType(G4String &interactionType)
void SetInteractivity(G4ITTrackingInteractivity *) override
G4double GetLimitingTimeStep() const override
void DoProcess()
void PrintWhyDoYouStop()
~G4Scheduler() override
void RegisterModel(G4VITStepModel *, G4double) override
static void DeleteInstance()
void Stop()
void Start()
virtual void UserPostTimeStepAction()
virtual void UserPreTimeStepAction()
#define DBL_MAX
Definition templates.hh:62
#define G4ThreadLocal
Definition tls.hh:77