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

#include <G4ParallelWorldScoringProcess.hh>

+ Inheritance diagram for G4ParallelWorldScoringProcess:

Public Member Functions

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

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- 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 52 of file G4ParallelWorldScoringProcess.hh.

Constructor & Destructor Documentation

◆ G4ParallelWorldScoringProcess()

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

Definition at line 50 of file G4ParallelWorldScoringProcess.cc.

52 : G4VProcess(processName, theType), fFieldTrack('0')
53{
54 pParticleChange = &aDummyParticleChange;
55
56 fGhostStep = new G4Step();
57 fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
58 fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
59
61 fPathFinder = G4PathFinder::GetInstance();
62
63 fGhostWorld = nullptr;
64 fGhostSafety = 0.;
65 fOnBoundary = false;
66
67 if (verboseLevel > 0) {
68 G4cout << GetProcessName() << " is created " << G4endl;
69 }
70}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4PathFinder * GetInstance()
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
static G4TransportationManager * GetTransportationManager()
G4int verboseLevel
G4VParticleChange * pParticleChange
const G4String & GetProcessName() const

◆ ~G4ParallelWorldScoringProcess()

G4ParallelWorldScoringProcess::~G4ParallelWorldScoringProcess ( )
override

Definition at line 75 of file G4ParallelWorldScoringProcess.cc.

76{
77 delete fGhostStep;
78}

Member Function Documentation

◆ AlongStepDoIt()

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

Implements G4VProcess.

Definition at line 335 of file G4ParallelWorldScoringProcess.cc.

336{
337 // Dummy ParticleChange ie: does nothing
338 // Expecting G4Transportation to move the track
340 return pParticleChange;
341}
virtual void Initialize(const G4Track &)

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 275 of file G4ParallelWorldScoringProcess.cc.

278{
279 static G4ThreadLocal G4FieldTrack* endTrack_G4MT_TLS_ = nullptr;
280 if (endTrack_G4MT_TLS_ == nullptr) endTrack_G4MT_TLS_ = new G4FieldTrack('0');
281 G4FieldTrack& endTrack = *endTrack_G4MT_TLS_;
282 static G4ThreadLocal ELimited* eLimited_G4MT_TLS_ = nullptr;
283 if (eLimited_G4MT_TLS_ == nullptr) eLimited_G4MT_TLS_ = new ELimited;
284 ELimited& eLimited = *eLimited_G4MT_TLS_;
285
286 *selection = NotCandidateForSelection;
287 G4double returnedStep = DBL_MAX;
288
289 if (previousStepSize > 0.) {
290 fGhostSafety -= previousStepSize;
291 }
292 if (fGhostSafety < 0.) fGhostSafety = 0.0;
293
294 // Determination of the proposed STEP LENGTH:
295 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.) {
296 // I have no chance to limit
297 returnedStep = currentMinimumStep;
298 fOnBoundary = false;
299 proposedSafety = fGhostSafety - currentMinimumStep;
300 }
301 else // (currentMinimumStep > fGhostSafety: I may limit the Step)
302 {
303 G4FieldTrackUpdator::Update(&fFieldTrack, &track);
304 // ComputeStep
305 returnedStep = fPathFinder->ComputeStep(fFieldTrack, currentMinimumStep, fNavigatorID,
306 track.GetCurrentStepNumber(), fGhostSafety, eLimited,
307 endTrack, track.GetVolume());
308 if (eLimited == kDoNot) {
309 // Track is not on the boundary
310 fOnBoundary = false;
311 fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
312 }
313 else {
314 // Track is on the boundary
315 fOnBoundary = true;
316 }
317 proposedSafety = fGhostSafety;
318 if (eLimited == kUnique || eLimited == kSharedOther) {
319 *selection = CandidateForSelection;
320 }
321 else if (eLimited == kSharedTransport) {
322 returnedStep *= (1.0 + 1.0e-9);
323 // Expand to disable its selection in Step Manager comparison
324 }
325 }
326
327 // ----------------------------------------------
328 // Returns the fGhostSafety as the proposedSafety
329 // The SteppingManager will take care of keeping
330 // the smallest one.
331 // ----------------------------------------------
332 return returnedStep;
333}
@ CandidateForSelection
@ NotCandidateForSelection
@ kDoNot
@ 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 * G4ParallelWorldScoringProcess::AtRestDoIt ( const G4Track & track,
const G4Step & step )
overridevirtual

Implements G4VProcess.

Definition at line 175 of file G4ParallelWorldScoringProcess.cc.

177{
178 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
179 G4VSensitiveDetector* aSD = nullptr;
180 if (fOldGhostTouchable->GetVolume() != nullptr) {
181 aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector();
182 }
183 fOnBoundary = false;
184 CopyStep(step);
185 fGhostPreStepPoint->SetSensitiveDetector(aSD);
186
187 fNewGhostTouchable = fOldGhostTouchable;
188
189 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
190 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
191 if (fNewGhostTouchable->GetVolume() != nullptr) {
192 fGhostPostStepPoint->SetSensitiveDetector(
193 fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
194 }
195 else {
196 fGhostPostStepPoint->SetSensitiveDetector(nullptr);
197 }
198
199 if (verboseLevel > 1) Verbose(step);
200
201 G4VSensitiveDetector* sd = fGhostPreStepPoint->GetSensitiveDetector();
202 if (sd != nullptr) {
203 sd->Hit(fGhostStep);
204 }
205
207 return pParticleChange;
208}
G4VSensitiveDetector * GetSensitiveDetector() const
void SetSensitiveDetector(G4VSensitiveDetector *)
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4TouchableHandle & GetTouchableHandle() const
G4VSensitiveDetector * GetSensitiveDetector() const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
G4LogicalVolume * GetLogicalVolume() const
G4bool Hit(G4Step *aStep)

◆ AtRestGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 163 of file G4ParallelWorldScoringProcess.cc.

165{
166 *condition = Forced;
167 return DBL_MAX;
168}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced

◆ IsAtRestRequired()

G4bool G4ParallelWorldScoringProcess::IsAtRestRequired ( G4ParticleDefinition * partDef)

Definition at line 107 of file G4ParallelWorldScoringProcess.cc.

108{
109 G4int pdgCode = partDef->GetPDGEncoding();
110 if (pdgCode == 0) {
111 G4String partName = partDef->GetParticleName();
112 if (partName == "geantino") return false;
113 if (partName == "chargedgeantino") return false;
114 }
115 else {
116 if (pdgCode == 11 || pdgCode == 2212) return false; // electrons and proton
117 pdgCode = std::abs(pdgCode);
118 if (pdgCode == 22) return false; // gamma and optical photons
119 if (pdgCode == 12 || pdgCode == 14 || pdgCode == 16) return false; // all neutronos
120 }
121 return true;
122}
int G4int
Definition G4Types.hh:85
const G4String & GetParticleName() const

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 228 of file G4ParallelWorldScoringProcess.cc.

230{
231 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
232 G4VSensitiveDetector* aSD = nullptr;
233 if (fOldGhostTouchable->GetVolume() != nullptr) {
234 aSD = fOldGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector();
235 }
236 CopyStep(step);
237 fGhostPreStepPoint->SetSensitiveDetector(aSD);
238
239 if (fOnBoundary) {
240 // Locate the point and get new touchable
241 fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
242 }
243 else {
244 // reuse the touchable
245 fNewGhostTouchable = fOldGhostTouchable;
246 }
247
248 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
249 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
250
251 if (fNewGhostTouchable->GetVolume() != nullptr) {
252 fGhostPostStepPoint->SetSensitiveDetector(
253 fNewGhostTouchable->GetVolume()->GetLogicalVolume()->GetSensitiveDetector());
254 }
255 else {
256 fGhostPostStepPoint->SetSensitiveDetector(nullptr);
257 }
258
259 if (verboseLevel > 1) Verbose(step);
260
261 G4VSensitiveDetector* sd = fGhostPreStepPoint->GetSensitiveDetector();
262 if (sd != nullptr) {
263 sd->Hit(fGhostStep);
264 }
265
266 pParticleChange->Initialize(track); // Does not change the track properties
267 return pParticleChange;
268}
G4TouchableHandle CreateTouchableHandle(G4int navId) const

◆ PostStepGetPhysicalInteractionLength()

G4double G4ParallelWorldScoringProcess::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 215 of file G4ParallelWorldScoringProcess.cc.

217{
218 // I must be invoked anyway to score the hit.
220 return DBL_MAX;
221}
@ StronglyForced

◆ SetParallelWorld() [1/2]

void G4ParallelWorldScoringProcess::SetParallelWorld ( G4String parallelWorldName)

Definition at line 85 of file G4ParallelWorldScoringProcess.cc.

86{
87 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
88 // Get pointers of the parallel world and its navigator
89 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
90 fGhostWorldName = parallelWorldName;
91 fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
92 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
93 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
94}
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
G4Navigator * GetNavigator(const G4String &worldName)

◆ SetParallelWorld() [2/2]

void G4ParallelWorldScoringProcess::SetParallelWorld ( G4VPhysicalVolume * parallelWorld)

Definition at line 96 of file G4ParallelWorldScoringProcess.cc.

97{
98 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
99 // Get pointer of navigator
100 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
101 fGhostWorldName = parallelWorld->GetName();
102 fGhostWorld = parallelWorld;
103 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
104 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
105}
const G4String & GetName() const

◆ StartTracking()

void G4ParallelWorldScoringProcess::StartTracking ( G4Track * trk)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 129 of file G4ParallelWorldScoringProcess.cc.

130{
131 // Activate navigator and get the navigator ID
132 if (fGhostNavigator != nullptr) {
133 fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator);
134 }
135 else {
136 G4Exception("G4ParallelWorldScoringProcess::StartTracking", "ProcParaWorld000", FatalException,
137 "G4ParallelWorldScoringProcess is used for tracking without having a parallel "
138 "world assigned");
139 }
140
141 // Let PathFinder initialize
142 fPathFinder->PrepareNewTrack(trk->GetPosition(), trk->GetMomentumDirection());
143
144 // Setup initial touchables for the first step
145 fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
146 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
147 fNewGhostTouchable = fOldGhostTouchable;
148 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
149
150 // Initialize
151 fGhostSafety = -1.;
152 fOnBoundary = false;
153 fGhostPreStepPoint->SetStepStatus(fUndefined);
154 fGhostPostStepPoint->SetStepStatus(fUndefined);
155}
@ 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 SetStepStatus(const G4StepStatus aValue)
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4int ActivateNavigator(G4Navigator *aNavigator)

◆ Verbose()

void G4ParallelWorldScoringProcess::Verbose ( const G4Step & step) const

Definition at line 366 of file G4ParallelWorldScoringProcess.cc.

367{
368 G4cout << "In mass geometry ------------------------------------------------" << G4endl;
369 G4cout << " StepLength : " << step.GetStepLength() / mm
370 << " TotalEnergyDeposit : " << step.GetTotalEnergyDeposit() / MeV << G4endl;
371 G4cout << " PreStepPoint : " << step.GetPreStepPoint()->GetPhysicalVolume()->GetName() << " - ";
372 if (step.GetPreStepPoint()->GetProcessDefinedStep() != nullptr) {
374 }
375 else {
376 G4cout << "NoProcessAssigned";
377 }
378 G4cout << G4endl;
379 G4cout << " " << step.GetPreStepPoint()->GetPosition() << G4endl;
380 G4cout << " PostStepPoint : ";
381 if (step.GetPostStepPoint()->GetPhysicalVolume() != nullptr) {
383 }
384 else {
385 G4cout << "OutOfWorld";
386 }
387 G4cout << " - ";
388 if (step.GetPostStepPoint()->GetProcessDefinedStep() != nullptr) {
390 }
391 else {
392 G4cout << "NoProcessAssigned";
393 }
394 G4cout << G4endl;
395 G4cout << " " << step.GetPostStepPoint()->GetPosition() << G4endl;
396
397 G4cout << "In ghost geometry ------------------------------------------------" << G4endl;
398 G4cout << " StepLength : " << fGhostStep->GetStepLength() / mm
399 << " TotalEnergyDeposit : " << fGhostStep->GetTotalEnergyDeposit() / MeV << G4endl;
400 G4cout << " PreStepPoint : " << fGhostStep->GetPreStepPoint()->GetPhysicalVolume()->GetName()
401 << " [" << fGhostStep->GetPreStepPoint()->GetTouchable()->GetReplicaNumber() << " ]"
402 << " - ";
403 if (fGhostStep->GetPreStepPoint()->GetProcessDefinedStep() != nullptr) {
405 }
406 else {
407 G4cout << "NoProcessAssigned";
408 }
409 G4cout << G4endl;
410 G4cout << " " << fGhostStep->GetPreStepPoint()->GetPosition() << G4endl;
411 G4cout << " PostStepPoint : ";
412 if (fGhostStep->GetPostStepPoint()->GetPhysicalVolume() != nullptr) {
413 G4cout << fGhostStep->GetPostStepPoint()->GetPhysicalVolume()->GetName() << " ["
414 << fGhostStep->GetPostStepPoint()->GetTouchable()->GetReplicaNumber() << " ]";
415 }
416 else {
417 G4cout << "OutOfWorld";
418 }
419 G4cout << " - ";
420 if (fGhostStep->GetPostStepPoint()->GetProcessDefinedStep() != nullptr) {
422 }
423 else {
424 G4cout << "NoProcessAssigned";
425 }
426 G4cout << G4endl;
427 G4cout << " " << fGhostStep->GetPostStepPoint()->GetPosition()
428 << " == " << fGhostStep->GetTrack()->GetMomentumDirection() << G4endl;
429}
const G4VTouchable * GetTouchable() const
const G4VProcess * GetProcessDefinedStep() const
const G4ThreeVector & GetPosition() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4Track * GetTrack() const
G4double GetStepLength() const
G4double GetTotalEnergyDeposit() const
virtual G4int GetReplicaNumber(G4int depth=0) const

Referenced by AtRestDoIt(), and PostStepDoIt().


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