Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ImportanceProcess.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4ImportanceProcess
27//
28// Author: Michael Dressel, 2002
29// --------------------------------------------------------------------
30
33#include "G4GeometryCell.hh"
35#include "G4VTrackTerminator.hh"
36#include "G4VIStore.hh"
37
38#include "G4Step.hh"
39#include "G4Navigator.hh"
40#include "G4VTouchable.hh"
41#include "G4VPhysicalVolume.hh"
42#include "G4ParticleChange.hh"
43#include "G4PathFinder.hh"
45#include "G4StepPoint.hh"
47
48
50G4ImportanceProcess(const G4VImportanceAlgorithm &aImportanceAlgorithm,
51 const G4VIStore &aIstore,
52 const G4VTrackTerminator *TrackTerminator,
53 const G4String &aName, G4bool para)
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}
93
95{
96 delete fPostStepAction;
97 delete fParticleChange;
98}
99
100//------------------------------------------------------
101//
102// SetParallelWorld
103//
104//------------------------------------------------------
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}
117
118//------------------------------------------------------
119//
120// StartTracking
121//
122//------------------------------------------------------
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}
160
168
171 const G4Step &aStep)
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}
254
256{
257 fParticleChange->ProposeTrackStatus(fStopAndKill);
258}
259
261{
262 return theProcessName;
263}
264
267 const G4Track& track, G4double previousStepSize,
268 G4double currentMinimumStep,
269 G4double& proposedSafety, G4GPILSelection* selection)
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}
338
345
347AtRestDoIt(const G4Track&, const G4Step&)
348{
349 return nullptr;
350}
351
353AlongStepDoIt(const G4Track& aTrack, const G4Step& )
354{
355 // Dummy ParticleChange ie: does nothing
356 // Expecting G4Transportation to move the track
357 pParticleChange->Initialize(aTrack);
358 return pParticleChange;
359}
360
361void G4ImportanceProcess::CopyStep(const G4Step& step)
362{
363 fGhostStep->SetTrack(step.GetTrack());
364 fGhostStep->SetStepLength(step.GetStepLength());
366 fGhostStep->SetControlFlag(step.GetControlFlag());
367
368 *fGhostPreStepPoint = *(step.GetPreStepPoint());
369 *fGhostPostStepPoint = *(step.GetPostStepPoint());
370
371//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
372// Set StepStatus for ghost world
373//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
374 if(fOnBoundary)
375 { fGhostPostStepPoint->SetStepStatus(fGeomBoundary); }
376 else if(fGhostPostStepPoint->GetStepStatus()==fGeomBoundary)
377 { fGhostPostStepPoint->SetStepStatus(fPostStepDoItProc); }
378//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
379}
G4double condition(const G4ErrorSymMatrix &m)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4ForceCondition
@ Forced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
@ kDoNot
@ kUnique
@ kSharedOther
@ kSharedTransport
@ fParallel
@ fGeomBoundary
@ fPostStepDoItProc
@ fStopAndKill
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static void Update(G4FieldTrack *, const G4Track *)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
void SetParallelWorld(const G4String &parallelWorldName)
virtual void KillTrack() const
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
G4ImportanceProcess(const G4VImportanceAlgorithm &aImportanceAlgorithm, const G4VIStore &aIstore, const G4VTrackTerminator *TrackTerminator, const G4String &aName="ImportanceProcess", G4bool para=false)
void StartTracking(G4Track *)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
virtual const G4String & GetName() const
static G4PathFinder * GetInstance()
G4StepStatus GetStepStatus() const
const G4VTouchable * GetTouchable() const
void SetStepStatus(const G4StepStatus aValue)
G4VPhysicalVolume * GetPhysicalVolume() const
G4SteppingControl GetControlFlag() const
G4Track * GetTrack() const
void SetStepLength(G4double value)
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4double GetTotalEnergyDeposit() const
void SetControlFlag(G4SteppingControl StepControlFlag)
void SetTotalEnergyDeposit(G4double value)
G4StepPoint * GetPostStepPoint() const
void SetTrack(G4Track *value)
virtual G4int GetReplicaNumber(G4int depth=0) const
G4TrackStatus GetTrackStatus() const
G4VPhysicalVolume * GetVolume() const
G4VPhysicalVolume * GetNextVolume() const
G4double GetWeight() const
const G4ThreeVector & GetPosition() const
G4int GetCurrentStepNumber() const
const G4ThreeVector & GetMomentumDirection() const
static G4TransportationManager * GetTransportationManager()
virtual void Initialize(const G4Track &)
G4VProcess(const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
Definition G4VProcess.cc:46
G4int verboseLevel
G4String theProcessName
G4VParticleChange * pParticleChange
const G4String & GetProcessName() const
#define DBL_MAX
Definition templates.hh:62