Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ITStepProcessor2.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: G4ITStepProcessor2.cc 94218 2015-11-09 08:24:48Z 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 "G4ITStepProcessor.hh"
37 #include "G4LossTableManager.hh"
38 #include "G4EnergyLossTables.hh"
39 #include "G4ProductionCuts.hh"
40 #include "G4ProductionCutsTable.hh"
41 #include "G4VITProcess.hh"
42 #include "G4TrackingInformation.hh"
43 #include "G4IT.hh"
44 #include "G4ITTrackingManager.hh"
45 #include "G4ITTransportation.hh"
46 
47 #include "G4ITNavigator.hh" // Include from 'geometry'
48 
49 #include "G4ITSteppingVerbose.hh"
50 #include "G4VITSteppingVerbose.hh"
51 
52 #include "G4ITTrackHolder.hh"
53 #include "G4ITReaction.hh"
54 
55 #include "G4TrackingInformation.hh"
56 
57 //#define DEBUG_MEM 1
58 
59 #ifdef DEBUG_MEM
60 #include "G4MemStat.hh"
61 using namespace G4MemStat;
62 using G4MemStat::MemStat;
63 #endif
64 
66 {
67  // Now Store the secondaries from ParticleChange to SecondaryList
68  G4Track* tempSecondaryTrack;
69 
70  for(G4int DSecLoop = 0; DSecLoop < fpParticleChange->GetNumberOfSecondaries();
71  DSecLoop++)
72  {
73  tempSecondaryTrack = fpParticleChange->GetSecondary(DSecLoop);
74 
75  if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
76  {
77  ApplyProductionCut(tempSecondaryTrack);
78  }
79 
80  // Set parentID
81  tempSecondaryTrack->SetParentID(fpTrack->GetTrackID());
82 
83  // Set the process pointer which created this track
84  tempSecondaryTrack->SetCreatorProcess(fpCurrentProcess);
85 
86  // If this 2ndry particle has 'zero' kinetic energy, make sure
87  // it invokes a rest process at the beginning of the tracking
88  if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN)
89  {
90  G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
91  if (pm->GetAtRestProcessVector()->entries()>0)
92  {
93  tempSecondaryTrack->SetTrackStatus( fStopButAlive );
94  fpSecondary->push_back( tempSecondaryTrack );
95  fN2ndariesAtRestDoIt++;
96  }
97  else
98  {
99  delete tempSecondaryTrack;
100  }
101  }
102  else
103  {
104  fpSecondary->push_back( tempSecondaryTrack );
105  counter++;
106  }
107  } //end of loop on secondary
108 }
109 
110 //_________________________________________________________________________
111 
112 void G4ITStepProcessor::DoIt(double timeStep)
113 
114 // Call the process having the min step length or just propagate the track on the given time step
115 
116 // If the track is "leading the step" (ie one of its process has been selected
117 // as the one having the minimum time step over all tracks and processes),
118 // it will undergo its selected processes. Otherwise, it will just propagate the track
119 // on the given time step.
120 
121 {
122  if(fpVerbose) fpVerbose->DoItStarted();
123 
124  G4TrackManyList* mainList = fpTrackContainer->GetMainList();
125  G4TrackManyList::iterator it = mainList->end();
126  it--;
127  size_t initialSize = mainList->size();
128 
129  // G4cout << "initialSize = " << initialSize << G4endl;
130 
131  for(size_t i = 0 ; i < initialSize ; ++i)
132  {
133 
134  // G4cout << "i = " << i << G4endl;
135 
136  G4Track* track = *it;
137  if (!track)
138  {
139  G4ExceptionDescription exceptionDescription;
140  exceptionDescription << "No track was pop back the main track list.";
141  G4Exception("G4ITStepProcessor::DoIt", "NO_TRACK",
142  FatalException, exceptionDescription);
143  }
144  // G4TrackManyList::iterator next_it (it);
145  // next_it--;
146  // it = next_it;
147 
148  it--;
149  // Must be called before EndTracking(track)
150  // Otherwise the iterator will point to the list of killed tracks
151 
152  if(track->GetTrackStatus() == fStopAndKill)
153  {
154  fpTrackingManager->EndTracking(track);
155 // G4cout << GetIT(track)->GetName() << G4endl;
156 // G4cout << " ************************ CONTINUE ********************" << G4endl;
157  continue;
158  }
159 
160 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DOIT)
161  MemStat mem_first, mem_second, mem_diff;
162  mem_first = MemoryUsage();
163 #endif
164 
165  Stepping(track, timeStep);
166 
167 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DOIT)
168  MemStat mem_intermediaire = MemoryUsage();
169  mem_diff = mem_intermediaire-mem_first;
170  G4cout << "\t\t >> || MEM || In DoIT with track "
171  << track->GetTrackID() << ", diff is : " << mem_diff << G4endl;
172 #endif
173 
174  ExtractDoItData();
175 
176 #if defined (DEBUG_MEM) && defined (DEBUG_MEM_DOIT)
177  mem_second = MemoryUsage();
178  mem_diff = mem_second-mem_first;
179  G4cout << "\t >> || MEM || In DoIT with track "
180  << track->GetTrackID()
181  << ", diff is : " << mem_diff << G4endl;
182 #endif
183  }
184 
185 
186  fpTrackContainer->MergeSecondariesWithMainList();
187  fpTrackContainer->KillTracks(); // (18-06-15 : MK) Remove it ?
188  fLeadingTracks.Reset();
189 }
190 
191 //_________________________________________________________________________
192 
194 {
195  if (!fpTrack)
196  {
197  CleanProcessor();
198  return;
199  }
200 
201  G4TrackStatus status = fpTrack->GetTrackStatus();
202 
203  switch (status)
204  {
205  case fAlive:
206  case fStopButAlive:
207  case fSuspend:
209  default:
210  PushSecondaries();
211  break;
212 
213  case fStopAndKill:
215  PushSecondaries();
216 // G4TrackList::Pop(fpTrack);
217  fpTrackingManager->EndTracking(fpTrack);
218 // fTrackContainer->PushToKill(fpTrack);
219  break;
220 
223  if (fpSecondary)
224  {
225  for (size_t i = 0; i < fpSecondary->size(); ++i)
226  {
227  delete (*fpSecondary)[i];
228  }
229  fpSecondary->clear();
230  }
231 // G4TrackList::Pop(fpTrack);
232  fpTrackingManager->EndTracking(fpTrack);
233 // fTrackContainer->PushToKill(fpTrack);
234  break;
235  }
236 
237  CleanProcessor();
238 }
239 
240 //_________________________________________________________________________
241 
243 {
244  if (!fpSecondary || fpSecondary->empty())
245  {
246  // DEBUG
247  // G4cout << "NO SECONDARIES !!! " << G4endl;
248  return;
249  }
250 
251  // DEBUG
252  // G4cout << "There are secondaries : "<< secondaries -> size() << G4endl ;
253 
254  G4TrackVector::iterator secondaries_i = fpSecondary->begin();
255 
256  for (; secondaries_i != fpSecondary->end(); ++secondaries_i)
257  {
258  G4Track* secondary = *secondaries_i;
259  fpTrackContainer->_PushTrack(secondary);
260  }
261 }
262 
263 //______________________________________________________________________________
264 
265 void G4ITStepProcessor::Stepping(G4Track* track, const double & timeStep)
266 {
267 
268 #ifdef DEBUG_MEM
269  MemStat mem_first, mem_second, mem_diff;
270 #endif
271 
272 #ifdef DEBUG_MEM
273  mem_first = MemoryUsage();
274 #endif
275 
276  CleanProcessor();
277 
278 #ifdef DEBUG_MEM
279  MemStat mem_intermediaire = MemoryUsage();
280  mem_diff = mem_intermediaire-mem_first;
281  G4cout << "\t\t\t >> || MEM || After CleanProcessor " << track->GetTrackID() << ", diff is : " << mem_diff << G4endl;
282 #endif
283 
284  if(track == 0) return; // maybe put an exception here
285  fTimeStep = timeStep;
286  SetTrack(track);
287  DoStepping();
288 }
289 //______________________________________________________________________________
290 
291 // ************************************************************************
292 // Stepping
293 // ************************************************************************
295 {
296  SetupMembers();
297 
298 #ifdef DEBUG_MEM
299  MemStat mem_first, mem_second, mem_diff;
300 #endif
301 
302 #ifdef DEBUG_MEM
303  mem_first = MemoryUsage();
304 #endif
305 
306 #ifdef G4VERBOSE
307  if(fpVerbose) fpVerbose->PreStepVerbose(fpTrack);
308 #endif
309 
310  if(!fpProcessInfo)
311  {
312  G4ExceptionDescription exceptionDescription;
313  exceptionDescription << "No process info found for particle :"
314  << fpTrack->GetDefinition()->GetParticleName();
315  G4Exception("G4ITStepProcessor::DoStepping",
316  "ITStepProcessor0012",
318  exceptionDescription);
319  return;
320  }
321 // else if(fpTrack->GetTrackStatus() == fStopAndKill)
322 // {
323 // fpState->fStepStatus = fUndefined;
324 // return;
325 // }
326 
327  if(fpProcessInfo->MAXofPostStepLoops == 0 &&
328  fpProcessInfo->MAXofAlongStepLoops == 0
329  && fpProcessInfo->MAXofAtRestLoops == 0)
330  {
331  G4ExceptionDescription exceptionDescription;
332  exceptionDescription << "No process was found for particle :"
333  << fpTrack->GetDefinition()->GetParticleName();
334  G4Exception("G4ITStepProcessor::DoStepping",
335  "ITStepProcessorNoProcess",
336  JustWarning,
337  exceptionDescription);
338 
339  fpTrack->SetTrackStatus(fStopAndKill);
340  fpState->fStepStatus = fUndefined;
341  return;
342  }
343 
344  //--------
345  // Prelude
346  //--------
347 #ifdef G4VERBOSE
348  // !!!!! Verbose
349  if(fpVerbose) fpVerbose->NewStep();
350 #endif
351 
352  //---------------------------------
353  // AtRestStep, AlongStep and PostStep Processes
354  //---------------------------------
355 
356  fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()->GetNavigatorState());
357 // fpNavigator->ResetHierarchyAndLocate( fpTrack->GetPosition(),
358 // fpTrack->GetMomentumDirection(),
359 // *((G4TouchableHistory*)fpTrack->GetTouchableHandle()()) );
360 // fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()->GetNavigatorState());
361  // We reset the navigator state before checking for AtRest
362  // in case a AtRest processe would use a navigator info
363 
364 #ifdef DEBUG_MEM
365  MemStat mem_intermediaire = MemoryUsage();
366  mem_diff = mem_intermediaire-mem_first;
367  G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After dealing with navigator with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
368 #endif
369 
370  if(fpTrack->GetTrackStatus() == fStopButAlive)
371  {
372  if(fpProcessInfo->MAXofAtRestLoops > 0 && fpProcessInfo->fpAtRestDoItVector
373  != 0) // second condition to make coverity happy
374  {
375  //-----------------
376  // AtRestStepDoIt
377  //-----------------
378  InvokeAtRestDoItProcs();
379  fpState->fStepStatus = fAtRestDoItProc;
380  fpStep->GetPostStepPoint()->SetStepStatus(fpState->fStepStatus);
381 
382 #ifdef G4VERBOSE
383  // !!!!! Verbose
384  if(fpVerbose) fpVerbose->AtRestDoItInvoked();
385 #endif
386 
387  }
388  // Make sure the track is killed
389  // fpTrack->SetTrackStatus(fStopAndKill);
390  }
391  else // if(fTimeStep > 0.) // Bye, because PostStepIL can return 0 => time =0
392  {
393  if(fpITrack == 0)
394  {
395  G4ExceptionDescription exceptionDescription;
396  exceptionDescription << " !!! TrackID : " << fpTrack->GetTrackID()
397  << G4endl<< " !!! Track status : "<< fpTrack->GetTrackStatus() << G4endl
398  << " !!! Particle Name : "<< fpTrack -> GetDefinition() -> GetParticleName() << G4endl
399  << "No G4ITStepProcessor::fpITrack found" << G4endl;
400 
401  G4Exception("G4ITStepProcessor::DoStepping",
402  "ITStepProcessor0013",
404  exceptionDescription);
405  return; // to make coverity happy
406  }
407 
408  if(fpITrack->GetTrackingInfo()->IsLeadingStep() == false)
409  {
410  // In case the track has NOT the minimum step length
411  // Given the final step time, the transportation
412  // will compute the final position of the particle
413  fpState->fStepStatus = fPostStepDoItProc;
414  fpStep->GetPostStepPoint()->SetProcessDefinedStep(fpTransportation);
415  FindTransportationStep();
416  }
417 
418 #ifdef DEBUG_MEM
419  mem_intermediaire = MemoryUsage();
420  mem_diff = mem_intermediaire-mem_first;
421  G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After FindTransportationStep() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
422 #endif
423 
424  // Store the Step length (geometrical length) to G4Step and G4Track
425  fpTrack->SetStepLength(fpState->fPhysicalStep);
426  fpStep->SetStepLength(fpState->fPhysicalStep);
427 
428  G4double GeomStepLength = fpState->fPhysicalStep;
429 
430  // Store StepStatus to PostStepPoint
431  fpStep->GetPostStepPoint()->SetStepStatus(fpState->fStepStatus);
432 
433  // Invoke AlongStepDoIt
434  InvokeAlongStepDoItProcs();
435 
436 #ifdef DEBUG_MEM
437  mem_intermediaire = MemoryUsage();
438  mem_diff = mem_intermediaire-mem_first;
439  G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After InvokeAlongStepDoItProcs() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
440 #endif
441 
442 #ifdef G4VERBOSE
443  // !!!!! Verbose
444  if(fpVerbose) fpVerbose->AlongStepDoItAllDone();
445 #endif
446 
447  // Update track by taking into account all changes by AlongStepDoIt
448  // fpStep->UpdateTrack(); // done in InvokeAlongStepDoItProcs
449 
450  // Update safety after invocation of all AlongStepDoIts
451  fpState->fEndpointSafOrigin = fpPostStepPoint->GetPosition();
452 
453  fpState->fEndpointSafety =
454  std::max(fpState->fProposedSafety - GeomStepLength, kCarTolerance);
455 
456  fpStep->GetPostStepPoint()->SetSafety(fpState->fEndpointSafety);
457 
458  if(GetIT(fpTrack)->GetTrackingInfo()->IsLeadingStep())
459  {
460  // Invoke PostStepDoIt including G4ITTransportation::PSDI
461  InvokePostStepDoItProcs();
462 
463 #ifdef DEBUG_MEM
464  mem_intermediaire = MemoryUsage();
465  mem_diff = mem_intermediaire-mem_first;
466  G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After InvokePostStepDoItProcs() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
467 #endif
468 #ifdef G4VERBOSE
469  // !!!!! Verbose
470  if(fpVerbose) fpVerbose->StepInfoForLeadingTrack();
471 #endif
472  }
473  else
474  {
475  // Only invoke transportation and all other forced processes
476  InvokeTransportationProc();
477  fpStep->GetPostStepPoint()->SetProcessDefinedStep(fpTransportation);
478 
479 #ifdef DEBUG_MEM
480  mem_intermediaire = MemoryUsage();
481  mem_diff = mem_intermediaire-mem_first;
482  G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After InvokeTransportationProc() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
483 #endif
484  }
485 
486 #ifdef G4VERBOSE
487  // !!!!! Verbose
488  if(fpVerbose) fpVerbose->PostStepDoItAllDone();
489 #endif
490  }
491 
492  fpNavigator->ResetNavigatorState();
493 
494 #ifdef DEBUG_MEM
495  mem_intermediaire = MemoryUsage();
496  mem_diff = mem_intermediaire-mem_first;
497  G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After fpNavigator->SetNavigatorState with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
498 #endif
499 
500  //-------
501  // Finale
502  //-------
503 
504  // Update 'TrackLength' and remeber the Step length of the current Step
505  fpTrack->AddTrackLength(fpStep->GetStepLength());
506  fpTrack->IncrementCurrentStepNumber();
507 
508 //#ifdef G4VERBOSE
509 // // !!!!! Verbose
510 // if(fpVerbose) fpVerbose->StepInfo();
511 //#endif
512 
513 #ifdef G4VERBOSE
514  if(fpVerbose) fpVerbose->PostStepVerbose(fpTrack);
515 #endif
516 
517 // G4cout << " G4ITStepProcessor::DoStepping -- " <<fpTrack->GetTrackID() << " tps = " << fpTrack->GetGlobalTime() << G4endl;
518 
519  // Send G4Step information to Hit/Dig if the volume is sensitive
520  /***
521  fpCurrentVolume = fpStep->GetPreStepPoint()->GetPhysicalVolume();
522  StepControlFlag = fpStep->GetControlFlag();
523 
524  if( fpCurrentVolume != 0 && StepControlFlag != AvoidHitInvocation)
525  {
526  fpSensitive = fpStep->GetPreStepPoint()->
527  GetSensitiveDetector();
528  if( fpSensitive != 0 )
529  {
530  fpSensitive->Hit(fpStep);
531  }
532  }
533 
534  User intervention process.
535  if( fpUserSteppingAction != 0 )
536  {
537  fpUserSteppingAction->UserSteppingAction(fpStep);
538  }
539  G4UserSteppingAction* regionalAction
540  = fpStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetRegion()
541  ->GetRegionalSteppingAction();
542  if( regionalAction ) regionalAction->UserSteppingAction(fpStep);
543  ***/
544  fpTrackingManager->AppendStep(fpTrack, fpStep);
545  // Stepping process finish. Return the value of the StepStatus.
546 
547 #ifdef DEBUG_MEM
548  MemStat mem_intermediaire = MemoryUsage();
549  mem_diff = mem_intermediaire-mem_first;
550  G4cout << "\t\t\t >> || MEM || End of DoStepping() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
551 #endif
552 
553  // return fpState->fStepStatus;
554 }
555 
556 //______________________________________________________________________________
557 
558 // ************************************************************************
559 // AtRestDoIt
560 // ************************************************************************
561 
563 {
564  fpStep->SetStepLength(0.); //the particle has stopped
565  fpTrack->SetStepLength(0.);
566 
567  G4SelectedAtRestDoItVector& selectedAtRestDoItVector =
568  fpState->fSelectedAtRestDoItVector;
569 
570  // invoke selected process
571  for(size_t np = 0; np < fpProcessInfo->MAXofAtRestLoops; np++)
572  {
573  //
574  // Note: DoItVector has inverse order against GetPhysIntVector
575  // and SelectedAtRestDoItVector.
576  //
577  if(selectedAtRestDoItVector[fpProcessInfo->MAXofAtRestLoops - np - 1] != InActivated)
578  {
579  fpCurrentProcess =
580  (G4VITProcess*) (*fpProcessInfo->fpAtRestDoItVector)[np];
581 
582 // G4cout << " Invoke : "
583 // << fpCurrentProcess->GetProcessName()
584 // << G4endl;
585 
586  // if(fpVerbose)
587  // {
588  // fpVerbose->AtRestDoItOneByOne();
589  // }
590 
591  fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess
592  ->GetProcessID()));
593  fpParticleChange = fpCurrentProcess->AtRestDoIt(*fpTrack, *fpStep);
594  fpCurrentProcess->ResetProcessState();
595 
596  // Set the current process as a process which defined this Step length
597  fpStep->GetPostStepPoint()->SetProcessDefinedStep(fpCurrentProcess);
598 
599  // Update Step
600  fpParticleChange->UpdateStepForAtRest(fpStep);
601 
602  // Now Store the secondaries from ParticleChange to SecondaryList
603  DealWithSecondaries(fN2ndariesAtRestDoIt);
604 
605  // Set the track status according to what the process defined
606  // if kinetic energy >0, otherwise set fStopButAlive
607  fpTrack->SetTrackStatus(fpParticleChange->GetTrackStatus());
608 
609  // clear ParticleChange
610  fpParticleChange->Clear();
611 
612  } //if(fSelectedAtRestDoItVector[np] != InActivated){
613  } //for(size_t np=0; np < MAXofAtRestLoops; np++){
614  fpStep->UpdateTrack();
615 
616  // Modification par rapport au transport standard :
617  // fStopAndKill doit etre propose par le modele
618  // sinon d autres processus AtRest seront appeles
619  // au pas suivant
620  // fpTrack->SetTrackStatus(fStopAndKill);
621 }
622 
623 //______________________________________________________________________________
624 
625 // ************************************************************************
626 // AlongStepDoIt
627 // ************************************************************************
628 
630 {
631 
632 #ifdef DEBUG_MEM
633  MemStat mem_first, mem_second, mem_diff;
634 #endif
635 
636 #ifdef DEBUG_MEM
637  mem_first = MemoryUsage();
638 #endif
639 
640  // If the current Step is defined by a 'ExclusivelyForced'
641  // PostStepDoIt, then don't invoke any AlongStepDoIt
642  if(fpState->fStepStatus == fExclusivelyForcedProc)
643  {
644  return; // Take note 'return' at here !!!
645  }
646 
647  // Invoke the all active continuous processes
648  for(size_t ci = 0; ci < fpProcessInfo->MAXofAlongStepLoops; ci++)
649  {
650  fpCurrentProcess =
651  (G4VITProcess*) (*fpProcessInfo->fpAlongStepDoItVector)[ci];
652  if(fpCurrentProcess == 0) continue;
653  // NULL means the process is inactivated by a user on fly.
654 
655  fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess
656  ->GetProcessID()));
657  fpParticleChange = fpCurrentProcess->AlongStepDoIt(*fpTrack, *fpStep);
658 
659 #ifdef DEBUG_MEM
660  MemStat mem_intermediaire = MemoryUsage();
661  mem_diff = mem_intermediaire-mem_first;
662  G4cout << "\t\t\t >> || MEM || After calling AlongStepDoIt for " << fpCurrentProcess->GetProcessName() << " and track "<< fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
663 #endif
664 
665 // fpCurrentProcess->SetProcessState(0);
666  fpCurrentProcess->ResetProcessState();
667  // Update the PostStepPoint of Step according to ParticleChange
668 
669  fpParticleChange->UpdateStepForAlongStep(fpStep);
670 
671 #ifdef G4VERBOSE
672  // !!!!! Verbose
673  if(fpVerbose) fpVerbose->AlongStepDoItOneByOne();
674 #endif
675 
676  // Now Store the secondaries from ParticleChange to SecondaryList
677  DealWithSecondaries(fN2ndariesAlongStepDoIt);
678 
679  // Set the track status according to what the process defined
680  // if kinetic energy >0, otherwise set fStopButAlive
681  fpTrack->SetTrackStatus(fpParticleChange->GetTrackStatus());
682 
683  // clear ParticleChange
684  fpParticleChange->Clear();
685  }
686 
687 #ifdef DEBUG_MEM
688  MemStat mem_intermediaire = MemoryUsage();
689  mem_diff = mem_intermediaire-mem_first;
690  G4cout << "\t\t\t >> || MEM || After looping on processes with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
691 #endif
692 
693  fpStep->UpdateTrack();
694 
695  G4TrackStatus fNewStatus = fpTrack->GetTrackStatus();
696 
697  if(fNewStatus == fAlive && fpTrack->GetKineticEnergy() <= DBL_MIN)
698  {
699  // G4cout << "G4ITStepProcessor::InvokeAlongStepDoItProcs : Track will be killed" << G4endl;
700  if(fpProcessInfo->MAXofAtRestLoops>0) fNewStatus = fStopButAlive;
701  else fNewStatus = fStopAndKill;
702  fpTrack->SetTrackStatus( fNewStatus );
703  }
704 
705 }
706 
707 //______________________________________________________________________________
708 
709 // ************************************************************************
710 // PostStepDoIt
711 // ************************************************************************
712 
714 {
715  size_t _MAXofPostStepLoops = fpProcessInfo->MAXofPostStepLoops;
716  G4SelectedPostStepDoItVector& selectedPostStepDoItVector = fpState
717  ->fSelectedPostStepDoItVector;
718  G4StepStatus& stepStatus = fpState->fStepStatus;
719 
720  // Invoke the specified discrete processes
721  for(size_t np = 0; np < _MAXofPostStepLoops; np++)
722  {
723  //
724  // Note: DoItVector has inverse order against GetPhysIntVector
725  // and SelectedPostStepDoItVector.
726  //
727  G4int Cond = selectedPostStepDoItVector[_MAXofPostStepLoops
728  - np - 1];
729  if(Cond != InActivated)
730  {
731  if(((Cond == NotForced) && (stepStatus == fPostStepDoItProc)) ||
732  ((Cond == Forced) && (stepStatus != fExclusivelyForcedProc))
733  ||
734  // ((Cond == Conditionally) && (stepStatus == fAlongStepDoItProc)) ||
735  ((Cond == ExclusivelyForced) && (stepStatus == fExclusivelyForcedProc))
736  || ((Cond == StronglyForced)))
737  {
738 
739  InvokePSDIP(np);
740  }
741  } //if(*fSelectedPostStepDoItVector(np)........
742 
743  // Exit from PostStepLoop if the track has been killed,
744  // but extra treatment for processes with Strongly Forced flag
745  if(fpTrack->GetTrackStatus() == fStopAndKill)
746  {
747  for(size_t np1 = np + 1; np1 < _MAXofPostStepLoops; np1++)
748  {
749  G4int Cond2 = selectedPostStepDoItVector[_MAXofPostStepLoops
750  - np1 - 1];
751  if(Cond2 == StronglyForced)
752  {
753  InvokePSDIP(np1);
754  }
755  }
756  break;
757  }
758  } //for(size_t np=0; np < MAXofPostStepLoops; np++){
759 }
760 
761 //______________________________________________________________________________
762 
764 {
765  fpCurrentProcess = (G4VITProcess*) (*fpProcessInfo->fpPostStepDoItVector)[np];
766 
767  fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess
768  ->GetProcessID()));
769  fpParticleChange = fpCurrentProcess->PostStepDoIt(*fpTrack, *fpStep);
770 // fpCurrentProcess->SetProcessState(0);
771  fpCurrentProcess->ResetProcessState();
772 
773  // Update PostStepPoint of Step according to ParticleChange
774  fpParticleChange->UpdateStepForPostStep(fpStep);
775 
776 #ifdef G4VERBOSE
777  // !!!!! Verbose
778  if(fpVerbose) fpVerbose->PostStepDoItOneByOne();
779 #endif
780 
781  // Update G4Track according to ParticleChange after each PostStepDoIt
782  fpStep->UpdateTrack();
783 
784  // Update safety after each invocation of PostStepDoIts
785  fpStep->GetPostStepPoint()->SetSafety(CalculateSafety());
786 
787  // Now Store the secondaries from ParticleChange to SecondaryList
788  DealWithSecondaries(fN2ndariesPostStepDoIt);
789 
790  // Set the track status according to what the process defined
791  fpTrack->SetTrackStatus(fpParticleChange->GetTrackStatus());
792 
793  // clear ParticleChange
794  fpParticleChange->Clear();
795 }
796 
797 //______________________________________________________________________________
798 
799 // ************************************************************************
800 // Transport on a given time
801 // ************************************************************************
802 
804 {
805  double physicalStep(0.);
806 
807  fpTransportation = fpProcessInfo->fpTransportation;
808  // dynamic_cast<G4ITTransportation*>((fpProcessInfo->fpAlongStepGetPhysIntVector)[MAXofAlongStepLoops-1]);
809 
810  if(!fpTrack)
811  {
812  G4ExceptionDescription exceptionDescription;
813  exceptionDescription << "No G4ITStepProcessor::fpTrack found";
814  G4Exception("G4ITStepProcessor::FindTransportationStep",
815  "ITStepProcessor0013",
817  exceptionDescription);
818  return;
819 
820  }
821  if(!fpITrack)
822  {
823  G4ExceptionDescription exceptionDescription;
824  exceptionDescription << "No G4ITStepProcessor::fITrack";
825  G4Exception("G4ITStepProcessor::FindTransportationStep",
826  "ITStepProcessor0014",
828  exceptionDescription);
829  return;
830  }
831  if(!(fpITrack->GetTrack()))
832  {
833  G4ExceptionDescription exceptionDescription;
834  exceptionDescription << "No G4ITStepProcessor::fITrack->GetTrack()";
835  G4Exception("G4ITStepProcessor::FindTransportationStep",
836  "ITStepProcessor0015",
838  exceptionDescription);
839  return;
840  }
841 
842  if(fpTransportation)
843  {
844  fpTransportation->SetProcessState(fpTrackingInfo->GetProcessState(fpTransportation
845  ->GetProcessID()));
846  fpTransportation->ComputeStep(*fpTrack, *fpStep, fTimeStep, physicalStep);
847 
848 // fpTransportation->SetProcessState(0);
849  fpTransportation->ResetProcessState();
850  }
851 
852  if(physicalStep >= DBL_MAX)
853  {
854 // G4cout << "---- 2) Setting stop and kill for " << GetIT(fpTrack)->GetName() << G4endl;
855  fpTrack -> SetTrackStatus(fStopAndKill);
856  return;
857  }
858 
859  fpState->fPhysicalStep = physicalStep;
860 }
861 
862 //______________________________________________________________________________
863 
865 {
866  size_t _MAXofPostStepLoops = fpProcessInfo->MAXofPostStepLoops;
867  G4SelectedPostStepDoItVector& selectedPostStepDoItVector = fpState
868  ->fSelectedPostStepDoItVector;
869  G4StepStatus& stepStatus = fpState->fStepStatus;
870 
871  // Invoke the specified discrete processes
872  for(size_t np = 0; np < _MAXofPostStepLoops; np++)
873  {
874  //
875  // Note: DoItVector has inverse order against GetPhysIntVector
876  // and SelectedPostStepDoItVector.
877  //
878  G4int Cond = selectedPostStepDoItVector[_MAXofPostStepLoops - np - 1];
879  if(Cond != InActivated)
880  {
881  if(((Cond == Forced) && (stepStatus != fExclusivelyForcedProc)) ||
882  // ((Cond == Conditionally) && (stepStatus == fAlongStepDoItProc)) ||
883  ((Cond == ExclusivelyForced) && (stepStatus == fExclusivelyForcedProc))
884  || ((Cond == StronglyForced)))
885  {
886 
887  InvokePSDIP(np);
888  }
889  } //if(Cond != InActivated)
890 
891  // Exit from PostStepLoop if the track has been killed,
892  // but extra treatment for processes with Strongly Forced flag
893  if(fpTrack->GetTrackStatus() == fStopAndKill)
894  {
895  for(size_t np1 = np + 1; np1 < _MAXofPostStepLoops; np1++)
896  {
897  G4int Cond2 = selectedPostStepDoItVector[_MAXofPostStepLoops - np1 - 1];
898  if(Cond2 == StronglyForced)
899  {
900  InvokePSDIP(np1);
901  }
902  }
903  break;
904  }
905  }
906 }
907 
908 //______________________________________________________________________________
909 
910 // ************************************************************************
911 // Apply cuts
912 // ************************************************************************
913 
915 {
916  G4bool tBelowCutEnergyAndSafety = false;
917  G4int tPtclIdx = G4ProductionCuts::GetIndex(aSecondary->GetDefinition());
918  if(tPtclIdx < 0)
919  {
920  return;
921  }
922  G4ProductionCutsTable* tCutsTbl =
924  G4int tCoupleIdx = tCutsTbl->GetCoupleIndex(fpPreStepPoint
925  ->GetMaterialCutsCouple());
926  G4double tProdThreshold =
927  (*(tCutsTbl->GetEnergyCutsVector(tPtclIdx)))[tCoupleIdx];
928  if(aSecondary->GetKineticEnergy() < tProdThreshold)
929  {
930  tBelowCutEnergyAndSafety = true;
931  if(std::abs(aSecondary->GetDynamicParticle()->GetCharge()) > DBL_MIN)
932  {
933  G4double currentRange
935  aSecondary->GetKineticEnergy(),
936  fpPreStepPoint->GetMaterialCutsCouple());
937  tBelowCutEnergyAndSafety = (currentRange < CalculateSafety() );
938  }
939  }
940 
941  if(tBelowCutEnergyAndSafety)
942  {
943  if(!(aSecondary->IsGoodForTracking()))
944  {
945  // Add kinetic energy to the total energy deposit
946  fpStep->AddTotalEnergyDeposit(aSecondary->GetKineticEnergy());
947  aSecondary->SetKineticEnergy(0.0);
948  }
949  }
950 }
void SetTrackStatus(const G4TrackStatus aTrackStatus)
G4ParticleDefinition * GetDefinition() const
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
static G4LossTableManager * Instance()
static G4int GetIndex(const G4String &name)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4DynamicParticle * GetDynamicParticle() const
class std::vector< int, std::allocator< int > > G4SelectedAtRestDoItVector
void RemoveReactionSet(G4Track *track)
G4TrackStatus GetTrackStatus() const
size_t size() const
int G4int
Definition: G4Types.hh:78
MemStat MemoryUsage()
Definition: G4MemStat.cc:55
G4StepStatus
Definition: G4StepStatus.hh:49
void SetCreatorProcess(const G4VProcess *aValue)
static G4ITReactionSet * Instance()
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:49
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4double GetKineticEnergy() const
G4int entries() const
class std::vector< int, std::allocator< int > > G4SelectedPostStepDoItVector
G4GLOB_DLL std::ostream G4cout
void DoIt(double timeStep)
void ApplyProductionCut(G4Track *)
bool G4bool
Definition: G4Types.hh:79
G4int GetCoupleIndex(const G4MaterialCutsCouple *aCouple) const
G4double GetCharge() const
void DealWithSecondaries(G4int &)
const G4double kCarTolerance
G4int GetTrackID() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4ProductionCutsTable * GetProductionCutsTable()
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4ProcessManager * GetProcessManager() const
void Stepping(G4Track *, const double &)
void SetParentID(const G4int aValue)
#define DBL_MIN
Definition: templates.hh:75
G4bool IsGoodForTracking() const
#define G4endl
Definition: G4ios.hh:61
G4bool GetApplyCutsFlag() const
void SetKineticEnergy(const G4double aValue)
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4TrackStatus