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