Geant4-11
|
#include <G4ITModelProcessor.hh>
Public Member Functions | |
G4double | CalculateMinTimeStep (G4double currentGlobalTime, G4double definedMinTimeStep) |
void | CleanProcessor () |
void | ComputeTrackReaction (G4ITStepStatus fITStepStatus, G4double fGlobalTime, G4double currentTimeStep, G4double previousTimeStep, G4bool reachedUserTimeLimit, G4double fTimeTolerance, G4UserTimeStepAction *fpUserTimeStepAction, G4int fVerbose) |
G4ITModelProcessor () | |
G4ITModelProcessor (const G4ITModelProcessor &other)=delete | |
bool | GetComputeTimeStep () const |
const G4Track * | GetTrack () const |
void | Initialize () |
void | InitializeStepper (G4double currentGlobalTime, G4double userMinTime) |
G4ITModelProcessor & | operator= (const G4ITModelProcessor &other)=delete |
void | RegisterModel (double time, G4VITStepModel *) |
void | SetModelHandler (G4ITModelHandler *) |
void | SetTrackingManager (G4ITTrackingManager *trackingManager) |
virtual | ~G4ITModelProcessor () |
Protected Member Functions | |
void | SetTrack (const G4Track *) |
Protected Attributes | |
std::vector< G4VITStepModel * > | fActiveModels |
bool | fComputeReaction |
bool | fComputeTimeStep |
G4bool | fInitialized |
G4VITStepModel * | fpActiveModelWithMinTimeStep |
G4ITModelHandler * | fpModelHandler |
const G4Track * | fpTrack |
G4ITTrackHolder * | fpTrackContainer |
G4ITTrackingManager * | fpTrackingManager |
std::vector< std::unique_ptr< G4ITReactionChange > > | fReactionInfo |
G4ITReactionSet * | fReactionSet |
G4double | fTSTimeStep |
G4double | fUserMinTimeStep |
The G4ITModelProcessor will call the two processes defined in G4VITModel. This processes act at the beginning and end of each step. The first one, the TimeStepper will calculate a time step to propagate all the track and eventually it can return some tracks that can likely react at the end of the step. The second one, the ReactionProcess will make the tracks reacting.
Definition at line 71 of file G4ITModelProcessor.hh.
G4ITModelProcessor::G4ITModelProcessor | ( | ) |
Definition at line 55 of file G4ITModelProcessor.cc.
References DBL_MAX, fComputeReaction, fComputeTimeStep, fInitialized, fpActiveModelWithMinTimeStep, fpModelHandler, fpTrack, fpTrackContainer, fpTrackingManager, fReactionSet, fTSTimeStep, and fUserMinTimeStep.
|
delete |
|
virtualdefault |
G4double G4ITModelProcessor::CalculateMinTimeStep | ( | G4double | currentGlobalTime, |
G4double | definedMinTimeStep | ||
) |
Definition at line 95 of file G4ITModelProcessor.cc.
References DBL_MAX, G4ITReactionSet::Empty(), fActiveModels, fpActiveModelWithMinTimeStep, fReactionSet, fTSTimeStep, G4cout, G4endl, G4VITStepModel::GetReactionProcess(), G4ITReactionSet::GetReactionsPerTime(), G4VITReactionProcess::Initialize(), InitializeStepper(), and G4MemStat::MemoryUsage().
Referenced by G4Scheduler::Stepping().
void G4ITModelProcessor::CleanProcessor | ( | ) |
Restore the original state. This method should be called only by G4Scheduler
Definition at line 375 of file G4ITModelProcessor.cc.
References fpTrack.
void G4ITModelProcessor::ComputeTrackReaction | ( | G4ITStepStatus | fITStepStatus, |
G4double | fGlobalTime, | ||
G4double | currentTimeStep, | ||
G4double | previousTimeStep, | ||
G4bool | reachedUserTimeLimit, | ||
G4double | fTimeTolerance, | ||
G4UserTimeStepAction * | fpUserTimeStepAction, | ||
G4int | fVerbose | ||
) |
Definition at line 173 of file G4ITModelProcessor.cc.
References eCollisionBetweenTracks, G4ITReactionSet::Empty(), G4ITTrackingManager::EndTracking(), FatalErrorInArgument, G4VITReactionProcess::FindReaction(), fpActiveModelWithMinTimeStep, fpTrackContainer, fpTrackingManager, fReactionInfo, fReactionSet, fStopAndKill, G4BestUnit, G4cout, G4endl, G4Exception(), G4Track::GetGlobalTime(), GetIT(), G4IT::GetName(), G4Track::GetParticleDefinition(), G4ParticleDefinition::GetParticleName(), G4VITStepModel::GetReactionProcess(), G4Track::GetStep(), G4Track::GetTrackID(), G4ITTrackHolder::KillTracks(), G4ITTrackHolder::MergeSecondariesWithMainList(), G4IT::SetParentID(), and G4UserTimeStepAction::UserReactionAction().
Referenced by G4Scheduler::Stepping().
bool G4ITModelProcessor::GetComputeTimeStep | ( | ) | const |
Definition at line 380 of file G4ITModelProcessor.cc.
References fComputeTimeStep.
Referenced by G4Scheduler::Stepping().
const G4Track * G4ITModelProcessor::GetTrack | ( | ) | const |
Definition at line 385 of file G4ITModelProcessor.cc.
References fpTrack.
void G4ITModelProcessor::Initialize | ( | ) |
Definition at line 77 of file G4ITModelProcessor.cc.
References fComputeReaction, fComputeTimeStep, fInitialized, fpModelHandler, fpTrackContainer, fReactionSet, G4ITModelHandler::GetReactionProcessFlag(), G4ITModelHandler::GetTimeStepComputerFlag(), G4ITModelHandler::Initialize(), G4ITReactionSet::Instance(), and G4ITTrackHolder::Instance().
Referenced by G4Scheduler::Process().
Definition at line 146 of file G4ITModelProcessor.cc.
References fActiveModels, fpModelHandler, G4cout, G4endl, G4ITModelHandler::GetActiveModels(), G4MemStat::MemoryUsage(), and G4VITTimeStepComputer::SetTimes().
Referenced by CalculateMinTimeStep().
|
delete |
void G4ITModelProcessor::RegisterModel | ( | double | time, |
G4VITStepModel * | model | ||
) |
Definition at line 72 of file G4ITModelProcessor.cc.
References fpModelHandler, and G4ITModelHandler::RegisterModel().
void G4ITModelProcessor::SetModelHandler | ( | G4ITModelHandler * | pModelHandler | ) |
Definition at line 362 of file G4ITModelProcessor.cc.
References FatalErrorInArgument, fInitialized, fpModelHandler, and G4Exception().
Referenced by G4Scheduler::Initialize().
|
protected |
Definition at line 357 of file G4ITModelProcessor.cc.
References fpTrack.
void G4ITModelProcessor::SetTrackingManager | ( | G4ITTrackingManager * | trackingManager | ) |
Definition at line 390 of file G4ITModelProcessor.cc.
References fpTrackingManager.
Referenced by G4Scheduler::Initialize().
|
protected |
Definition at line 123 of file G4ITModelProcessor.hh.
Referenced by CalculateMinTimeStep(), and InitializeStepper().
|
protected |
Definition at line 129 of file G4ITModelProcessor.hh.
Referenced by G4ITModelProcessor(), and Initialize().
|
protected |
Definition at line 128 of file G4ITModelProcessor.hh.
Referenced by G4ITModelProcessor(), GetComputeTimeStep(), and Initialize().
|
protected |
Definition at line 117 of file G4ITModelProcessor.hh.
Referenced by G4ITModelProcessor(), Initialize(), and SetModelHandler().
|
protected |
Definition at line 124 of file G4ITModelProcessor.hh.
Referenced by CalculateMinTimeStep(), ComputeTrackReaction(), and G4ITModelProcessor().
|
protected |
Definition at line 118 of file G4ITModelProcessor.hh.
Referenced by G4ITModelProcessor(), Initialize(), InitializeStepper(), RegisterModel(), and SetModelHandler().
|
protected |
Definition at line 120 of file G4ITModelProcessor.hh.
Referenced by CleanProcessor(), G4ITModelProcessor(), GetTrack(), and SetTrack().
|
protected |
Definition at line 115 of file G4ITModelProcessor.hh.
Referenced by ComputeTrackReaction(), G4ITModelProcessor(), and Initialize().
|
protected |
Definition at line 114 of file G4ITModelProcessor.hh.
Referenced by ComputeTrackReaction(), G4ITModelProcessor(), and SetTrackingManager().
|
protected |
Definition at line 126 of file G4ITModelProcessor.hh.
Referenced by ComputeTrackReaction().
|
protected |
Definition at line 113 of file G4ITModelProcessor.hh.
Referenced by CalculateMinTimeStep(), ComputeTrackReaction(), G4ITModelProcessor(), and Initialize().
|
protected |
Definition at line 112 of file G4ITModelProcessor.hh.
Referenced by CalculateMinTimeStep(), and G4ITModelProcessor().
|
protected |
Definition at line 121 of file G4ITModelProcessor.hh.
Referenced by G4ITModelProcessor().