Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParallelWorldScoringProcess Class Reference

#include <G4ParallelWorldScoringProcess.hh>

Inheritance diagram for G4ParallelWorldScoringProcess:
Collaboration diagram for G4ParallelWorldScoringProcess:

Public Member Functions

 G4ParallelWorldScoringProcess (const G4String &processName="ParaWorldScore", G4ProcessType theType=fParameterisation)
 
virtual ~G4ParallelWorldScoringProcess ()
 
void SetParallelWorld (G4String parallelWorldName)
 
void SetParallelWorld (G4VPhysicalVolume *parallelWorld)
 
G4bool IsAtRestRequired (G4ParticleDefinition *partDef)
 
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 &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
void Verbose (const G4Step &) const
 
- 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 &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- 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 68 of file G4ParallelWorldScoringProcess.hh.

Constructor & Destructor Documentation

G4ParallelWorldScoringProcess::G4ParallelWorldScoringProcess ( const G4String processName = "ParaWorldScore",
G4ProcessType  theType = fParameterisation 
)

Definition at line 54 of file G4ParallelWorldScoringProcess.cc.

55 :G4VProcess(processName,theType), fGhostNavigator(0), fNavigatorID(-1), fFieldTrack('0')
56 {
57  pParticleChange = &aDummyParticleChange;
58 
59  fGhostStep = new G4Step();
60  fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
61  fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
62 
63  fTransportationManager = G4TransportationManager::GetTransportationManager();
64  fPathFinder = G4PathFinder::GetInstance();
65 
66  fGhostWorld = 0;
67  fGhostSafety = 0.;
68  fOnBoundary = false;
69 
70  if (verboseLevel>0)
71  {
72  G4cout << GetProcessName() << " is created " << G4endl;
73  }
74 }
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
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
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4ParallelWorldScoringProcess::~G4ParallelWorldScoringProcess ( )
virtual

Definition at line 79 of file G4ParallelWorldScoringProcess.cc.

80 {
81  delete fGhostStep;
82 }

Member Function Documentation

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

Implements G4VProcess.

Definition at line 394 of file G4ParallelWorldScoringProcess.cc.

396 {
397  // Dummy ParticleChange ie: does nothing
398  // Expecting G4Transportation to move the track
399  pParticleChange->Initialize(track);
400  return pParticleChange;
401 }
virtual void Initialize(const G4Track &)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

Here is the call graph for this function:

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

Implements G4VProcess.

Definition at line 328 of file G4ParallelWorldScoringProcess.cc.

331 {
332  static G4ThreadLocal G4FieldTrack *endTrack_G4MT_TLS_ = 0 ; if (!endTrack_G4MT_TLS_) endTrack_G4MT_TLS_ = new G4FieldTrack ('0') ; G4FieldTrack &endTrack = *endTrack_G4MT_TLS_;
333  static G4ThreadLocal ELimited *eLimited_G4MT_TLS_ = 0 ; if (!eLimited_G4MT_TLS_) eLimited_G4MT_TLS_ = new ELimited ; ELimited &eLimited = *eLimited_G4MT_TLS_;
334 
335  *selection = NotCandidateForSelection;
336  G4double returnedStep = DBL_MAX;
337 
338  if (previousStepSize > 0.)
339  { fGhostSafety -= previousStepSize; }
340 // else
341 // { fGhostSafety = -1.; }
342  if (fGhostSafety < 0.) fGhostSafety = 0.0;
343 
344  // ------------------------------------------
345  // Determination of the proposed STEP LENGTH:
346  // ------------------------------------------
347  if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
348  {
349  // I have no chance to limit
350  returnedStep = currentMinimumStep;
351  fOnBoundary = false;
352  proposedSafety = fGhostSafety - currentMinimumStep;
353  }
354  else // (currentMinimumStep > fGhostSafety: I may limit the Step)
355  {
356  G4FieldTrackUpdator::Update(&fFieldTrack,&track);
357 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
358 // ComputeStep
359 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
360  returnedStep
361  = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
362  track.GetCurrentStepNumber(),fGhostSafety,eLimited,
363  endTrack,track.GetVolume());
364 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
365  if(eLimited == kDoNot)
366  {
367  // Track is not on the boundary
368  fOnBoundary = false;
369  fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
370  }
371  else
372  {
373  // Track is on the boundary
374  fOnBoundary = true;
375  // fEndGhostSafety = 0.0; // Will apply at the end of the step ...
376  }
377  proposedSafety = fGhostSafety;
378  if(eLimited == kUnique || eLimited == kSharedOther) {
379  *selection = CandidateForSelection;
380  }else if (eLimited == kSharedTransport) {
381  returnedStep *= (1.0 + 1.0e-9);
382  // Expand to disable its selection in Step Manager comparison
383  }
384  }
385 
386  // ----------------------------------------------
387  // Returns the fGhostSafety as the proposedSafety
388  // The SteppingManager will take care of keeping
389  // the smallest one.
390  // ----------------------------------------------
391  return returnedStep;
392 }
ELimited
#define G4ThreadLocal
Definition: tls.hh:89
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 *)

Here is the call graph for this function:

G4VParticleChange * G4ParallelWorldScoringProcess::AtRestDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VProcess.

Definition at line 202 of file G4ParallelWorldScoringProcess.cc.

205 {
206  fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
207  G4VSensitiveDetector* aSD = 0;
208  if(fOldGhostTouchable->GetVolume())
209  { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
210  fOnBoundary = false;
211  CopyStep(step);
212  fGhostPreStepPoint->SetSensitiveDetector(aSD);
213 
214  fNewGhostTouchable = fOldGhostTouchable;
215 
216  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
217  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
218  if(fNewGhostTouchable->GetVolume())
219  {
220  fGhostPostStepPoint->SetSensitiveDetector(
221  fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
222  }
223  else
224  { fGhostPostStepPoint->SetSensitiveDetector(0); }
225 
226  if (verboseLevel>1) Verbose(step);
227 
228  G4VSensitiveDetector* sd = fGhostPreStepPoint->GetSensitiveDetector();
229  if(sd)
230  {
231  sd->Hit(fGhostStep);
232  }
233 
234  pParticleChange->Initialize(track);
235  return pParticleChange;
236 }
virtual void Initialize(const G4Track &)
G4int verboseLevel
Definition: G4VProcess.hh:368
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
G4VSensitiveDetector * GetSensitiveDetector() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4VSensitiveDetector * GetSensitiveDetector() const
const G4TouchableHandle & GetTouchableHandle() const

Here is the call graph for this function:

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

Implements G4VProcess.

Definition at line 189 of file G4ParallelWorldScoringProcess.cc.

192 {
193  *condition = Forced;
194  return DBL_MAX;
195 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
G4bool G4ParallelWorldScoringProcess::IsAtRestRequired ( G4ParticleDefinition partDef)

Definition at line 114 of file G4ParallelWorldScoringProcess.cc.

115 {
116  G4int pdgCode = partDef->GetPDGEncoding();
117  if(pdgCode==0)
118  {
119  G4String partName = partDef->GetParticleName();
120  if(partName=="opticalphoton") return false;
121  if(partName=="geantino") return false;
122  if(partName=="chargedgeantino") return false;
123  }
124  else
125  {
126  if(pdgCode==22) return false; // gamma
127  if(pdgCode==11) return false; // electron
128  if(pdgCode==2212) return false; // proton
129  if(pdgCode==-12) return false; // anti_nu_e
130  if(pdgCode==12) return false; // nu_e
131  if(pdgCode==-14) return false; // anti_nu_mu
132  if(pdgCode==14) return false; // nu_mu
133  if(pdgCode==-16) return false; // anti_nu_tau
134  if(pdgCode==16) return false; // nu_tau
135  }
136  return true;
137 }
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const

Here is the call graph for this function:

G4VParticleChange * G4ParallelWorldScoringProcess::PostStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4VProcess.

Definition at line 259 of file G4ParallelWorldScoringProcess.cc.

262 {
263  fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
264  G4VSensitiveDetector* aSD = 0;
265  if(fOldGhostTouchable->GetVolume())
266  { aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector(); }
267  CopyStep(step);
268  fGhostPreStepPoint->SetSensitiveDetector(aSD);
269 
270  // fPathFinder->Locate( track.GetPosition(),
271  // track.GetMomentumDirection(),
272  // true);
273 
274  // fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
275  // step.GetPostStepPoint()->GetMomentumDirection());
276 
277  if(fOnBoundary)
278  {
279 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
280 // Locate the point and get new touchable
281 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
282  //?? fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
283  //?? step.GetPostStepPoint()->GetMomentumDirection());
284  fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
285 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
286  }
287  else
288  {
289 // Do I need this ??????????????????????????????????????????????????????????
290 // fGhostNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
291 // ?????????????????????????????????????????????????????????????????????????
292 
293  // fPathFinder->ReLocate(track.GetPosition());
294 
295  // reuse the touchable
296  fNewGhostTouchable = fOldGhostTouchable;
297  }
298 
299  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
300  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
301 
302  if(fNewGhostTouchable->GetVolume())
303  {
304  fGhostPostStepPoint->SetSensitiveDetector(
305  fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
306  }
307  else
308  { fGhostPostStepPoint->SetSensitiveDetector(0); }
309 
310  if (verboseLevel>1) Verbose(step);
311 
312  G4VSensitiveDetector* sd = fGhostPreStepPoint->GetSensitiveDetector();
313  if(sd)
314  {
315  sd->Hit(fGhostStep);
316  }
317 
318  pParticleChange->Initialize(track); // Does not change the track properties
319  return pParticleChange;
320 }
virtual void Initialize(const G4Track &)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4TouchableHandle CreateTouchableHandle(G4int navId) const
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
G4VSensitiveDetector * GetSensitiveDetector() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4VSensitiveDetector * GetSensitiveDetector() const
const G4TouchableHandle & GetTouchableHandle() const

Here is the call graph for this function:

G4double G4ParallelWorldScoringProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 244 of file G4ParallelWorldScoringProcess.cc.

248 {
249  // I must be invoked anyway to score the hit.
251  return DBL_MAX;
252 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
void G4ParallelWorldScoringProcess::SetParallelWorld ( G4String  parallelWorldName)

Definition at line 90 of file G4ParallelWorldScoringProcess.cc.

91 {
92 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
93 // Get pointers of the parallel world and its navigator
94 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
95  fGhostWorldName = parallelWorldName;
96  fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
97  fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
98 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
99 }
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
G4Navigator * GetNavigator(const G4String &worldName)

Here is the call graph for this function:

void G4ParallelWorldScoringProcess::SetParallelWorld ( G4VPhysicalVolume parallelWorld)

Definition at line 102 of file G4ParallelWorldScoringProcess.cc.

103 {
104 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
105 // Get pointer of navigator
106 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
107  fGhostWorldName = parallelWorld->GetName();
108  fGhostWorld = parallelWorld;
109  fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
110 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
111 }
const G4String & GetName() const
G4Navigator * GetNavigator(const G4String &worldName)

Here is the call graph for this function:

void G4ParallelWorldScoringProcess::StartTracking ( G4Track trk)
virtual

Reimplemented from G4VProcess.

Definition at line 145 of file G4ParallelWorldScoringProcess.cc.

146 {
147 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
148 // Activate navigator and get the navigator ID
149 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
150 // G4cout << " G4ParallelWorldScoringProcess::StartTracking" << G4endl;
151  if(fGhostNavigator)
152  { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
153  else
154  {
155  G4Exception("G4ParallelWorldScoringProcess::StartTracking",
156  "ProcParaWorld000",FatalException,
157  "G4ParallelWorldScoringProcess is used for tracking without having a parallel world assigned");
158  }
159 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
160 
161 // G4cout << "G4ParallelWorldScoringProcess::StartTracking <<<<<<<<<<<<<<<<<< " << G4endl;
162 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
163 // Let PathFinder initialize
164 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
165  fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
166 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
167 
168 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
169 // Setup initial touchables for the first step
170 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
171  fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
172  fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
173  fNewGhostTouchable = fOldGhostTouchable;
174  fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
175 
176  // Initialize
177  fGhostSafety = -1.;
178  fOnBoundary = false;
179  fGhostPreStepPoint->SetStepStatus(fUndefined);
180  fGhostPostStepPoint->SetStepStatus(fUndefined);
181 }
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=0)
const G4ThreeVector & GetPosition() const
G4TouchableHandle CreateTouchableHandle(G4int navId) const
void SetStepStatus(const G4StepStatus aValue)
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
void SetTouchableHandle(const G4TouchableHandle &apValue)

Here is the call graph for this function:

void G4ParallelWorldScoringProcess::Verbose ( const G4Step step) const

Definition at line 428 of file G4ParallelWorldScoringProcess.cc.

429 {
430  G4cout << "In mass geometry ------------------------------------------------" << G4endl;
431  G4cout << " StepLength : " << step.GetStepLength()/mm << " TotalEnergyDeposit : "
432  << step.GetTotalEnergyDeposit()/MeV << G4endl;
433  G4cout << " PreStepPoint : "
434  << step.GetPreStepPoint()->GetPhysicalVolume()->GetName() << " - ";
437  else
438  { G4cout << "NoProcessAssigned"; }
439  G4cout << G4endl;
440  G4cout << " " << step.GetPreStepPoint()->GetPosition() << G4endl;
441  G4cout << " PostStepPoint : ";
442  if(step.GetPostStepPoint()->GetPhysicalVolume())
443  { G4cout << step.GetPostStepPoint()->GetPhysicalVolume()->GetName(); }
444  else
445  { G4cout << "OutOfWorld"; }
446  G4cout << " - ";
449  else
450  { G4cout << "NoProcessAssigned"; }
451  G4cout << G4endl;
452  G4cout << " " << step.GetPostStepPoint()->GetPosition() << G4endl;
453 
454  G4cout << "In ghost geometry ------------------------------------------------" << G4endl;
455  G4cout << " StepLength : " << fGhostStep->GetStepLength()/mm
456  << " TotalEnergyDeposit : "
457  << fGhostStep->GetTotalEnergyDeposit()/MeV << G4endl;
458  G4cout << " PreStepPoint : "
459  << fGhostStep->GetPreStepPoint()->GetPhysicalVolume()->GetName() << " ["
460  << fGhostStep->GetPreStepPoint()->GetTouchable()->GetReplicaNumber()
461  << " ]" << " - ";
462  if(fGhostStep->GetPreStepPoint()->GetProcessDefinedStep())
463  { G4cout << fGhostStep->GetPreStepPoint()->GetProcessDefinedStep()->GetProcessName(); }
464  else
465  { G4cout << "NoProcessAssigned"; }
466  G4cout << G4endl;
467  G4cout << " " << fGhostStep->GetPreStepPoint()->GetPosition() << G4endl;
468  G4cout << " PostStepPoint : ";
469  if(fGhostStep->GetPostStepPoint()->GetPhysicalVolume())
470  {
471  G4cout << fGhostStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() << " ["
472  << fGhostStep->GetPostStepPoint()->GetTouchable()->GetReplicaNumber()
473  << " ]";
474  }
475  else
476  { G4cout << "OutOfWorld"; }
477  G4cout << " - ";
478  if(fGhostStep->GetPostStepPoint()->GetProcessDefinedStep())
479  { G4cout << fGhostStep->GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName(); }
480  else
481  { G4cout << "NoProcessAssigned"; }
482  G4cout << G4endl;
483  G4cout << " " << fGhostStep->GetPostStepPoint()->GetPosition() << " == "
484  << fGhostStep->GetTrack()->GetMomentumDirection()
485  << G4endl;
486 
487 }
static constexpr double mm
Definition: G4SIunits.hh:115
G4double GetStepLength() const
const G4VTouchable * GetTouchable() const
G4StepPoint * GetPreStepPoint() const
G4GLOB_DLL std::ostream G4cout
G4VPhysicalVolume * GetPhysicalVolume() const
const G4String & GetName() const
const G4ThreeVector & GetPosition() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4double GetTotalEnergyDeposit() const
const G4VProcess * GetProcessDefinedStep() const
const G4ThreeVector & GetMomentumDirection() const
G4StepPoint * GetPostStepPoint() const
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:58
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4Track * GetTrack() const

Here is the call graph for this function:

Here is the caller graph for this function:


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