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

#include <G4ImportanceProcess.hh>

+ Inheritance diagram for G4ImportanceProcess:

Public Member Functions

 G4ImportanceProcess (const G4VImportanceAlgorithm &aImportanceAlgorithm, const G4VIStore &aIstore, const G4VTrackTerminator *TrackTerminator, const G4String &aName="ImportanceProcess", G4bool para=false)
 
virtual ~G4ImportanceProcess ()
 
 G4ImportanceProcess (const G4ImportanceProcess &)=delete
 
G4ImportanceProcessoperator= (const G4ImportanceProcess &)=delete
 
void SetParallelWorld (const G4String &parallelWorldName)
 
void StartTracking (G4Track *)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual void KillTrack () const
 
virtual const G4StringGetName () const
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- 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 &)
 
- Public Member Functions inherited from G4VTrackTerminator
 G4VTrackTerminator ()
 
virtual ~G4VTrackTerminator ()
 

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
 
- Protected Attributes inherited from G4VTrackTerminator
G4double kCarTolerance
 

Detailed Description

Definition at line 57 of file G4ImportanceProcess.hh.

Constructor & Destructor Documentation

◆ G4ImportanceProcess() [1/2]

G4ImportanceProcess::G4ImportanceProcess ( const G4VImportanceAlgorithm & aImportanceAlgorithm,
const G4VIStore & aIstore,
const G4VTrackTerminator * TrackTerminator,
const G4String & aName = "ImportanceProcess",
G4bool para = false )

Definition at line 49 of file G4ImportanceProcess.cc.

54 : G4VProcess(aName, fParallel),
55 fParticleChange(new G4ParticleChange),
56 fImportanceAlgorithm(aImportanceAlgorithm),
57 fIStore(aIstore),
58 fParaflag(para)
59{
60 G4cout << "### G4ImportanceProcess:: Creating " << G4endl;
61 if (TrackTerminator != nullptr)
62 {
63 fPostStepAction = new G4SamplingPostStepAction(*TrackTerminator);
64 }
65 else
66 {
67 fPostStepAction = new G4SamplingPostStepAction(*this);
68 }
69 if (fParticleChange == nullptr)
70 {
71 G4Exception("G4ImportanceProcess::G4ImportanceProcess()",
72 "FatalError", FatalException,
73 "Failed allocation of G4ParticleChange !");
74 }
75 G4VProcess::pParticleChange = fParticleChange;
76
77
78 fGhostStep = new G4Step();
79 fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
80 fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
81
83 fPathFinder = G4PathFinder::GetInstance();
84
85 if (verboseLevel>0)
86 {
87 G4cout << GetProcessName() << " is created " << G4endl;
88 }
89
90 G4cout << "G4ImportanceProcess:: importance process paraflag is: " << fParaflag << G4endl;
91
92}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
@ fParallel
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4PathFinder * GetInstance()
static G4TransportationManager * GetTransportationManager()
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition G4VProcess.cc:46
G4int verboseLevel
G4VParticleChange * pParticleChange
const G4String & GetProcessName() const

Referenced by G4ImportanceProcess(), and operator=().

◆ ~G4ImportanceProcess()

G4ImportanceProcess::~G4ImportanceProcess ( )
virtual

Definition at line 94 of file G4ImportanceProcess.cc.

95{
96 delete fPostStepAction;
97 delete fParticleChange;
98}

◆ G4ImportanceProcess() [2/2]

G4ImportanceProcess::G4ImportanceProcess ( const G4ImportanceProcess & )
delete

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4ImportanceProcess::AlongStepDoIt ( const G4Track & aTrack,
const G4Step &  )
virtual

Implements G4VProcess.

Definition at line 352 of file G4ImportanceProcess.cc.

354{
355 // Dummy ParticleChange ie: does nothing
356 // Expecting G4Transportation to move the track
357 pParticleChange->Initialize(aTrack);
358 return pParticleChange;
359}

◆ AlongStepGetPhysicalInteractionLength()

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

Implements G4VProcess.

Definition at line 265 of file G4ImportanceProcess.cc.

270{
271 if(fParaflag)
272 {
273 *selection = NotCandidateForSelection;
274 G4double returnedStep = DBL_MAX;
275
276 if (previousStepSize > 0.)
277 { fGhostSafety -= previousStepSize; }
278 // else
279 // { fGhostSafety = -1.; }
280 if (fGhostSafety < 0.) { fGhostSafety = 0.0; }
281
282 // ------------------------------------------
283 // Determination of the proposed STEP LENGTH:
284 // ------------------------------------------
285 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
286 {
287 // I have no chance to limit
288 returnedStep = currentMinimumStep;
289 fOnBoundary = false;
290 proposedSafety = fGhostSafety - currentMinimumStep;
291 }
292 else // (currentMinimumStep > fGhostSafety: I may limit the Step)
293 {
294 G4FieldTrackUpdator::Update(&fFieldTrack,&track);
295 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
296 // ComputeStep
297 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
298 returnedStep
299 = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
300 track.GetCurrentStepNumber(),fGhostSafety,feLimited,
301 fEndTrack,track.GetVolume());
302 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
303 if(feLimited == kDoNot)
304 {
305 // Track is not on the boundary
306 fOnBoundary = false;
307 fGhostSafety = fGhostNavigator->ComputeSafety(fEndTrack.GetPosition());
308 }
309 else
310 {
311 // Track is on the boundary
312 fOnBoundary = true;
313 }
314 proposedSafety = fGhostSafety;
315 if(feLimited == kUnique || feLimited == kSharedOther)
316 {
317 *selection = CandidateForSelection;
318 }
319 else if (feLimited == kSharedTransport)
320 {
321 returnedStep *= (1.0 + 1.0e-9);
322 // Expand to disable its selection in Step Manager comparison
323 }
324 }
325
326 // ----------------------------------------------
327 // Returns the fGhostSafety as the proposedSafety
328 // The SteppingManager will take care of keeping
329 // the smallest one.
330 // ----------------------------------------------
331 return returnedStep;
332 }
333 else
334 {
335 return DBL_MAX;
336 }
337}
@ CandidateForSelection
@ NotCandidateForSelection
@ kDoNot
@ kUnique
@ kSharedOther
@ kSharedTransport
double G4double
Definition G4Types.hh:83
static void Update(G4FieldTrack *, const G4Track *)
G4VPhysicalVolume * GetVolume() const
G4int GetCurrentStepNumber() const
#define DBL_MAX
Definition templates.hh:62

◆ AtRestDoIt()

G4VParticleChange * G4ImportanceProcess::AtRestDoIt ( const G4Track & ,
const G4Step &  )
virtual

Implements G4VProcess.

Definition at line 346 of file G4ImportanceProcess.cc.

348{
349 return nullptr;
350}

◆ AtRestGetPhysicalInteractionLength()

G4double G4ImportanceProcess::AtRestGetPhysicalInteractionLength ( const G4Track & ,
G4ForceCondition *  )
virtual

Implements G4VProcess.

Definition at line 339 of file G4ImportanceProcess.cc.

342{
343 return -1.0;
344}

◆ GetName()

const G4String & G4ImportanceProcess::GetName ( ) const
virtual

Implements G4VTrackTerminator.

Definition at line 260 of file G4ImportanceProcess.cc.

261{
262 return theProcessName;
263}
G4String theProcessName

◆ KillTrack()

void G4ImportanceProcess::KillTrack ( ) const
virtual

Implements G4VTrackTerminator.

Definition at line 255 of file G4ImportanceProcess.cc.

256{
257 fParticleChange->ProposeTrackStatus(fStopAndKill);
258}
@ fStopAndKill

◆ operator=()

G4ImportanceProcess & G4ImportanceProcess::operator= ( const G4ImportanceProcess & )
delete

◆ PostStepDoIt()

G4VParticleChange * G4ImportanceProcess::PostStepDoIt ( const G4Track & aTrack,
const G4Step & aStep )
virtual

Implements G4VProcess.

Definition at line 170 of file G4ImportanceProcess.cc.

172{
173 fParticleChange->Initialize(aTrack);
174
175 if(aTrack.GetNextVolume() == nullptr)
176 {
177 return fParticleChange;
178 }
179
180 if(fParaflag)
181 {
182 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
183 //xbug? fOnBoundary = false;
184 CopyStep(aStep);
185
186 if(fOnBoundary)
187 {
188//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
189// Locate the point and get new touchable
190//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
191 //?? fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
192 //?? step.GetPostStepPoint()->GetMomentumDirection());
193 fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
194 }
195 else
196 {
197 // reuse the touchable
198 fNewGhostTouchable = fOldGhostTouchable;
199 }
200
201 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
202 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
203
204 if ( (fGhostPostStepPoint->GetStepStatus() == fGeomBoundary)
205 && (aStep.GetStepLength() > kCarTolerance) )
206 {
207 if (aTrack.GetTrackStatus()==fStopAndKill)
208 {
209 G4cout << "WARNING - G4ImportanceProcess::PostStepDoIt()"
210 << " StopAndKill track. on boundary" << G4endl;
211 }
212
213 G4GeometryCell prekey(*(fGhostPreStepPoint->GetPhysicalVolume()),
214 fGhostPreStepPoint->GetTouchable()->GetReplicaNumber());
215 G4GeometryCell postkey(*(fGhostPostStepPoint->GetPhysicalVolume()),
216 fGhostPostStepPoint->GetTouchable()->GetReplicaNumber());
217
218 G4Nsplit_Weight nw =
219 fImportanceAlgorithm.Calculate(fIStore.GetImportance(prekey),
220 fIStore.GetImportance(postkey),
221 aTrack.GetWeight());
222 fPostStepAction->DoIt(aTrack, fParticleChange, nw);
223 }
224 }
225 else
226 {
227 if ( (aStep.GetPostStepPoint()->GetStepStatus() == fGeomBoundary)
228 && (aStep.GetStepLength() > kCarTolerance) )
229 {
230 if (aTrack.GetTrackStatus()==fStopAndKill)
231 {
232 G4cout << "WARNING - G4ImportanceProcess::PostStepDoIt()"
233 << " StopAndKill track. on boundary non-parallel"
234 << G4endl;
235 }
236
237 G4StepPoint *prepoint = aStep.GetPreStepPoint();
238 G4StepPoint *postpoint = aStep.GetPostStepPoint();
239
240 G4GeometryCell prekey(*(prepoint->GetPhysicalVolume()),
241 prepoint->GetTouchable()->GetReplicaNumber());
242 G4GeometryCell postkey(*(postpoint->GetPhysicalVolume()),
243 postpoint->GetTouchable()->GetReplicaNumber());
244
245 G4Nsplit_Weight nw =
246 fImportanceAlgorithm.Calculate(fIStore.GetImportance(prekey),
247 fIStore.GetImportance(postkey),
248 aTrack.GetWeight());
249 fPostStepAction->DoIt(aTrack, fParticleChange, nw);
250 }
251 }
252 return fParticleChange;
253}
@ fGeomBoundary
G4StepStatus GetStepStatus() const
const G4VTouchable * GetTouchable() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
virtual G4int GetReplicaNumber(G4int depth=0) const
G4TrackStatus GetTrackStatus() const
G4VPhysicalVolume * GetNextVolume() const
G4double GetWeight() const

◆ PostStepGetPhysicalInteractionLength()

G4double G4ImportanceProcess::PostStepGetPhysicalInteractionLength ( const G4Track & aTrack,
G4double previousStepSize,
G4ForceCondition * condition )
virtual

Implements G4VProcess.

Definition at line 161 of file G4ImportanceProcess.cc.

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

◆ SetParallelWorld()

void G4ImportanceProcess::SetParallelWorld ( const G4String & parallelWorldName)

Definition at line 105 of file G4ImportanceProcess.cc.

106{
107 G4cout << G4endl << G4endl << G4endl;
108 G4cout << "G4ImportanceProcess:: SetParallelWorld name = " << parallelWorldName << G4endl;
109//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
110// Get pointers of the parallel world and its navigator
111//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
112 fGhostWorldName = parallelWorldName;
113 fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
114 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
115//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
116}

◆ StartTracking()

void G4ImportanceProcess::StartTracking ( G4Track * trk)
virtual

Reimplemented from G4VProcess.

Definition at line 123 of file G4ImportanceProcess.cc.

124{
125//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
126// Activate navigator and get the navigator ID
127//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
128
129 if(fParaflag)
130 {
131 if(fGhostNavigator != nullptr)
132 {
133 fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator);
134 }
135 else
136 {
137 G4Exception("G4ImportanceProcess::StartTracking",
138 "ProcParaWorld000",FatalException,
139 "G4ImportanceProcess is used for tracking without having a parallel world assigned");
140 }
141
142//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
143// Let PathFinder initialize
144//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
145 fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
146
147//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
148// Setup initial touchables for the first step
149//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
150 fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
151 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
152 fNewGhostTouchable = fOldGhostTouchable;
153 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
154
155 // Initialize
156 fGhostSafety = -1.;
157 fOnBoundary = false;
158 }
159}
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const

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