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