G4ITStepProcessor.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id: G4ITStepProcessor.cc 65022 2012-11-12 16:43:12Z gcosmo $
00027 //
00028 // Author: Mathieu Karamitros (kara (AT) cenbg . in2p3 . fr) 
00029 //
00030 // History:
00031 // -----------
00032 // 10 Oct 2011 M.Karamitros created
00033 //
00034 // -------------------------------------------------------------------
00035 
00036 #include "G4ITStepProcessor.hh"
00037 #include "G4UImanager.hh"
00038 #include "G4ForceCondition.hh"
00039 #include "G4GPILSelection.hh"
00040 #include "G4ITTransportationManager.hh"
00041 // #include "G4VSensitiveDetector.hh"    // Include from 'hits/digi'
00042 #include "G4GeometryTolerance.hh"
00043 #include "G4ParticleTable.hh"
00044 #include "G4ITTrackingManager.hh"
00045 #include "G4TrackingInformation.hh"
00046 #include "G4IT.hh"
00047 #include "G4ITNavigator.hh"             // Include from 'geometry'
00048 
00049 #include "G4VITProcess.hh"
00050 #include "G4VProcess.hh"
00051 #include "G4ITTransportation.hh"
00052 
00053 #include <iomanip>              // Include from 'system'
00054 #include <vector>               // Include from 'system'
00055 
00056 using namespace std;
00057 
00058 static const size_t SizeOfSelectedDoItVector=100;
00059 static const size_t& gMaxNProcesses(G4VITProcess::GetMaxProcessIndex());
00060 
00061 //____________________________________________________________________________________
00062 
00063 G4ITStepProcessor::G4ITStepProcessor()
00064 {
00065     verboseLevel = 0 ;
00066     //    fpUserSteppingAction = 0 ;
00067     fStoreTrajectory = 0;
00068     fpTrackingManager = 0;
00069     fpNavigator = 0;
00070     kCarTolerance = -1.;
00071     fInitialized = false;
00072     fPreviousTimeStep = DBL_MAX;
00073     CleanProcessor();
00074     ResetSecondaries();
00075 }
00076 
00077 G4ITStepProcessor::G4ITStepProcessorState::G4ITStepProcessorState() :
00078     G4ITStepProcessorState_Lock(),
00079     fSelectedAtRestDoItVector (gMaxNProcesses,0),
00080     fSelectedPostStepDoItVector (gMaxNProcesses,0)
00081 {
00082     fPhysicalStep = -1.;
00083     fPreviousStepSize = -1.;
00084 
00085     fSafety = -1.;
00086     proposedSafety = -1.;
00087     endpointSafety = -1;
00088 
00089     fStepStatus = fUndefined;
00090 
00091     fTouchableHandle = 0;
00092 }
00093 
00094 // should not be used
00095 G4ITStepProcessor::G4ITStepProcessorState::G4ITStepProcessorState(const G4ITStepProcessorState& ) :
00096     G4ITStepProcessorState_Lock(),
00097     fSelectedAtRestDoItVector (gMaxNProcesses,0),
00098     fSelectedPostStepDoItVector (gMaxNProcesses,0)
00099 {
00100     fPhysicalStep = -1.;
00101     fPreviousStepSize = -1.;
00102 
00103     fSafety = -1.;
00104     proposedSafety = -1.;
00105     endpointSafety = -1;
00106 
00107     fStepStatus = fUndefined;
00108 
00109     fTouchableHandle = 0;
00110 }
00111 
00112 // should not be used
00113 G4ITStepProcessor::G4ITStepProcessorState&  G4ITStepProcessor::G4ITStepProcessorState::operator=(const G4ITStepProcessorState& rhs)
00114 {
00115     if(this == &rhs) return *this;
00116 
00117     fSelectedAtRestDoItVector.clear();
00118     fSelectedAtRestDoItVector.resize(gMaxNProcesses,0);
00119     fSelectedPostStepDoItVector.clear();
00120     fSelectedPostStepDoItVector.resize(gMaxNProcesses,0);
00121 
00122     fPhysicalStep = -1.;
00123     fPreviousStepSize = -1.;
00124 
00125     fSafety = -1.;
00126     proposedSafety = -1.;
00127     endpointSafety = -1;
00128 
00129     fStepStatus = fUndefined;
00130 
00131     fTouchableHandle = 0;
00132     return *this;
00133 }
00134 //____________________________________________________________________________________
00135 
00136 G4ITStepProcessor::G4ITStepProcessorState::~G4ITStepProcessorState()
00137 {;}
00138 //____________________________________________________________________________________
00139 
00140 void G4ITStepProcessor::ClearProcessInfo()
00141 {
00142     std::map<const G4ParticleDefinition*, ProcessGeneralInfo*> ::iterator it;
00143 
00144     for(it = fProcessGeneralInfoMap.begin();it != fProcessGeneralInfoMap.end();it++)
00145     {
00146         if(it->second)
00147         {
00148             delete it->second;
00149             it->second = 0;
00150         }
00151     }
00152 
00153     fProcessGeneralInfoMap.clear();
00154 }
00155 
00156 //____________________________________________________________________________________
00157 
00158 void G4ITStepProcessor::ForceReInitialization()
00159 {
00160     fInitialized = false;
00161     ClearProcessInfo();
00162     Initialize();
00163 }
00164 
00165 //____________________________________________________________________________________
00166 
00167 void G4ITStepProcessor::Initialize()
00168 {
00169     CleanProcessor();
00170     if(fInitialized) return;
00171     //    ActiveOnlyITProcess();
00172 
00173     SetNavigator(G4ITTransportationManager::GetTransportationManager()
00174                  ->GetNavigatorForTracking());
00175 
00176     fPhysIntLength = DBL_MAX;
00177     kCarTolerance = 0.5*G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
00178 
00179     fInitialized = true;
00180 }
00181 //______________________________________________________________________________
00182 
00183 G4ITStepProcessor::~G4ITStepProcessor()
00184 {
00185     if(fpStep)
00186     {
00187         fpStep->DeleteSecondaryVector();
00188         delete fpStep;
00189     }
00190 
00191     if(fpSecondary)                      delete fpSecondary;
00192     ClearProcessInfo();
00193     G4ITTransportationManager::DeleteInstance();
00194 
00195     //    if(fpUserSteppingAction)             delete fpUserSteppingAction;
00196 }
00197 //______________________________________________________________________________
00198 // should not be used
00199 G4ITStepProcessor::G4ITStepProcessor(const G4ITStepProcessor& rhs)
00200 {
00201     verboseLevel = rhs.verboseLevel ;
00202     fStoreTrajectory = rhs.fStoreTrajectory ;
00203 
00204     //    fpUserSteppingAction = 0 ;
00205     fpTrackingManager = 0;
00206     fpNavigator = 0;
00207     fInitialized = false;
00208 
00209     kCarTolerance = rhs.kCarTolerance;
00210     fInitialized = false;
00211     fPreviousTimeStep = DBL_MAX;
00212 
00213     CleanProcessor();
00214     ResetSecondaries();
00215 }
00216 //______________________________________________________________________________
00217 
00218 G4ITStepProcessor& G4ITStepProcessor::operator=(const G4ITStepProcessor& rhs)
00219 {
00220     if (this == &rhs) return *this; // handle self assignment
00221     //assignment operator
00222     return *this;
00223 }
00224 // ******************************************************************
00225 
00226 void G4ITStepProcessor::ActiveOnlyITProcess()
00227 {
00228     // Method not used for the time being
00229 #ifdef debug
00230     G4cout<<"G4ITStepProcessor::CloneProcesses: is called"<<G4endl;
00231 #endif
00232 
00233     G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
00234     G4ParticleTable::G4PTblDicIterator* theParticleIterator = theParticleTable->GetIterator();
00235 
00236     theParticleIterator->reset();
00237     // TODO : Ne faire la boucle que sur les IT **** !!!
00238     while( (*theParticleIterator)() )
00239     {
00240         G4ParticleDefinition* particle = theParticleIterator->value();
00241         G4ProcessManager* pm= particle->GetProcessManager();
00242 
00243         if(!pm)
00244         {
00245             G4cerr << "ERROR - G4ITStepProcessor::GetProcessNumber()" << G4endl
00246                    << "        ProcessManager is NULL for particle = "
00247                    << particle->GetParticleName() << ", PDG_code = "
00248                    << particle->GetPDGEncoding() << G4endl;
00249             G4Exception("G4ITStepProcessor::GetProcessNumber()", "ITStepProcessor0001",
00250                         FatalException, "Process Manager is not found.");
00251             return;
00252         }
00253 
00254         ActiveOnlyITProcess(pm);
00255     }
00256 }
00257 // ******************************************************************
00258 
00259 void G4ITStepProcessor::ActiveOnlyITProcess(G4ProcessManager* processManager)
00260 {
00261     // Method not used for the time being
00262     G4ProcessVector* processVector = processManager->GetProcessList();
00263 
00264     G4VITProcess* itProcess = 0 ;
00265     for(int i = 0 ; i < processVector->size() ; i++)
00266     {
00267         G4VProcess* base_process = (*processVector)[i];
00268         itProcess = dynamic_cast<G4VITProcess*>(base_process);
00269 
00270         if(!itProcess)
00271         {
00272             processManager->SetProcessActivation(base_process, false);
00273         }
00274     }
00275 }
00276 // ******************************************************************
00277 void G4ITStepProcessor::SetupGeneralProcessInfo(G4ParticleDefinition* particle,
00278                                                 G4ProcessManager* pm)
00279 {
00280 
00281 #ifdef debug
00282     G4cout<<"G4ITStepProcessor::GetProcessNumber: is called track"<<G4endl;
00283 #endif
00284     if(!pm)
00285     {
00286         G4cerr << "ERROR - G4SteppingManager::GetProcessNumber()" << G4endl
00287                << "        ProcessManager is NULL for particle = "
00288                << particle->GetParticleName() << ", PDG_code = "
00289                << particle->GetPDGEncoding() << G4endl;
00290         G4Exception("G4SteppingManager::GetProcessNumber()", "ITStepProcessor0002",
00291                     FatalException, "Process Manager is not found.");
00292         return;
00293     }
00294 
00295     std::map<const G4ParticleDefinition*, ProcessGeneralInfo*>::iterator it = fProcessGeneralInfoMap.find(particle);
00296     if(it != fProcessGeneralInfoMap.end())
00297     {
00298         G4Exception("G4SteppingManager::SetupGeneralProcessInfo()", "ITStepProcessor0003",
00299                     FatalException, "Process info already registered.");
00300         return;
00301     }
00302 
00303     // here used as temporary
00304     fpProcessInfo = new ProcessGeneralInfo();
00305 
00306     // AtRestDoits
00307     fpProcessInfo->MAXofAtRestLoops =        pm->GetAtRestProcessVector()->entries();
00308     fpProcessInfo->fpAtRestDoItVector =       pm->GetAtRestProcessVector(typeDoIt);
00309     fpProcessInfo->fpAtRestGetPhysIntVector = pm->GetAtRestProcessVector(typeGPIL);
00310 #ifdef debug
00311     G4cout << "G4ITStepProcessor::GetProcessNumber: #ofAtRest="
00312            << fpProcessInfo->MAXofAtRestLoops << G4endl;
00313 #endif
00314 
00315     // AlongStepDoits
00316     fpProcessInfo->MAXofAlongStepLoops = pm->GetAlongStepProcessVector()->entries();
00317     fpProcessInfo->fpAlongStepDoItVector = pm->GetAlongStepProcessVector(typeDoIt);
00318     fpProcessInfo->fpAlongStepGetPhysIntVector = pm->GetAlongStepProcessVector(typeGPIL);
00319 #ifdef debug
00320     G4cout << "G4ITStepProcessor::GetProcessNumber:#ofAlongStp="
00321            << fpProcessInfo->MAXofAlongStepLoops << G4endl;
00322 #endif
00323 
00324     // PostStepDoits
00325     fpProcessInfo->MAXofPostStepLoops = pm->GetPostStepProcessVector()->entries();
00326     fpProcessInfo->fpPostStepDoItVector = pm->GetPostStepProcessVector(typeDoIt);
00327     fpProcessInfo->fpPostStepGetPhysIntVector = pm->GetPostStepProcessVector(typeGPIL);
00328 #ifdef debug
00329     G4cout << "G4ITStepProcessor::GetProcessNumber: #ofPostStep="
00330            << fpProcessInfo->MAXofPostStepLoops << G4endl;
00331 #endif
00332 
00333     if (SizeOfSelectedDoItVector<fpProcessInfo->MAXofAtRestLoops    ||
00334             SizeOfSelectedDoItVector<fpProcessInfo->MAXofAlongStepLoops ||
00335             SizeOfSelectedDoItVector<fpProcessInfo->MAXofPostStepLoops  )
00336     {
00337         G4cerr << "ERROR - G4ITStepProcessor::GetProcessNumber()" << G4endl
00338                << "        SizeOfSelectedDoItVector= " << SizeOfSelectedDoItVector
00339                << " ; is smaller then one of MAXofAtRestLoops= "
00340                << fpProcessInfo->MAXofAtRestLoops << G4endl
00341                << "        or MAXofAlongStepLoops= " << fpProcessInfo->MAXofAlongStepLoops
00342                << " or MAXofPostStepLoops= " << fpProcessInfo->MAXofPostStepLoops << G4endl;
00343         G4Exception("G4ITStepProcessor::GetProcessNumber()",
00344                     "ITStepProcessor0004", FatalException,
00345                     "The array size is smaller than the actual No of processes.");
00346     }
00347 
00348     if(!fpProcessInfo->fpAtRestDoItVector       &&
00349             !fpProcessInfo->fpAlongStepDoItVector    &&
00350             !fpProcessInfo->fpPostStepDoItVector)
00351     {
00352         G4ExceptionDescription exceptionDescription ;
00353         exceptionDescription << "No DoIt process found " ;
00354         G4Exception("G4ITStepProcessor::DoStepping","ITStepProcessor0005",
00355                     FatalErrorInArgument,exceptionDescription);
00356         return ;
00357     }
00358 
00359     if(fpProcessInfo->fpAlongStepGetPhysIntVector && fpProcessInfo->MAXofAlongStepLoops>0)
00360     {
00361         fpProcessInfo->fpTransportation = dynamic_cast<G4ITTransportation*>
00362                 ((*fpProcessInfo->fpAlongStepGetPhysIntVector)[fpProcessInfo->MAXofAlongStepLoops-1]);
00363 
00364         if(fpProcessInfo->fpTransportation == 0)
00365         {
00366             G4ExceptionDescription exceptionDescription ;
00367             exceptionDescription << "No transportation process found " ;
00368             G4Exception("G4ITStepProcessor::SetupGeneralProcessInfo","ITStepProcessor0006",
00369                         FatalErrorInArgument,exceptionDescription);
00370         }
00371     }
00372     fProcessGeneralInfoMap[particle] = fpProcessInfo;
00373     //    fpProcessInfo = 0;
00374 }
00375 
00376 // ******************************************************************
00377 
00378 void G4ITStepProcessor::SetTrack(G4Track* track)
00379 {
00380     fpTrack = track ;
00381     if(fpTrack)
00382     {
00383         fpITrack = GetIT(fpTrack) ;
00384         fpStep = const_cast<G4Step*>(fpTrack -> GetStep());
00385 
00386         if(fpITrack)
00387         {
00388             fpTrackingInfo = fpITrack->GetTrackingInfo() ;
00389         }
00390         else
00391         {
00392             fpTrackingInfo = 0;
00393             G4cerr << "Track ID : " << fpTrack->GetTrackID() << G4endl;
00394 
00395             G4ExceptionDescription exceptionDescription ("No IT pointer was attached to the track you try to process.");
00396             G4Exception("G4ITStepProcessor::SetTrack","ITStepProcessor0007",
00397                         FatalErrorInArgument,exceptionDescription);
00398         }
00399     }
00400     else
00401     {
00402         fpITrack = 0;
00403         fpStep = 0 ;
00404     }
00405 }
00406 //______________________________________________________________________________
00407 
00408 void G4ITStepProcessor::GetProcessInfo()
00409 {
00410     G4ParticleDefinition* particle = fpTrack->GetDefinition();
00411     std::map<const G4ParticleDefinition*, ProcessGeneralInfo*>::iterator it = fProcessGeneralInfoMap.find(particle);
00412 
00413     if(it == fProcessGeneralInfoMap.end())
00414     {
00415         SetupGeneralProcessInfo(particle,fpTrack->GetDefinition()->GetProcessManager());
00416         if(fpProcessInfo == 0)
00417         {
00418             G4ExceptionDescription exceptionDescription ("...");
00419             G4Exception("G4ITStepProcessor::GetProcessNumber","ITStepProcessor0008",
00420                         FatalErrorInArgument,exceptionDescription);
00421             return;
00422         }
00423     }
00424     else
00425     {
00426         fpProcessInfo = it->second;
00427     }
00428 }
00429 //______________________________________________________________________________
00430 
00431 void G4ITStepProcessor::SetupMembers()
00432 {
00433     fpSecondary      = fpStep->GetfSecondary();
00434     fpPreStepPoint   = fpStep->GetPreStepPoint();
00435     fpPostStepPoint  = fpStep->GetPostStepPoint();
00436 
00437     fpState = (G4ITStepProcessorState*) fpITrack->GetTrackingInfo()->GetStepProcessorState();
00438 
00439     GetProcessInfo();
00440     ResetSecondaries();
00441 }
00442 //______________________________________________________________________________
00443 
00444 void G4ITStepProcessor::ResetSecondaries()
00445 {
00446     // Reset the secondary particles
00447     fN2ndariesAtRestDoIt    = 0;
00448     fN2ndariesAlongStepDoIt = 0;
00449     fN2ndariesPostStepDoIt  = 0;
00450 }
00451 //______________________________________________________________________________
00452 
00453 void G4ITStepProcessor::GetAtRestIL()
00454 {
00455     // Select the rest process which has the shortest time before
00456     // it is invoked. In rest processes, GPIL()
00457     // returns the time before a process occurs.
00458     G4double lifeTime (DBL_MAX), shortestLifeTime (DBL_MAX);
00459 
00460     fAtRestDoItProcTriggered = 0;
00461     shortestLifeTime = DBL_MAX;
00462 
00463     unsigned int NofInactiveProc=0;
00464 
00465     for( size_t ri=0 ; ri < fpProcessInfo->MAXofAtRestLoops ; ri++ )
00466     {
00467         fpCurrentProcess = (G4VITProcess*) (*fpProcessInfo->fpAtRestGetPhysIntVector)[ri];
00468         if (fpCurrentProcess== 0)
00469         {
00470             (fpState->fSelectedAtRestDoItVector)[ri] = InActivated;
00471             NofInactiveProc++;
00472             continue;
00473         }   // NULL means the process is inactivated by a user on fly.
00474 
00475         fCondition=NotForced;
00476         fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess->GetProcessID()));
00477         lifeTime = fpCurrentProcess->AtRestGPIL( *fpTrack, &fCondition );
00478         fpCurrentProcess->SetProcessState(0);
00479 
00480         if(fCondition==Forced)
00481         {
00482             (fpState->fSelectedAtRestDoItVector)[ri] = Forced;
00483         }
00484         else
00485         {
00486             (fpState->fSelectedAtRestDoItVector)[ri] = InActivated;
00487             if(lifeTime < shortestLifeTime )
00488             {
00489                 shortestLifeTime = lifeTime;
00490                 fAtRestDoItProcTriggered =  G4int(ri);
00491                 (fpState->fSelectedAtRestDoItVector)[fAtRestDoItProcTriggered] = NotForced;
00492             }
00493         }
00494     }
00495 
00496     fTimeStep = shortestLifeTime ;
00497 
00498     // at least one process is necessary to destroy the particle
00499     // exit with warning
00500     if(NofInactiveProc==fpProcessInfo->MAXofAtRestLoops)
00501     {
00502         G4cerr << "ERROR - G4ITStepProcessor::InvokeAtRestDoItProcs()" << G4endl
00503                << "        No AtRestDoIt process is active!" << G4endl;
00504     }
00505 }
00506 //___________________________________________________________________________
00507 
00508 void G4ITStepProcessor::DefinePhysicalStepLength(G4Track* track)
00509 {
00510     SetTrack(track);
00511     DoDefinePhysicalStepLength();
00512 }
00513 //______________________________________________________________________________
00514 
00515 
00516 void G4ITStepProcessor::SetInitialStep()
00517 {
00518     // DEBUG
00519     //    G4cout << "SetInitialStep for : " << fpITrack-> GetName() << G4endl;
00520     //________________________________________________________
00521     // Initialize geometry
00522 
00523 
00524     if ( ! fpTrack->GetTouchableHandle())
00525     {
00526         G4ThreeVector direction= fpTrack->GetMomentumDirection();
00527         fpNavigator->LocateGlobalPointAndSetup( fpTrack->GetPosition(),
00528                                                 &direction, false, false );
00529         fpState->fTouchableHandle = fpNavigator->CreateTouchableHistory();
00530 
00531         fpTrack->SetTouchableHandle( fpState->fTouchableHandle );
00532         fpTrack->SetNextTouchableHandle( fpState->fTouchableHandle );
00533     }
00534     else
00535     {
00536         fpState->fTouchableHandle = fpTrack->GetTouchableHandle();
00537         fpTrack->SetNextTouchableHandle( fpState->fTouchableHandle );
00538         G4VPhysicalVolume* oldTopVolume= fpTrack->GetTouchableHandle()->GetVolume();
00539         G4VPhysicalVolume* newTopVolume=
00540                 fpNavigator->ResetHierarchyAndLocate( fpTrack->GetPosition(),
00541                                                       fpTrack->GetMomentumDirection(),
00542                                                       *((G4TouchableHistory*)fpTrack->GetTouchableHandle()()) );
00543         if(newTopVolume != oldTopVolume || oldTopVolume->GetRegularStructureId() == 1 )
00544         {
00545             fpState->fTouchableHandle = fpNavigator->CreateTouchableHistory();
00546             fpTrack->SetTouchableHandle( fpState->fTouchableHandle );
00547             fpTrack->SetNextTouchableHandle( fpState->fTouchableHandle );
00548         }
00549     }
00550 
00551     fpCurrentVolume = fpState->fTouchableHandle->GetVolume();
00552 
00553     //________________________________________________________
00554     // If the primary track has 'Suspend' or 'PostponeToNextEvent' state,
00555     // set the track state to 'Alive'.
00556     if( (fpTrack->GetTrackStatus()==fSuspend) ||
00557             (fpTrack->GetTrackStatus()==fPostponeToNextEvent) )
00558     {
00559         fpTrack->SetTrackStatus(fAlive);
00560     }
00561 
00562     // If the primary track has 'zero' kinetic energy, set the track
00563     // state to 'StopButAlive'.
00564     if(fpTrack->GetKineticEnergy() <= 0.0)
00565     {
00566         fpTrack->SetTrackStatus( fStopButAlive );
00567     }
00568     //________________________________________________________
00569     // Set vertex information of G4Track at here
00570     if ( fpTrack->GetCurrentStepNumber() == 0 )
00571     {
00572         fpTrack->SetVertexPosition( fpTrack->GetPosition() );
00573         fpTrack->SetVertexMomentumDirection( fpTrack->GetMomentumDirection() );
00574         fpTrack->SetVertexKineticEnergy( fpTrack->GetKineticEnergy() );
00575         fpTrack->SetLogicalVolumeAtVertex( fpTrack->GetVolume()->GetLogicalVolume() );
00576     }
00577     //________________________________________________________
00578     // If track is already outside the world boundary, kill it
00579     if( fpCurrentVolume==0 )
00580     {
00581         // If the track is a primary, stop processing
00582         if(fpTrack->GetParentID()==0)
00583         {
00584             G4cerr << "ERROR - G4ITStepProcessor::SetInitialStep()" << G4endl
00585                    << "        Primary particle starting at - "
00586                    << fpTrack->GetPosition()
00587                    << " - is outside of the world volume." << G4endl;
00588             G4Exception("G4ITStepProcessor::SetInitialStep()", "ITStepProcessor0011",
00589                         FatalException, "Primary vertex outside of the world!");
00590         }
00591 
00592         fpTrack->SetTrackStatus( fStopAndKill );
00593         G4cout << "WARNING - G4ITStepProcessor::SetInitialStep()" << G4endl
00594                << "          Initial track position is outside world! - "
00595                << fpTrack->GetPosition() << G4endl;
00596     }
00597     else{
00598         // Initial set up for attribues of 'Step'
00599         fpStep->InitializeStep( fpTrack );
00600     }
00601 
00602 
00603     if( fpTrack->GetTrackStatus() == fStopAndKill ) return ;
00604 
00605     fpTrackingManager->StartTracking(fpTrack);
00606 
00607     fpState->fStepStatus = fUndefined;
00608 }
00609 //______________________________________________________________________________
00610 
00611 void G4ITStepProcessor::InitDefineStep()
00612 {
00613 
00614     if(!fpStep)
00615     {
00616 
00617         // Create new Step and give it to the track
00618         fpStep = new G4Step();
00619         fpTrack->SetStep(fpStep);
00620         fpSecondary = fpStep->NewSecondaryVector();
00621 
00622         // Create new state and set it in the trackingInfo
00623         fpState = new G4ITStepProcessorState();
00624         fpITrack->GetTrackingInfo()->SetStepProcessorState((G4ITStepProcessorState_Lock*)fpState);
00625 
00626         SetupMembers();
00627         fpNavigator->NewNavigatorState();
00628 
00629         SetInitialStep();
00630     }
00631     else
00632     {
00633         SetupMembers();
00634 
00635         fpState->fPreviousStepSize = fpTrack->GetStepLength();
00636 /***
00637         // Send G4Step information to Hit/Dig if the volume is sensitive
00638         fpCurrentVolume = fpStep->GetPreStepPoint()->GetPhysicalVolume();
00639         StepControlFlag =  fpStep->GetControlFlag();
00640         if( fpCurrentVolume != 0 && StepControlFlag != AvoidHitInvocation)
00641         {
00642             fpSensitive = fpStep->GetPreStepPoint()->
00643                     GetSensitiveDetector();
00644 
00645             //            if( fSensitive != 0 ) {
00646             //              fSensitive->Hit(fStep);
00647             //            }
00648         }
00649 ***/
00650         // Store last PostStepPoint to PreStepPoint, and swap current and next
00651         // volume information of G4Track. Reset total energy deposit in one Step.
00652         fpStep->CopyPostToPreStepPoint();
00653         fpStep->ResetTotalEnergyDeposit();
00654 
00655         //JA Set the volume before it is used (in DefineStepLength() for User Limit)
00656         fpCurrentVolume = fpStep->GetPreStepPoint()->GetPhysicalVolume();
00657 /*
00658         G4cout << G4endl;
00659         G4cout << "!!!!!!!!!!!!!!!!!!!!!!!!!!!" << G4endl;
00660         G4cout << "PreStepPoint Volume : " << fpCurrentVolume->GetName() << G4endl;
00661         G4cout << "Track Touchable : " << fpTrack->GetTouchableHandle()->GetVolume()->GetName() << G4endl;
00662         G4cout << "Track NextTouchable : " << fpTrack->GetNextTouchableHandle()->GetVolume()->GetName() << G4endl;
00663 */
00664         // Reset the step's auxiliary points vector pointer
00665         fpStep->SetPointerToVectorOfAuxiliaryPoints(0);
00666 
00667         // Switch next touchable in track to current one
00668         fpTrack->SetTouchableHandle(fpTrack->GetNextTouchableHandle());
00669         fpState->fTouchableHandle = fpTrack->GetTouchableHandle();
00670         fpTrack->SetNextTouchableHandle( fpState->fTouchableHandle );
00671         G4VPhysicalVolume* oldTopVolume= fpTrack->GetTouchableHandle()->GetVolume();
00672         fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()->GetNavigatorState());
00673 
00674         G4VPhysicalVolume* newTopVolume=
00675                 fpNavigator->ResetHierarchyAndLocate( fpTrack->GetPosition(),
00676                                                       fpTrack->GetMomentumDirection(),
00677                                                       *((G4TouchableHistory*)fpTrack->GetTouchableHandle()()) );
00678 
00679         //        G4cout << "New Top Volume : " << newTopVolume->GetName() << G4endl;
00680 
00681         if(newTopVolume != oldTopVolume || oldTopVolume->GetRegularStructureId() == 1 )
00682         {
00683             fpState->fTouchableHandle = fpNavigator->CreateTouchableHistory();
00684             fpTrack->SetTouchableHandle( fpState->fTouchableHandle );
00685             fpTrack->SetNextTouchableHandle( fpState->fTouchableHandle );
00686         }
00687 
00688         fpNavigator->SetNavigatorState(fpITrack->GetTrackingInfo()->GetNavigatorState());
00689     }
00690 }
00691 
00692 //______________________________________________________________________________
00693 
00694 
00695 // ************************************************************************
00696 //      Compute Interaction Length
00697 // ************************************************************************
00698 void G4ITStepProcessor::DoDefinePhysicalStepLength()
00699 {
00700 
00701     InitDefineStep();
00702 
00703     G4TrackStatus trackStatus = fpTrack -> GetTrackStatus()  ;
00704 
00705     if(trackStatus == fStopAndKill)
00706     {
00707         return ;
00708     }
00709 
00710     if(trackStatus == fStopButAlive)
00711     {
00712         fpITrack->GetTrackingInfo()->SetNavigatorState(fpNavigator->GetNavigatorState());
00713         fpNavigator->SetNavigatorState(0);
00714         return GetAtRestIL() ;
00715     }
00716 
00717 
00718     // Find minimum Step length and corresponding time
00719     // demanded by active disc./cont. processes
00720 
00721     // ReSet the counter etc.
00722     fpState->fPhysicalStep  = DBL_MAX;          // Initialize by a huge number
00723     fPhysIntLength = DBL_MAX;          // Initialize by a huge number
00724 
00725     double proposedTimeStep = DBL_MAX;
00726     G4VProcess* processWithPostStepGivenByTimeStep(0);
00727 
00728     // GPIL for PostStep
00729     fPostStepDoItProcTriggered = fpProcessInfo->MAXofPostStepLoops;
00730     fPostStepAtTimeDoItProcTriggered = fpProcessInfo->MAXofPostStepLoops;
00731 
00732     //    G4cout << "fpProcessInfo->MAXofPostStepLoops : " << fpProcessInfo->MAXofPostStepLoops
00733     //           << " mol : " << fpITrack -> GetName() << " id : " << fpTrack->GetTrackID()
00734     //           << G4endl;
00735 
00736     for(size_t np=0; np < fpProcessInfo->MAXofPostStepLoops; np++)
00737     {
00738         fpCurrentProcess = (G4VITProcess*) (*fpProcessInfo->fpPostStepGetPhysIntVector)[np];
00739         if (fpCurrentProcess== 0)
00740         {
00741             (fpState->fSelectedPostStepDoItVector)[np] = InActivated;
00742             continue;
00743         }   // NULL means the process is inactivated by a user on fly.
00744 
00745         fCondition=NotForced;
00746         fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess->GetProcessID()));
00747 
00748         //        G4cout << "Is going to call : " << fpCurrentProcess -> GetProcessName() << G4endl;
00749         fPhysIntLength = fpCurrentProcess->
00750                 PostStepGPIL( *fpTrack,
00751                               fpState->fPreviousStepSize,
00752                               &fCondition );
00753         fpCurrentProcess->SetProcessState(0);
00754 
00755         switch (fCondition)
00756         {
00757         case ExclusivelyForced: // Will need special treatment
00758             (fpState->fSelectedPostStepDoItVector)[np] = ExclusivelyForced;
00759             fpState->fStepStatus = fExclusivelyForcedProc;
00760             fpStep->GetPostStepPoint()
00761                     ->SetProcessDefinedStep(fpCurrentProcess);
00762             break;
00763 
00764         case Conditionally:
00765             //       (fpState->fSelectedPostStepDoItVector)[np] = Conditionally;
00766             G4Exception("G4ITStepProcessor::DefinePhysicalStepLength()", "ITStepProcessor0008",
00767                         FatalException, "This feature is no more supported");
00768             break;
00769 
00770         case Forced:
00771             (fpState->fSelectedPostStepDoItVector)[np] = Forced;
00772             break;
00773 
00774         case StronglyForced:
00775             (fpState->fSelectedPostStepDoItVector)[np] = StronglyForced;
00776             break;
00777 
00778         default:
00779             (fpState->fSelectedPostStepDoItVector)[np] = InActivated;
00780             break;
00781         }
00782 
00783         if (fCondition==ExclusivelyForced)
00784         {
00785             for(size_t nrest=np+1; nrest < fpProcessInfo->MAXofPostStepLoops; nrest++)
00786             {
00787                 (fpState->fSelectedPostStepDoItVector)[nrest] = InActivated;
00788             }
00789             return;  // Please note the 'return' at here !!!
00790         }
00791         else
00792         {
00793             if(fPhysIntLength < fpState->fPhysicalStep )
00794             {
00795                 // To avoid checking whether the process is actually
00796                 // proposing a time step, the returned time steps are
00797                 // negative (just for tagging)
00798                 if(fpCurrentProcess->ProposesTimeStep())
00799                 {
00800                     fPhysIntLength *= -1;
00801                     if(fPhysIntLength < proposedTimeStep)
00802                     {
00803                         proposedTimeStep = fPhysIntLength;
00804                         fPostStepAtTimeDoItProcTriggered = np;
00805                         processWithPostStepGivenByTimeStep = fpCurrentProcess;
00806                     }
00807                 }
00808                 else
00809                 {
00810                     fpState->fPhysicalStep = fPhysIntLength;
00811                     fpState->fStepStatus = fPostStepDoItProc;
00812                     fPostStepDoItProcTriggered = G4int(np);
00813                     fpStep->GetPostStepPoint()
00814                             ->SetProcessDefinedStep(fpCurrentProcess);
00815                 }
00816             }
00817         }
00818     }
00819 
00820     // GPIL for AlongStep
00821     fpState->proposedSafety = DBL_MAX;
00822     G4double safetyProposedToAndByProcess = fpState->proposedSafety;
00823 
00824     for(size_t kp=0; kp < fpProcessInfo->MAXofAlongStepLoops; kp++)
00825     {
00826         fpCurrentProcess = (G4VITProcess*) (*fpProcessInfo->fpAlongStepGetPhysIntVector)[kp];
00827         if (fpCurrentProcess== 0) continue;
00828         // NULL means the process is inactivated by a user on fly.
00829 
00830         fpCurrentProcess->SetProcessState(fpTrackingInfo->GetProcessState(fpCurrentProcess->GetProcessID()));
00831         fPhysIntLength = fpCurrentProcess-> AlongStepGPIL( *fpTrack,
00832                                                            fpState->fPreviousStepSize,
00833                                                            fpState->fPhysicalStep,
00834                                                            safetyProposedToAndByProcess,
00835                                                            &fGPILSelection );
00836 
00837         if(fPhysIntLength < fpState->fPhysicalStep)
00838         {
00839             fpState->fPhysicalStep      = fPhysIntLength;
00840             // Should save PS and TS in IT
00841 
00842             // Check if the process wants to be the GPIL winner. For example,
00843             // multi-scattering proposes Step limit, but won't be the winner.
00844             if(fGPILSelection==CandidateForSelection)
00845             {
00846                 fpState->fStepStatus = fAlongStepDoItProc;
00847                 fpStep->GetPostStepPoint()
00848                         ->SetProcessDefinedStep(fpCurrentProcess);
00849             }
00850 
00851             // Transportation is assumed to be the last process in the vector
00852             if(kp == fpProcessInfo->MAXofAlongStepLoops-1)
00853             {
00854                 fpTransportation = dynamic_cast<G4ITTransportation*>(fpCurrentProcess);
00855 
00856                 if(! fpTransportation)
00857                 {
00858                     G4ExceptionDescription exceptionDescription ;
00859                     exceptionDescription << "No transportation process found " ;
00860                     G4Exception("G4ITStepProcessor::DoDefinePhysicalStepLength","ITStepProcessor0009",
00861                                 FatalErrorInArgument,exceptionDescription);
00862                 }
00863 
00864                 fTimeStep               = fpTransportation->GetInteractionTimeLeft();
00865 
00866 
00867                 if (fpTrack->GetNextVolume() != 0)
00868                     fpState->fStepStatus = fGeomBoundary;
00869                 else
00870                     fpState->fStepStatus = fWorldBoundary;
00871             }
00872         }
00873         else
00874         {
00875             if(kp == fpProcessInfo->MAXofAlongStepLoops-1)
00876             {
00877                 fpTransportation = dynamic_cast<G4ITTransportation*>(fpCurrentProcess);
00878 
00879                 if(! fpTransportation)
00880                 {
00881                     G4ExceptionDescription exceptionDescription ;
00882                     exceptionDescription << "No transportation process found " ;
00883                     G4Exception("G4ITStepProcessor::DoDefinePhysicalStepLength","ITStepProcessor0010",
00884                                 FatalErrorInArgument,exceptionDescription);
00885                 }
00886 
00887                 fTimeStep               = fpTransportation->GetInteractionTimeLeft();
00888             }
00889         }
00890 
00891         if(proposedTimeStep < fTimeStep)
00892         {
00893             if(fPostStepAtTimeDoItProcTriggered<fpProcessInfo->MAXofPostStepLoops)
00894             {
00895                 if ((fpState->fSelectedPostStepDoItVector)[fPostStepAtTimeDoItProcTriggered] ==
00896                         InActivated)
00897                 {
00898                     (fpState->fSelectedPostStepDoItVector)[fPostStepAtTimeDoItProcTriggered] = NotForced;
00899                     //                    (fpState->fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] = InActivated;
00900 
00901                     fpState->fStepStatus = fPostStepDoItProc;
00902                     fpStep->GetPostStepPoint()->SetProcessDefinedStep(processWithPostStepGivenByTimeStep);
00903 
00904                     fTimeStep = proposedTimeStep;
00905 
00906                     fpTransportation->ComputeStep(*fpTrack,*fpStep,fTimeStep,fpState->fPhysicalStep);
00907                 }
00908             }
00909         }
00910         else
00911         {
00912             if (fPostStepDoItProcTriggered<fpProcessInfo->MAXofPostStepLoops)
00913             {
00914                 if ((fpState->fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] ==
00915                         InActivated)
00916                 {
00917                     (fpState->fSelectedPostStepDoItVector)[fPostStepDoItProcTriggered] =
00918                             NotForced;
00919                 }
00920             }
00921         }
00922 
00923         fpCurrentProcess->SetProcessState(0);
00924 
00925         // Make sure to check the safety, even if Step is not limited
00926         //  by this process.                      J. Apostolakis, June 20, 1998
00927         //
00928         if (safetyProposedToAndByProcess < fpState->proposedSafety)
00929             // proposedSafety keeps the smallest value:
00930             fpState->proposedSafety               = safetyProposedToAndByProcess;
00931         else
00932             // safetyProposedToAndByProcess always proposes a valid safety:
00933             safetyProposedToAndByProcess = fpState->proposedSafety;
00934 
00935     }
00936 
00937     fpITrack->GetTrackingInfo()->SetNavigatorState(fpNavigator->GetNavigatorState());
00938     fpNavigator->SetNavigatorState(0);
00939 }
00940 
00941 //______________________________________________________________________________

Generated on Mon May 27 17:48:42 2013 for Geant4 by  doxygen 1.4.7