Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Member Functions
G4ParallelWorldProcess Class Reference

#include <G4ParallelWorldProcess.hh>

Inheritance diagram for G4ParallelWorldProcess:
G4VProcess

Public Member Functions

 G4ParallelWorldProcess (const G4String &processName="ParaWorld", G4ProcessType theType=fParameterisation)
 
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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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 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
 
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)
 

Additional Inherited Members

- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 73 of file G4ParallelWorldProcess.hh.

Constructor & Destructor Documentation

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

Definition at line 61 of file G4ParallelWorldProcess.cc.

References G4cout, G4endl, G4PathFinder::GetInstance(), G4TransportationManager::GetNavigatorForTracking(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4VProcess::GetProcessName(), G4TransportationManager::GetTransportationManager(), G4VProcess::pParticleChange, G4Navigator::SetPushVerbosity(), and G4VProcess::verboseLevel.

62 :G4VProcess(processName,theType), fGhostNavigator(0), fNavigatorID(-1),
63  fFieldTrack('0'),layeredMaterialFlag(false)
64 {
65  if(!fpHyperStep) fpHyperStep = new G4Step();
66  iParallelWorld = ++nParallelWorlds;
67 
68  pParticleChange = &aDummyParticleChange;
69 
70  fGhostStep = new G4Step();
71  fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
72  fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
73 
74  fTransportationManager = G4TransportationManager::GetTransportationManager();
75  fTransportationManager->GetNavigatorForTracking()->SetPushVerbosity(false);
76  fPathFinder = G4PathFinder::GetInstance();
77 
78  if (verboseLevel>0)
79  {
80  G4cout << GetProcessName() << " is created " << G4endl;
81  }
82 }
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:57
G4int verboseLevel
Definition: G4VProcess.hh:368
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition: G4VProcess.cc:52
G4Navigator * GetNavigatorForTracking() const
G4StepPoint * GetPreStepPoint() const
G4GLOB_DLL std::ostream G4cout
Definition: G4Step.hh:76
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
static G4TransportationManager * GetTransportationManager()
G4StepPoint * GetPostStepPoint() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
void SetPushVerbosity(G4bool mode)
#define G4endl
Definition: G4ios.hh:61
G4ParallelWorldProcess::~G4ParallelWorldProcess ( )
virtual

Definition at line 84 of file G4ParallelWorldProcess.cc.

85 {
86  delete fGhostStep;
87  nParallelWorlds--;
88  if(nParallelWorlds==0)
89  {
90  delete fpHyperStep;
91  fpHyperStep = 0;
92  }
93 }

Member Function Documentation

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

Implements G4VProcess.

Definition at line 314 of file G4ParallelWorldProcess.cc.

References G4VParticleChange::Initialize(), and G4VProcess::pParticleChange.

316 {
317  pParticleChange->Initialize(track);
318  return pParticleChange;
319 }
virtual void Initialize(const G4Track &)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4double G4ParallelWorldProcess::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
)
virtual

Implements G4VProcess.

Definition at line 259 of file G4ParallelWorldProcess.cc.

References CandidateForSelection, G4Navigator::ComputeSafety(), G4PathFinder::ComputeStep(), DBL_MAX, G4ThreadLocal, G4Track::GetCurrentStepNumber(), G4Track::GetVolume(), kDoNot, kSharedOther, kSharedTransport, kUndefLimited, kUnique, NotCandidateForSelection, and G4FieldTrackUpdator::Update().

262 {
263  static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ; if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ; G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
264  //static ELimited eLimited;
265  ELimited eLimited;
266  ELimited eLim = kUndefLimited;
267 
268  *selection = NotCandidateForSelection;
269  G4double returnedStep = DBL_MAX;
270 
271  if (previousStepSize > 0.)
272  { fGhostSafety -= previousStepSize; }
273  if (fGhostSafety < 0.) fGhostSafety = 0.0;
274 
275  if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
276  {
277  // I have no chance to limit
278  returnedStep = currentMinimumStep;
279  fOnBoundary = false;
280  proposedSafety = fGhostSafety - currentMinimumStep;
281  eLim = kDoNot;
282  }
283  else
284  {
285  G4FieldTrackUpdator::Update(&fFieldTrack,&track);
286  returnedStep
287  = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
288  track.GetCurrentStepNumber(),fGhostSafety,eLimited,
289  endTrack,track.GetVolume());
290  if(eLimited == kDoNot)
291  {
292  fOnBoundary = false;
293  fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
294  }
295  else
296  {
297  fOnBoundary = true;
298  }
299  proposedSafety = fGhostSafety;
300  if(eLimited == kUnique || eLimited == kSharedOther) {
301  *selection = CandidateForSelection;
302  }
303  else if (eLimited == kSharedTransport) {
304  returnedStep *= (1.0 + 1.0e-9);
305  }
306  eLim = eLimited;
307  }
308 
309  if(iParallelWorld==nParallelWorlds) fNavIDHyp = 0;
310  if(eLim == kUnique || eLim == kSharedOther) fNavIDHyp = fNavigatorID;
311  return returnedStep;
312 }
ELimited
#define G4ThreadLocal
Definition: tls.hh:52
G4int GetCurrentStepNumber() const
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
G4VPhysicalVolume * GetVolume() const
double G4double
Definition: G4Types.hh:76
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
#define DBL_MAX
Definition: templates.hh:83
static void Update(G4FieldTrack *, const G4Track *)
G4VParticleChange * G4ParallelWorldProcess::AtRestDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VProcess.

Definition at line 165 of file G4ParallelWorldProcess.cc.

References G4VPhysicalVolume::GetLogicalVolume(), G4LogicalVolume::GetSensitiveDetector(), G4StepPoint::GetTouchableHandle(), G4VTouchable::GetVolume(), G4VSensitiveDetector::Hit(), G4VParticleChange::Initialize(), G4VProcess::pParticleChange, G4StepPoint::SetSensitiveDetector(), and G4StepPoint::SetTouchableHandle().

168 {
169 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
170 // At Rest must be registered ONLY for the particle which has other At Rest
171 // process(es).
172 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
173  fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
174  G4VSensitiveDetector* aSD = 0;
175  if(fOldGhostTouchable->GetVolume())
176  { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
177  fOnBoundary = false;
178  if(aSD)
179  {
180  CopyStep(step);
181  fGhostPreStepPoint->SetSensitiveDetector(aSD);
182 
183  fNewGhostTouchable = fOldGhostTouchable;
184 
185  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
186  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
187  if(fNewGhostTouchable->GetVolume())
188  {
189  fGhostPostStepPoint->SetSensitiveDetector(
190  fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
191  }
192  else
193  { fGhostPostStepPoint->SetSensitiveDetector(0); }
194 
195  aSD->Hit(fGhostStep);
196  }
197 
198  pParticleChange->Initialize(track);
199  return pParticleChange;
200 }
virtual void Initialize(const G4Track &)
void SetSensitiveDetector(G4VSensitiveDetector *)
G4bool Hit(G4Step *aStep)
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4VSensitiveDetector * GetSensitiveDetector() const
const G4TouchableHandle & GetTouchableHandle() const
G4double G4ParallelWorldProcess::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 153 of file G4ParallelWorldProcess.cc.

References DBL_MAX, and Forced.

156 {
157 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
158 // At Rest must be registered ONLY for the particle which has other At Rest
159 // process(es).
160 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
161  *condition = Forced;
162  return DBL_MAX;
163 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
const G4Step * G4ParallelWorldProcess::GetHyperStep ( )
static

Definition at line 55 of file G4ParallelWorldProcess.cc.

Referenced by G4OpBoundaryProcess::PostStepDoIt().

56 { return fpHyperStep; }
G4int G4ParallelWorldProcess::GetHypNavigatorID ( )
static

Definition at line 57 of file G4ParallelWorldProcess.cc.

Referenced by G4OpBoundaryProcess::PostStepDoIt().

58 { return fNavIDHyp; }
G4bool G4ParallelWorldProcess::GetLayeredMaterialFlag ( ) const
inline

Definition at line 127 of file G4ParallelWorldProcess.hh.

128  { return layeredMaterialFlag; }
G4bool G4ParallelWorldProcess::IsAtRestRequired ( G4ParticleDefinition partDef)

Definition at line 400 of file G4ParallelWorldProcess.cc.

References G4ParticleDefinition::GetParticleName(), and G4ParticleDefinition::GetPDGEncoding().

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

401 {
402  G4int pdgCode = partDef->GetPDGEncoding();
403  if(pdgCode==0)
404  {
405  G4String partName = partDef->GetParticleName();
406  if(partName=="opticalphoton") return false;
407  if(partName=="geantino") return false;
408  if(partName=="chargedgeantino") return false;
409  }
410  else
411  {
412  if(pdgCode==22) return false; // gamma
413  if(pdgCode==11) return false; // electron
414  if(pdgCode==2212) return false; // proton
415  if(pdgCode==-12) return false; // anti_nu_e
416  if(pdgCode==12) return false; // nu_e
417  if(pdgCode==-14) return false; // anti_nu_mu
418  if(pdgCode==14) return false; // nu_mu
419  if(pdgCode==-16) return false; // anti_nu_tau
420  if(pdgCode==16) return false; // nu_tau
421  }
422  return true;
423 }
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4VParticleChange * G4ParallelWorldProcess::PostStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VProcess.

Definition at line 212 of file G4ParallelWorldProcess.cc.

References G4PathFinder::CreateTouchableHandle(), G4VPhysicalVolume::GetLogicalVolume(), G4StepPoint::GetSensitiveDetector(), G4LogicalVolume::GetSensitiveDetector(), G4Track::GetStep(), G4StepPoint::GetTouchableHandle(), G4VTouchable::GetVolume(), G4VSensitiveDetector::Hit(), G4VParticleChange::Initialize(), G4VProcess::pParticleChange, G4StepPoint::SetSensitiveDetector(), and G4StepPoint::SetTouchableHandle().

215 {
216  fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
217  G4VSensitiveDetector* aSD = 0;
218  if(fOldGhostTouchable->GetVolume())
219  { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
220  CopyStep(step);
221  fGhostPreStepPoint->SetSensitiveDetector(aSD);
222 
223  if(fOnBoundary)
224  {
225  fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
226  }
227  else
228  {
229  fNewGhostTouchable = fOldGhostTouchable;
230  }
231 
232  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
233  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
234 
235  if(fNewGhostTouchable->GetVolume())
236  {
237  fGhostPostStepPoint->SetSensitiveDetector(
238  fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
239  }
240  else
241  { fGhostPostStepPoint->SetSensitiveDetector(0); }
242 
243  G4VSensitiveDetector* sd = fGhostPreStepPoint->GetSensitiveDetector();
244  if(sd)
245  {
246  sd->Hit(fGhostStep);
247  }
248 
249  pParticleChange->Initialize(track);
250  if(layeredMaterialFlag)
251  {
252  G4StepPoint* realWorldPostStepPoint =
253  ((G4Step*)(track.GetStep()))->GetPostStepPoint();
254  SwitchMaterial(realWorldPostStepPoint);
255  }
256  return pParticleChange;
257 }
virtual void Initialize(const G4Track &)
G4TouchableHandle CreateTouchableHandle(G4int navId) const
const G4Step * GetStep() const
void SetSensitiveDetector(G4VSensitiveDetector *)
G4bool Hit(G4Step *aStep)
Definition: G4Step.hh:76
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4VSensitiveDetector * GetSensitiveDetector() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4VSensitiveDetector * GetSensitiveDetector() const
const G4TouchableHandle & GetTouchableHandle() const
G4double G4ParallelWorldProcess::PostStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 203 of file G4ParallelWorldProcess.cc.

References DBL_MAX, and StronglyForced.

207 {
209  return DBL_MAX;
210 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
void G4ParallelWorldProcess::SetLayeredMaterialFlag ( G4bool  flg = true)
inline

Definition at line 125 of file G4ParallelWorldProcess.hh.

Referenced by G4ParallelWorldPhysics::ConstructProcess().

126  { layeredMaterialFlag = flg; }
void G4ParallelWorldProcess::SetParallelWorld ( G4String  parallelWorldName)

Definition at line 96 of file G4ParallelWorldProcess.cc.

References G4TransportationManager::GetNavigator(), G4TransportationManager::GetParallelWorld(), and G4Navigator::SetPushVerbosity().

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

97 {
98  fGhostWorldName = parallelWorldName;
99  fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
100  fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
101  fGhostNavigator->SetPushVerbosity(false);
102 }
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
G4Navigator * GetNavigator(const G4String &worldName)
void SetPushVerbosity(G4bool mode)
void G4ParallelWorldProcess::SetParallelWorld ( G4VPhysicalVolume parallelWorld)

Definition at line 105 of file G4ParallelWorldProcess.cc.

References G4VPhysicalVolume::GetName(), G4TransportationManager::GetNavigator(), and G4Navigator::SetPushVerbosity().

106 {
107  fGhostWorldName = parallelWorld->GetName();
108  fGhostWorld = parallelWorld;
109  fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
110  fGhostNavigator->SetPushVerbosity(false);
111 }
const G4String & GetName() const
G4Navigator * GetNavigator(const G4String &worldName)
void SetPushVerbosity(G4bool mode)
void G4ParallelWorldProcess::StartTracking ( G4Track trk)
virtual

Reimplemented from G4VProcess.

Definition at line 113 of file G4ParallelWorldProcess.cc.

References G4TransportationManager::ActivateNavigator(), G4PathFinder::CreateTouchableHandle(), FatalException, fUndefined, G4Exception(), G4Track::GetMomentumDirection(), G4Track::GetPosition(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4Track::GetStep(), G4PathFinder::PrepareNewTrack(), G4StepPoint::SetStepStatus(), and G4StepPoint::SetTouchableHandle().

114 {
115  if(fGhostNavigator)
116  { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
117  else
118  {
119  G4Exception("G4ParallelWorldProcess::StartTracking",
120  "ProcParaWorld000",FatalException,
121  "G4ParallelWorldProcess is used for tracking without having a parallel world assigned");
122  }
123  fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
124 
125  fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
126  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
127  fNewGhostTouchable = fOldGhostTouchable;
128  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
129 
130  fGhostSafety = -1.;
131  fOnBoundary = false;
132  fGhostPreStepPoint->SetStepStatus(fUndefined);
133  fGhostPostStepPoint->SetStepStatus(fUndefined);
134 
135 // G4VPhysicalVolume* thePhys = fNewGhostTouchable->GetVolume();
136 // if(thePhys)
137 // {
138 // G4Material* ghostMaterial = thePhys->GetLogicalVolume()->GetMaterial();
139 // if(ghostMaterial)
140 // { G4cout << " --- Material : " << ghostMaterial->GetName() << G4endl; }
141 // }
142 
143  *(fpHyperStep->GetPostStepPoint()) = *(trk->GetStep()->GetPostStepPoint());
144  if(layeredMaterialFlag)
145  {
146  G4StepPoint* realWorldPostStepPoint = trk->GetStep()->GetPostStepPoint();
147  SwitchMaterial(realWorldPostStepPoint);
148  }
149  *(fpHyperStep->GetPreStepPoint()) = *(fpHyperStep->GetPostStepPoint());
150 }
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=0)
const G4ThreeVector & GetPosition() const
G4TouchableHandle CreateTouchableHandle(G4int navId) const
const G4Step * GetStep() const
void SetStepStatus(const G4StepStatus aValue)
G4StepPoint * GetPreStepPoint() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4int ActivateNavigator(G4Navigator *aNavigator)
const G4ThreeVector & GetMomentumDirection() const
G4StepPoint * GetPostStepPoint() const
void SetTouchableHandle(const G4TouchableHandle &apValue)

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