Geant4  10.01.p01
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 <sstream>
57 
58 //#define DEBUG_MEM 1
59 #define DEBUG_MEM_STEPPING
60 #define DEBUG_MEM_DETAILED_STEPPING
61 //#define DEBUG_MEM_DOIT
62 
63 #ifdef DEBUG_MEM
64 #include "G4MemStat.hh"
65 using namespace G4MemStat;
66 using G4MemStat::MemStat;
67 #endif
68 
69 //COLOR FOR DEBUGING
70 //#define USE_COLOR 1
71 
72 #ifdef USE_COLOR
73 #define RED "\033[0;31m"
74 #define LIGHT_RED "\33[1;31m"
75 #define GREEN "\033[32;40m"
76 #define GREEN_ON_BLUE "\033[1;32;44m"
77 #define RESET_COLOR "\033[0m"
78 #else
79 #define RED ""
80 #define LIGHT_RED ""
81 #define GREEN ""
82 #define GREEN_ON_BLUE ""
83 #define RESET_COLOR ""
84 #endif
85 
86 using namespace std;
87 
89 
91 {
92 public:
93  explicit IosFlagSaver(std::ostream& _ios) :
94  ios(_ios), f(_ios.flags())
95  {
96  }
98  {
99  ios.flags(f);
100  }
101 
102  // IosFlagSaver(const IosFlagSaver &rhs) = delete;
103  // IosFlagSaver& operator= (const IosFlagSaver& rhs) = delete;
104 
105 private:
106  std::ostream& ios;
107  std::ios::fmtflags f;
108 };
109 
110 template<typename T>
111  inline bool IsInf(T value)
112  {
113  return std::numeric_limits<T>::has_infinity
114  && value == std::numeric_limits<T>::infinity();
115  }
116 //_________________________________________________________________________
117 
119 {
120  if (fgScheduler == 0) fgScheduler = new G4Scheduler();
121  return fgScheduler;
122 }
123 //_________________________________________________________________________
124 
126 {
127  if (requestedState == G4State_Quit)
128  {
129  if (fVerbose >= 4)
130  {
131  G4cout << "G4Scheduler received G4State_Quit" << G4endl;
132  }
133  Clear();
134  //DeleteInstance();
135  }
136  return true;
137 }
138 //_________________________________________________________________________
139 
141 {
142  if (fgScheduler)
143  {
144  delete fgScheduler;
145  }
146 }
147 //_________________________________________________________________________
148 
151  fTrackContainer((G4ITTrackHolder&)*G4ITTrackHolder::Instance())
152 {
153  Create();
154 }
155 
157 {
158  fUseDefaultTimeSteps = true;
159  fUserUpperTimeLimit = -1;
160  fpGun = 0;
161  fContinue = true;
162 // fpMainList = 0;
163 // fpWaitingList = 0;
165 
167 
168  fpUserTimeSteps = 0;
169 
170  fTimeStep = DBL_MAX;
174 
175  fZeroTimeCount = 0;
177 
178  fStartTime = 0;
179  fComputeTimeStep = false;
180  fComputeReaction = false;
182  fEndTime = 1 * microsecond;
183  fGlobalTime = -1;
184  fInteractionStep = true;
185  fUsePreDefinedTimeSteps = false;
186 
188 
189  fpStepProcessor = 0;
190  // fpMasterStepProcessor = 0 ;
191 
192  fpModelProcessor = 0;
193  // fpMasterModelProcessor= 0;
194 
195  fNbSteps = 0;
196  fMaxSteps = -1;
197 
198  fRunning = false;
199  fInitialized = false;
200 
202  fpTrackingManager = 0;
203 
204  fNbTracks = -1;
207 
208  fVerbose = 0;
209  fWhyDoYouStop = false;
210  fDefinedMinTimeStep = -1.;
211  fReachedUserTimeLimit = false;
212  fTmpEndTime = -1.;
213  fTmpGlobalTime = -1.;
214 
216 
218 }
219 
220 //_________________________________________________________________________
221 
223 {
224 
225  if (fSteppingMsg) // is used as a flag to know whether the manager was cleared
226  {
227  Clear();
228  }
229  fgScheduler = 0;
230 
231 // if (fVerbose >= 1)
232 // {
233 // G4cout << "G4Scheduler is being deleted. Bye :) !" << G4endl;
234 // }
235 }
236 
238 {
239 // if (fVerbose) G4cout << "*** G4Scheduler is being cleared ***" << G4endl;
240 
241  if(fSteppingMsg)
242  {
243  delete fSteppingMsg;
244  fSteppingMsg = 0;
245  }
246  if(fpStepProcessor)
247  {
248  delete fpStepProcessor;
249  fpStepProcessor = 0;
250  }
251  if(fpModelProcessor)
252  {
253  delete fpModelProcessor;
254  fpModelProcessor = 0;
255  }
256  if(fpModelHandler)
257  {
258  delete fpModelHandler;
259  fpModelHandler = 0;
260  }
261  // if(fpMasterStepProcessor) delete fpMasterStepProcessor ;
262  // if(fpMasterModelProcessor) delete fpMasterModelProcessor ;
263  G4ITTypeManager::Instance() -> ReleaseRessource();
264  ClearList();
266  {
267  delete fpTrackingManager;
268  fpTrackingManager = 0;
269  }
270 
271  /*
272  * DEBUG
273  * assert(G4StateManager::GetStateManager()->
274  * DeregisterDependent(this) == true);
275  */
276 }
277 
278 //_________________________________________________________________________
279 
281 {
282 // if (fNbTracks == 0) return;
283 
285 
286  fNbTracks = -1;
287 
289 }
290 
291 //_________________________________________________________________________
293 {
294  fpModelHandler->RegisterModel(model, time);
295 }
296 
297 //_________________________________________________________________________
298 
300 {
301  if (fpStepProcessor)
302  {
303  delete fpStepProcessor;
304  }
305  if (fpModelProcessor)
306  {
307  delete fpModelProcessor;
308  }
309  // if(fpMasterModelProcessor)
310  // {
311  // delete fpMasterModelProcessor;
312  // }
313 
314  //______________________________________________________________
315 
318 
319  // fpMasterModelProcessor = new G4ITModelProcessor();
320  // fpMasterModelProcessor->SetModelHandler(fpModelHandler);
321  // fpModelProcessor = fpMasterModelProcessor;
322 
323  //______________________________________________________________
324 
327 
329 
330  // fpMasterStepProcessor = new G4ITStepProcessor();
331  // fpMasterStepProcessor->SetTrackingManager(fpTrackingManager);
332  // fpStepProcessor = fpMasterStepProcessor ;
333 
336  //______________________________________________________________
337 
338  fInitialized = true;
339 }
340 //_________________________________________________________________________
341 
343 {
344  fStartTime = 0;
345  fUserUpperTimeLimit = -1;
346  fTimeStep = DBL_MAX;
350  fGlobalTime = -1;
351  fInteractionStep = true;
353  fZeroTimeCount = 0;
354 
355  fNbSteps = 0;
356  fContinue = true;
357  fReactingTracks.clear();
358  fLeadingTracks.clear();
359 }
360 //_________________________________________________________________________
361 
363 {
364 
365 #ifdef G4VERBOSE
366  if (fVerbose)
367  {
368  G4cout << "*** G4Scheduler starts processing " << G4endl;
369  if(fVerbose > 2)
370  G4cout << "___________________________________________"
371  "___________________________" << G4endl;
372  }
373 #endif
374 
375  if (fInitialized == false) Initialize();
376 
377  // fpTrackingManager->Initialize();
380 
381  // TODO
382  // fpMasterModelProcessor->Initialize();
383  // fpMasterStepProcessor->Initialize();
384 
385  if (fpGun) fpGun->DefineTracks();
386 
388 
389  // ___________________
390  fRunning = true;
391  Reset();
392 
394 
395 #ifdef G4VERBOSE
396  G4bool trackFound = false;
397 #endif
398 
399  //===========================================================================
400  // By default, before the G4Scheduler is launched, the tracks are pushed to
401  // the delayed lists
402  //===========================================================================
403 
405  {
407 #ifdef G4VERBOSE
408  trackFound = true;
409  G4Timer localtimer;
410  if(fVerbose>1)
411  {
412  localtimer.Start();
413  }
414 #endif
416 #ifdef G4VERBOSE
417  if(fVerbose>1)
418  {
419  localtimer.Stop();
420  G4cout << "G4Scheduler: process time= "<< localtimer << G4endl;
421  }
422 #endif
423  }
424 
425 // //---------------------------------
426 // // TODO: This could be removed ?
427 // if(fTrackContainer.MainListsNOTEmpty()
428 // && (fMaxSteps == -1 ? true : fNbSteps < fMaxSteps)
429 // && fGlobalTime < fEndTime
430 // && fContinue == true)
431 //{
432 //#ifdef G4VERBOSE
433 // trackFound = true;
434 //#endif
435 // DoProcess();
436 //}
437 // //---------------------------------
438 
439 #ifdef G4VERBOSE
440  if(fVerbose)
441  {
442  if(trackFound)
443  {
444  G4cout << "*** G4Scheduler ends at time : "
445  << G4BestUnit(fGlobalTime,"Time") << G4endl;
446  G4cout << "___________________________________" << G4endl;
447  }
448  else
449  {
450  G4cout << "*** G4Scheduler did not start because no "
451  "track was found to be processed"<< G4endl;
452  G4cout << "___________________________________" << G4endl;
453  }
454  }
455 #endif
456 
457  fRunning = false;
458 
460 
461  // ___________________
462  EndTracking();
463  ClearList();
464 
465  Reset();
466 
468 }
469 //_________________________________________________________________________
470 
472 {
473 // if(fTrackContainer.WaitingListsNOTEmpty())
474 // {
475 // G4ExceptionDescription exceptionDescription;
476 // exceptionDescription
477 // << "There is a waiting track list (fpWaitingList != 0).";
478 // exceptionDescription
479 // << " When G4Scheduler::SynchronizeTracks() is called, ";
480 // exceptionDescription
481 // << "no more tracks should remain in the fpWaitingList.";
482 // G4Exception("G4Scheduler::SynchronizeTracks", "ITScheduler002",
483 // FatalErrorInArgument, exceptionDescription);
484 // }
485 
486  // Backup main list and time feature
487  // Reminder : the main list here, should
488  // have the biggest global time
489  // fTrackContainer.MoveMainToWaitingList();
490  // TODO: not yet supported
491 
494 
496  G4double tmpGlobalTime = fGlobalTime;
497  while(fTrackContainer.MergeNextTimeToMainList(tmpGlobalTime))
498  {
499 // assert(tmpGlobalTime == fGlobalTime);
501  DoProcess();
502  }
503 }
504 //_________________________________________________________________________
505 
507 // We split it from the Process() method to avoid repeating code in SynchronizeTracks
508 {
510 
511 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_STEPPING)
512  MemStat mem_first, mem_second, mem_diff;
513 #endif
514 
515 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_STEPPING)
516  mem_first = MemoryUsage();
517 #endif
518 
519  while (fGlobalTime < fEndTime
521  && (fMaxSteps == -1 ? true : fNbSteps < fMaxSteps)
522  && fContinue == true)
523  {
524 // G4cout << "Mainlist size : " << fTrackContainer.GetMainList()->size()
525 // << G4endl;
526 
527  Stepping();
528 
529 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_STEPPING)
530  mem_second = MemoryUsage();
531  mem_diff = mem_second-mem_first;
532  G4cout << "\t || MEM || After step " << fNbSteps << ", diff is : "
533  << mem_diff << G4endl;
534 #endif
535  }
536 
537 #ifdef G4VERBOSE
538  if (fWhyDoYouStop)
539  {
540  G4cout << "G4Scheduler has reached a stage: it might be"
541  " a transition or the end"
542  << G4endl;
543 
544  G4bool normalStop = false;
545 
546  if(fGlobalTime >= fEndTime)
547  {
548  G4cout << "== G4Scheduler: I stop because I reached the final time : "
549  << G4BestUnit(fEndTime,"Time") << " =="<< G4endl;
550  normalStop = true;
551  }
552  if(fTrackContainer.MainListsNOTEmpty() == false) // is empty
553  {
554  G4cout << "G4Scheduler: I stop because the current main list of tracks"
555  "is empty"
556  << G4endl;
557  normalStop = true;
558  }
559  if(fMaxSteps == -1 ? false : fNbSteps >= fMaxSteps)
560  {
561  G4cout << "G4Scheduler: I stop because I reached the maximum allowed "
562  "number of steps=" << fMaxSteps
563  << G4endl;
564  normalStop = true;
565  }
566  if(fContinue && normalStop == false)
567  {
568  G4cout << "G4Scheduler: It might be that I stop because "
569  "I have been told so. You may check "
570  "member fContinue and usage of the method G4Scheduler::Stop()."
571  << G4endl;
572  }
573  }
574 #endif
575 
576 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_STEPPING)
577  mem_second = MemoryUsage();
578  mem_diff = mem_second-mem_first;
579  G4cout << "\t || MEM || After stepping, diff is : " << mem_diff << G4endl;
580 #endif
581 
582 #ifdef G4VERBOSE
583  if (fVerbose > 2) G4cout
584  << "*** G4Scheduler has finished processing a track list at time : "
585  << G4BestUnit(fGlobalTime, "Time") << G4endl;
586 #endif
587 }
588 //_________________________________________________________________________
589 
591 {
592  IosFlagSaver iosfs(G4cout);
593  fTimeStep = DBL_MAX;
594 
597 
598  fInteractionStep = false;
599  fReachedUserTimeLimit = false;
600 
602 
603  // Start of step
604 #ifdef G4VERBOSE
605  if (fVerbose > 2)
606  {
607 #ifdef USE_COLOR
608  G4cout << LIGHT_RED;
609 #endif
610  G4cout << "*** Start Of Step N°" << fNbSteps + 1 << " ***" << G4endl;
611  G4cout << "Current Global time : " << G4BestUnit(fGlobalTime, "Time")
612  <<G4endl;
613 #ifdef USE_COLOR
614  G4cout << RESET_COLOR;
615 #endif
616  }
617 #endif
618 
619 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
620  MemStat mem_first, mem_second, mem_diff;
621 #endif
622 
623 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
624  mem_first = MemoryUsage();
625 #endif
626 
628 
630  {
631  if (fpUserTimeSteps == 0) // Extra checking, is it necessary ?
632  {
633  G4ExceptionDescription exceptionDescription;
634  exceptionDescription
635  << "You are asking to use user defined steps but you did not give any.";
636  G4Exception("G4Scheduler::FindUserPreDefinedTimeStep",
637  "ITScheduler004", FatalErrorInArgument,
638  exceptionDescription);
639  return; // makes coverity happy
640  }
641 #ifdef G4VERBOSE
642  if (fVerbose > 2)
643  {
644 #ifdef USE_COLOR
645  G4cout << LIGHT_RED;
646 #endif
647  G4cout << "*** At time : " << G4BestUnit(fGlobalTime, "Time")
648  << " the chosen user time step is : "
649  << G4BestUnit(fDefinedMinTimeStep, "Time") << " ***" << G4endl;
650 #ifdef USE_COLOR
651  G4cout << RESET_COLOR;
652 #endif
653  }
654 #endif
655  }
656 
657  if (fComputeTimeStep)
658  {
659  CalculateMinTimeStep(); // => at least N (N = nb of tracks) loops
660  }
661  else if(fUseDefaultTimeSteps)
662  {
664  }
665 
666 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
667  mem_second = MemoryUsage();
668  mem_diff = mem_second-mem_first;
669  G4cout << "|| MEM || After computing TS, diff is : " << mem_diff << G4endl;
670 #endif
671 
672 #ifdef G4VERBOSE
673  if (fVerbose > 2)
674  {
675 #ifdef USE_COLOR
676  G4cout << LIGHT_RED;
677 #endif
678  G4cout << "*** Time stepper returned : " << G4BestUnit(fTSTimeStep, "Time")
679  << " ***" << G4endl;
680 #ifdef USE_COLOR
681  G4cout << RESET_COLOR;
682 #endif
683  }
684 #endif
685 
686 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
687  mem_first = MemoryUsage();
688 #endif
689 
690  // Call IL even if fTSTimeStep == 0
691  // if fILTimeStep == 0 give the priority to DoIt processes
692  ComputeInteractionLength(); // => at least N loops
693  // All process returns the physical step of interaction
694  // using the approximation : E = cste along the step
695  // Only the transportation calculates the corresponding
696  // time step
697 
698 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
699  mem_second = MemoryUsage();
700  mem_diff = mem_second-mem_first;
701  G4cout << "|| MEM || After IL, diff is : " << mem_diff << G4endl;
702 #endif
703 
704 #ifdef G4VERBOSE
705  if (fVerbose > 2)
706  {
707 #ifdef USE_COLOR
708  G4cout << LIGHT_RED;
709 #endif
710  G4cout << "*** The minimum time returned by the processes is : "
711  << G4BestUnit(fILTimeStep, "Time") << " ***" << G4endl;
712 #ifdef USE_COLOR
713  G4cout << RESET_COLOR;
714 #endif
715  }
716 #endif
717 
718 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
719  mem_first = MemoryUsage();
720 #endif
721 
722  if (fILTimeStep <= fTSTimeStep)
723  // Give the priority to the IL
724  {
725  fInteractionStep = true;
726  fReactingTracks.clear(); // Give the priority to the IL
729  }
730  else
731  {
732  fInteractionStep = false;
736  }
737 
739  // This check is done at very time step
740  {
742  fITStepStatus = eInteractionWithMedium; // ie: transportation
743  fInteractionStep = true;
744  fReactingTracks.clear();
746  }
747 
748  if (fTimeStep == 0) // < fTimeTolerance)
749  {
750  fZeroTimeCount++;
752  {
753  G4ExceptionDescription exceptionDescription;
754 
755  exceptionDescription << "Too many zero time steps were detected. ";
756  exceptionDescription << "The simulation is probably stuck. ";
757  exceptionDescription
758  << "The maximum number of zero time steps is currently : "
760  exceptionDescription << ".";
761 
762  G4Exception("G4Scheduler::Stepping", "ITSchedulerNullTimeSteps",
763  FatalErrorInArgument, exceptionDescription);
764  }
765  }
766  else
767  {
768  fZeroTimeCount = 0;
769  }
770 
772  ((fTimeStep <= fDefinedMinTimeStep) || ((fTimeStep > fDefinedMinTimeStep)
774  true : false;
775 
777  // TODO: pre/post
778 
779 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
780  mem_second = MemoryUsage();
781  mem_diff = mem_second-mem_first;
782  G4cout << "|| MEM || After LeadingTracks and UserPreTimeStepAction: "
783  << mem_diff << G4endl;
784 #endif
785 
786 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
787  mem_first = MemoryUsage();
788 #endif
789 
790  // if fTSTimeStep > 0 => needs to call the transportation process
791  // if fILTimeStep < fTSTimeStep => call DoIt processes
792  // if fILTimeStep == fTSTimeStep => give the priority to the DoIt processes
793  if (fTSTimeStep > 0 || fILTimeStep <= fTSTimeStep)
794  {
795  // G4cout << "Will call DoIT" << G4endl;
796  DoIt();
797 
799  KillTracks();
800  }
801  // else
802  // {
803  // G4cout << "fTSTimeStep : " << fTSTimeStep << G4endl;
804  // G4cout << "fILTimeStep : " << fILTimeStep << G4endl;
805  // }
806 
807 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
808  mem_second = MemoryUsage();
809  mem_diff = mem_second-mem_first;
810  G4cout << "|| MEM || After DoIT, diff is : " << mem_diff << G4endl;
811 #endif
812 
814 
815 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
816  mem_first = MemoryUsage();
817 #endif
818 
820 
822 
824 
825  fNbSteps++;
826 
828  {
830  }
831 
833 
834 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
835  mem_second = MemoryUsage();
836  mem_diff = mem_second-mem_first;
837  G4cout << "|| MEM || After computing reactions + UserPostTimeStepAction, "
838  "diff is : " << mem_diff << G4endl;
839 #endif
840 
841  // End of step
842 #ifdef G4VERBOSE
843  if (fVerbose >= 2)
844  {
845 #ifdef USE_COLOR
846  G4cout << LIGHT_RED;
847 #endif
848 
849  G4String interactionType;
850  GetCollisionType(interactionType);
851 
852  std::stringstream finalOutput;
853 
854  finalOutput << "*** End of step N°" << fNbSteps
855  << "\t T_i= " << G4BestUnit(fGlobalTime-fTimeStep, "Time")
856  << "\t dt= " << G4BestUnit(fTimeStep, "Time")
857  << "\t T_f= " << G4BestUnit(fGlobalTime, "Time")
858  << "\t " << interactionType
859  << G4endl;
860 
861  if(fVerbose>2)
862  {
864  {
865  finalOutput << "It has also reached the user time limit" << G4endl;
866  }
867  finalOutput << "_______________________________________________________________"
868  "_______"<< G4endl;
869  }
870 
871  G4cout << finalOutput.str();
872 
873 #ifdef USE_COLOR
874  G4cout << RESET_COLOR;
875 #endif
876  }
877 #endif
878 
879 }
880 //_________________________________________________________________________
881 
883 {
884  if (fpUserTimeSteps == 0) return fDefaultMinTimeStep;
886  return fDefinedMinTimeStep;
887 
888  map<double, double>::const_iterator it_fpUserTimeSteps_i = fpUserTimeSteps
889  ->upper_bound(fGlobalTime);
890  map<double, double>::const_iterator it_fpUserTimeSteps_low = fpUserTimeSteps
891  ->lower_bound(fGlobalTime);
892 
893  // DEBUG
894  // G4cout << "fGlobalTime : " << G4BestUnit(fGlobalTime,"Time")
895  // << G4endl;
896  // G4cout << "fpUserTimeSteps_i : "
897  // <<"<"<<G4BestUnit(it_fpUserTimeSteps->first,"Time")
898  // <<", "<< G4BestUnit(it_fpUserTimeSteps->second,"Time")<<">"
899  // << "\t fpUserTimeSteps_low : "
900  // <<"<"<<G4BestUnit(fpUserTimeSteps_low->first,"Time")<<", "*
901  // << G4BestUnit(fpUserTimeSteps_low->second,"Time")<<">"
902  // << G4endl;
903 
904  if (it_fpUserTimeSteps_i == fpUserTimeSteps->end())
905  {
906  it_fpUserTimeSteps_i--;
908  }
909  else if (fabs(fGlobalTime - it_fpUserTimeSteps_low->first) < fTimeTolerance)
910  {
911  // Case : fGlobalTime = X picosecond
912  // and fpUserTimeSteps_low->first = X picosecond
913  // but the precision is not good enough
914  it_fpUserTimeSteps_i = it_fpUserTimeSteps_low;
915  map<double, double>::const_iterator tmp_it = it_fpUserTimeSteps_low;
916  tmp_it++;
917  if (tmp_it == fpUserTimeSteps->end())
918  {
920  }
921  else
922  {
923  fUserUpperTimeLimit = tmp_it->second;
924  }
925  }
926  else if (it_fpUserTimeSteps_i == it_fpUserTimeSteps_low)
927  {
928  // "Normal" cases
929  fUserUpperTimeLimit = it_fpUserTimeSteps_i->second;
930  it_fpUserTimeSteps_i--;
931  }
932  else
933  {
934  fUserUpperTimeLimit = it_fpUserTimeSteps_i->second;
935  it_fpUserTimeSteps_i = it_fpUserTimeSteps_low;
936  }
937 
938  return it_fpUserTimeSteps_i->second;
939 }
940 
941 //_________________________________________________________________________
942 
944 {
945 
946  if (fpUserTimeSteps == 0)
947  {
948  G4ExceptionDescription exceptionDescription;
949  exceptionDescription
950  << "You are asking to use user defined steps but you did not give any.";
951  G4Exception("G4Scheduler::FindUserPreDefinedTimeStep",
952  "ITScheduler004", FatalErrorInArgument, exceptionDescription);
953  return; // makes coverity happy
954  }
955  map<double, double>::iterator fpUserTimeSteps_i =
956  fpUserTimeSteps->upper_bound(fGlobalTime);
957  map<double, double>::iterator fpUserTimeSteps_low = fpUserTimeSteps
958  ->lower_bound(fGlobalTime);
959 
960  // DEBUG
961  // G4cout << "fGlobalTime : " << G4BestUnit(fGlobalTime,"Time") << G4endl;
962  // G4cout << "fpUserTimeSteps_i : "
963  // <<"<"<<G4BestUnit(fpUserTimeSteps_i->first,"Time")<<", "
964  // << G4BestUnit(fpUserTimeSteps_i->second,"Time")<<">"
965  // << "\t fpUserTimeSteps_low : "
966  // <<"<"<<G4BestUnit(fpUserTimeSteps_low->first,"Time")<<", "
967  // << G4BestUnit(fpUserTimeSteps_low->second,"Time")<<">"
968  // << G4endl;
969 
970  if (fpUserTimeSteps_i == fpUserTimeSteps->end())
971  {
972  fpUserTimeSteps_i--;
973  }
974  else if (fabs(fGlobalTime - fpUserTimeSteps_low->first) < fTimeTolerance)
975  {
976  // Case : fGlobalTime = X picosecond
977  // and fpUserTimeSteps_low->first = X picosecond
978  // but the precision is not good enough
979  fpUserTimeSteps_i = fpUserTimeSteps_low;
980  }
981  else if (fpUserTimeSteps_i == fpUserTimeSteps_low)
982  {
983  // "Normal" cases
984  fpUserTimeSteps_i--;
985  }
986  else
987  {
988  fpUserTimeSteps_i = fpUserTimeSteps_low;
989  }
990 
991  fDefinedMinTimeStep = fpUserTimeSteps_i->second;
992 }
993 //_________________________________________________________________________
994 
996 {
997 
998  if (fpModelProcessor == 0)
999  // if(fpMasterModelProcessor == 0)
1000  {
1001  G4ExceptionDescription exceptionDescription;
1002  exceptionDescription
1003  << "There is no G4ITModelProcessor to hande IT reaction. ";
1004  exceptionDescription
1005  << "You probably did not initialize the G4Scheduler. ";
1006  exceptionDescription
1007  << "Just do G4Scheduler::Instance()->Initialize(); ";
1008  exceptionDescription << " but only after initializing the run manager.";
1009  G4Exception("G4Scheduler::CalculateMinStep", "ITScheduler005",
1010  FatalErrorInArgument, exceptionDescription);
1011  return; // makes coverity happy
1012  }
1013 
1014 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
1015  MemStat mem_first, mem_second, mem_diff;
1016 #endif
1017 
1018 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
1019  mem_first = MemoryUsage();
1020 #endif
1021 
1023  // TODO
1024  // fpMasterModelProcessor -> InitializeStepper(fGlobalTime, fDefinedMinTimeStep) ;
1025 
1026 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
1027  mem_second = MemoryUsage();
1028  mem_diff = mem_second-mem_first;
1029  G4cout << "\t || MEM || G4Scheduler::CalculateMinTimeStep || After "
1030  "computing fpModelProcessor -> InitializeStepper, diff is : "
1031  << mem_diff
1032  << G4endl;
1033 #endif
1034 
1035 // G4TrackList::iterator fpMainList_i = fpMainList->begin();
1036 
1038  G4TrackManyList::iterator it = mainList->begin();
1039  G4TrackManyList::iterator end = mainList->end();
1040 
1041  for (; it != end; it++)
1042  {
1043  G4Track * track = *it;
1044 
1045  if (track == 0)
1046  {
1047  G4ExceptionDescription exceptionDescription;
1048  exceptionDescription << "No track found.";
1049  G4Exception("G4Scheduler::CalculateMinStep", "ITScheduler006",
1050  FatalErrorInArgument, exceptionDescription);
1051  return; // makes coverity happy
1052  }
1053 
1054 #ifdef DEBUG
1055  G4cout << "*_* " << GetIT(track)->GetName()
1056  << " ID: " << track->GetTrackID()
1057  << " at time : " << track->GetGlobalTime()
1058  << G4endl;
1059 #endif
1060 
1061  G4TrackStatus trackStatus = track->GetTrackStatus();
1062  if (trackStatus == fStopAndKill || trackStatus == fStopButAlive)
1063  {
1064  continue;
1065  }
1066 
1067  // This Extract... was thought for MT mode at the track level
1068  //ExtractTimeStepperData(fpModelProcessor) ;
1069 
1071 
1072  // if MT mode at track level, this command should be displaced
1074  }
1075 
1076 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
1077  mem_second = MemoryUsage();
1078  mem_diff = mem_second-mem_first;
1079  G4cout << "\t || MEM || G4Scheduler::CalculateMinTimeStep || "
1080  "After looping on tracks, diff is : " << mem_diff << G4endl;
1081 #endif
1082 }
1083 //_________________________________________________________________________
1084 
1086 {
1087  G4Track* track = (G4Track*) MP->GetTrack();
1088  if (track == 0)
1089  {
1090  MP->CleanProcessor();
1091  return;
1092  }
1093 
1094  const std::vector<std::vector<G4VITStepModel*> >* model =
1095  MP->GetCurrentModel();
1096 
1097  for (unsigned i = 0; i < model->size(); i++)
1098  {
1099  for (unsigned j = 0; j < (*model)[i].size(); j++)
1100  {
1101  G4VITStepModel* mod = (*model)[i][j];
1102 
1103  if (mod == 0)
1104  {
1105  continue;
1106  }
1107 
1108  G4VITTimeStepComputer* stepper(mod->GetTimeStepper());
1109 
1110  G4double sampledMinTimeStep(stepper->GetSampledMinTimeStep());
1111  G4TrackVectorHandle reactants(stepper->GetReactants());
1112 
1113  if (sampledMinTimeStep < fTSTimeStep)
1114  {
1115  /*
1116  // DEBUG SPECIAL CASE
1117  if(!reactants)
1118  {
1119  G4ExceptionDescription exceptionDescription ;
1120  exceptionDescription << "No reactants were found by the time stepper.";
1121  G4Exception("G4Scheduler::ExtractTimeStepperData","ITScheduler007",
1122  FatalErrorInArgument,exceptionDescription);
1123  continue ;
1124  }
1125  */
1126 
1127  fTSTimeStep = sampledMinTimeStep;
1128  fReactingTracks.clear();
1129  if (bool(reactants))
1130  {
1131  // G4cout << "*** (1) G4Scheduler::ExtractTimeStepperData insert
1132  // reactants for " << GetIT(track)->GetName() << G4endl;
1133  // G4cout << "bool(reactants) = " << bool(reactants) << G4endl;
1134  // G4cout << reactants->size() << G4endl;
1135  // G4cout << GetIT(reactants->operator[](0))->GetName() << G4endl;
1136 
1137  fReactingTracks.insert(make_pair(track, reactants));
1138  stepper->ResetReactants();
1139  }
1140  }
1141  else if (fTSTimeStep == sampledMinTimeStep)
1142  {
1143  /*
1144  // DEBUG SPECIAL CASE
1145  if(!reactants)
1146  {
1147  G4ExceptionDescription exceptionDescription ;
1148  exceptionDescription << "No reactants were found by the time stepper.";
1149  G4Exception("G4Scheduler::ExtractTimeStepperData","ITScheduler008",
1150  FatalErrorInArgument,exceptionDescription);
1151  continue ;
1152  }
1153  */
1154  if (bool(reactants))
1155  {
1156  // G4cout << "*** (2) G4Scheduler::ExtractTimeStepperData insert
1157  // reactants for " << GetIT(track)->GetName() << G4endl;
1158  // G4cout << "bool(reactants) = " << bool(reactants) << G4endl;
1159  // G4cout << "trackA : " << GetIT(track)->GetName()
1160  // << " ("<< track->GetTrackID() << ")" << G4endl;
1161  // G4cout << reactants->size() << G4endl;
1162  // G4cout << GetIT(reactants->operator[](0))->GetName() << G4endl;
1163 
1164  fReactingTracks.insert(make_pair(track, reactants));
1165  stepper->ResetReactants();
1166  }
1167  }
1168  else
1169  {
1170  if (bool(reactants))
1171  {
1172  stepper->ResetReactants();
1173  }
1174  }
1175  }
1176  }
1177 
1178  MP->CleanProcessor();
1179 }
1180 //_________________________________________________________________________
1181 
1183 {
1185  G4TrackManyList::iterator it = mainList ->begin();
1186  G4TrackManyList::iterator end = mainList ->end();
1187 
1189 
1190  for (; it != end; it++)
1191  {
1192  G4Track * track = *it;
1193 
1194 #ifdef DEBUG
1195  G4cout << "*CIL* " << GetIT(track)->GetName()
1196  << " ID: " << track->GetTrackID()
1197  << " at time : " << track->GetGlobalTime()
1198  << G4endl;
1199 #endif
1200 
1202 
1204  }
1205 }
1206 //_________________________________________________________________________
1207 
1209 {
1210  G4Track* track = SP->GetTrack();
1211  if (!track)
1212  {
1213  SP->CleanProcessor();
1214  return;
1215  }
1216  if (track->GetTrackStatus() == fStopAndKill)
1217  {
1218  fTrackContainer.GetMainList()->pop(track);
1219  EndTracking(track);
1220  //return;
1221  }
1222 
1223  if (IsInf(SP->GetInteractionTime()))
1224  {
1225  // G4cout << "!!!!!!!!!!!! IS INF " << track->GetTrackID() << G4endl;
1226  SP->CleanProcessor();
1227  return;
1228  }
1229  else if (SP->GetInteractionTime() < fILTimeStep - DBL_EPSILON)
1230  {
1231  // G4cout << "!!!!!!!!!!!! TEMPS DIFFERENTS "
1232  // << track->GetTrackID() << G4endl;
1233 
1235 
1236  fILTimeStep = SP -> GetInteractionTime();
1237 
1238 // G4cout << "Will set leading step to true for time :"
1239 // << SP -> GetInteractionTime() << " against fTimeStep : "
1240 // << G4BestUnit(fILTimeStep, "Time") << " the trackID is : " << track->GetTrackID()
1241 // << G4endl;
1242 
1243  GetIT(track)->GetTrackingInfo()->SetLeadingStep(true);
1244  fLeadingTracks.push_back(track);
1245  }
1246  else if(fabs(fILTimeStep - SP -> GetInteractionTime()) < DBL_EPSILON )
1247  {
1248 
1249  // G4cout << "!!!!!!!!!!!! MEME TEMPS " << track->GetTrackID() << G4endl;
1250  // G4cout << "Will set leading step to true for time :"
1251  // << SP -> GetInteractionTime() << " against fTimeStep : "
1252  // << fTimeStep << " the trackID is : " << track->GetTrackID()<< G4endl;
1253  GetIT(track)->GetTrackingInfo()->SetLeadingStep(true);
1254  fLeadingTracks.push_back(track);
1255  }
1256  // else
1257  // {
1258  // G4cout << "!!!! Bigger time : " << "currentTime : "<<fILTimeStep
1259  // << " proposedTime : " << SP -> GetInteractionTime() << G4endl;
1260  // }
1261 
1262  SP->CleanProcessor();
1263 }
1264 //_________________________________________________________________________
1265 
1267 
1268 // Call the process having the min step length or just propagate the track on the given time step
1269 
1270 // If the track is "leading the step" (ie one of its process has been selected
1271 // as the one having the minimum time step over all tracks and processes),
1272 // it will undergo its selected processes. Otherwise, it will just propagate the track
1273 // on the given time step.
1274 
1275 {
1276 
1277 #ifdef G4VERBOSE
1278  if (fVerbose > 3)
1279  {
1280 #ifdef USE_COLOR
1281  G4cout << LIGHT_RED;
1282 #endif
1283  G4cout << "*** G4Scheduler::DoIt ***" << G4endl;
1284  G4cout << setw(18) << left << "#Name"
1285  << setw(15) << "trackID"
1286  << setw(35) << "Position"
1287  << setw(25) << "Pre step volume"
1288  << setw(25) << "Post step volume"
1289  << setw(22) << "Process"
1290  << G4endl;
1291 #ifdef USE_COLOR
1292  G4cout << RESET_COLOR;
1293 #endif
1294  }
1295 #endif
1296 
1298  G4TrackManyList::iterator it = mainList->end();
1299  it--;
1300  size_t initialSize = mainList->size();
1301 
1302 // G4cout << "initialSize = " << initialSize << G4endl;
1303 
1304  for(size_t i = 0 ; i < initialSize ; i++)
1305  {
1306 
1307 // G4cout << "i = " << i << G4endl;
1308 
1309  G4Track* track = *it;
1310  if (!track)
1311  {
1312  G4ExceptionDescription exceptionDescription;
1313  exceptionDescription << "No track was pop back the main track list.";
1314  G4Exception("G4Scheduler::DoIt", "ITScheduler009",
1315  FatalErrorInArgument, exceptionDescription);
1316  }
1317 
1318 #ifdef DEBUG
1319 #ifdef USE_COLOR
1320  G4cout << LIGHT_RED;
1321 #endif
1322  G4cout << "*DoIt* " << GetIT(track)->GetName()
1323  << " ID: " << track->GetTrackID()
1324  << " at time : " << track->GetGlobalTime()
1325  << G4endl;
1326 #ifdef USE_COLOR
1327  G4cout << RESET_COLOR;
1328 #endif
1329 #endif
1330 // G4TrackManyList::iterator next_it (it);
1331 // next_it--;
1332 // it = next_it;
1333 
1334  it--;
1335 
1336 #ifdef G4VERBOSE
1337  // PRE STEP VERBOSE
1339  if(fVerbose > 3)
1340  {
1341  G4String volumeName;
1342 
1343  G4TouchableHandle nextTouchable = track->GetNextTouchableHandle();
1344  G4VPhysicalVolume* volume (0);
1345 
1346  if(nextTouchable
1347  && (volume = nextTouchable->GetVolume()))
1348  {
1349  volumeName = volume->GetName();
1350 
1351  if(volume->IsParameterised() || volume->IsReplicated())
1352  {
1353  volumeName+=" ";
1354  volumeName+=nextTouchable->GetReplicaNumber();
1355  }
1356  }
1357  else
1358  {
1359  volumeName = "OutOfWorld";
1360  }
1361 
1362  G4cout << setw(18) << left << GetIT(track)->GetName()
1363  << setw(15) << track->GetTrackID()
1364  << std::setprecision(3)
1365  << setw(35) << G4String(G4BestUnit(track->GetPosition(), "Length"))
1366  << setw(25) << volumeName
1367  << setw(25) << "---"
1368  << G4endl;
1369  }
1370 #endif
1371 
1372 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DOIT)
1373  MemStat mem_first, mem_second, mem_diff;
1374 #endif
1375 
1376 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DOIT)
1377  mem_first = MemoryUsage();
1378 #endif
1379 
1380  // G4cout << " G4Scheduler::DoIt -- fTimeStep="
1381 // << fTimeStep << G4endl;
1383 
1384 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DOIT)
1385  MemStat mem_intermediaire = MemoryUsage();
1386  mem_diff = mem_intermediaire-mem_first;
1387  G4cout << "\t\t >> || MEM || In DoIT with track "
1388  << track->GetTrackID() << ", diff is : " << mem_diff << G4endl;
1389 #endif
1390 
1392 
1393 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DOIT)
1394  mem_second = MemoryUsage();
1395  mem_diff = mem_second-mem_first;
1396  G4cout << "\t >> || MEM || In DoIT with track "
1397  << track->GetTrackID()
1398  << ", diff is : " << mem_diff << G4endl;
1399 #endif
1400 
1401 #ifdef G4VERBOSE
1402  // POST STEP VERBOSE
1404  if(fVerbose > 3)
1405  {
1406  G4cout << setw(18) << left << GetIT(track)->GetName()
1407  << setw(15) << track->GetTrackID()
1408  << std::setprecision(3)
1409  << setw(35) << G4String(G4BestUnit(track->GetPosition(), "Length" ))
1410  << setw(25) << "---";
1411 
1412  G4TouchableHandle nextTouchable = track->GetNextTouchableHandle();
1413  G4VPhysicalVolume* volume (0);
1414 
1415  if(nextTouchable
1416  && (volume = nextTouchable->GetVolume()))
1417  {
1418  G4String volumeName = volume->GetName();
1419 
1420  if(volume->IsParameterised() || volume->IsReplicated())
1421  {
1422  volumeName+=" ";
1423  volumeName+=nextTouchable->GetReplicaNumber();
1424  }
1425 
1426  G4cout << setw(25) << volumeName;
1427  }
1428  else
1429  {
1430  G4cout << setw(25) << "OutOfWorld";
1431  }
1433  {
1434  G4cout << setw(22)
1436  ->GetProcessName();
1437  }
1438  else
1439  {
1440  G4cout << "---";
1441  }
1442  G4cout << G4endl;
1443 
1444  if(fVerbose > 4)
1445  {
1446  const G4TrackVector* secondaries = 0;
1447  if((secondaries = track->GetStep()->GetSecondary()))
1448  {
1449  if(secondaries->empty() == false)
1450  {
1451  G4cout << "\t\t ---->";
1452  for(size_t j = 0; j < secondaries->size(); j++)
1453  {
1454  G4cout << GetIT((*secondaries)[j])->GetName()
1455  << "("<<(*secondaries)[j]->GetTrackID() << ")"<< " ";
1456  }
1457  G4cout << G4endl;
1458  }
1459  }
1460  }
1461  }
1462 #endif
1463  }
1464 }
1465 //_________________________________________________________________________
1466 
1468 {
1469 
1470  G4Track* track = SP->GetTrack();
1471  if (!track)
1472  {
1473  SP->CleanProcessor();
1474  return;
1475  }
1476 
1477  G4TrackStatus status = track->GetTrackStatus();
1478 
1479  switch (status)
1480  {
1481  case fAlive:
1482  case fStopButAlive:
1483  case fSuspend:
1484  case fPostponeToNextEvent:
1485  default:
1486  PushSecondaries(SP);
1487  break;
1488 
1489  case fStopAndKill:
1490  PushSecondaries(SP);
1491  G4TrackList::Pop(track);
1492  EndTracking(track);
1493  break;
1494 
1496  G4TrackVector* secondaries = SP->GetSecondaries();
1497  if (secondaries)
1498  {
1499  for (size_t i = 0; i < secondaries->size(); i++)
1500  {
1501  delete (*secondaries)[i];
1502  }
1503  secondaries->clear();
1504  }
1505  G4TrackList::Pop(track);
1506  EndTracking(track);
1507  break;
1508  }
1509  SP->CleanProcessor();
1510 }
1511 //_________________________________________________________________________
1512 
1514 {
1515  G4TrackVector* secondaries = SP->GetSecondaries();
1516  if (!secondaries || secondaries->empty())
1517  {
1518  // DEBUG
1519  // G4cout << "NO SECONDARIES !!! " << G4endl;
1520  return;
1521  }
1522 
1523  // DEBUG
1524  // G4cout << "There are secondaries : "<< secondaries -> size() << G4endl ;
1525 
1526  G4TrackVector::iterator secondaries_i = secondaries->begin();
1527 
1528  for (; secondaries_i != secondaries->end(); secondaries_i++)
1529  {
1530  G4Track* secondary = *secondaries_i;
1531  fTrackContainer._PushTrack(secondary);
1532  }
1533 }
1534 
1535 //_________________________________________________________________________
1536 
1538 {
1539  if (fReactingTracks.empty())
1540  {
1541  return;
1542  }
1543 
1545  // if(fInteractionStep == false)
1546  {
1549  // TODO
1550  // A ne faire uniquement si le temps choisis est celui calculé par le time stepper
1551  // Sinon utiliser quelque chose comme : fModelProcessor->FindReaction(&fMainList);
1552 
1553  std::vector<G4ITReactionChange*> * reactionInfo_v = fpModelProcessor
1554  ->GetReactionInfo();
1555  std::vector<G4ITReactionChange*>::iterator reactionInfo_i = reactionInfo_v
1556  ->begin();
1557 
1558  for (; reactionInfo_i != reactionInfo_v->end(); reactionInfo_i++)
1559  {
1560  G4ITReactionChange* changes = (*reactionInfo_i);
1561  G4Track* trackA = const_cast<G4Track*>(changes->GetTrackA());
1562  G4Track* trackB = const_cast<G4Track*>(changes->GetTrackB());
1563 
1564  if (trackA == 0 || trackB == 0 || trackA->GetTrackStatus() == fStopAndKill
1565  || trackB->GetTrackStatus() == fStopAndKill) continue;
1566 
1567  G4int nbSecondaries = changes->GetNumberOfSecondaries();
1568 
1570  {
1571  const vector<G4Track*>* productsVector = changes->GetfSecondary();
1572  fpUserTimeStepAction->UserReactionAction(*trackA, *trackB,
1573  *productsVector);
1574  }
1575 
1576 #ifdef G4VERBOSE
1577  if (fVerbose)
1578  {
1579  G4cout << "At time : " << setw(7) << G4BestUnit(fGlobalTime, "Time")
1580  << " Reaction : " << GetIT(trackA)->GetName() << " ("
1581  << trackA->GetTrackID() << ") + " << GetIT(trackB)->GetName()
1582  << " (" << trackB->GetTrackID() << ") -> ";
1583  }
1584 #endif
1585 
1586  if (nbSecondaries > 0)
1587  {
1588  for (int i = 0; i < nbSecondaries; i++)
1589  {
1590 #ifdef G4VERBOSE
1591  if (fVerbose && i != 0) G4cout << " + ";
1592 #endif
1593 
1594  G4Track* secondary = changes->GetSecondary(i);
1595  fTrackContainer._PushTrack(secondary);
1596  GetIT(secondary)->SetParentID(trackA->GetTrackID(),
1597  trackB->GetTrackID());
1598 
1599  if (secondary->GetGlobalTime() - fGlobalTime > fTimeTolerance)
1600  {
1601  G4ExceptionDescription exceptionDescription;
1602  exceptionDescription
1603  << "The time of the secondary should not be bigger than the"
1604  " current global time."
1605  << " This may cause synchronization problem. If the process you"
1606  " are using required "
1607  << "such feature please contact the developpers."
1608  << G4endl
1609  << "The global time in the step manager : "
1610  << G4BestUnit(fGlobalTime,"Time")
1611  << G4endl
1612  << "The global time of the track : "
1613  << G4BestUnit(secondary->GetGlobalTime(),"Time")
1614  << G4endl;
1615 
1616  G4Exception("G4Scheduler::ComputeInteractionBetweenTracks",
1617  "ITScheduler010", FatalErrorInArgument,
1618  exceptionDescription);
1619  }
1620 
1621 #ifdef G4VERBOSE
1622  if (fVerbose) G4cout << GetIT(secondary)->GetName() << " ("
1623  << secondary->GetTrackID() << ")";
1624 #endif
1625  }
1626  }
1627  else
1628  {
1629 #ifdef G4VERBOSE
1630  if (fVerbose) G4cout << "No product";
1631 #endif
1632  }
1633 #ifdef G4VERBOSE
1634  if (fVerbose) G4cout << G4endl;
1635 #endif
1636  if (trackA->GetTrackID() == 0 || trackB->GetTrackID() == 0)
1637  {
1638  G4Track* track = 0;
1639  if (trackA->GetTrackID() == 0) track = trackA;
1640  else track = trackB;
1641 
1642  G4ExceptionDescription exceptionDescription;
1643  exceptionDescription
1644  << "The problem was found for the reaction between tracks :"
1645  << trackA->GetParticleDefinition()->GetParticleName() << " ("
1646  << trackA->GetTrackID() << ") & "
1647  << trackB->GetParticleDefinition()->GetParticleName() << " ("
1648  << trackB->GetTrackID() << "). \n";
1649 
1650  if (track->GetStep() == 0)
1651  {
1652  exceptionDescription << "Also no step was found"
1653  << " ie track->GetStep() == 0 \n";
1654  }
1655 
1656  exceptionDescription << "Parent ID of trackA : "
1657  << trackA->GetParentID() << "\n";
1658  exceptionDescription << "Parent ID of trackB : "
1659  << trackB->GetParentID() << "\n";
1660 
1661  exceptionDescription
1662  << "The ID of one of the reaction track was not setup.";
1663  G4Exception("G4Scheduler::ComputeInteractionBetweenTracks",
1664  "ITScheduler011", FatalErrorInArgument,
1665  exceptionDescription);
1666  }
1667 
1668  if (changes->WereParentsKilled())
1669  {
1670  // DEBUG
1671  // G4cout << "Erasing tracks : "
1672  // << "trackA at time : " << G4BestUnit(trackA->GetGlobalTime() , "Time")
1673  // << "\t trackB at time : "<< G4BestUnit(trackB->GetGlobalTime(), "Time")
1674  // << "\t GlobalTime : " << G4BestUnit(fGlobalTime, "Time")
1675  // << G4endl;
1676 
1677  trackA->SetTrackStatus(fStopAndKill);
1678  trackB->SetTrackStatus(fStopAndKill);
1679 
1680  G4TrackList::Pop(trackA);
1681  G4TrackList::Pop(trackB);
1682  EndTracking(trackA);
1683  EndTracking(trackB);
1684  }
1685 
1686  delete changes;
1687  }
1688 
1689  reactionInfo_v->clear(); // never null because point to a non-pointer member of ModelProcessor
1690  }
1691  // DEBUG
1692  //else
1693  //{
1694  // G4cout << "fInteractionStep == true" << G4endl ;
1695  //}
1696 
1697  fReactingTracks.clear();
1698 }
1699 
1700 //_________________________________________________________________________
1701 
1703 {
1704  if (fNbTracks == 0) fNbTracks = -1;
1705  track->SetTrackID(fNbTracks);
1706  fNbTracks--;
1707 }
1708 
1709 //_________________________________________________________________________
1710 
1712 {
1713  if (fRunning)
1714  {
1715  G4ExceptionDescription exceptionDescription;
1716  exceptionDescription
1717  << "G4Scheduler::PushTrack : You are trying to push tracks while the "
1718  "ITScheduler is running";
1719  G4Exception("G4Scheduler::PushTrack", "ITScheduler012",
1720  FatalErrorInArgument, exceptionDescription);
1721  }
1722  fTrackContainer._PushTrack(track);
1723 }
1724 //_________________________________________________________________________
1725 
1727 {
1729 }
1730 //_________________________________________________________________________
1731 
1733 {
1734  if (fLeadingTracks.empty() == false)
1735  {
1736  std::vector<G4Track*>::iterator fLeadingTracks_i = fLeadingTracks.begin();
1737 
1738  while (fLeadingTracks_i != fLeadingTracks.end())
1739  {
1740  G4Track* track = *fLeadingTracks_i;
1741  if (track)
1742  {
1743  G4IT* ITrack = GetIT(*fLeadingTracks_i);
1744  if (ITrack)
1745  {
1746  GetIT(*fLeadingTracks_i)->GetTrackingInfo()->SetLeadingStep(false);
1747  }
1748  }
1749 
1750  fLeadingTracks_i++;
1751  continue;
1752  }
1753 
1754  fLeadingTracks.clear();
1755  }
1756 }
1757 
1758 //_________________________________________________________________________
1759 
1760 void G4Scheduler::EndTracking(G4Track* trackToBeKilled)
1761 {
1762  fpTrackingManager->EndTracking(trackToBeKilled);
1763  fTrackContainer.PushToKill(trackToBeKilled);
1764 }
1765 //_________________________________________________________________________
1766 
1768 {
1769  if (fRunning)
1770  {
1771  G4ExceptionDescription exceptionDescription;
1772  exceptionDescription
1773  << "End tracking is called while G4Scheduler is still running."
1774  << G4endl;
1775 
1776  G4Exception("G4Scheduler::EndTracking", "ITScheduler017",
1777  FatalErrorInArgument, exceptionDescription);
1778  }
1779 
1781 
1783  {
1785  G4TrackManyList::iterator it = mainList->begin();
1786  G4TrackManyList::iterator end = mainList->end();
1787  for (; it != end; it++)
1788  {
1790  }
1791  }
1792 
1793  if (fTrackContainer.SecondaryListsNOTEmpty()) // should be empty
1794  {
1796  G4TrackManyList::iterator it = secondaries->begin();
1797  G4TrackManyList::iterator end = secondaries->end();
1798 
1799  for (; it != end; it++)
1800  {
1802  }
1803  }
1804 }
1805 
1806 //_________________________________________________________________________
1808 {
1809  fpTrackingInteractivity = interactivity;
1812 }
1813 
1814 //_________________________________________________________________________
1816 {
1817  fInitialized = false;
1818  Initialize();
1819 }
1820 
1821 //_________________________________________________________________________
1824 fTrackContainer((G4ITTrackHolder&)*G4ITTrackHolder::Instance())
1825 
1826 {
1827  Create();
1828 }
1829 
1830 //_________________________________________________________________________
1832 {
1833  if (this != &right)
1834  {
1835  Create();
1836  }
1837  return *this;
1838 }
1839 
1841 {
1842  return fTrackContainer.GetNTracks();
1843 }
1844 //_________________________________________________________________________
1845 
1847 {
1848  switch(fITStepStatus)
1849  {
1851  interactionType = "eInteractionWithMedium";
1852  break;
1854  interactionType = "eCollisionBetweenTracks";
1855  break;
1856  default:
1857  interactionType = "eCollisionBetweenTracks";
1858  break;
1859  }
1860 }
1861 
G4ITModelHandler holds for two IT types the corresponding model manager.
void RegisterModel(G4VITStepModel *aModel, const G4double globalTime)
G4bool fUseDefaultTimeSteps
Definition: G4Scheduler.hh:264
void Process()
Definition: G4Scheduler.cc:362
void SetTrackStatus(const G4TrackStatus aTrackStatus)
std::map< double, double > * fpUserTimeSteps
Definition: G4Scheduler.hh:231
Its role is the same as G4StepManager :
G4int GetParentID() const
bool fInteractionStep
Definition: G4Scheduler.hh:247
virtual void DefineTracks()
Definition: G4ITGun.hh:56
void MergeSecondariesWithMainList()
virtual size_t GetNTracks()
Similar to G4ParticleChange, but deal with two tracks rather than one.
void FindUserPreDefinedTimeStep()
Definition: G4Scheduler.cc:943
const G4Track * GetTrackA()
The G4ITModelProcessor will call the two processes defined in G4VITModel.
void EndTracking()
virtual G4bool Notify(G4ApplicationState requestedState)
Definition: G4Scheduler.cc:125
G4TrackManyList * GetSecondariesList()
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void SetTrackingManager(G4ITTrackingManager *trackMan)
void SetParentID(int, int)
Definition: G4IT.hh:250
virtual G4bool IsReplicated() const =0
void PushSecondaries(G4ITStepProcessor *)
bool fReachedUserTimeLimit
Definition: G4Scheduler.hh:236
bool fComputeTimeStep
Definition: G4Scheduler.hh:223
virtual void UserReactionAction(const G4Track &, const G4Track &, const std::vector< G4Track * > &)
Inform about a reaction.
double fTSTimeStep
Definition: G4Scheduler.hh:238
double fTmpEndTime
Definition: G4Scheduler.hh:216
void Clear(Node *)
Define what to do before stepping and after stepping.
const G4ThreeVector & GetPosition() const
G4TrackStatus GetTrackStatus() const
virtual ~G4Scheduler()
Definition: G4Scheduler.cc:222
void CalculateTimeStep(const G4Track *, const G4double)
G4int GetNumberOfSecondaries() const
static void Pop(G4Track *)
void SynchronizeTracks()
Definition: G4Scheduler.cc:471
double fPreviousTimeStep
Definition: G4Scheduler.hh:218
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
virtual const G4String & GetName() const =0
void SetPreviousStepTime(G4double)
G4ITStepStatus fITStepStatus
Definition: G4Scheduler.hh:208
#define G4ThreadLocal
Definition: tls.hh:89
bool IsInf(T value)
Definition: G4Scheduler.cc:111
int fMaxNZeroTimeStepsAllowed
Definition: G4Scheduler.hh:220
const G4Step * GetStep() const
void ComputeTrackReaction()
int G4int
Definition: G4Types.hh:78
double fILTimeStep
Definition: G4Scheduler.hh:240
static G4Scheduler * Instance()
Definition: G4Scheduler.cc:118
MemStat MemoryUsage()
Definition: G4MemStat.cc:55
iterator pop(OBJECT *)
virtual void PushTrack(G4Track *)
const G4String & GetParticleName() const
CLHEP::shared_ptr< std::vector< G4Track * > > G4TrackVectorHandle
void KillTracks()
static G4ITTypeManager * Instance()
Definition: G4ITType.cc:58
static G4ThreadLocal G4Scheduler * fgScheduler
Definition: G4Scheduler.hh:198
void DefinePhysicalStepLength(G4Track *)
G4ITGun * fpGun
Definition: G4Scheduler.hh:262
static const double microsecond
Definition: G4SIunits.hh:141
virtual void EndTracking(G4Track *)
bool fUsePreDefinedTimeSteps
Definition: G4Scheduler.hh:229
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:49
virtual void UserPostTimeStepAction()
bool fInitialized
Definition: G4Scheduler.hh:201
G4GLOB_DLL std::ostream G4cout
double fUserUpperTimeLimit
Definition: G4Scheduler.hh:233
static void DeleteInstance()
const G4String & GetName() const
#define RESET_COLOR
Definition: G4Scheduler.cc:83
void Clear()
Definition: G4Scheduler.cc:237
virtual void UserPreTimeStepAction()
In this method, the user can use : G4Scheduler::Instance()->GetGlobalTime(), to know the current simu...
void _PushTrack(G4Track *track)
bool G4bool
Definition: G4Types.hh:79
G4ITStepManager enables to synchronize in time the step of tracks.
Definition: G4Scheduler.hh:75
G4ITTrackingManager * fpTrackingManager
Definition: G4Scheduler.hh:258
const std::vector< std::vector< G4VITStepModel * > > * GetCurrentModel()
void ExtractILData(G4ITStepProcessor *)
double GetLimitingTimeStep() const
Definition: G4Scheduler.cc:882
void ForceReinitialization()
void ExtractDoItData(G4ITStepProcessor *)
void ClearList()
Definition: G4Scheduler.cc:280
void AddTrackID(G4Track *)
double fEndTime
Definition: G4Scheduler.hh:215
bool MergeNextTimeToMainList(double &time)
G4ITModelProcessor * fpModelProcessor
Definition: G4Scheduler.hh:254
const G4TouchableHandle & GetNextTouchableHandle() const
const G4ParticleDefinition * GetParticleDefinition() const
G4TrackVector * GetSecondaries()
#define DBL_EPSILON
Definition: templates.hh:87
bool SecondaryListsNOTEmpty()
G4SchedulerMessenger * fSteppingMsg
Definition: G4Scheduler.hh:196
double fTimeTolerance
Definition: G4Scheduler.hh:211
void PushToKill(G4Track *track)
void Create()
Definition: G4Scheduler.cc:156
G4Track * GetSecondary(G4int) const
G4int GetTrackID() const
static void DeleteInstance()
DeleteInstance should be used instead of the destructor.
Definition: G4Scheduler.cc:140
G4double GetGlobalTime() const
const G4Track * GetTrack() const
void ReserveRessource()
Definition: G4ITType.cc:77
G4ITModelHandler * fpModelHandler
Definition: G4Scheduler.hh:256
IosFlagSaver(std::ostream &_ios)
Definition: G4Scheduler.cc:93
virtual size_t GetNTracks()
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4VITTimeStepComputer * GetTimeStepper()
G4ITStepProcessor * fpStepProcessor
Definition: G4Scheduler.hh:253
double fTimeStep
Definition: G4Scheduler.hh:226
const G4Track * GetTrackB()
std::vector< G4ITReactionChange * > * GetReactionInfo()
std::vector< G4Track * > fLeadingTracks
Definition: G4Scheduler.hh:245
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4bool IsParameterised() const =0
virtual void Initialize()
G4TrackingInformation * GetTrackingInfo()
Definition: G4IT.hh:152
double fStartTime
Definition: G4Scheduler.hh:214
const G4VProcess * GetProcessDefinedStep() const
G4ITTrackingInteractivity * fpTrackingInteractivity
Definition: G4Scheduler.hh:260
bool fWhyDoYouStop
Definition: G4Scheduler.hh:200
G4UserTimeStepAction * fpUserTimeStepAction
Definition: G4Scheduler.hh:259
std::vector< G4Track * > G4TrackVector
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
void CalculateMinTimeStep()
Definition: G4Scheduler.cc:995
void Stepping(G4Track *, const double &)
void Stop()
G4StepPoint * GetPostStepPoint() const
double fDefinedMinTimeStep
Definition: G4Scheduler.hh:234
G4bool fContinue
Definition: G4Scheduler.hh:263
size_t size() const
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:58
void FindReaction(std::map< G4Track *, G4TrackVectorHandle > *, const double currentStepTime, const double previousStepTime, const bool reachedUserStepTimeLimit)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
void InitializeStepper(const G4double &currentGlobalTime, const G4double &userMinTime)
void Reset()
Definition: G4Scheduler.cc:342
G4bool WereParentsKilled() const
G4Scheduler & operator=(const G4Scheduler &)
virtual void RegisterModel(G4VITStepModel *, double)
Definition: G4Scheduler.cc:292
#define G4endl
Definition: G4ios.hh:61
double fDefaultMinTimeStep
Definition: G4Scheduler.hh:230
std::map< G4Track *, G4TrackVectorHandle > fReactingTracks
Definition: G4Scheduler.hh:244
void SetModelHandler(G4ITModelHandler *)
double fGlobalTime
Definition: G4Scheduler.hh:212
static const double picosecond
Definition: G4SIunits.hh:142
void Start()
G4TrackList * GetMainList(Key)
std::ostream & ios
Definition: G4Scheduler.cc:106
double G4double
Definition: G4Types.hh:76
std::vector< G4Track * > * GetfSecondary()
void DoProcess()
Definition: G4Scheduler.cc:506
void SetInteractivity(G4ITTrackingInteractivity *)
void GetCollisionType(G4String &interactionType)
int fZeroTimeCount
Definition: G4Scheduler.hh:219
Before stepping all tracks G4Scheduler calls all the G4VITModel which may contain a G4VITTimeStepper ...
#define DBL_MAX
Definition: templates.hh:83
void ComputeInteractionLength()
G4TrackStatus
{ Class description:
bool GetTimeStepComputerFlag()
std::ios::fmtflags f
Definition: G4Scheduler.cc:107
bool fComputeReaction
Definition: G4Scheduler.hh:224
void Initialize()
Definition: G4Scheduler.cc:299
void Stepping()
Definition: G4Scheduler.cc:590
G4ITTrackHolder & fTrackContainer
Definition: G4Scheduler.hh:251
G4ApplicationState
void SetInteractivity(G4ITTrackingInteractivity *)
void SetTrackID(const G4int aValue)
double fTmpGlobalTime
Definition: G4Scheduler.hh:213
void ExtractTimeStepperData(G4ITModelProcessor *)
const G4TrackVector * GetSecondary() const
void ResetLeadingTracks()
#define LIGHT_RED
Definition: G4Scheduler.cc:80
virtual void StartProcessing()
void CleanProcessor()
Restaure original state of the modelProcessor.