Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ITModelProcessor.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: G4ITModelProcessor.cc 94010 2015-11-05 10:08:33Z 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 "G4ITModelProcessor.hh"
37 #include "G4VITTimeStepComputer.hh"
38 #include "G4VITReactionProcess.hh"
39 //#include "G4ITTimeStepper.hh"
40 #include "G4ITReaction.hh"
41 #include "G4ITTrackHolder.hh"
42 #include "G4ITTrackingManager.hh"
43 #include "G4UserTimeStepAction.hh"
44 #include "G4UnitsTable.hh"
45 #include <vector>
46 
47 using namespace std;
48 
49 //#define DEBUG_MEM
50 
51 #ifdef DEBUG_MEM
52 #include "G4MemStat.hh"
53 using namespace G4MemStat;
54 #endif
55 
57 {
58  fpTrack = 0;
59  fpModel = 0;
60  fInitialized = false;
61  fpModelManager = 0;
62  fCurrentModel.assign(G4ITType::size(), std::vector<G4VITStepModel*>());
63 
64  for(int i = 0; i < (int) G4ITType::size(); i++)
65  {
66  fCurrentModel[i].assign(G4ITType::size(), 0);
67  }
68  fUserMinTimeStep = -1.;
69  fTSTimeStep = DBL_MAX;
70  fpTrackingManager = 0;
71  fReactionSet = 0;
72  fpTrackContainer = 0;
73  fpModelHandler = 0;
74  fComputeTimeStep = false;
75  fComputeReaction = false;
76 }
77 
79 {
80  //dtor
81  // if(fpModelHandler) delete fpModelHandler; deleted by G4Scheduler
82  fCurrentModel.clear();
83  fReactionInfo.clear();
84 }
85 
86 // Should not be used
88 {
89  //copy ctorr
90  fpTrack = 0;
91  fpModelHandler = 0;
92  fpModel = 0;
93  fInitialized = false;
94  fpModelManager = 0;
95  fUserMinTimeStep = -1.;
96  fTSTimeStep = DBL_MAX;
97  fpTrackingManager = 0;
98  fReactionSet = 0;
99  fpTrackContainer = 0;
100  fComputeTimeStep = false;
101  fComputeReaction = false;
102 }
103 
105 {
106  fpModelHandler->RegisterModel(model, time);
107 }
108 
109 // Should not be used
111 {
112  if(this == &rhs) return *this; // handle self assignment
113  //assignment operator
114  return *this;
115 }
116 //______________________________________________________________________________
118 {
119  fpModelHandler->Initialize();
120  fReactionSet = G4ITReactionSet::Instance();
121  fpTrackContainer = G4ITTrackHolder::Instance();
122  fInitialized = true;
123  fComputeTimeStep = false;
124  fComputeReaction = false;
125  if(fpModelHandler->GetTimeStepComputerFlag()) fComputeTimeStep = true;
126  if(fpModelHandler->GetReactionProcessFlag()) fComputeReaction = true;
127 }
128 
129 //_________________________________________________________________________
130 
132  G4double definedMinTimeStep)
133 {
134 
135 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
136  MemStat mem_first, mem_second, mem_diff;
137 #endif
138 
139 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
140  mem_first = MemoryUsage();
141 #endif
142 
143  fTSTimeStep = DBL_MAX;
144 
145  InitializeStepper(currentGlobalTime, definedMinTimeStep);
146  // TODO
147  // fpMasterModelProcessor -> InitializeStepper(fGlobalTime, fDefinedMinTimeStep) ;
148 
149 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
150  mem_second = MemoryUsage();
151  mem_diff = mem_second-mem_first;
152  G4cout << "\t || MEM || G4Scheduler::CalculateMinTimeStep || After "
153  "computing fpModelProcessor -> InitializeStepper, diff is : "
154  << mem_diff
155  << G4endl;
156 #endif
157 
158 // G4TrackList::iterator fpMainList_i = fpMainList->begin();
159 
160  G4TrackManyList* mainList = fpTrackContainer->GetMainList();
161  G4TrackManyList::iterator it = mainList->begin();
162  G4TrackManyList::iterator end = mainList->end();
163 
164  for (; it != end; ++it)
165  {
166  G4Track * track = *it;
167 
168  if (track == 0)
169  {
170  G4ExceptionDescription exceptionDescription;
171  exceptionDescription << "No track found.";
172  G4Exception("G4Scheduler::CalculateMinStep", "ITScheduler006",
173  FatalErrorInArgument, exceptionDescription);
174  return 0; // makes coverity happy
175  }
176 
177 #ifdef DEBUG
178  G4cout << "*_* " << GetIT(track)->GetName()
179  << " ID: " << track->GetTrackID()
180  << " at time : " << track->GetGlobalTime()
181  << G4endl;
182 #endif
183 
184  G4TrackStatus trackStatus = track->GetTrackStatus();
185  if (trackStatus == fStopAndKill || trackStatus == fStopButAlive)
186  {
187  continue;
188  }
189 
190  // This Extract... was thought for MT mode at the track level
191  //ExtractTimeStepperData(fpModelProcessor) ;
192 
193  CalculateTimeStep(track, definedMinTimeStep);
194 
195  // if MT mode at track level, this command should be displaced
196  ExtractTimeStepperData();
197  }
198 
199 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DETAILED_STEPPING)
200  mem_second = MemoryUsage();
201  mem_diff = mem_second-mem_first;
202  G4cout << "\t || MEM || G4Scheduler::CalculateMinTimeStep || "
203  "After looping on tracks, diff is : " << mem_diff << G4endl;
204 #endif
205  return fTSTimeStep;
206 }
207 
208 //_________________________________________________________________________
209 
211 {
212  if(fpTrack == 0)
213  {
214  CleanProcessor();
215  return;
216  }
217 
218  const std::vector<std::vector<G4VITStepModel*> >* model = GetCurrentModel();
219 
220  for(unsigned i = 0; i < model->size(); ++i)
221  {
222  for(unsigned j = 0; j < (*model)[i].size(); ++j)
223  {
224  G4VITStepModel* mod = (*model)[i][j];
225 
226  if(mod == 0)
227  {
228  continue;
229  }
230 
231  G4VITTimeStepComputer* stepper(mod->GetTimeStepper());
232 
233  G4double sampledMinTimeStep(stepper->GetSampledMinTimeStep());
234  G4TrackVectorHandle reactants(stepper->GetReactants());
235 
236  if(sampledMinTimeStep < fTSTimeStep)
237  {
238  /*
239  // DEBUG SPECIAL CASE
240  if(!reactants)
241  {
242  G4ExceptionDescription exceptionDescription ;
243  CalculateMinTimeStep(); // => at least N (N = nb of tracks) loops exceptionDescription << "No reactants were found by the time stepper.";
244  G4Exception("G4Scheduler::ExtractTimeStepperData","ITScheduler007",
245  FatalErrorInArgument,exceptionDescription);
246  continue ;
247  }
248  */
249 
250  fTSTimeStep = sampledMinTimeStep;
251  //fReactingTracks.clear();
252 
253  fReactionSet->CleanAllReaction();
254  if(bool(reactants))
255  {
256  // G4cout << "*** (1) G4Scheduler::ExtractTimeStepperData insert
257  // reactants for " << GetIT(track)->GetName() << G4endl;
258  // G4cout << "bool(reactants) = " << bool(reactants) << G4endl;
259  // G4cout << reactants->size() << G4endl;
260  // G4cout << GetIT(reactants->operator[](0))->GetName() << G4endl;
261 
262  // fReactingTracks.insert(make_pair(track, reactants));
263  fReactionSet->AddReactions(fTSTimeStep,
264  const_cast<G4Track*>(fpTrack),
265  reactants);
266  stepper->ResetReactants();
267  }
268  }
269  else if(fTSTimeStep == sampledMinTimeStep)
270  {
271  /*
272  // DEBUG SPECIAL CASE
273  if(!reactants)
274  {
275  G4ExceptionDescription exceptionDescription ;
276  exceptionDescription << "No reactants were found by the time stepper.";
277  G4Exception("G4Scheduler::ExtractTimeStepperData","ITScheduler008",
278  FatalErrorInArgument,exceptionDescription);
279  continue ;
280  }
281  */
282  if(bool(reactants))
283  {
284  // G4cout << "*** (2) G4Scheduler::ExtractTimeStepperData insert
285  // reactants for " << GetIT(track)->GetName() << G4endl;
286  // G4cout << "bool(reactants) = " << bool(reactants) << G4endl;
287  // G4cout << "trackA : " << GetIT(track)->GetName()
288  // << " ("<< track->GetTrackID() << ")" << G4endl;
289  // G4cout << reactants->size() << G4endl;
290  // G4cout << GetIT(reactants->operator[](0))->GetName() << G4endl;
291 
292  // fReactingTracks.insert(make_pair(track, reactants));
293 
294  fReactionSet->AddReactions(fTSTimeStep,
295  const_cast<G4Track*>(fpTrack),
296  reactants);
297  stepper->ResetReactants();
298  }
299  }
300  else
301  {
302  if(bool(reactants))
303  {
304  stepper->ResetReactants();
305  }
306  }
307  }
308  }
309 
310  CleanProcessor();
311 }
312 
313 //______________________________________________________________________________
314 
316  G4double userMinTime)
317 {
318  // G4cout << "G4ITModelProcessor::InitializeStepper" << G4endl;
319  if(fpModelHandler == 0)
320  {
321  G4ExceptionDescription exceptionDescription;
322  exceptionDescription
323  << "No G4ITModelHandler was passed to the modelProcessor.";
324  G4Exception("G4ITModelProcessor::InitializeStepper",
325  "ITModelProcessor002",
327  exceptionDescription);
328  }
329  const std::vector<std::vector<G4ITModelManager*> >* modelManager =
330  fpModelHandler->GetAllModelManager();
331 
332  if(modelManager == 0)
333  {
334  G4ExceptionDescription exceptionDescription;
335  exceptionDescription
336  << "No G4ITModelManager was register to G4ITModelHandler.";
337  G4Exception("G4ITModelProcessor::InitializeStepper",
338  "ITModelProcessor003",
340  exceptionDescription);
341  }
342 
343  int nbModels1 = modelManager->size();
344 
345  G4VITTimeStepComputer::SetTimes(currentGlobalTime, userMinTime);
346 
347  // TODO !!!
348  // if( nbModels1 != 1 || (nbModels1 == 1 && !fpModelManager) )
349  {
350  int nbModels2 = -1;
351  G4VITStepModel* model = 0;
352  G4ITModelManager* modman = 0;
353 
354  for(int i = 0; i < nbModels1; i++)
355  {
356  nbModels2 = (*modelManager)[i].size();
357 
358  for(int j = 0; j <= i; j++)
359  {
360  modman = (*modelManager)[i][j];
361 
362  if(modman == 0) continue;
363 
364  model = modman->GetModel(currentGlobalTime);
365  G4VITTimeStepComputer* stepper = model->GetTimeStepper();
366 
367 #if defined (DEBUG_MEM)
368  MemStat mem_first, mem_second, mem_diff;
369 #endif
370 
371 #if defined (DEBUG_MEM)
372  mem_first = MemoryUsage();
373 #endif
374 
375  // stepper->PrepareForAllProcessors() ;
376  stepper->Prepare();
377 
378 #if defined (DEBUG_MEM)
379  mem_second = MemoryUsage();
380  mem_diff = mem_second-mem_first;
381  G4cout << "\t || MEM || G4ITModelProcessor::InitializeStepper || After computing stepper -> Prepare(), diff is : " << mem_diff << G4endl;
382 #endif
383  fCurrentModel[i][j] = model;
384  }
385  }
386 
387  if(nbModels1 == 1 && nbModels2 == 1)
388  {
389  fpModelManager = modman;
390  fpModel = model;
391  }
392  else fpModel = 0;
393  }
394 }
395 
396 //______________________________________________________________________________
398  const G4double userMinTimeStep)
399 {
400  // G4cout << "G4ITModelProcessor::CalculateStep" << G4endl;
401  CleanProcessor();
402  if(track == 0)
403  {
404  G4ExceptionDescription exceptionDescription;
405  exceptionDescription << "No track was passed to the method (track == 0).";
406  G4Exception("G4ITModelProcessor::CalculateStep",
407  "ITModelProcessor004",
409  exceptionDescription);
410  }
411  SetTrack(track);
412  fUserMinTimeStep = userMinTimeStep;
413 
414  DoCalculateStep();
415 }
416 
417 //______________________________________________________________________________
419 {
420  if(fpModel) // ie only one model has been declared and will be used
421  {
422  fpModel->GetTimeStepper()->CalculateStep(*fpTrack, fUserMinTimeStep);
423  }
424  else // ie many models have been declared and will be used
425  {
426  std::vector<G4VITStepModel*>& model = fCurrentModel[GetIT(fpTrack)
427  ->GetITType()];
428 
429  for(int i = 0; i < (int) model.size(); i++)
430  {
431  if(model[i] == 0) continue;
432  model[i]->GetTimeStepper()->CalculateStep(*fpTrack, fUserMinTimeStep);
433  }
434  }
435 }
436 
437 //_________________________________________________________________________
438 
440  G4double fGlobalTime,
441  G4double currentTimeStep,
442  G4double previousTimeStep,
443  G4bool reachedUserTimeLimit,
444  G4double fTimeTolerance,
445  G4UserTimeStepAction* fpUserTimeStepAction,
446  G4int
447 #ifdef G4VERBOSE
448  fVerbose
449 #endif
450  )
451 {
452 // if (fReactingTracks.empty())
453  if(fReactionSet->Empty())
454  {
455  return;
456  }
457 
458  if(fITStepStatus == eCollisionBetweenTracks)
459  // if(fInteractionStep == false)
460  {
461  // TODO
462  FindReaction(fReactionSet,
463  currentTimeStep,
464  previousTimeStep,
465  reachedUserTimeLimit);
466  // TODO
467  // A ne faire uniquement si le temps choisis est celui calculé par le time stepper
468  // Sinon utiliser quelque chose comme : fModelProcessor->FindReaction(&fMainList);
469 
470  std::vector<G4ITReactionChange*> * reactionInfo_v = GetReactionInfo();
471  std::vector<G4ITReactionChange*>::iterator reactionInfo_i = reactionInfo_v
472  ->begin();
473 
474  for(; reactionInfo_i != reactionInfo_v->end(); ++reactionInfo_i)
475  {
476  G4ITReactionChange* changes = (*reactionInfo_i);
477  G4Track* trackA = const_cast<G4Track*>(changes->GetTrackA());
478  G4Track* trackB = const_cast<G4Track*>(changes->GetTrackB());
479 
480  if(trackA == 0 || trackB == 0 || trackA->GetTrackStatus() == fStopAndKill
481  || trackB->GetTrackStatus() == fStopAndKill) continue;
482 
483  G4int nbSecondaries = changes->GetNumberOfSecondaries();
484  const vector<G4Track*>* productsVector = changes->GetfSecondary();
485 
486  if(fpUserTimeStepAction)
487  {
488  fpUserTimeStepAction->UserReactionAction(*trackA,
489  *trackB,
490  productsVector);
491  }
492 
493 #ifdef G4VERBOSE
494  if(fVerbose)
495  {
496  G4cout << "At time : " << setw(7) << G4BestUnit(fGlobalTime, "Time")
497  << " Reaction : " << GetIT(trackA)->GetName() << " ("
498  << trackA->GetTrackID() << ") + " << GetIT(trackB)->GetName() << " ("
499  << trackB->GetTrackID() << ") -> ";
500  }
501 #endif
502 
503  if(nbSecondaries > 0)
504  {
505  for(int i = 0; i < nbSecondaries; ++i)
506  {
507 #ifdef G4VERBOSE
508  if(fVerbose && i != 0) G4cout << " + ";
509 #endif
510 
511  G4Track* secondary = (*productsVector)[i]; //changes->GetSecondary(i);
512  fpTrackContainer->_PushTrack(secondary);
513  GetIT(secondary)->SetParentID(trackA->GetTrackID(),
514  trackB->GetTrackID());
515 
516  if(secondary->GetGlobalTime() - fGlobalTime > fTimeTolerance)
517  {
518  G4ExceptionDescription exceptionDescription;
519  exceptionDescription << "The time of the secondary should not be bigger than the"
520  " current global time."
521  << " This may cause synchronization problem. If the process you"
522  " are using required "
523  << "such feature please contact the developpers." << G4endl<< "The global time in the step manager : "
524  << G4BestUnit(fGlobalTime,"Time")
525  << G4endl
526  << "The global time of the track : "
527  << G4BestUnit(secondary->GetGlobalTime(),"Time")
528  << G4endl;
529 
530  G4Exception("G4Scheduler::ComputeInteractionBetweenTracks",
531  "ITScheduler010",
533  exceptionDescription);
534  }
535 
536 #ifdef G4VERBOSE
537  if(fVerbose)
538  G4cout << GetIT(secondary)->GetName() << " ("
539  << secondary->GetTrackID() << ")";
540 #endif
541  }
542  }
543  else
544  {
545 #ifdef G4VERBOSE
546  if(fVerbose) G4cout << "No product";
547 #endif
548  }
549 #ifdef G4VERBOSE
550  if(fVerbose) G4cout << G4endl;
551 #endif
552  if(trackA->GetTrackID() == 0 || trackB->GetTrackID() == 0)
553  {
554  G4Track* track = 0;
555  if(trackA->GetTrackID() == 0) track = trackA;
556  else track = trackB;
557 
558  G4ExceptionDescription exceptionDescription;
559  exceptionDescription
560  << "The problem was found for the reaction between tracks :"
561  << trackA->GetParticleDefinition()->GetParticleName() << " ("
562  << trackA->GetTrackID() << ") & "
563  << trackB->GetParticleDefinition()->GetParticleName() << " ("
564  << trackB->GetTrackID() << "). \n";
565 
566  if(track->GetStep() == 0)
567  {
568  exceptionDescription << "Also no step was found"
569  << " ie track->GetStep() == 0 \n";
570  }
571 
572  exceptionDescription << "Parent ID of trackA : "
573  << trackA->GetParentID() << "\n";
574  exceptionDescription << "Parent ID of trackB : "
575  << trackB->GetParentID() << "\n";
576 
577  exceptionDescription
578  << "The ID of one of the reaction track was not setup.";
579  G4Exception("G4Scheduler::ComputeInteractionBetweenTracks",
580  "ITScheduler011",
582  exceptionDescription);
583  }
584 
585  if(changes->WereParentsKilled())
586  {
587  // DEBUG
588  // G4cout << "Erasing tracks : "
589  // << "trackA at time : " << G4BestUnit(trackA->GetGlobalTime() , "Time")
590  // << "\t trackB at time : "<< G4BestUnit(trackB->GetGlobalTime(), "Time")
591  // << "\t GlobalTime : " << G4BestUnit(fGlobalTime, "Time")
592  // << G4endl;
593 
594  trackA->SetTrackStatus(fStopAndKill);
595  trackB->SetTrackStatus(fStopAndKill);
596 
597 // G4TrackList::Pop(trackA);
598 // G4TrackList::Pop(trackB);
599  fpTrackingManager->EndTracking(trackA);
600  fpTrackingManager->EndTracking(trackB);
601  }
602 
603  delete changes;
604  }
605 
606  reactionInfo_v->clear(); // never null because point to a non-pointer member of ModelProcessor
607  }
608  // DEBUG
609  //else
610  //{
611  // G4cout << "fInteractionStep == true" << G4endl ;
612  //}
613 
614 // fReactingTracks.clear();
615  fReactionSet->CleanAllReaction();
616 
617  fpTrackContainer->MergeSecondariesWithMainList();
618  fpTrackContainer->KillTracks();
619 }
620 
621 //______________________________________________________________________________
623  const double currentStepTime,
624  const double previousStepTime,
625  const bool reachedUserStepTimeLimit)
626 {
627  // DEBUG
628  // G4cout << "G4ITReactionManager::FindReaction" << G4endl;
629  //if (tracks == 0) return;
630  if(reactionSet == 0) return;
631  if(fpModelHandler->GetAllModelManager()->empty()) return;
632 
633  G4ITReactionPerTrackMap& reactionPerTrackMap = reactionSet->GetReactionMap();
634 
635  std::map<G4Track*, G4ITReactionPerTrackPtr, compTrackPerID>::iterator tracks_i =
636  reactionPerTrackMap.begin();
637 
638  //std::map<G4Track*, G4TrackVectorHandle, compTrackPerID>::iterator tracks_i = tracks->begin();
639 
640  // G4cout << "G4ITModelProcessor::FindReaction at step :" << G4ITTimeStepper::Instance()->GetNbSteps() << G4endl;
641 
642  // for (tracks_i = tracks->begin(); tracks_i != tracks->end(); tracks_i++)
643  for(tracks_i = reactionPerTrackMap.begin();
644  tracks_i != reactionPerTrackMap.end();
645  tracks_i = reactionPerTrackMap.begin())
646  {
647  //G4cout << "here" << G4endl;
648  G4Track* trackA = tracks_i->first;
649  if(trackA->GetTrackStatus() == fStopAndKill)
650  {
651  //G4cout << "continue 1" << G4endl;
652  continue;
653  }
654 
655  G4ITReactionPerTrackPtr reactionPerTrack = tracks_i->second;
656  G4ITReactionList& reactionList = reactionPerTrack->GetReactionList();
657 
658  G4IT* ITA = GetIT(trackA);
659  G4ITType ITypeA = ITA->GetITType();
660 
661  const std::vector<G4VITStepModel*> model = fCurrentModel[ITypeA];
662 
663  G4Track* trackB = 0;
664  G4ITType ITypeB(-1);
665  G4VITReactionProcess* process = 0;
666  G4ITReactionChange* changes = 0;
667 
668  assert(reactionList.begin() != reactionList.end());
669 
670  for(G4ITReactionList::iterator it = reactionList.begin();
671  it != reactionList.end(); it = reactionList.begin())
672  {
673  G4ITReactionPtr reaction(*it);
674  trackB = reaction->GetReactant(trackA);
675  if(trackB->GetTrackStatus() == fStopAndKill)
676  {
677  //G4cout << "continue 2" << G4endl;
678  continue;
679  }
680 
681  if(trackB == trackA)
682  {
683  G4ExceptionDescription exceptionDescription;
684  exceptionDescription
685  << "The IT reaction process sent back a reaction between trackA and trackB. ";
686  exceptionDescription << "The problem is trackA == trackB";
687  G4Exception("G4ITModelProcessor::FindReaction",
688  "ITModelProcessor005",
690  exceptionDescription);
691  }
692 
693  G4IT* ITB = GetIT(trackB);
694  G4ITType ITypeBtmp = ITB->GetITType();
695 
696  if(ITypeB != ITypeBtmp)
697  {
698  ITypeB = ITypeBtmp;
699 
700  if(model[ITypeB]) process = model[ITypeB]->GetReactionProcess();
701  }
702 
703  reactionSet->SelectThisReaction(reaction);
704 
705  if(process && process->TestReactibility(*trackA,
706  *trackB,
707  currentStepTime,
708  previousStepTime,
709  reachedUserStepTimeLimit))
710  {
711  changes = process->MakeReaction(*trackA, *trackB);
712  }
713 
714  if(changes)
715  {
716  fReactionInfo.push_back(changes);
717 
718  // G4cout << "pushing reaction for trackA (" << trackA->GetTrackID() << ") and trackB ("
719  // << trackB->GetTrackID() << ")" << G4endl;
720  //
721  // G4cout << "nb of secondaries : " << changes->GetNumberOfSecondaries() << G4endl;
722  // G4cout << "with track 0 = " << changes->GetSecondary(0) << G4endl;
723 
724  process->ResetChanges();
725  changes = 0;
726 
727  break;
728  }
729  }
730  }
731 
732  //assert(G4ITReaction::gAll->empty() == true);
733 }
void SetTrackStatus(const G4TrackStatus aTrackStatus)
virtual G4bool TestReactibility(const G4Track &, const G4Track &, const double, const double, bool)=0
G4int GetParentID() const
Definition: G4IT.hh:88
const G4Track * GetTrackA()
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void SetParentID(int, int)
Definition: G4IT.hh:229
G4shared_ptr< std::vector< G4Track * > > G4TrackVectorHandle
Definition: G4ITReaction.hh:43
G4TrackStatus GetTrackStatus() const
void CalculateTimeStep(const G4Track *, const G4double)
G4ITStepStatus
G4shared_ptr< G4ITReaction > G4ITReactionPtr
Definition: G4ITReaction.hh:59
G4int GetNumberOfSecondaries() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
virtual const G4String & GetName() const =0
const G4Step * GetStep() const
int G4int
Definition: G4Types.hh:78
MemStat MemoryUsage()
Definition: G4MemStat.cc:55
std::map< G4Track *, G4ITReactionPerTrackPtr, compTrackPerID > G4ITReactionPerTrackMap
Definition: G4ITReaction.hh:66
const G4String & GetParticleName() const
G4ITReactionPerTrackMap & GetReactionMap()
static void SetTimes(const G4double &, const G4double &)
G4double CalculateMinTimeStep(G4double currentGlobalTime, G4double definedMinTimeStep)
static G4ITReactionSet * Instance()
void SelectThisReaction(G4ITReactionPtr reaction)
void FindReaction(G4ITReactionSet *reactionSet, const double currentStepTime, const double previousStepTime, const bool reachedUserStepTimeLimit)
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:49
std::list< G4ITReactionPtr > G4ITReactionList
Definition: G4ITReaction.hh:63
G4GLOB_DLL std::ostream G4cout
static size_t size()
Definition: G4ITType.cc:46
bool G4bool
Definition: G4Types.hh:79
typedef int(XMLCALL *XML_NotStandaloneHandler)(void *userData)
const G4ParticleDefinition * GetParticleDefinition() const
G4int GetTrackID() const
G4double GetGlobalTime() const
G4VITTimeStepComputer * GetTimeStepper()
const G4Track * GetTrackB()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4ITModelProcessor & operator=(const G4ITModelProcessor &other)
G4shared_ptr< G4ITReactionPerTrack > G4ITReactionPerTrackPtr
Definition: G4ITReaction.hh:61
void RegisterModel(double time, G4VITStepModel *)
virtual void UserReactionAction(const G4Track &, const G4Track &, const std::vector< G4Track * > *)
virtual const G4ITType GetITType() const =0
void ComputeTrackReaction(G4ITStepStatus fITStepStatus, G4double fGlobalTime, G4double currentTimeStep, G4double previousTimeStep, G4bool reachedUserTimeLimit, G4double fTimeTolerance, G4UserTimeStepAction *fpUserTimeStepAction, G4int fVerbose)
void InitializeStepper(G4double currentGlobalTime, G4double userMinTime)
G4bool WereParentsKilled() const
static G4ITTrackHolder * Instance()
#define G4endl
Definition: G4ios.hh:61
G4VITStepModel * GetModel(const G4double globalTime)
double G4double
Definition: G4Types.hh:76
std::vector< G4Track * > * GetfSecondary()
virtual G4ITReactionChange * MakeReaction(const G4Track &, const G4Track &)=0
#define DBL_MAX
Definition: templates.hh:83
G4TrackStatus
const XML_Char XML_Content * model
Definition: expat.h:151