Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4WeightWindowProcess.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// G4WeightWindowProcess
27//
28// Author: Michael Dressel, 2002
29// --------------------------------------------------------------------
30
33#include "G4GeometryCell.hh"
35#include "G4VTrackTerminator.hh"
36#include "G4PlaceOfAction.hh"
38
39#include "G4Step.hh"
40#include "G4Navigator.hh"
41#include "G4VTouchable.hh"
42#include "G4VPhysicalVolume.hh"
43#include "G4ParticleChange.hh"
44#include "G4PathFinder.hh"
46#include "G4StepPoint.hh"
48
49
51 const G4VWeightWindowAlgorithm &aWeightWindowAlgorithm,
52 const G4VWeightWindowStore &aWWStore,
53 const G4VTrackTerminator *TrackTerminator,
54 G4PlaceOfAction placeOfAction,
55 const G4String &aName, G4bool para)
56 : G4VProcess(aName),
57 fParticleChange(new G4ParticleChange),
58 fWeightWindowAlgorithm(aWeightWindowAlgorithm),
59 fWeightWindowStore(aWWStore),
60 fPlaceOfAction(placeOfAction)
61{
62 if (TrackTerminator != nullptr)
63 {
64 fPostStepAction = new G4SamplingPostStepAction(*TrackTerminator);
65 }
66 else
67 {
68 fPostStepAction = new G4SamplingPostStepAction(*this);
69 }
70 if (fParticleChange == nullptr)
71 {
72 G4Exception("G4WeightWindowProcess::G4WeightWindowProcess()",
73 "FatalError", FatalException,
74 "Failed allocation of G4ParticleChange !");
75 }
76 G4VProcess::pParticleChange = fParticleChange;
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 fParaflag = para;
91
92}
93
95{
96
97 delete fPostStepAction;
98 delete fParticleChange;
99 // delete fGhostStep;
100
101}
102
103
104//------------------------------------------------------
105//
106// SetParallelWorld
107//
108//------------------------------------------------------
110SetParallelWorld(const G4String &parallelWorldName)
111{
112//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
113// Get pointers of the parallel world and its navigator
114//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
115 fGhostWorldName = parallelWorldName;
116 fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
117 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
118//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
119}
120
123{
124//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
125// Get pointer of navigator
126//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
127 fGhostWorldName = parallelWorld->GetName();
128 fGhostWorld = parallelWorld;
129 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
130//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
131}
132
133//------------------------------------------------------
134//
135// StartTracking
136//
137//------------------------------------------------------
139{
140//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
141// Activate navigator and get the navigator ID
142//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
143// G4cout << " G4ParallelWorldScoringProcess::StartTracking" << G4endl;
144
145 if(fParaflag) {
146 if(fGhostNavigator != nullptr)
147 { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
148 else
149 {
150 G4Exception("G4WeightWindowProcess::StartTracking",
151 "ProcParaWorld000",FatalException,
152 "G4WeightWindowProcess is used for tracking without having a parallel world assigned");
153 }
154//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
155
156// G4cout << "G4ParallelWorldScoringProcess::StartTracking <<<<<<<<<<<<<<<<<< " << G4endl;
157//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
158// Let PathFinder initialize
159//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
160 fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
161//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
162
163//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
164// Setup initial touchables for the first step
165//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
166 fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
167 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
168 fNewGhostTouchable = fOldGhostTouchable;
169 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
170
171 // Initialize
172 fGhostSafety = -1.;
173 fOnBoundary = false;
174 }
175}
176
177
180 G4double ,
182{
183// *condition = Forced;
184// return kInfinity;
185
186// *condition = StronglyForced;
187 *condition = Forced;
188 return DBL_MAX;
189}
190
193 const G4Step &aStep)
194{
195
196 fParticleChange->Initialize(aTrack);
197
198 if(fParaflag) {
199 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
200 //xbug? fOnBoundary = false;
201 CopyStep(aStep);
202
203 if(fOnBoundary)
204 {
205//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
206// Locate the point and get new touchable
207//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
208 //?? fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
209 //?? step.GetPostStepPoint()->GetMomentumDirection());
210 fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
211//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
212 }
213 else
214 {
215 // Do I need this ??????????????????????????????????????????????????????????
216 // fGhostNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
217 // ?????????????????????????????????????????????????????????????????????????
218
219 // fPathFinder->ReLocate(track.GetPosition());
220
221 // reuse the touchable
222 fNewGhostTouchable = fOldGhostTouchable;
223 }
224
225 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
226 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
227
228 }
229
230 if (aStep.GetStepLength() > kCarTolerance)
231 {
232// if ( ( fPlaceOfAction == onBoundaryAndCollision)
233// || ( (fPlaceOfAction == onBoundary) &&
234// (aStep.GetPostStepPoint()->GetStepStatus() == fGeomBoundary) )
235// || ( (fPlaceOfAction == onCollision) &&
236// (aStep.GetPostStepPoint()->GetStepStatus() != fGeomBoundary) ) )
237 if(fParaflag) {
238 if ( ( fPlaceOfAction == onBoundaryAndCollision)
239 || ( (fPlaceOfAction == onBoundary) &&
240 (fGhostPostStepPoint->GetStepStatus() == fGeomBoundary) )
241 || ( (fPlaceOfAction == onCollision) &&
242 (fGhostPostStepPoint->GetStepStatus() != fGeomBoundary) ) )
243 {
244
245// G4StepPoint *postpoint = aStep.GetPostStepPoint();
246
247// G4GeometryCell postCell(*(postpoint->GetPhysicalVolume()),
248// postpoint->GetTouchable()->GetReplicaNumber());
249
250 G4GeometryCell postCell(*(fGhostPostStepPoint->GetPhysicalVolume()),
251 fGhostPostStepPoint->GetTouchable()->GetReplicaNumber());
252 G4Nsplit_Weight nw =
253 fWeightWindowAlgorithm.Calculate(aTrack.GetWeight(),
254 fWeightWindowStore.GetLowerWeight(postCell,
255 aTrack.GetKineticEnergy()));
256 fPostStepAction->DoIt(aTrack, fParticleChange, nw);
257 }
258 } else {
259 if ( ( fPlaceOfAction == onBoundaryAndCollision)
260 || ( (fPlaceOfAction == onBoundary) &&
262 || ( (fPlaceOfAction == onCollision) &&
264 {
265
266 G4StepPoint *postpoint = aStep.GetPostStepPoint();
267
268 G4GeometryCell postCell(*(postpoint->GetPhysicalVolume()),
269 postpoint->GetTouchable()->GetReplicaNumber());
270
271 G4Nsplit_Weight nw =
272 fWeightWindowAlgorithm.Calculate(aTrack.GetWeight(),
273 fWeightWindowStore.GetLowerWeight(postCell,
274 aTrack.GetKineticEnergy()));
275 fPostStepAction->DoIt(aTrack, fParticleChange, nw);
276 }
277 }
278 }
279 return fParticleChange;
280}
281
283{
284 fParticleChange->ProposeTrackStatus(fStopAndKill);
285}
286
288{
290}
291
294 const G4Track& track, G4double previousStepSize, G4double currentMinimumStep,
295 G4double& proposedSafety, G4GPILSelection* selection)
296{
297 if(fParaflag)
298 {
299 *selection = NotCandidateForSelection;
300 G4double returnedStep = DBL_MAX;
301
302 if (previousStepSize > 0.)
303 { fGhostSafety -= previousStepSize; }
304 // else
305 // { fGhostSafety = -1.; }
306 if (fGhostSafety < 0.) fGhostSafety = 0.0;
307
308 // ------------------------------------------
309 // Determination of the proposed STEP LENGTH:
310 // ------------------------------------------
311 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
312 {
313 // I have no chance to limit
314 returnedStep = currentMinimumStep;
315 fOnBoundary = false;
316 proposedSafety = fGhostSafety - currentMinimumStep;
317 }
318 else // (currentMinimumStep > fGhostSafety: I may limit the Step)
319 {
320 G4FieldTrackUpdator::Update(&fFieldTrack,&track);
321 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
322 // ComputeStep
323 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
324 returnedStep
325 = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
326 track.GetCurrentStepNumber(),fGhostSafety,feLimited,
327 fEndTrack,track.GetVolume());
328 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
329 if(feLimited == kDoNot)
330 {
331 // Track is not on the boundary
332 fOnBoundary = false;
333 fGhostSafety = fGhostNavigator->ComputeSafety(fEndTrack.GetPosition());
334 // fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition(), DBL_MAX, true);
335 }
336 else
337 {
338 // Track is on the boundary
339 fOnBoundary = true;
340 proposedSafety = fGhostSafety;
341 }
342 //xbug? proposedSafety = fGhostSafety;
343 if(feLimited == kUnique || feLimited == kSharedOther) {
344 *selection = CandidateForSelection;
345 }else if (feLimited == kSharedTransport) {
346 returnedStep *= (1.0 + 1.0e-9);
347 // Expand to disable its selection in Step Manager comparison
348 }
349
350 }
351
352 // ----------------------------------------------
353 // Returns the fGhostSafety as the proposedSafety
354 // The SteppingManager will take care of keeping
355 // the smallest one.
356 // ----------------------------------------------
357 return returnedStep;
358
359 } else {
360 return DBL_MAX;
361 // not sensible - time goes backwards! return -1.0;
362 }
363
364}
365
369{
370 return -1.0;
371}
372
374AtRestDoIt(const G4Track&, const G4Step&)
375{
376 return nullptr;
377}
378
380AlongStepDoIt(const G4Track& track, const G4Step&)
381{
382 // Dummy ParticleChange ie: does nothing
383 // Expecting G4Transportation to move the track
385 return pParticleChange;
386}
387
388void G4WeightWindowProcess::CopyStep(const G4Step & step)
389{
390 fGhostStep->SetTrack(step.GetTrack());
391 fGhostStep->SetStepLength(step.GetStepLength());
393 fGhostStep->SetControlFlag(step.GetControlFlag());
394
395 *fGhostPreStepPoint = *(step.GetPreStepPoint());
396 *fGhostPostStepPoint = *(step.GetPostStepPoint());
397
398//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
399// Set StepStatus for ghost world
400//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
401 if(fOnBoundary)
402 { fGhostPostStepPoint->SetStepStatus(fGeomBoundary); }
403 else if(fGhostPostStepPoint->GetStepStatus()==fGeomBoundary)
404 { fGhostPostStepPoint->SetStepStatus(fPostStepDoItProc); }
405//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
406}
G4double condition(const G4ErrorSymMatrix &m)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
G4ForceCondition
@ Forced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
@ kDoNot
@ kUnique
@ kSharedOther
@ kSharedTransport
G4PlaceOfAction
@ onCollision
@ onBoundaryAndCollision
@ onBoundary
@ fGeomBoundary
Definition: G4StepStatus.hh:43
@ fPostStepDoItProc
Definition: G4StepStatus.hh:49
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static void Update(G4FieldTrack *, const G4Track *)
G4ThreeVector GetPosition() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
void Initialize(const G4Track &) override
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
void PrepareNewTrack(const G4ThreeVector &position, const G4ThreeVector &direction, G4VPhysicalVolume *massStartVol=nullptr)
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:52
G4TouchableHandle CreateTouchableHandle(G4int navId) const
void DoIt(const G4Track &aTrack, G4ParticleChange *aParticleChange, const G4Nsplit_Weight &nw)
G4StepStatus GetStepStatus() const
const G4VTouchable * GetTouchable() const
void SetStepStatus(const G4StepStatus aValue)
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4TouchableHandle & GetTouchableHandle() const
G4VPhysicalVolume * GetPhysicalVolume() const
Definition: G4Step.hh:62
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)
G4VPhysicalVolume * GetVolume() const
G4double GetWeight() const
const G4ThreeVector & GetPosition() const
G4int GetCurrentStepNumber() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
static G4TransportationManager * GetTransportationManager()
G4int ActivateNavigator(G4Navigator *aNavigator)
G4Navigator * GetNavigator(const G4String &worldName)
void ProposeTrackStatus(G4TrackStatus status)
virtual void Initialize(const G4Track &)
const G4String & GetName() const
G4int verboseLevel
Definition: G4VProcess.hh:360
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:50
virtual G4Nsplit_Weight Calculate(G4double init_w, G4double lowerWeightBound) const =0
virtual G4double GetLowerWeight(const G4GeometryCell &gCell, G4double partEnergy) const =0
G4WeightWindowProcess(const G4VWeightWindowAlgorithm &aWeightWindowAlgorithm, const G4VWeightWindowStore &aWWStore, const G4VTrackTerminator *TrackTerminator, G4PlaceOfAction placeOfAction, const G4String &aName="WeightWindowProcess", G4bool para=false)
void SetParallelWorld(const G4String &parallelWorldName)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
virtual void KillTrack() const
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
virtual const G4String & GetName() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
#define DBL_MAX
Definition: templates.hh:62