Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ITStepProcessor.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: G4ITStepProcessor.cc 65022 2012-11-12 16:43:12Z 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 "G4UImanager.hh"
38 #include "G4ForceCondition.hh"
39 #include "G4GPILSelection.hh"
41 // #include "G4VSensitiveDetector.hh" // Include from 'hits/digi'
42 #include "G4GeometryTolerance.hh"
43 #include "G4ParticleTable.hh"
44 #include "G4ITTrackingManager.hh"
45 #include "G4TrackingInformation.hh"
46 #include "G4IT.hh"
47 #include "G4ITNavigator.hh" // Include from 'geometry'
48 
49 #include "G4VITProcess.hh"
50 #include "G4VProcess.hh"
51 #include "G4ITTransportation.hh"
52 
53 #include <iomanip> // Include from 'system'
54 #include <vector> // Include from 'system'
55 
56 using namespace std;
57 
58 static const size_t SizeOfSelectedDoItVector=100;
59 static const size_t& gMaxNProcesses(G4VITProcess::GetMaxProcessIndex());
60 
61 //____________________________________________________________________________________
62 
64 {
65  verboseLevel = 0 ;
66  // fpUserSteppingAction = 0 ;
67  fStoreTrajectory = 0;
68  fpTrackingManager = 0;
69  fpNavigator = 0;
70  kCarTolerance = -1.;
71  fInitialized = false;
72  fPreviousTimeStep = DBL_MAX;
73  CleanProcessor();
74  ResetSecondaries();
75 }
76 
77 G4ITStepProcessor::G4ITStepProcessorState::G4ITStepProcessorState() :
79  fSelectedAtRestDoItVector (gMaxNProcesses,0),
80  fSelectedPostStepDoItVector (gMaxNProcesses,0)
81 {
82  fPhysicalStep = -1.;
83  fPreviousStepSize = -1.;
84 
85  fSafety = -1.;
86  proposedSafety = -1.;
87  endpointSafety = -1;
88 
89  fStepStatus = fUndefined;
90 
91  fTouchableHandle = 0;
92 }
93 
94 // should not be used
95 G4ITStepProcessor::G4ITStepProcessorState::G4ITStepProcessorState(const G4ITStepProcessorState& ) :
97  fSelectedAtRestDoItVector (gMaxNProcesses,0),
98  fSelectedPostStepDoItVector (gMaxNProcesses,0)
99 {
100  fPhysicalStep = -1.;
101  fPreviousStepSize = -1.;
102 
103  fSafety = -1.;
104  proposedSafety = -1.;
105  endpointSafety = -1;
106 
107  fStepStatus = fUndefined;
108 
109  fTouchableHandle = 0;
110 }
111 
112 // should not be used
113 G4ITStepProcessor::G4ITStepProcessorState& G4ITStepProcessor::G4ITStepProcessorState::operator=(const G4ITStepProcessorState& rhs)
114 {
115  if(this == &rhs) return *this;
116 
117  fSelectedAtRestDoItVector.clear();
118  fSelectedAtRestDoItVector.resize(gMaxNProcesses,0);
119  fSelectedPostStepDoItVector.clear();
120  fSelectedPostStepDoItVector.resize(gMaxNProcesses,0);
121 
122  fPhysicalStep = -1.;
123  fPreviousStepSize = -1.;
124 
125  fSafety = -1.;
126  proposedSafety = -1.;
127  endpointSafety = -1;
128 
129  fStepStatus = fUndefined;
130 
131  fTouchableHandle = 0;
132  return *this;
133 }
134 //____________________________________________________________________________________
135 
136 G4ITStepProcessor::G4ITStepProcessorState::~G4ITStepProcessorState()
137 {;}
138 //____________________________________________________________________________________
139 
141 {
142  std::map<const G4ParticleDefinition*, ProcessGeneralInfo*> ::iterator it;
143 
144  for(it = fProcessGeneralInfoMap.begin();it != fProcessGeneralInfoMap.end();it++)
145  {
146  if(it->second)
147  {
148  delete it->second;
149  it->second = 0;
150  }
151  }
152 
153  fProcessGeneralInfoMap.clear();
154 }
155 
156 //____________________________________________________________________________________
157 
159 {
160  fInitialized = false;
162  Initialize();
163 }
164 
165 //____________________________________________________________________________________
166 
168 {
169  CleanProcessor();
170  if(fInitialized) return;
171  // ActiveOnlyITProcess();
172 
174  ->GetNavigatorForTracking());
175 
176  fPhysIntLength = DBL_MAX;
177  kCarTolerance = 0.5*G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
178 
179  fInitialized = true;
180 }
181 //______________________________________________________________________________
182 
184 {
185  if(fpStep)
186  {
187  fpStep->DeleteSecondaryVector();
188  delete fpStep;
189  }
190 
191  if(fpSecondary) delete fpSecondary;
194 
195  // if(fpUserSteppingAction) delete fpUserSteppingAction;
196 }
197 //______________________________________________________________________________
198 // should not be used
200 {
201  verboseLevel = rhs.verboseLevel ;
202  fStoreTrajectory = rhs.fStoreTrajectory ;
203 
204  // fpUserSteppingAction = 0 ;
205  fpTrackingManager = 0;
206  fpNavigator = 0;
207  fInitialized = false;
208 
209  kCarTolerance = rhs.kCarTolerance;
210  fInitialized = false;
211  fPreviousTimeStep = DBL_MAX;
212 
213  CleanProcessor();
215 }
216 //______________________________________________________________________________
217 
219 {
220  if (this == &rhs) return *this; // handle self assignment
221  //assignment operator
222  return *this;
223 }
224 // ******************************************************************
225 
227 {
228  // Method not used for the time being
229 #ifdef debug
230  G4cout<<"G4ITStepProcessor::CloneProcesses: is called"<<G4endl;
231 #endif
232 
234  G4ParticleTable::G4PTblDicIterator* theParticleIterator = theParticleTable->GetIterator();
235 
236  theParticleIterator->reset();
237  // TODO : Ne faire la boucle que sur les IT **** !!!
238  while( (*theParticleIterator)() )
239  {
240  G4ParticleDefinition* particle = theParticleIterator->value();
241  G4ProcessManager* pm= particle->GetProcessManager();
242 
243  if(!pm)
244  {
245  G4cerr << "ERROR - G4ITStepProcessor::GetProcessNumber()" << G4endl
246  << " ProcessManager is NULL for particle = "
247  << particle->GetParticleName() << ", PDG_code = "
248  << particle->GetPDGEncoding() << G4endl;
249  G4Exception("G4ITStepProcessor::GetProcessNumber()", "ITStepProcessor0001",
250  FatalException, "Process Manager is not found.");
251  return;
252  }
253 
255  }
256 }
257 // ******************************************************************
258 
260 {
261  // Method not used for the time being
262  G4ProcessVector* processVector = processManager->GetProcessList();
263 
264  G4VITProcess* itProcess = 0 ;
265  for(int i = 0 ; i < processVector->size() ; i++)
266  {
267  G4VProcess* base_process = (*processVector)[i];
268  itProcess = dynamic_cast<G4VITProcess*>(base_process);
269 
270  if(!itProcess)
271  {
272  processManager->SetProcessActivation(base_process, false);
273  }
274  }
275 }
276 // ******************************************************************
278  G4ProcessManager* pm)
279 {
280 
281 #ifdef debug
282  G4cout<<"G4ITStepProcessor::GetProcessNumber: is called track"<<G4endl;
283 #endif
284  if(!pm)
285  {
286  G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
287  << " ProcessManager is NULL for particle = "
288  << particle->GetParticleName() << ", PDG_code = "
289  << particle->GetPDGEncoding() << G4endl;
290  G4Exception("G4SteppingManager::GetProcessNumber()", "ITStepProcessor0002",
291  FatalException, "Process Manager is not found.");
292  return;
293  }
294 
295  std::map<const G4ParticleDefinition*, ProcessGeneralInfo*>::iterator it = fProcessGeneralInfoMap.find(particle);
296  if(it != fProcessGeneralInfoMap.end())
297  {
298  G4Exception("G4SteppingManager::SetupGeneralProcessInfo()", "ITStepProcessor0003",
299  FatalException, "Process info already registered.");
300  return;
301  }
302 
303  // here used as temporary
304  fpProcessInfo = new ProcessGeneralInfo();
305 
306  // AtRestDoits
307  fpProcessInfo->MAXofAtRestLoops = pm->GetAtRestProcessVector()->entries();
308  fpProcessInfo->fpAtRestDoItVector = pm->GetAtRestProcessVector(typeDoIt);
309  fpProcessInfo->fpAtRestGetPhysIntVector = pm->GetAtRestProcessVector(typeGPIL);
310 #ifdef debug
311  G4cout << "G4ITStepProcessor::GetProcessNumber: #ofAtRest="
312  << fpProcessInfo->MAXofAtRestLoops << G4endl;
313 #endif
314 
315  // AlongStepDoits
316  fpProcessInfo->MAXofAlongStepLoops = pm->GetAlongStepProcessVector()->entries();
317  fpProcessInfo->fpAlongStepDoItVector = pm->GetAlongStepProcessVector(typeDoIt);
318  fpProcessInfo->fpAlongStepGetPhysIntVector = pm->GetAlongStepProcessVector(typeGPIL);
319 #ifdef debug
320  G4cout << "G4ITStepProcessor::GetProcessNumber:#ofAlongStp="
321  << fpProcessInfo->MAXofAlongStepLoops << G4endl;
322 #endif
323 
324  // PostStepDoits
325  fpProcessInfo->MAXofPostStepLoops = pm->GetPostStepProcessVector()->entries();
326  fpProcessInfo->fpPostStepDoItVector = pm->GetPostStepProcessVector(typeDoIt);
327  fpProcessInfo->fpPostStepGetPhysIntVector = pm->GetPostStepProcessVector(typeGPIL);
328 #ifdef debug
329  G4cout << "G4ITStepProcessor::GetProcessNumber: #ofPostStep="
330  << fpProcessInfo->MAXofPostStepLoops << G4endl;
331 #endif
332 
333  if (SizeOfSelectedDoItVector<fpProcessInfo->MAXofAtRestLoops ||
334  SizeOfSelectedDoItVector<fpProcessInfo->MAXofAlongStepLoops ||
335  SizeOfSelectedDoItVector<fpProcessInfo->MAXofPostStepLoops )
336  {
337  G4cerr << "ERROR - G4ITStepProcessor::GetProcessNumber()" << G4endl
338  << " SizeOfSelectedDoItVector= " << SizeOfSelectedDoItVector
339  << " ; is smaller then one of MAXofAtRestLoops= "
340  << fpProcessInfo->MAXofAtRestLoops << G4endl
341  << " or MAXofAlongStepLoops= " << fpProcessInfo->MAXofAlongStepLoops
342  << " or MAXofPostStepLoops= " << fpProcessInfo->MAXofPostStepLoops << G4endl;
343  G4Exception("G4ITStepProcessor::GetProcessNumber()",
344  "ITStepProcessor0004", FatalException,
345  "The array size is smaller than the actual No of processes.");
346  }
347 
348  if(!fpProcessInfo->fpAtRestDoItVector &&
349  !fpProcessInfo->fpAlongStepDoItVector &&
350  !fpProcessInfo->fpPostStepDoItVector)
351  {
352  G4ExceptionDescription exceptionDescription ;
353  exceptionDescription << "No DoIt process found " ;
354  G4Exception("G4ITStepProcessor::DoStepping","ITStepProcessor0005",
355  FatalErrorInArgument,exceptionDescription);
356  return ;
357  }
358 
359  if(fpProcessInfo->fpAlongStepGetPhysIntVector && fpProcessInfo->MAXofAlongStepLoops>0)
360  {
361  fpProcessInfo->fpTransportation = dynamic_cast<G4ITTransportation*>
362  ((*fpProcessInfo->fpAlongStepGetPhysIntVector)[fpProcessInfo->MAXofAlongStepLoops-1]);
363 
364  if(fpProcessInfo->fpTransportation == 0)
365  {
366  G4ExceptionDescription exceptionDescription ;
367  exceptionDescription << "No transportation process found " ;
368  G4Exception("G4ITStepProcessor::SetupGeneralProcessInfo","ITStepProcessor0006",
369  FatalErrorInArgument,exceptionDescription);
370  }
371  }
372  fProcessGeneralInfoMap[particle] = fpProcessInfo;
373  // fpProcessInfo = 0;
374 }
375 
376 // ******************************************************************
377 
379 {
380  fpTrack = track ;
381  if(fpTrack)
382  {
383  fpITrack = GetIT(fpTrack) ;
384  fpStep = const_cast<G4Step*>(fpTrack -> GetStep());
385 
386  if(fpITrack)
387  {
388  fpTrackingInfo = fpITrack->GetTrackingInfo() ;
389  }
390  else
391  {
392  fpTrackingInfo = 0;
393  G4cerr << "Track ID : " << fpTrack->GetTrackID() << G4endl;
394 
395  G4ExceptionDescription exceptionDescription ("No IT pointer was attached to the track you try to process.");
396  G4Exception("G4ITStepProcessor::SetTrack","ITStepProcessor0007",
397  FatalErrorInArgument,exceptionDescription);
398  }
399  }
400  else
401  {
402  fpITrack = 0;
403  fpStep = 0 ;
404  }
405 }
406 //______________________________________________________________________________
407 
409 {
410  G4ParticleDefinition* particle = fpTrack->GetDefinition();
411  std::map<const G4ParticleDefinition*, ProcessGeneralInfo*>::iterator it = fProcessGeneralInfoMap.find(particle);
412 
413  if(it == fProcessGeneralInfoMap.end())
414  {
416  if(fpProcessInfo == 0)
417  {
418  G4ExceptionDescription exceptionDescription ("...");
419  G4Exception("G4ITStepProcessor::GetProcessNumber","ITStepProcessor0008",
420  FatalErrorInArgument,exceptionDescription);
421  return;
422  }
423  }
424  else
425  {
426  fpProcessInfo = it->second;
427  }
428 }
429 //______________________________________________________________________________
430 
432 {
433  fpSecondary = fpStep->GetfSecondary();
434  fpPreStepPoint = fpStep->GetPreStepPoint();
435  fpPostStepPoint = fpStep->GetPostStepPoint();
436 
437  fpState = (G4ITStepProcessorState*) fpITrack->GetTrackingInfo()->GetStepProcessorState();
438 
439  GetProcessInfo();
441 }
442 //______________________________________________________________________________
443 
445 {
446  // Reset the secondary particles
447  fN2ndariesAtRestDoIt = 0;
448  fN2ndariesAlongStepDoIt = 0;
449  fN2ndariesPostStepDoIt = 0;
450 }
451 //______________________________________________________________________________
452 
454 {
455  // Select the rest process which has the shortest time before
456  // it is invoked. In rest processes, GPIL()
457  // returns the time before a process occurs.
458  G4double lifeTime (DBL_MAX), shortestLifeTime (DBL_MAX);
459 
460  fAtRestDoItProcTriggered = 0;
461  shortestLifeTime = DBL_MAX;
462 
463  unsigned int NofInactiveProc=0;
464 
465  for( size_t ri=0 ; ri < fpProcessInfo->MAXofAtRestLoops ; ri++ )
466  {
467  fpCurrentProcess = (G4VITProcess*) (*fpProcessInfo->fpAtRestGetPhysIntVector)[ri];
468  if (fpCurrentProcess== 0)
469  {
470  (fpState->fSelectedAtRestDoItVector)[ri] = InActivated;
471  NofInactiveProc++;
472  continue;
473  } // NULL means the process is inactivated by a user on fly.
474 
475  fCondition=NotForced;
476  fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess->GetProcessID()));
477  lifeTime = fpCurrentProcess->AtRestGPIL( *fpTrack, &fCondition );
478  fpCurrentProcess->SetProcessState(0);
479 
480  if(fCondition==Forced)
481  {
482  (fpState->fSelectedAtRestDoItVector)[ri] = Forced;
483  }
484  else
485  {
486  (fpState->fSelectedAtRestDoItVector)[ri] = InActivated;
487  if(lifeTime < shortestLifeTime )
488  {
489  shortestLifeTime = lifeTime;
490  fAtRestDoItProcTriggered = G4int(ri);
491  (fpState->fSelectedAtRestDoItVector)[fAtRestDoItProcTriggered] = NotForced;
492  }
493  }
494  }
495 
496  fTimeStep = shortestLifeTime ;
497 
498  // at least one process is necessary to destroy the particle
499  // exit with warning
500  if(NofInactiveProc==fpProcessInfo->MAXofAtRestLoops)
501  {
502  G4cerr << "ERROR - G4ITStepProcessor::InvokeAtRestDoItProcs()" << G4endl
503  << " No AtRestDoIt process is active!" << G4endl;
504  }
505 }
506 //___________________________________________________________________________
507 
509 {
510  SetTrack(track);
512 }
513 //______________________________________________________________________________
514 
515 
517 {
518  // DEBUG
519  // G4cout << "SetInitialStep for : " << fpITrack-> GetName() << G4endl;
520  //________________________________________________________
521  // Initialize geometry
522 
523 
524  if ( ! fpTrack->GetTouchableHandle())
525  {
526  G4ThreeVector direction= fpTrack->GetMomentumDirection();
527  fpNavigator->LocateGlobalPointAndSetup( fpTrack->GetPosition(),
528  &direction, false, false );
529  fpState->fTouchableHandle = fpNavigator->CreateTouchableHistory();
530 
531  fpTrack->SetTouchableHandle( fpState->fTouchableHandle );
532  fpTrack->SetNextTouchableHandle( fpState->fTouchableHandle );
533  }
534  else
535  {
536  fpState->fTouchableHandle = fpTrack->GetTouchableHandle();
537  fpTrack->SetNextTouchableHandle( fpState->fTouchableHandle );
538  G4VPhysicalVolume* oldTopVolume= fpTrack->GetTouchableHandle()->GetVolume();
539  G4VPhysicalVolume* newTopVolume=
540  fpNavigator->ResetHierarchyAndLocate( fpTrack->GetPosition(),
541  fpTrack->GetMomentumDirection(),
542  *((G4TouchableHistory*)fpTrack->GetTouchableHandle()()) );
543  if(newTopVolume != oldTopVolume || oldTopVolume->GetRegularStructureId() == 1 )
544  {
545  fpState->fTouchableHandle = fpNavigator->CreateTouchableHistory();
546  fpTrack->SetTouchableHandle( fpState->fTouchableHandle );
547  fpTrack->SetNextTouchableHandle( fpState->fTouchableHandle );
548  }
549  }
550 
551  fpCurrentVolume = fpState->fTouchableHandle->GetVolume();
552 
553  //________________________________________________________
554  // If the primary track has 'Suspend' or 'PostponeToNextEvent' state,
555  // set the track state to 'Alive'.
556  if( (fpTrack->GetTrackStatus()==fSuspend) ||
557  (fpTrack->GetTrackStatus()==fPostponeToNextEvent) )
558  {
559  fpTrack->SetTrackStatus(fAlive);
560  }
561 
562  // If the primary track has 'zero' kinetic energy, set the track
563  // state to 'StopButAlive'.
564  if(fpTrack->GetKineticEnergy() <= 0.0)
565  {
566  fpTrack->SetTrackStatus( fStopButAlive );
567  }
568  //________________________________________________________
569  // Set vertex information of G4Track at here
570  if ( fpTrack->GetCurrentStepNumber() == 0 )
571  {
572  fpTrack->SetVertexPosition( fpTrack->GetPosition() );
573  fpTrack->SetVertexMomentumDirection( fpTrack->GetMomentumDirection() );
574  fpTrack->SetVertexKineticEnergy( fpTrack->GetKineticEnergy() );
575  fpTrack->SetLogicalVolumeAtVertex( fpTrack->GetVolume()->GetLogicalVolume() );
576  }
577  //________________________________________________________
578  // If track is already outside the world boundary, kill it
579  if( fpCurrentVolume==0 )
580  {
581  // If the track is a primary, stop processing
582  if(fpTrack->GetParentID()==0)
583  {
584  G4cerr << "ERROR - G4ITStepProcessor::SetInitialStep()" << G4endl
585  << " Primary particle starting at - "
586  << fpTrack->GetPosition()
587  << " - is outside of the world volume." << G4endl;
588  G4Exception("G4ITStepProcessor::SetInitialStep()", "ITStepProcessor0011",
589  FatalException, "Primary vertex outside of the world!");
590  }
591 
592  fpTrack->SetTrackStatus( fStopAndKill );
593  G4cout << "WARNING - G4ITStepProcessor::SetInitialStep()" << G4endl
594  << " Initial track position is outside world! - "
595  << fpTrack->GetPosition() << G4endl;
596  }
597  else{
598  // Initial set up for attribues of 'Step'
599  fpStep->InitializeStep( fpTrack );
600  }
601 
602 
603  if( fpTrack->GetTrackStatus() == fStopAndKill ) return ;
604 
605  fpTrackingManager->StartTracking(fpTrack);
606 
607  fpState->fStepStatus = fUndefined;
608 }
609 //______________________________________________________________________________
610 
612 {
613 
614  if(!fpStep)
615  {
616 
617  // Create new Step and give it to the track
618  fpStep = new G4Step();
619  fpTrack->SetStep(fpStep);
620  fpSecondary = fpStep->NewSecondaryVector();
621 
622  // Create new state and set it in the trackingInfo
623  fpState = new G4ITStepProcessorState();
625 
626  SetupMembers();
627  fpNavigator->NewNavigatorState();
628 
629  SetInitialStep();
630  }
631  else
632  {
633  SetupMembers();
634 
635  fpState->fPreviousStepSize = fpTrack->GetStepLength();
636 /***
637  // Send G4Step information to Hit/Dig if the volume is sensitive
638  fpCurrentVolume = fpStep->GetPreStepPoint()->GetPhysicalVolume();
639  StepControlFlag = fpStep->GetControlFlag();
640  if( fpCurrentVolume != 0 && StepControlFlag != AvoidHitInvocation)
641  {
642  fpSensitive = fpStep->GetPreStepPoint()->
643  GetSensitiveDetector();
644 
645  // if( fSensitive != 0 ) {
646  // fSensitive->Hit(fStep);
647  // }
648  }
649 ***/
650  // Store last PostStepPoint to PreStepPoint, and swap current and next
651  // volume information of G4Track. Reset total energy deposit in one Step.
652  fpStep->CopyPostToPreStepPoint();
653  fpStep->ResetTotalEnergyDeposit();
654 
655  //JA Set the volume before it is used (in DefineStepLength() for User Limit)
656  fpCurrentVolume = fpStep->GetPreStepPoint()->GetPhysicalVolume();
657 /*
658  G4cout << G4endl;
659  G4cout << "!!!!!!!!!!!!!!!!!!!!!!!!!!!" << G4endl;
660  G4cout << "PreStepPoint Volume : " << fpCurrentVolume->GetName() << G4endl;
661  G4cout << "Track Touchable : " << fpTrack->GetTouchableHandle()->GetVolume()->GetName() << G4endl;
662  G4cout << "Track NextTouchable : " << fpTrack->GetNextTouchableHandle()->GetVolume()->GetName() << G4endl;
663 */
664  // Reset the step's auxiliary points vector pointer
666 
667  // Switch next touchable in track to current one
668  fpTrack->SetTouchableHandle(fpTrack->GetNextTouchableHandle());
669  fpState->fTouchableHandle = fpTrack->GetTouchableHandle();
670  fpTrack->SetNextTouchableHandle( fpState->fTouchableHandle );
671  G4VPhysicalVolume* oldTopVolume= fpTrack->GetTouchableHandle()->GetVolume();
672  fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()->GetNavigatorState());
673 
674  G4VPhysicalVolume* newTopVolume=
675  fpNavigator->ResetHierarchyAndLocate( fpTrack->GetPosition(),
676  fpTrack->GetMomentumDirection(),
677  *((G4TouchableHistory*)fpTrack->GetTouchableHandle()()) );
678 
679  // G4cout << "New Top Volume : " << newTopVolume->GetName() << G4endl;
680 
681  if(newTopVolume != oldTopVolume || oldTopVolume->GetRegularStructureId() == 1 )
682  {
683  fpState->fTouchableHandle = fpNavigator->CreateTouchableHistory();
684  fpTrack->SetTouchableHandle( fpState->fTouchableHandle );
685  fpTrack->SetNextTouchableHandle( fpState->fTouchableHandle );
686  }
687 
688  fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()->GetNavigatorState());
689  }
690 }
691 
692 //______________________________________________________________________________
693 
694 
695 // ************************************************************************
696 // Compute Interaction Length
697 // ************************************************************************
699 {
700 
701  InitDefineStep();
702 
703  G4TrackStatus trackStatus = fpTrack -> GetTrackStatus() ;
704 
705  if(trackStatus == fStopAndKill)
706  {
707  return ;
708  }
709 
710  if(trackStatus == fStopButAlive)
711  {
712  fpITrack->GetTrackingInfo()->SetNavigatorState(fpNavigator->GetNavigatorState());
713  fpNavigator->SetNavigatorState(0);
714  return GetAtRestIL() ;
715  }
716 
717 
718  // Find minimum Step length and corresponding time
719  // demanded by active disc./cont. processes
720 
721  // ReSet the counter etc.
722  fpState->fPhysicalStep = DBL_MAX; // Initialize by a huge number
723  fPhysIntLength = DBL_MAX; // Initialize by a huge number
724 
725  double proposedTimeStep = DBL_MAX;
726  G4VProcess* processWithPostStepGivenByTimeStep(0);
727 
728  // GPIL for PostStep
729  fPostStepDoItProcTriggered = fpProcessInfo->MAXofPostStepLoops;
730  fPostStepAtTimeDoItProcTriggered = fpProcessInfo->MAXofPostStepLoops;
731 
732  // G4cout << "fpProcessInfo->MAXofPostStepLoops : " << fpProcessInfo->MAXofPostStepLoops
733  // << " mol : " << fpITrack -> GetName() << " id : " << fpTrack->GetTrackID()
734  // << G4endl;
735 
736  for(size_t np=0; np < fpProcessInfo->MAXofPostStepLoops; np++)
737  {
738  fpCurrentProcess = (G4VITProcess*) (*fpProcessInfo->fpPostStepGetPhysIntVector)[np];
739  if (fpCurrentProcess== 0)
740  {
741  (fpState->fSelectedPostStepDoItVector)[np] = InActivated;
742  continue;
743  } // NULL means the process is inactivated by a user on fly.
744 
745  fCondition=NotForced;
746  fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess->GetProcessID()));
747 
748  // G4cout << "Is going to call : " << fpCurrentProcess -> GetProcessName() << G4endl;
749  fPhysIntLength = fpCurrentProcess->
750  PostStepGPIL( *fpTrack,
751  fpState->fPreviousStepSize,
752  &fCondition );
753  fpCurrentProcess->SetProcessState(0);
754 
755  switch (fCondition)
756  {
757  case ExclusivelyForced: // Will need special treatment
758  (fpState->fSelectedPostStepDoItVector)[np] = ExclusivelyForced;
759  fpState->fStepStatus = fExclusivelyForcedProc;
760  fpStep->GetPostStepPoint()
761  ->SetProcessDefinedStep(fpCurrentProcess);
762  break;
763 
764  case Conditionally:
765  // (fpState->fSelectedPostStepDoItVector)[np] = Conditionally;
766  G4Exception("G4ITStepProcessor::DefinePhysicalStepLength()", "ITStepProcessor0008",
767  FatalException, "This feature is no more supported");
768  break;
769 
770  case Forced:
771  (fpState->fSelectedPostStepDoItVector)[np] = Forced;
772  break;
773 
774  case StronglyForced:
775  (fpState->fSelectedPostStepDoItVector)[np] = StronglyForced;
776  break;
777 
778  default:
779  (fpState->fSelectedPostStepDoItVector)[np] = InActivated;
780  break;
781  }
782 
783  if (fCondition==ExclusivelyForced)
784  {
785  for(size_t nrest=np+1; nrest < fpProcessInfo->MAXofPostStepLoops; nrest++)
786  {
787  (fpState->fSelectedPostStepDoItVector)[nrest] = InActivated;
788  }
789  return; // Please note the 'return' at here !!!
790  }
791  else
792  {
793  if(fPhysIntLength < fpState->fPhysicalStep )
794  {
795  // To avoid checking whether the process is actually
796  // proposing a time step, the returned time steps are
797  // negative (just for tagging)
798  if(fpCurrentProcess->ProposesTimeStep())
799  {
800  fPhysIntLength *= -1;
801  if(fPhysIntLength < proposedTimeStep)
802  {
803  proposedTimeStep = fPhysIntLength;
804  fPostStepAtTimeDoItProcTriggered = np;
805  processWithPostStepGivenByTimeStep = fpCurrentProcess;
806  }
807  }
808  else
809  {
810  fpState->fPhysicalStep = fPhysIntLength;
811  fpState->fStepStatus = fPostStepDoItProc;
812  fPostStepDoItProcTriggered = G4int(np);
813  fpStep->GetPostStepPoint()
814  ->SetProcessDefinedStep(fpCurrentProcess);
815  }
816  }
817  }
818  }
819 
820  // GPIL for AlongStep
821  fpState->proposedSafety = DBL_MAX;
822  G4double safetyProposedToAndByProcess = fpState->proposedSafety;
823 
824  for(size_t kp=0; kp < fpProcessInfo->MAXofAlongStepLoops; kp++)
825  {
826  fpCurrentProcess = (G4VITProcess*) (*fpProcessInfo->fpAlongStepGetPhysIntVector)[kp];
827  if (fpCurrentProcess== 0) continue;
828  // NULL means the process is inactivated by a user on fly.
829 
830  fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess->GetProcessID()));
831  fPhysIntLength = fpCurrentProcess-> AlongStepGPIL( *fpTrack,
832  fpState->fPreviousStepSize,
833  fpState->fPhysicalStep,
834  safetyProposedToAndByProcess,
835  &fGPILSelection );
836 
837  if(fPhysIntLength < fpState->fPhysicalStep)
838  {
839  fpState->fPhysicalStep = fPhysIntLength;
840  // Should save PS and TS in IT
841 
842  // Check if the process wants to be the GPIL winner. For example,
843  // multi-scattering proposes Step limit, but won't be the winner.
844  if(fGPILSelection==CandidateForSelection)
845  {
846  fpState->fStepStatus = fAlongStepDoItProc;
847  fpStep->GetPostStepPoint()
848  ->SetProcessDefinedStep(fpCurrentProcess);
849  }
850 
851  // Transportation is assumed to be the last process in the vector
852  if(kp == fpProcessInfo->MAXofAlongStepLoops-1)
853  {
854  fpTransportation = dynamic_cast<G4ITTransportation*>(fpCurrentProcess);
855 
856  if(! fpTransportation)
857  {
858  G4ExceptionDescription exceptionDescription ;
859  exceptionDescription << "No transportation process found " ;
860  G4Exception("G4ITStepProcessor::DoDefinePhysicalStepLength","ITStepProcessor0009",
861  FatalErrorInArgument,exceptionDescription);
862  }
863 
864  fTimeStep = fpTransportation->GetInteractionTimeLeft();
865 
866 
867  if (fpTrack->GetNextVolume() != 0)
868  fpState->fStepStatus = fGeomBoundary;
869  else
870  fpState->fStepStatus = fWorldBoundary;
871  }
872  }
873  else
874  {
875  if(kp == fpProcessInfo->MAXofAlongStepLoops-1)
876  {
877  fpTransportation = dynamic_cast<G4ITTransportation*>(fpCurrentProcess);
878 
879  if(! fpTransportation)
880  {
881  G4ExceptionDescription exceptionDescription ;
882  exceptionDescription << "No transportation process found " ;
883  G4Exception("G4ITStepProcessor::DoDefinePhysicalStepLength","ITStepProcessor0010",
884  FatalErrorInArgument,exceptionDescription);
885  }
886 
887  fTimeStep = fpTransportation->GetInteractionTimeLeft();
888  }
889  }
890 
891  if(proposedTimeStep < fTimeStep)
892  {
893  if(fPostStepAtTimeDoItProcTriggered<fpProcessInfo->MAXofPostStepLoops)
894  {
895  if ((fpState->fSelectedPostStepDoItVector)[fPostStepAtTimeDoItProcTriggered] ==
896  InActivated)
897  {
898  (fpState->fSelectedPostStepDoItVector)[fPostStepAtTimeDoItProcTriggered] = NotForced;
899  // (fpState->fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] = InActivated;
900 
901  fpState->fStepStatus = fPostStepDoItProc;
902  fpStep->GetPostStepPoint()->SetProcessDefinedStep(processWithPostStepGivenByTimeStep);
903 
904  fTimeStep = proposedTimeStep;
905 
906  fpTransportation->ComputeStep(*fpTrack,*fpStep,fTimeStep,fpState->fPhysicalStep);
907  }
908  }
909  }
910  else
911  {
912  if (fPostStepDoItProcTriggered<fpProcessInfo->MAXofPostStepLoops)
913  {
914  if ((fpState->fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] ==
915  InActivated)
916  {
917  (fpState->fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] =
918  NotForced;
919  }
920  }
921  }
922 
923  fpCurrentProcess->SetProcessState(0);
924 
925  // Make sure to check the safety, even if Step is not limited
926  // by this process. J. Apostolakis, June 20, 1998
927  //
928  if (safetyProposedToAndByProcess < fpState->proposedSafety)
929  // proposedSafety keeps the smallest value:
930  fpState->proposedSafety = safetyProposedToAndByProcess;
931  else
932  // safetyProposedToAndByProcess always proposes a valid safety:
933  safetyProposedToAndByProcess = fpState->proposedSafety;
934 
935  }
936 
937  fpITrack->GetTrackingInfo()->SetNavigatorState(fpNavigator->GetNavigatorState());
938  fpNavigator->SetNavigatorState(0);
939 }
940 
941 //______________________________________________________________________________