Geant4  10.01.p01
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 87375 2014-12-02 08:17:28Z 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 //#define DEBUG_MEM 1
50 
51 #ifdef DEBUG_MEM
52 #include "G4MemStat.hh"
53 using namespace G4MemStat;
54 using G4MemStat::MemStat;
55 #endif
56 
58 {
59  // Now Store the secondaries from ParticleChange to SecondaryList
60  G4Track* tempSecondaryTrack;
61 
62  for(G4int DSecLoop=0 ;
63  DSecLoop<fpParticleChange->GetNumberOfSecondaries() ;
64  DSecLoop++)
65  {
66  tempSecondaryTrack = fpParticleChange->GetSecondary(DSecLoop);
67 
68  if(tempSecondaryTrack->GetDefinition()->GetApplyCutsFlag())
69  {
70  ApplyProductionCut(tempSecondaryTrack);
71  }
72 
73  // Set parentID
74  tempSecondaryTrack->SetParentID( fpTrack->GetTrackID() );
75 
76  // Set the process pointer which created this track
77  tempSecondaryTrack->SetCreatorProcess( fpCurrentProcess );
78 
79  // If this 2ndry particle has 'zero' kinetic energy, make sure
80  // it invokes a rest process at the beginning of the tracking
81  if(tempSecondaryTrack->GetKineticEnergy() <= DBL_MIN)
82  {
83  G4ProcessManager* pm = tempSecondaryTrack->GetDefinition()->GetProcessManager();
84  if (pm->GetAtRestProcessVector()->entries()>0){
85  tempSecondaryTrack->SetTrackStatus( fStopButAlive );
86  fpSecondary->push_back( tempSecondaryTrack );
87  fN2ndariesAtRestDoIt++;
88  } else {
89  delete tempSecondaryTrack;
90  }
91  }
92  else
93  {
94  fpSecondary->push_back( tempSecondaryTrack );
95  counter++;
96  }
97  } //end of loop on secondary
98 }
99 
100 //______________________________________________________________________________
101 
102 void G4ITStepProcessor::Stepping(G4Track* track, const double & timeStep)
103 {
104 
105 #ifdef DEBUG_MEM
106  MemStat mem_first, mem_second, mem_diff;
107 #endif
108 
109 #ifdef DEBUG_MEM
110  mem_first = MemoryUsage();
111 #endif
112 
113  CleanProcessor();
114 
115 #ifdef DEBUG_MEM
116  MemStat mem_intermediaire = MemoryUsage();
117  mem_diff = mem_intermediaire-mem_first;
118  G4cout << "\t\t\t >> || MEM || After CleanProcessor " << track->GetTrackID() << ", diff is : " << mem_diff << G4endl;
119 #endif
120 
121  if(track == 0) return ; // maybe put an exception here
122  fTimeStep = timeStep ;
123  SetTrack(track);
124  DoStepping();
125 }
126 //______________________________________________________________________________
127 
128 // ************************************************************************
129 // Stepping
130 // ************************************************************************
132 {
133  SetupMembers() ;
134 
135 #ifdef DEBUG_MEM
136  MemStat mem_first, mem_second, mem_diff;
137 #endif
138 
139 #ifdef DEBUG_MEM
140  mem_first = MemoryUsage();
141 #endif
142 
143  if(!fpProcessInfo)
144  {
145  G4ExceptionDescription exceptionDescription ;
146  exceptionDescription << "No process info found for particle :"
147  << fpTrack->GetDefinition()->GetParticleName();
148  G4Exception("G4ITStepProcessor::DoStepping","ITStepProcessor0012",
149  FatalErrorInArgument,exceptionDescription);
150  return ;
151  }
152  else if(fpTrack->GetTrackStatus() == fStopAndKill )
153  {
154  fpState->fStepStatus = fUndefined;
155  return ;
156  }
157 
158  if(fpProcessInfo->MAXofPostStepLoops == 0
159  && fpProcessInfo->MAXofAlongStepLoops == 0
160  && fpProcessInfo->MAXofAtRestLoops == 0)
161  {
162  G4ExceptionDescription exceptionDescription ;
163  exceptionDescription << "No process was found for particle :"
164  << fpTrack->GetDefinition()->GetParticleName();
165  G4Exception("G4ITStepProcessor::DoStepping","ITStepProcessorNoProcess",
166  JustWarning,exceptionDescription);
167 
168  fpTrack -> SetTrackStatus(fStopAndKill) ;
169  fpState->fStepStatus = fUndefined;
170  return ;
171  }
172  //---------------------------------
173  // AtRestStep, AlongStep and PostStep Processes
174  //---------------------------------
175  else
176  {
177  fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()->GetNavigatorState());
178 // fpNavigator->ResetHierarchyAndLocate( fpTrack->GetPosition(),
179 // fpTrack->GetMomentumDirection(),
180 // *((G4TouchableHistory*)fpTrack->GetTouchableHandle()()) );
181 // fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()->GetNavigatorState());
182  // We reset the navigator state before checking for AtRest
183  // in case a AtRest processe would use a navigator info
184 
185 #ifdef DEBUG_MEM
186  MemStat mem_intermediaire = MemoryUsage();
187  mem_diff = mem_intermediaire-mem_first;
188  G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After dealing with navigator with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
189 #endif
190 
191  if( fpTrack->GetTrackStatus() == fStopButAlive )
192  {
193  if( fpProcessInfo->MAXofAtRestLoops>0 &&
194  fpProcessInfo->fpAtRestDoItVector != 0) // second condition to make coverity happy
195  {
196  //-----------------
197  // AtRestStepDoIt
198  //-----------------
199  InvokeAtRestDoItProcs();
200  fpState->fStepStatus = fAtRestDoItProc;
201  fpStep->GetPostStepPoint()->SetStepStatus( fpState->fStepStatus );
202 
203  }
204  // Make sure the track is killed
205  fpTrack->SetTrackStatus( fStopAndKill );
206  }
207  else // if(fTimeStep > 0.) // Bye, because PostStepIL can return 0 => time =0
208  {
209  if(fpITrack == 0)
210  {
211  G4ExceptionDescription exceptionDescription ;
212  exceptionDescription
213  << " !!! TrackID : "<< fpTrack->GetTrackID() << G4endl
214  << " !!! Track status : "<< fpTrack->GetTrackStatus() << G4endl
215  << " !!! Particle Name : "<< fpTrack -> GetDefinition() -> GetParticleName() << G4endl
216  << "No G4ITStepProcessor::fpITrack found" << G4endl;
217 
218  G4Exception("G4ITStepProcessor::DoStepping","ITStepProcessor0013",
219  FatalErrorInArgument,exceptionDescription);
220  return ; // to make coverity happy
221  }
222 
223  if(fpITrack->GetTrackingInfo()->IsLeadingStep() == false)
224  {
225  // In case the track has NOT the minimum step length
226  // Given the final step time, the transportation
227  // will compute the final position of the particle
228  fpState->fStepStatus = fPostStepDoItProc;
229  fpStep->GetPostStepPoint()
230  ->SetProcessDefinedStep(fpTransportation);
231  FindTransportationStep();
232  }
233 
234 #ifdef DEBUG_MEM
235  mem_intermediaire = MemoryUsage();
236  mem_diff = mem_intermediaire-mem_first;
237  G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After FindTransportationStep() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
238 #endif
239 
240  // Store the Step length (geometrical length) to G4Step and G4Track
241  fpTrack->SetStepLength( fpState->fPhysicalStep );
242  fpStep->SetStepLength( fpState->fPhysicalStep );
243 
244  G4double GeomStepLength = fpState->fPhysicalStep;
245 
246  // Store StepStatus to PostStepPoint
247  fpStep->GetPostStepPoint()->SetStepStatus( fpState->fStepStatus );
248 
249  // Invoke AlongStepDoIt
250  InvokeAlongStepDoItProcs();
251 
252 #ifdef DEBUG_MEM
253  mem_intermediaire = MemoryUsage();
254  mem_diff = mem_intermediaire-mem_first;
255  G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After InvokeAlongStepDoItProcs() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
256 #endif
257 
258  // Update track by taking into account all changes by AlongStepDoIt
259  // fpStep->UpdateTrack(); // done in InvokeAlongStepDoItProcs
260 
261  // Update safety after invocation of all AlongStepDoIts
262  fpState->endpointSafOrigin= fpPostStepPoint->GetPosition();
263 
264  fpState->endpointSafety= std::max( fpState->proposedSafety - GeomStepLength, kCarTolerance);
265 
266  fpStep->GetPostStepPoint()->SetSafety( fpState->endpointSafety );
267 
268  if(GetIT(fpTrack)->GetTrackingInfo()->IsLeadingStep())
269  {
270  // Invoke PostStepDoIt including G4ITTransportation::PSDI
271  InvokePostStepDoItProcs();
272 
273 #ifdef DEBUG_MEM
274  mem_intermediaire = MemoryUsage();
275  mem_diff = mem_intermediaire-mem_first;
276  G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After InvokePostStepDoItProcs() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
277 #endif
278  }
279  else
280  {
281  // Only invoke transportation
282  InvokeTransportationProc();
283 
284 #ifdef DEBUG_MEM
285  mem_intermediaire = MemoryUsage();
286  mem_diff = mem_intermediaire-mem_first;
287  G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After InvokeTransportationProc() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
288 #endif
289  }
290  }
291 
292  fpNavigator->ResetNavigatorState();
293 
294 #ifdef DEBUG_MEM
295  mem_intermediaire = MemoryUsage();
296  mem_diff = mem_intermediaire-mem_first;
297  G4cout << "\t\t\t >> || MEM || G4ITStepProcessor::DoStepping || After fpNavigator->SetNavigatorState with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
298 #endif
299  }
300  //-------
301  // Finale
302  //-------
303 
304  // Update 'TrackLength' and remeber the Step length of the current Step
305  fpTrack->AddTrackLength(fpStep->GetStepLength());
306  fpTrack->IncrementCurrentStepNumber();
307 
308 // G4cout << " G4ITStepProcessor::DoStepping -- " <<fpTrack->GetTrackID() << " tps = " << fpTrack->GetGlobalTime() << G4endl;
309 
310  // Send G4Step information to Hit/Dig if the volume is sensitive
311 /***
312  fpCurrentVolume = fpStep->GetPreStepPoint()->GetPhysicalVolume();
313  StepControlFlag = fpStep->GetControlFlag();
314 
315  if( fpCurrentVolume != 0 && StepControlFlag != AvoidHitInvocation)
316  {
317  fpSensitive = fpStep->GetPreStepPoint()->
318  GetSensitiveDetector();
319  if( fpSensitive != 0 )
320  {
321  fpSensitive->Hit(fpStep);
322  }
323  }
324 
325  User intervention process.
326  if( fpUserSteppingAction != 0 )
327  {
328  fpUserSteppingAction->UserSteppingAction(fpStep);
329  }
330  G4UserSteppingAction* regionalAction
331  = fpStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetRegion()
332  ->GetRegionalSteppingAction();
333  if( regionalAction ) regionalAction->UserSteppingAction(fpStep);
334 ***/
335  fpTrackingManager->AppendStep(fpTrack,fpStep);
336  // Stepping process finish. Return the value of the StepStatus.
337 
338 
339 #ifdef DEBUG_MEM
340  MemStat mem_intermediaire = MemoryUsage();
341  mem_diff = mem_intermediaire-mem_first;
342  G4cout << "\t\t\t >> || MEM || End of DoStepping() with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
343 #endif
344 
345 
346  // return fpState->fStepStatus;
347 }
348 
349 //______________________________________________________________________________
350 
351 // ************************************************************************
352 // AtRestDoIt
353 // ************************************************************************
354 
356 {
357  fpStep->SetStepLength( 0. ); //the particle has stopped
358  fpTrack->SetStepLength( 0. );
359 
360  G4SelectedAtRestDoItVector& selectedAtRestDoItVector = fpState->fSelectedAtRestDoItVector;
361 
362  // invoke selected process
363  for(size_t np=0; np < fpProcessInfo->MAXofAtRestLoops; np++)
364  {
365  //
366  // Note: DoItVector has inverse order against GetPhysIntVector
367  // and SelectedAtRestDoItVector.
368  //
369  if( selectedAtRestDoItVector[fpProcessInfo->MAXofAtRestLoops-np-1] != InActivated)
370  {
371  fpCurrentProcess = (G4VITProcess*) (*fpProcessInfo->fpAtRestDoItVector)[np];
372 
373  fpCurrentProcess->SetProcessState(
374  fpTrackingInfo->GetProcessState(fpCurrentProcess->GetProcessID()));
375  fpParticleChange
376  = fpCurrentProcess->AtRestDoIt( *fpTrack, *fpStep);
377 // fpCurrentProcess->SetProcessState(0);
378  fpCurrentProcess->ResetProcessState();
379 
380  // Set the current process as a process which defined this Step length
381  fpStep->GetPostStepPoint()
382  ->SetProcessDefinedStep(fpCurrentProcess);
383 
384  // Update Step
385  fpParticleChange->UpdateStepForAtRest(fpStep);
386 
387  // Now Store the secondaries from ParticleChange to SecondaryList
388  DealWithSecondaries(fN2ndariesAtRestDoIt) ;
389 
390  // clear ParticleChange
391  fpParticleChange->Clear();
392 
393  } //if(fSelectedAtRestDoItVector[np] != InActivated){
394  } //for(size_t np=0; np < MAXofAtRestLoops; np++){
395  fpStep->UpdateTrack();
396 
397  fpTrack->SetTrackStatus( fStopAndKill );
398 }
399 
400 //______________________________________________________________________________
401 
402 // ************************************************************************
403 // AlongStepDoIt
404 // ************************************************************************
405 
407 {
408 
409 #ifdef DEBUG_MEM
410  MemStat mem_first, mem_second, mem_diff;
411 #endif
412 
413 #ifdef DEBUG_MEM
414  mem_first = MemoryUsage();
415 #endif
416 
417  // If the current Step is defined by a 'ExclusivelyForced'
418  // PostStepDoIt, then don't invoke any AlongStepDoIt
419  if(fpState->fStepStatus == fExclusivelyForcedProc)
420  {
421  return; // Take note 'return' at here !!!
422  }
423 
424  // Invoke the all active continuous processes
425  for( size_t ci=0 ; ci<fpProcessInfo->MAXofAlongStepLoops ; ci++ )
426  {
427  fpCurrentProcess = (G4VITProcess*) (*fpProcessInfo->fpAlongStepDoItVector)[ci];
428  if (fpCurrentProcess== 0) continue;
429  // NULL means the process is inactivated by a user on fly.
430 
431  fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess->GetProcessID()));
432  fpParticleChange = fpCurrentProcess->AlongStepDoIt( *fpTrack, *fpStep );
433 
434 #ifdef DEBUG_MEM
435  MemStat mem_intermediaire = MemoryUsage();
436  mem_diff = mem_intermediaire-mem_first;
437  G4cout << "\t\t\t >> || MEM || After calling AlongStepDoIt for " << fpCurrentProcess->GetProcessName() << " and track "<< fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
438 #endif
439 
440 // fpCurrentProcess->SetProcessState(0);
441  fpCurrentProcess->ResetProcessState();
442  // Update the PostStepPoint of Step according to ParticleChange
443 
444  fpParticleChange->UpdateStepForAlongStep(fpStep);
445 
446  // Now Store the secondaries from ParticleChange to SecondaryList
447  DealWithSecondaries(fN2ndariesAlongStepDoIt) ;
448 
449  // Set the track status according to what the process defined
450  // if kinetic energy >0, otherwise set fStopButAlive
451  fpTrack->SetTrackStatus( fpParticleChange->GetTrackStatus() );
452 
453  // clear ParticleChange
454  fpParticleChange->Clear();
455  }
456 
457 #ifdef DEBUG_MEM
458  MemStat mem_intermediaire = MemoryUsage();
459  mem_diff = mem_intermediaire-mem_first;
460  G4cout << "\t\t\t >> || MEM || After looping on processes with " << fpTrack->GetTrackID() << ", diff is : " << mem_diff << G4endl;
461 #endif
462 
463  fpStep->UpdateTrack();
464 
465  G4TrackStatus fNewStatus = fpTrack->GetTrackStatus();
466 
467  if ( fNewStatus == fAlive && fpTrack->GetKineticEnergy() <= DBL_MIN )
468  {
469  // G4cout << "G4ITStepProcessor::InvokeAlongStepDoItProcs : Track will be killed" << G4endl;
470  if(fpProcessInfo->MAXofAtRestLoops>0) fNewStatus = fStopButAlive;
471  else fNewStatus = fStopAndKill;
472  fpTrack->SetTrackStatus( fNewStatus );
473  }
474 
475 }
476 
477 //______________________________________________________________________________
478 
479 
480 // ************************************************************************
481 // PostStepDoIt
482 // ************************************************************************
483 
484 
486 {
487 
488  G4SelectedPostStepDoItVector& selectedPostStepDoItVector = fpState->fSelectedPostStepDoItVector;
489  G4StepStatus& stepStatus = fpState->fStepStatus;
490 
491  // Invoke the specified discrete processes
492  for(size_t np=0; np < fpProcessInfo->MAXofPostStepLoops; np++)
493  {
494  //
495  // Note: DoItVector has inverse order against GetPhysIntVector
496  // and SelectedPostStepDoItVector.
497  //
498  G4int Cond = selectedPostStepDoItVector[fpProcessInfo->MAXofPostStepLoops-np-1];
499  if(Cond != InActivated)
500  {
501  if( ((Cond == NotForced) && (stepStatus == fPostStepDoItProc)) ||
502  ((Cond == Forced) && (stepStatus != fExclusivelyForcedProc)) ||
503  // ((Cond == Conditionally) && (stepStatus == fAlongStepDoItProc)) ||
504  ((Cond == ExclusivelyForced) && (stepStatus == fExclusivelyForcedProc)) ||
505  ((Cond == StronglyForced) )
506  )
507  {
508 
509  InvokePSDIP(np);
510  }
511  } //if(*fSelectedPostStepDoItVector(np)........
512 
513  // Exit from PostStepLoop if the track has been killed,
514  // but extra treatment for processes with Strongly Forced flag
515  if(fpTrack->GetTrackStatus() == fStopAndKill)
516  {
517  for(size_t np1=np+1; np1 < fpProcessInfo->MAXofPostStepLoops; np1++)
518  {
519  G4int Cond2 = selectedPostStepDoItVector[fpProcessInfo->MAXofPostStepLoops-np1-1];
520  if (Cond2 == StronglyForced)
521  {
522  InvokePSDIP(np1);
523  }
524  }
525  break;
526  }
527  } //for(size_t np=0; np < MAXofPostStepLoops; np++){
528 }
529 
530 //______________________________________________________________________________
531 
533 {
534  fpCurrentProcess = (G4VITProcess*) (*fpProcessInfo->fpPostStepDoItVector)[np];
535 
536  fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess->GetProcessID()));
537  fpParticleChange
538  = fpCurrentProcess->PostStepDoIt( *fpTrack, *fpStep);
539 // fpCurrentProcess->SetProcessState(0);
540  fpCurrentProcess->ResetProcessState();
541 
542  // Update PostStepPoint of Step according to ParticleChange
543  fpParticleChange->UpdateStepForPostStep(fpStep);
544 
545  // Update G4Track according to ParticleChange after each PostStepDoIt
546  fpStep->UpdateTrack();
547 
548  // Update safety after each invocation of PostStepDoIts
549  fpStep->GetPostStepPoint()->SetSafety( CalculateSafety() );
550 
551  // Now Store the secondaries from ParticleChange to SecondaryList
552  DealWithSecondaries(fN2ndariesPostStepDoIt) ;
553 
554  // Set the track status according to what the process defined
555  fpTrack->SetTrackStatus( fpParticleChange->GetTrackStatus() );
556 
557  // clear ParticleChange
558  fpParticleChange->Clear();
559 }
560 
561 //______________________________________________________________________________
562 
563 // ************************************************************************
564 // Transport on a given time
565 // ************************************************************************
566 
568 {
569  double physicalStep(0.) ;
570 
571  fpTransportation = fpProcessInfo->fpTransportation;
572  // dynamic_cast<G4ITTransportation*>((fpProcessInfo->fpAlongStepGetPhysIntVector)[MAXofAlongStepLoops-1]);
573 
574  if(!fpTrack)
575  {
576  G4ExceptionDescription exceptionDescription ;
577  exceptionDescription
578  << "No G4ITStepProcessor::fpTrack found";
579  G4Exception("G4ITStepProcessor::FindTransportationStep","ITStepProcessor0013",
580  FatalErrorInArgument,exceptionDescription);
581  return;
582 
583  }
584  if(!fpITrack)
585  {
586  G4ExceptionDescription exceptionDescription ;
587  exceptionDescription
588  << "No G4ITStepProcessor::fITrack" ;
589  G4Exception("G4ITStepProcessor::FindTransportationStep","ITStepProcessor0014",
590  FatalErrorInArgument,exceptionDescription);
591  return;
592  }
593  if(!(fpITrack->GetTrack()))
594  {
595  G4ExceptionDescription exceptionDescription ;
596  exceptionDescription
597  << "No G4ITStepProcessor::fITrack->GetTrack()" ;
598  G4Exception("G4ITStepProcessor::FindTransportationStep","ITStepProcessor0015",
599  FatalErrorInArgument,exceptionDescription);
600  return;
601  }
602 
603  if(fpTransportation)
604  {
605  fpTransportation->SetProcessState(fpTrackingInfo->GetProcessState(fpTransportation->GetProcessID()));
606  fpTransportation->ComputeStep(*fpTrack, *fpStep, fTimeStep, physicalStep) ;
607 // fpTransportation->SetProcessState(0);
608  fpTransportation->ResetProcessState();
609  }
610 
611  if(physicalStep >= DBL_MAX)
612  {
613  G4cout << "---- 2) Setting stop and kill for " << GetIT(fpTrack)->GetName() << G4endl;
614 
615  fpTrack -> SetTrackStatus(fStopAndKill) ;
616  return ;
617  }
618 
619  fpState->fPhysicalStep = physicalStep ;
620 }
621 
622 //______________________________________________________________________________
623 
625 {
626  size_t _MAXofPostStepLoops = fpProcessInfo->MAXofPostStepLoops;
627  G4SelectedPostStepDoItVector& selectedPostStepDoItVector = fpState->fSelectedPostStepDoItVector;
628  G4StepStatus& stepStatus = fpState->fStepStatus;
629 
630  // Invoke the specified discrete processes
631  for(size_t np=0; np < _MAXofPostStepLoops; np++)
632  {
633  //
634  // Note: DoItVector has inverse order against GetPhysIntVector
635  // and SelectedPostStepDoItVector.
636  //
637  G4int Cond = selectedPostStepDoItVector[_MAXofPostStepLoops-np-1];
638  if(Cond != InActivated)
639  {
640  if(
641  ((Cond == Forced) && (stepStatus != fExclusivelyForcedProc)) ||
642  // ((Cond == Conditionally) && (stepStatus == fAlongStepDoItProc)) ||
643  ((Cond == ExclusivelyForced) && (stepStatus == fExclusivelyForcedProc)) ||
644  ((Cond == StronglyForced) )
645  )
646  {
647 
648  InvokePSDIP(np);
649  }
650  } //if(*fSelectedPostStepDoItVector(np)........
651 
652  // Exit from PostStepLoop if the track has been killed,
653  // but extra treatment for processes with Strongly Forced flag
654  if(fpTrack->GetTrackStatus() == fStopAndKill)
655  {
656  for(size_t np1=np+1; np1 < _MAXofPostStepLoops; np1++)
657  {
658  G4int Cond2 = selectedPostStepDoItVector[_MAXofPostStepLoops-np1-1];
659  if (Cond2 == StronglyForced)
660  {
661  InvokePSDIP(np1);
662  }
663  }
664  break;
665  }
666  }
667 }
668 
669 //______________________________________________________________________________
670 
671 // ************************************************************************
672 // Apply cuts
673 // ************************************************************************
674 
675 
677 {
678  G4bool tBelowCutEnergyAndSafety = false;
679  G4int tPtclIdx
680  = G4ProductionCuts::GetIndex(aSecondary->GetDefinition());
681  if (tPtclIdx<0)
682  {
683  return;
684  }
685  G4ProductionCutsTable* tCutsTbl
687  G4int tCoupleIdx
688  = tCutsTbl->GetCoupleIndex(fpPreStepPoint->GetMaterialCutsCouple());
689  G4double tProdThreshold
690  = (*(tCutsTbl->GetEnergyCutsVector(tPtclIdx)))[tCoupleIdx];
691  if( aSecondary->GetKineticEnergy()<tProdThreshold )
692  {
693  tBelowCutEnergyAndSafety = true;
694  if(std::abs(aSecondary->GetDynamicParticle()->GetCharge()) > DBL_MIN)
695  {
696  G4double currentRange
698  aSecondary->GetKineticEnergy(),
699  fpPreStepPoint->GetMaterialCutsCouple());
700  tBelowCutEnergyAndSafety = (currentRange < CalculateSafety() );
701  }
702  }
703 
704  if( tBelowCutEnergyAndSafety )
705  {
706  if( !(aSecondary->IsGoodForTracking()) )
707  {
708  // Add kinetic energy to the total energy deposit
709  fpStep->AddTotalEnergyDeposit(
710  aSecondary->GetKineticEnergy() );
711  aSecondary->SetKineticEnergy(0.0);
712  }
713  }
714 }
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
virtual const G4String & GetName() const =0
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
MemStat MemoryUsage()
Definition: G4MemStat.cc:55
G4StepStatus
Definition: G4StepStatus.hh:49
void SetCreatorProcess(const G4VProcess *aValue)
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 ApplyProductionCut(G4Track *)
bool G4bool
Definition: G4Types.hh:79
G4int GetCoupleIndex(const G4MaterialCutsCouple *aCouple) const
G4double GetCharge() const
void DealWithSecondaries(G4int &)
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
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
G4VITProcess inherits from G4VProcess.
Definition: G4VITProcess.hh:99
#define DBL_MAX
Definition: templates.hh:83
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4TrackStatus