Geant4 11.1.1
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)
 
virtual ~G4ParallelWorldProcess ()
 
void SetParallelWorld (G4String parallelWorldName)
 
void SetParallelWorld (G4VPhysicalVolume *parallelWorld)
 
void StartTracking (G4Track *)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 StartTracking (G4Track *)
 
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
 
G4NavigatorfGhostNavigator
 
G4int fNavigatorID
 
G4TouchableHandle fOldGhostTouchable
 
G4TouchableHandle fNewGhostTouchable
 
G4FieldTrack fFieldTrack
 
G4double fGhostSafety
 
G4bool fOnBoundary
 
G4bool layeredMaterialFlag
 
- 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 71 of file G4ParallelWorldProcess.hh.

Constructor & Destructor Documentation

◆ G4ParallelWorldProcess()

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

Definition at line 60 of file G4ParallelWorldProcess.cc.

62:G4VProcess(processName,theType),fGhostWorld(nullptr),fGhostNavigator(nullptr),
65{
67 if(!fpHyperStep) fpHyperStep = new G4Step();
68 iParallelWorld = ++nParallelWorlds;
69
71
72 fGhostStep = new G4Step();
75
79
80 fGhostWorldName = "** NotDefined **";
82
83 if (verboseLevel>0)
84 {
85 G4cout << GetProcessName() << " is created " << G4endl;
86 }
87}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetPushVerbosity(G4bool mode)
static G4ParallelWorldProcessStore * GetInstance()
void SetParallelWorld(G4ParallelWorldProcess *proc, G4String parallelWorldName)
G4VParticleChange aDummyParticleChange
G4TransportationManager * fTransportationManager
G4VPhysicalVolume * fGhostWorld
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:52
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
static G4TransportationManager * GetTransportationManager()
G4Navigator * GetNavigatorForTracking() const
G4int verboseLevel
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

◆ ~G4ParallelWorldProcess()

G4ParallelWorldProcess::~G4ParallelWorldProcess ( )
virtual

Definition at line 89 of file G4ParallelWorldProcess.cc.

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

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 343 of file G4ParallelWorldProcess.cc.

345{
347 return pParticleChange;
348}
virtual void Initialize(const G4Track &)

◆ AlongStepGetPhysicalInteractionLength()

G4double G4ParallelWorldProcess::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
)
virtual

Implements G4VProcess.

Definition at line 270 of file G4ParallelWorldProcess.cc.

273{
274 static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ; if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ; G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
275 //static ELimited eLimited;
276 ELimited eLimited;
277 ELimited eLim = kUndefLimited;
278
279 *selection = NotCandidateForSelection;
280 G4double returnedStep = DBL_MAX;
281
282 if (previousStepSize > 0.)
283 { fGhostSafety -= previousStepSize; }
284 if (fGhostSafety < 0.) fGhostSafety = 0.0;
285
286 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
287 {
288 // I have no chance to limit
289 returnedStep = currentMinimumStep;
290 fOnBoundary = false;
291 proposedSafety = fGhostSafety - currentMinimumStep;
292 eLim = kDoNot;
293 }
294 else
295 {
297
298#ifdef G4DEBUG_PARALLEL_WORLD_PROCESS
299 if( verboseLevel > 0 ){
300 int localVerb = verboseLevel-1;
301
302 if( localVerb == 1 ) {
303 G4cout << " Pll Wrl proc::AlongStepGPIL " << this->GetProcessName() << G4endl;
304 }else if( localVerb > 1 ) {
305 G4cout << "----------------------------------------------" << G4endl;
306 G4cout << " ParallelWorldProcess: field Track set to : " << G4endl;
307 G4cout << "----------------------------------------------" << G4endl;
309 G4cout << "----------------------------------------------" << G4endl;
310 }
311 }
312#endif
313
314 returnedStep
315 = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
316 track.GetCurrentStepNumber(),fGhostSafety,eLimited,
317 endTrack,track.GetVolume());
318 if(eLimited == kDoNot)
319 {
320 fOnBoundary = false;
322 }
323 else
324 {
325 fOnBoundary = true;
326 // fGhostSafetyEnd = 0.0; // At end-point of expected step only
327 }
328 proposedSafety = fGhostSafety;
329 if(eLimited == kUnique || eLimited == kSharedOther) {
330 *selection = CandidateForSelection;
331 }
332 else if (eLimited == kSharedTransport) {
333 returnedStep *= (1.0 + 1.0e-9);
334 }
335 eLim = eLimited;
336 }
337
338 if(iParallelWorld==nParallelWorlds) fNavIDHyp = 0;
339 if(eLim == kUnique || eLim == kSharedOther) fNavIDHyp = fNavigatorID;
340 return returnedStep;
341}
@ CandidateForSelection
@ NotCandidateForSelection
ELimited
@ 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 
)
virtual

Implements G4VProcess.

Definition at line 176 of file G4ParallelWorldProcess.cc.

179{
180//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
181// At Rest must be registered ONLY for the particle which has other At Rest
182// process(es).
183//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
185 G4VSensitiveDetector* aSD = 0;
188 fOnBoundary = false;
189 if(aSD)
190 {
191 CopyStep(step);
193
195
199 {
202 }
203 else
205
206 aSD->Hit(fGhostStep);
207 }
208
210 return pParticleChange;
211}
G4VSensitiveDetector * GetSensitiveDetector() const
G4TouchableHandle fNewGhostTouchable
void CopyStep(const G4Step &step)
G4TouchableHandle fOldGhostTouchable
void SetSensitiveDetector(G4VSensitiveDetector *)
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4TouchableHandle & GetTouchableHandle() const
G4LogicalVolume * GetLogicalVolume() const
G4bool Hit(G4Step *aStep)
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34

◆ AtRestGetPhysicalInteractionLength()

G4double G4ParallelWorldProcess::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 164 of file G4ParallelWorldProcess.cc.

167{
168//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
169// At Rest must be registered ONLY for the particle which has other At Rest
170// process(es).
171//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
172 *condition = Forced;
173 return DBL_MAX;
174}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced

◆ CopyStep()

void G4ParallelWorldProcess::CopyStep ( const G4Step step)
protected

Definition at line 350 of file G4ParallelWorldProcess.cc.

351{
353
359 fGhostStep->SetSecondary((const_cast<G4Step&>(step)).GetfSecondary());
360
363
365 if(fOnBoundary)
369
370 if(iParallelWorld==1)
371 {
372 G4StepStatus prevStatHyp = fpHyperStep->GetPostStepPoint()->GetStepStatus();
373
374 fpHyperStep->SetTrack(step.GetTrack());
375 fpHyperStep->SetStepLength(step.GetStepLength());
376 fpHyperStep->SetTotalEnergyDeposit(step.GetTotalEnergyDeposit());
378 fpHyperStep->SetControlFlag(step.GetControlFlag());
379
380 *(fpHyperStep->GetPreStepPoint()) = *(fpHyperStep->GetPostStepPoint());
381 *(fpHyperStep->GetPostStepPoint()) = *(step.GetPostStepPoint());
382
383 fpHyperStep->GetPreStepPoint()->SetStepStatus(prevStatHyp);
384 }
385
386 if(fOnBoundary)
387 { fpHyperStep->GetPostStepPoint()->SetStepStatus(fGeomBoundary); }
388}
G4StepStatus
Definition: G4StepStatus.hh:40
@ fGeomBoundary
Definition: G4StepStatus.hh:43
@ fPostStepDoItProc
Definition: G4StepStatus.hh:49
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

◆ GetHypNavigatorID()

G4int G4ParallelWorldProcess::GetHypNavigatorID ( )
static

◆ GetLayeredMaterialFlag()

G4bool G4ParallelWorldProcess::GetLayeredMaterialFlag ( ) const
inline

Definition at line 125 of file G4ParallelWorldProcess.hh.

126 { return layeredMaterialFlag; }

◆ IsAtRestRequired()

G4bool G4ParallelWorldProcess::IsAtRestRequired ( G4ParticleDefinition partDef)

Definition at line 430 of file G4ParallelWorldProcess.cc.

431{
432 G4int pdgCode = partDef->GetPDGEncoding();
433 if(pdgCode==0)
434 {
435 G4String partName = partDef->GetParticleName();
436 if(partName=="geantino") return false;
437 if(partName=="chargedgeantino") return false;
438 }
439 else
440 {
441 if(pdgCode==11 || pdgCode==2212) return false; // electrons and proton
442 pdgCode = std::abs(pdgCode);
443 if(pdgCode==22) return false; // gamma and optical photons
444 if(pdgCode==12 || pdgCode==14 || pdgCode==16) return false; // all neutronos
445 }
446 return true;
447}
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 
)
virtual

Implements G4VProcess.

Definition at line 223 of file G4ParallelWorldProcess.cc.

226{
228 G4VSensitiveDetector* aSD = 0;
231 CopyStep(step);
233
234 if(fOnBoundary)
235 {
237 }
238 else
239 {
241 }
242
245
247 {
250 }
251 else
253
255 if(sd)
256 {
257 sd->Hit(fGhostStep);
258 }
259
262 {
263 G4StepPoint* realWorldPostStepPoint =
264 ((G4Step*)(track.GetStep()))->GetPostStepPoint();
265 SwitchMaterial(realWorldPostStepPoint);
266 }
267 return pParticleChange;
268}
void SwitchMaterial(G4StepPoint *)
G4TouchableHandle CreateTouchableHandle(G4int navId) const
G4VSensitiveDetector * GetSensitiveDetector() const
const G4Step * GetStep() const

◆ PostStepGetPhysicalInteractionLength()

G4double G4ParallelWorldProcess::PostStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 214 of file G4ParallelWorldProcess.cc.

218{
220 return DBL_MAX;
221}
@ 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 109 of file G4ParallelWorldProcess.cc.

111{
112 fGhostWorldName = parallelWorld->GetName();
113 fGhostWorld = parallelWorld;
116}
const G4String & GetName() const

◆ StartTracking()

void G4ParallelWorldProcess::StartTracking ( G4Track trk)
virtual

Reimplemented from G4VProcess.

Definition at line 118 of file G4ParallelWorldProcess.cc.

119{
122 else
123 {
124 G4Exception("G4ParallelWorldProcess::StartTracking",
125 "ProcParaWorld000",FatalException,
126 "G4ParallelWorldProcess is used for tracking without having a parallel world assigned");
127 }
129
134
135 fGhostSafety = -1.;
136 fOnBoundary = false;
139
140// G4VPhysicalVolume* thePhys = fNewGhostTouchable->GetVolume();
141// if(thePhys)
142// {
143// G4Material* ghostMaterial = thePhys->GetLogicalVolume()->GetMaterial();
144// if(ghostMaterial)
145// { G4cout << " --- Material : " << ghostMaterial->GetName() << G4endl; }
146// }
147
148 *(fpHyperStep->GetPostStepPoint()) = *(trk->GetStep()->GetPostStepPoint());
150 {
151 G4StepPoint* realWorldPostStepPoint = trk->GetStep()->GetPostStepPoint();
152 SwitchMaterial(realWorldPostStepPoint);
153 G4StepPoint *realWorldPreStepPoint = trk->GetStep()->GetPreStepPoint();
154 SwitchMaterial(realWorldPreStepPoint);
155 G4double velocity = trk->CalculateVelocity();
156 realWorldPostStepPoint->SetVelocity(velocity);
157 realWorldPreStepPoint->SetVelocity(velocity);
158 trk->SetVelocity(velocity);
159 }
160 *(fpHyperStep->GetPreStepPoint()) = *(fpHyperStep->GetPostStepPoint());
161}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
@ fUndefined
Definition: G4StepStatus.hh:55
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 390 of file G4ParallelWorldProcess.cc.

391{
392 if(realWorldStepPoint->GetStepStatus()==fWorldBoundary) return;
394 if(thePhys)
395 {
396 G4Material* ghostMaterial = thePhys->GetLogicalVolume()->GetMaterial();
397 if(ghostMaterial)
398 {
399 G4Region* ghostRegion = thePhys->GetLogicalVolume()->GetRegion();
400 G4ProductionCuts* prodCuts =
401 realWorldStepPoint->GetMaterialCutsCouple()->GetProductionCuts();
402 if(ghostRegion)
403 {
404 G4ProductionCuts* ghostProdCuts = ghostRegion->GetProductionCuts();
405 if(ghostProdCuts) prodCuts = ghostProdCuts;
406 }
407 const G4MaterialCutsCouple* ghostMCCouple =
409 ->GetMaterialCutsCouple(ghostMaterial,prodCuts);
410 if(ghostMCCouple)
411 {
412 realWorldStepPoint->SetMaterial(ghostMaterial);
413 realWorldStepPoint->SetMaterialCutsCouple(ghostMCCouple);
414 *(fpHyperStep->GetPostStepPoint()) = *(fGhostPostStepPoint);
415 fpHyperStep->GetPostStepPoint()->SetMaterial(ghostMaterial);
416 fpHyperStep->GetPostStepPoint()->SetMaterialCutsCouple(ghostMCCouple);
417 }
418 else
419 {
420 G4cout << "!!! MaterialCutsCouple is not found for "
421 << ghostMaterial->GetName() << "." << G4endl
422 << " Material in real world ("
423 << realWorldStepPoint->GetMaterial()->GetName()
424 << ") is used." << G4endl;
425 }
426 }
427 }
428}
@ fWorldBoundary
Definition: G4StepStatus.hh:41
G4Region * GetRegion() const
G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
const G4String & GetName() const
Definition: G4Material.hh:172
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 143 of file G4ParallelWorldProcess.hh.

Referenced by G4ParallelWorldProcess().

◆ fFieldTrack

G4FieldTrack G4ParallelWorldProcess::fFieldTrack
protected

Definition at line 158 of file G4ParallelWorldProcess.hh.

Referenced by AlongStepGetPhysicalInteractionLength().

◆ fGhostNavigator

G4Navigator* G4ParallelWorldProcess::fGhostNavigator
protected

◆ fGhostPostStepPoint

G4StepPoint* G4ParallelWorldProcess::fGhostPostStepPoint
protected

◆ fGhostPreStepPoint

G4StepPoint* G4ParallelWorldProcess::fGhostPreStepPoint
protected

◆ fGhostSafety

G4double G4ParallelWorldProcess::fGhostSafety
protected

◆ fGhostStep

G4Step* G4ParallelWorldProcess::fGhostStep
protected

◆ fGhostWorld

G4VPhysicalVolume* G4ParallelWorldProcess::fGhostWorld
protected

Definition at line 153 of file G4ParallelWorldProcess.hh.

Referenced by SetParallelWorld().

◆ fGhostWorldName

G4String G4ParallelWorldProcess::fGhostWorldName
protected

Definition at line 152 of file G4ParallelWorldProcess.hh.

Referenced by G4ParallelWorldProcess(), and SetParallelWorld().

◆ fNavigatorID

G4int G4ParallelWorldProcess::fNavigatorID
protected

◆ fNewGhostTouchable

G4TouchableHandle G4ParallelWorldProcess::fNewGhostTouchable
protected

Definition at line 157 of file G4ParallelWorldProcess.hh.

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

◆ fOldGhostTouchable

G4TouchableHandle G4ParallelWorldProcess::fOldGhostTouchable
protected

Definition at line 156 of file G4ParallelWorldProcess.hh.

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

◆ fOnBoundary

G4bool G4ParallelWorldProcess::fOnBoundary
protected

◆ fPathFinder

G4PathFinder* G4ParallelWorldProcess::fPathFinder
protected

◆ fTransportationManager

G4TransportationManager* G4ParallelWorldProcess::fTransportationManager
protected

◆ layeredMaterialFlag

G4bool G4ParallelWorldProcess::layeredMaterialFlag
protected

◆ xParticleChange

G4ParticleChange G4ParallelWorldProcess::xParticleChange
protected

Definition at line 144 of file G4ParallelWorldProcess.hh.


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