Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParallelWorldProcess Class Reference

#include <G4ParallelWorldProcess.hh>

+ Inheritance diagram for G4ParallelWorldProcess:

Public Member Functions

 G4ParallelWorldProcess (const G4String &processName="ParaWorld", G4ProcessType theType=fParallel)
 
 ~G4ParallelWorldProcess () override
 
void SetParallelWorld (G4String parallelWorldName)
 
void SetParallelWorld (G4VPhysicalVolume *parallelWorld)
 
void StartTracking (G4Track *) override
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *) override
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
void SetLayeredMaterialFlag (G4bool flg=true)
 
G4bool GetLayeredMaterialFlag () const
 
G4bool IsAtRestRequired (G4ParticleDefinition *)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Static Public Member Functions

static const G4StepGetHyperStep ()
 
static G4int GetHypNavigatorID ()
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Protected Member Functions

void CopyStep (const G4Step &step)
 
void SwitchMaterial (G4StepPoint *)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4StepfGhostStep
 
G4StepPointfGhostPreStepPoint
 
G4StepPointfGhostPostStepPoint
 
G4VParticleChange aDummyParticleChange
 
G4ParticleChange xParticleChange
 
G4TransportationManagerfTransportationManager
 
G4PathFinderfPathFinder
 
G4String fGhostWorldName
 
G4VPhysicalVolumefGhostWorld {nullptr}
 
G4NavigatorfGhostNavigator {nullptr}
 
G4int fNavigatorID {-1}
 
G4TouchableHandle fOldGhostTouchable
 
G4TouchableHandle fNewGhostTouchable
 
G4FieldTrack fFieldTrack
 
G4double fGhostSafety {0.}
 
G4bool fOnBoundary {false}
 
G4bool layeredMaterialFlag {false}
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 56 of file G4ParallelWorldProcess.hh.

Constructor & Destructor Documentation

◆ G4ParallelWorldProcess()

G4ParallelWorldProcess::G4ParallelWorldProcess ( const G4String & processName = "ParaWorld",
G4ProcessType theType = fParallel )

Definition at line 64 of file G4ParallelWorldProcess.cc.

65 : G4VProcess(processName, theType), fFieldTrack('0')
66{
68 if (fpHyperStep == nullptr) fpHyperStep = new G4Step();
69 iParallelWorld = ++nParallelWorlds;
70
72
73 fGhostStep = new G4Step();
76
80
81 fGhostWorldName = "** NotDefined **";
83
84 if (verboseLevel > 0) {
85 G4cout << GetProcessName() << " is created " << G4endl;
86 }
87}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetPushVerbosity(G4bool mode)
static G4ParallelWorldProcessStore * GetInstance()
void SetParallelWorld(G4ParallelWorldProcess *proc, const G4String &parallelWorldName)
G4VParticleChange aDummyParticleChange
G4TransportationManager * fTransportationManager
static G4PathFinder * GetInstance()
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4int verboseLevel
void SetProcessSubType(G4int)
G4VParticleChange * pParticleChange
const G4String & GetProcessName() const

◆ ~G4ParallelWorldProcess()

G4ParallelWorldProcess::~G4ParallelWorldProcess ( )
override

Definition at line 89 of file G4ParallelWorldProcess.cc.

90{
91 delete fGhostStep;
92 nParallelWorlds--;
93 if (nParallelWorlds == 0) {
94 delete fpHyperStep;
95 fpHyperStep = nullptr;
96 }
97}

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4ParallelWorldProcess::AlongStepDoIt ( const G4Track & track,
const G4Step &  )
overridevirtual

Implements G4VProcess.

Definition at line 320 of file G4ParallelWorldProcess.cc.

321{
323 return pParticleChange;
324}
virtual void Initialize(const G4Track &)

◆ AlongStepGetPhysicalInteractionLength()

G4double G4ParallelWorldProcess::AlongStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4double currentMinimumStep,
G4double & proposedSafety,
G4GPILSelection * selection )
overridevirtual

Implements G4VProcess.

Definition at line 246 of file G4ParallelWorldProcess.cc.

251{
252 static G4ThreadLocal G4FieldTrack* endTrack_G4MT_TLS_ = nullptr;
253 if (endTrack_G4MT_TLS_ == nullptr) endTrack_G4MT_TLS_ = new G4FieldTrack('0');
254 G4FieldTrack& endTrack = *endTrack_G4MT_TLS_;
255 // static ELimited eLimited;
256 ELimited eLimited;
257 ELimited eLim = kUndefLimited;
258
259 *selection = NotCandidateForSelection;
260 G4double returnedStep = DBL_MAX;
261
262 if (previousStepSize > 0.) {
263 fGhostSafety -= previousStepSize;
264 }
265 if (fGhostSafety < 0.) fGhostSafety = 0.0;
266
267 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.) {
268 // I have no chance to limit
269 returnedStep = currentMinimumStep;
270 fOnBoundary = false;
271 proposedSafety = fGhostSafety - currentMinimumStep;
272 eLim = kDoNot;
273 }
274 else {
276
277#ifdef G4DEBUG_PARALLEL_WORLD_PROCESS
278 if (verboseLevel > 0) {
279 int localVerb = verboseLevel - 1;
280
281 if (localVerb == 1) {
282 G4cout << " Pll Wrl proc::AlongStepGPIL " << this->GetProcessName() << G4endl;
283 }
284 else if (localVerb > 1) {
285 G4cout << "----------------------------------------------" << G4endl;
286 G4cout << " ParallelWorldProcess: field Track set to : " << G4endl;
287 G4cout << "----------------------------------------------" << G4endl;
289 G4cout << "----------------------------------------------" << G4endl;
290 }
291 }
292#endif
293
294 returnedStep = fPathFinder->ComputeStep(fFieldTrack, currentMinimumStep, fNavigatorID,
295 track.GetCurrentStepNumber(), fGhostSafety, eLimited,
296 endTrack, track.GetVolume());
297 if (eLimited == kDoNot) {
298 fOnBoundary = false;
300 }
301 else {
302 fOnBoundary = true;
303 // fGhostSafetyEnd = 0.0; // At end-point of expected step only
304 }
305 proposedSafety = fGhostSafety;
306 if (eLimited == kUnique || eLimited == kSharedOther) {
307 *selection = CandidateForSelection;
308 }
309 else if (eLimited == kSharedTransport) {
310 returnedStep *= (1.0 + 1.0e-9);
311 }
312 eLim = eLimited;
313 }
314
315 if (iParallelWorld == nParallelWorlds) fNavIDHyp = 0;
316 if (eLim == kUnique || eLim == kSharedOther) fNavIDHyp = fNavigatorID;
317 return returnedStep;
318}
@ CandidateForSelection
@ NotCandidateForSelection
@ kDoNot
@ kUndefLimited
@ kUnique
@ kSharedOther
@ kSharedTransport
double G4double
Definition G4Types.hh:83
static void Update(G4FieldTrack *, const G4Track *)
G4ThreeVector GetPosition() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
G4VPhysicalVolume * GetVolume() const
G4int GetCurrentStepNumber() const
#define DBL_MAX
Definition templates.hh:62
#define G4ThreadLocal
Definition tls.hh:77

◆ AtRestDoIt()

G4VParticleChange * G4ParallelWorldProcess::AtRestDoIt ( const G4Track & track,
const G4Step & step )
overridevirtual

Implements G4VProcess.

Definition at line 162 of file G4ParallelWorldProcess.cc.

163{
164 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
165 // At Rest must be registered ONLY for the particle which has other At Rest
166 // process(es).
167 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
169 G4VSensitiveDetector* aSD = nullptr;
170 if (fOldGhostTouchable->GetVolume() != nullptr) {
172 }
173 fOnBoundary = false;
174 if (aSD != nullptr) {
175 CopyStep(step);
177
179
182 if (fNewGhostTouchable->GetVolume() != nullptr) {
185 }
186 else {
188 }
189
190 aSD->Hit(fGhostStep);
191 }
192
194 return pParticleChange;
195}
G4VSensitiveDetector * GetSensitiveDetector() const
void CopyStep(const G4Step &step)
void SetSensitiveDetector(G4VSensitiveDetector *)
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4TouchableHandle & GetTouchableHandle() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
G4LogicalVolume * GetLogicalVolume() const
G4bool Hit(G4Step *aStep)

◆ AtRestGetPhysicalInteractionLength()

G4double G4ParallelWorldProcess::AtRestGetPhysicalInteractionLength ( const G4Track & ,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 151 of file G4ParallelWorldProcess.cc.

153{
154 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
155 // At Rest must be registered ONLY for the particle which has other At Rest
156 // process(es).
157 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
158 *condition = Forced;
159 return DBL_MAX;
160}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced

◆ CopyStep()

void G4ParallelWorldProcess::CopyStep ( const G4Step & step)
protected

Definition at line 326 of file G4ParallelWorldProcess.cc.

327{
329
335 fGhostStep->SetSecondary((const_cast<G4Step&>(step)).GetfSecondary());
336
339
341 if (fOnBoundary) {
343 }
346 }
347
348 if (iParallelWorld == 1) {
349 G4StepStatus prevStatHyp = fpHyperStep->GetPostStepPoint()->GetStepStatus();
350
351 fpHyperStep->SetTrack(step.GetTrack());
352 fpHyperStep->SetStepLength(step.GetStepLength());
353 fpHyperStep->SetTotalEnergyDeposit(step.GetTotalEnergyDeposit());
355 fpHyperStep->SetControlFlag(step.GetControlFlag());
356
357 *(fpHyperStep->GetPreStepPoint()) = *(fpHyperStep->GetPostStepPoint());
358 *(fpHyperStep->GetPostStepPoint()) = *(step.GetPostStepPoint());
359
360 fpHyperStep->GetPreStepPoint()->SetStepStatus(prevStatHyp);
361 }
362
363 if (fOnBoundary) {
365 }
366}
G4StepStatus
@ fGeomBoundary
@ fPostStepDoItProc
G4StepStatus GetStepStatus() const
void SetStepStatus(const G4StepStatus aValue)
G4SteppingControl GetControlFlag() const
G4Track * GetTrack() const
void SetStepLength(G4double value)
void SetNonIonizingEnergyDeposit(G4double value)
G4double GetNonIonizingEnergyDeposit() const
G4double GetStepLength() const
G4double GetTotalEnergyDeposit() const
void SetSecondary(G4TrackVector *value)
void SetControlFlag(G4SteppingControl StepControlFlag)
void SetTotalEnergyDeposit(G4double value)
void SetTrack(G4Track *value)

Referenced by AtRestDoIt(), and PostStepDoIt().

◆ GetHyperStep()

const G4Step * G4ParallelWorldProcess::GetHyperStep ( )
static

Definition at line 54 of file G4ParallelWorldProcess.cc.

55{
56 return fpHyperStep;
57}

Referenced by G4RayTrajectory::AppendStep(), G4OpBoundaryProcess::PostStepDoIt(), and G4UCNBoundaryProcess::PostStepDoIt().

◆ GetHypNavigatorID()

G4int G4ParallelWorldProcess::GetHypNavigatorID ( )
static

◆ GetLayeredMaterialFlag()

G4bool G4ParallelWorldProcess::GetLayeredMaterialFlag ( ) const
inline

Definition at line 105 of file G4ParallelWorldProcess.hh.

◆ IsAtRestRequired()

G4bool G4ParallelWorldProcess::IsAtRestRequired ( G4ParticleDefinition * partDef)

Definition at line 400 of file G4ParallelWorldProcess.cc.

401{
402 G4int pdgCode = partDef->GetPDGEncoding();
403 if (pdgCode == 0) {
404 G4String partName = partDef->GetParticleName();
405 if (partName == "geantino") return false;
406 if (partName == "chargedgeantino") return false;
407 }
408 else {
409 if (pdgCode == 11 || pdgCode == 2212) return false; // electrons and proton
410 pdgCode = std::abs(pdgCode);
411 if (pdgCode == 22) return false; // gamma and optical photons
412 if (pdgCode == 12 || pdgCode == 14 || pdgCode == 16) return false; // all neutronos
413 }
414 return true;
415}
int G4int
Definition G4Types.hh:85
const G4String & GetParticleName() const

Referenced by G4ParallelWorldPhysics::ConstructProcess(), G4RunManager::ConstructScoringWorlds(), and G4WorkerRunManager::ConstructScoringWorlds().

◆ PostStepDoIt()

G4VParticleChange * G4ParallelWorldProcess::PostStepDoIt ( const G4Track & track,
const G4Step & step )
overridevirtual

Implements G4VProcess.

Definition at line 205 of file G4ParallelWorldProcess.cc.

206{
208 G4VSensitiveDetector* aSD = nullptr;
209 if (fOldGhostTouchable->GetVolume() != nullptr) {
211 }
212 CopyStep(step);
214
215 if (fOnBoundary) {
217 }
218 else {
220 }
221
224
225 if (fNewGhostTouchable->GetVolume() != nullptr) {
228 }
229 else {
231 }
232
234 if (sd != nullptr) {
235 sd->Hit(fGhostStep);
236 }
237
240 G4StepPoint* realWorldPostStepPoint = ((G4Step*)(track.GetStep()))->GetPostStepPoint();
241 SwitchMaterial(realWorldPostStepPoint);
242 }
243 return pParticleChange;
244}
void SwitchMaterial(G4StepPoint *)
G4TouchableHandle CreateTouchableHandle(G4int navId) const
G4VSensitiveDetector * GetSensitiveDetector() const
const G4Step * GetStep() const

◆ PostStepGetPhysicalInteractionLength()

G4double G4ParallelWorldProcess::PostStepGetPhysicalInteractionLength ( const G4Track & ,
G4double ,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 197 of file G4ParallelWorldProcess.cc.

200{
202 return DBL_MAX;
203}
@ StronglyForced

◆ SetLayeredMaterialFlag()

void G4ParallelWorldProcess::SetLayeredMaterialFlag ( G4bool flg = true)
inline

◆ SetParallelWorld() [1/2]

void G4ParallelWorldProcess::SetParallelWorld ( G4String parallelWorldName)

◆ SetParallelWorld() [2/2]

void G4ParallelWorldProcess::SetParallelWorld ( G4VPhysicalVolume * parallelWorld)

Definition at line 107 of file G4ParallelWorldProcess.cc.

108{
109 fGhostWorldName = parallelWorld->GetName();
110 fGhostWorld = parallelWorld;
113}
const G4String & GetName() const

◆ StartTracking()

void G4ParallelWorldProcess::StartTracking ( G4Track * trk)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 115 of file G4ParallelWorldProcess.cc.

116{
117 if (fGhostNavigator != nullptr) {
119 }
120 else {
122 "G4ParallelWorldProcess::StartTracking", "ProcParaWorld000", FatalException,
123 "G4ParallelWorldProcess is used for tracking without having a parallel world assigned");
124 }
126
131
132 fGhostSafety = -1.;
133 fOnBoundary = false;
136
137 *(fpHyperStep->GetPostStepPoint()) = *(trk->GetStep()->GetPostStepPoint());
139 G4StepPoint* realWorldPostStepPoint = trk->GetStep()->GetPostStepPoint();
140 SwitchMaterial(realWorldPostStepPoint);
141 G4StepPoint* realWorldPreStepPoint = trk->GetStep()->GetPreStepPoint();
142 SwitchMaterial(realWorldPreStepPoint);
143 G4double velocity = trk->CalculateVelocity();
144 realWorldPostStepPoint->SetVelocity(velocity);
145 realWorldPreStepPoint->SetVelocity(velocity);
146 trk->SetVelocity(velocity);
147 }
148 *(fpHyperStep->GetPreStepPoint()) = *(fpHyperStep->GetPostStepPoint());
149}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
@ fUndefined
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=nullptr)
void SetVelocity(G4double v)
void SetVelocity(G4double val)
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4double CalculateVelocity() const
G4int ActivateNavigator(G4Navigator *aNavigator)

◆ SwitchMaterial()

void G4ParallelWorldProcess::SwitchMaterial ( G4StepPoint * realWorldStepPoint)
protected

Definition at line 368 of file G4ParallelWorldProcess.cc.

369{
370 if (realWorldStepPoint->GetStepStatus() == fWorldBoundary) return;
372 if (thePhys != nullptr) {
373 G4Material* ghostMaterial = thePhys->GetLogicalVolume()->GetMaterial();
374 if (ghostMaterial != nullptr) {
375 G4Region* ghostRegion = thePhys->GetLogicalVolume()->GetRegion();
376 G4ProductionCuts* prodCuts = realWorldStepPoint->GetMaterialCutsCouple()->GetProductionCuts();
377 if (ghostRegion != nullptr) {
378 G4ProductionCuts* ghostProdCuts = ghostRegion->GetProductionCuts();
379 if (ghostProdCuts != nullptr) prodCuts = ghostProdCuts;
380 }
381 const G4MaterialCutsCouple* ghostMCCouple =
383 prodCuts);
384 if (ghostMCCouple != nullptr) {
385 realWorldStepPoint->SetMaterial(ghostMaterial);
386 realWorldStepPoint->SetMaterialCutsCouple(ghostMCCouple);
387 *(fpHyperStep->GetPostStepPoint()) = *(fGhostPostStepPoint);
388 fpHyperStep->GetPostStepPoint()->SetMaterial(ghostMaterial);
389 fpHyperStep->GetPostStepPoint()->SetMaterialCutsCouple(ghostMCCouple);
390 }
391 else {
392 G4cout << "!!! MaterialCutsCouple is not found for " << ghostMaterial->GetName() << "."
393 << G4endl << " Material in real world ("
394 << realWorldStepPoint->GetMaterial()->GetName() << ") is used." << G4endl;
395 }
396 }
397 }
398}
@ fWorldBoundary
G4Region * GetRegion() const
G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
const G4String & GetName() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ProductionCuts * GetProductionCuts() const
void SetMaterial(G4Material *)
G4Material * GetMaterial() const
void SetMaterialCutsCouple(const G4MaterialCutsCouple *)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const

Referenced by PostStepDoIt(), and StartTracking().

Member Data Documentation

◆ aDummyParticleChange

G4VParticleChange G4ParallelWorldProcess::aDummyParticleChange
protected

Definition at line 129 of file G4ParallelWorldProcess.hh.

Referenced by G4ParallelWorldProcess().

◆ fFieldTrack

G4FieldTrack G4ParallelWorldProcess::fFieldTrack
protected

Definition at line 144 of file G4ParallelWorldProcess.hh.

Referenced by AlongStepGetPhysicalInteractionLength().

◆ fGhostNavigator

G4Navigator* G4ParallelWorldProcess::fGhostNavigator {nullptr}
protected

◆ fGhostPostStepPoint

G4StepPoint* G4ParallelWorldProcess::fGhostPostStepPoint
protected

◆ fGhostPreStepPoint

G4StepPoint* G4ParallelWorldProcess::fGhostPreStepPoint
protected

◆ fGhostSafety

G4double G4ParallelWorldProcess::fGhostSafety {0.}
protected

Definition at line 145 of file G4ParallelWorldProcess.hh.

145{0.};

Referenced by AlongStepGetPhysicalInteractionLength(), and StartTracking().

◆ fGhostStep

G4Step* G4ParallelWorldProcess::fGhostStep
protected

◆ fGhostWorld

G4VPhysicalVolume* G4ParallelWorldProcess::fGhostWorld {nullptr}
protected

Definition at line 139 of file G4ParallelWorldProcess.hh.

139{nullptr};

Referenced by SetParallelWorld(), and SetParallelWorld().

◆ fGhostWorldName

G4String G4ParallelWorldProcess::fGhostWorldName
protected

◆ fNavigatorID

G4int G4ParallelWorldProcess::fNavigatorID {-1}
protected

Definition at line 141 of file G4ParallelWorldProcess.hh.

141{-1};

Referenced by AlongStepGetPhysicalInteractionLength(), PostStepDoIt(), and StartTracking().

◆ fNewGhostTouchable

G4TouchableHandle G4ParallelWorldProcess::fNewGhostTouchable
protected

Definition at line 143 of file G4ParallelWorldProcess.hh.

Referenced by AtRestDoIt(), PostStepDoIt(), StartTracking(), and SwitchMaterial().

◆ fOldGhostTouchable

G4TouchableHandle G4ParallelWorldProcess::fOldGhostTouchable
protected

Definition at line 142 of file G4ParallelWorldProcess.hh.

Referenced by AtRestDoIt(), PostStepDoIt(), and StartTracking().

◆ fOnBoundary

G4bool G4ParallelWorldProcess::fOnBoundary {false}
protected

◆ fPathFinder

G4PathFinder* G4ParallelWorldProcess::fPathFinder
protected

◆ fTransportationManager

G4TransportationManager* G4ParallelWorldProcess::fTransportationManager
protected

◆ layeredMaterialFlag

G4bool G4ParallelWorldProcess::layeredMaterialFlag {false}
protected

◆ xParticleChange

G4ParticleChange G4ParallelWorldProcess::xParticleChange
protected

Definition at line 130 of file G4ParallelWorldProcess.hh.


The documentation for this class was generated from the following files: