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
G4WeightCutOffProcess.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// G4WeightCutOffProcess.cc
33//
34// ----------------------------------------------------------------------
35
37//#include "G4VScorer.hh"
38#include "G4GeometryCellStep.hh"
39//#include "G4GCellFinder.hh"
40#include "G4TouchableHandle.hh"
41#include "G4VIStore.hh"
42
43#include "G4Step.hh"
44#include "G4Navigator.hh"
45#include "G4VTouchable.hh"
46#include "G4VPhysicalVolume.hh"
47#include "G4ParticleChange.hh"
48#include "G4PathFinder.hh"
50#include "G4StepPoint.hh"
52
53
56 G4double wlimit,
57 G4double isource,
58 G4VIStore *istore,
59 // const G4VGCellFinder &aGCellFinder,
60 const G4String &aName, G4bool para)
61 : G4VProcess(aName),
62 fParticleChange(new G4ParticleChange),
63 fWeightSurvival(wsurvival),
64 fWeightLimit(wlimit),
65 fSourceImportance(isource),
66 fIStore(istore),
67 // fGCellFinder(aGCellFinder),
68 fGhostNavigator(0), fNavigatorID(-1), fFieldTrack('0')
69{
70 if (!fParticleChange)
71 {
72 G4Exception("G4WeightCutOffProcess::G4WeightCutOffProcess()",
73 "FatalError", FatalException,
74 "Failed to allocate G4ParticleChange !");
75 }
76
77 G4VProcess::pParticleChange = fParticleChange;
78
79 fGhostStep = new G4Step();
80 fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
81 fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
82
84 fPathFinder = G4PathFinder::GetInstance();
85
86 if (verboseLevel>0)
87 {
88 G4cout << GetProcessName() << " is created " << G4endl;
89 }
90
91 paraflag = para;
92
93}
94
96{
97 delete fParticleChange;
98 delete fGhostStep;
99}
100
101
102//------------------------------------------------------
103//
104// SetParallelWorld
105//
106//------------------------------------------------------
108SetParallelWorld(G4String parallelWorldName)
109{
110//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
111// Get pointers of the parallel world and its navigator
112//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
113 fGhostWorldName = parallelWorldName;
114 fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
115 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
116//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
117}
118
121{
122//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
123// Get pointer of navigator
124//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
125 fGhostWorldName = parallelWorld->GetName();
126 fGhostWorld = parallelWorld;
127 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
128//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
129}
130
131//------------------------------------------------------
132//
133// StartTracking
134//
135//------------------------------------------------------
137{
138//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
139// Activate navigator and get the navigator ID
140//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
141// G4cout << " G4ParallelWorldScoringProcess::StartTracking" << G4endl;
142
143 if(paraflag) {
144 if(fGhostNavigator)
145 { fNavigatorID = fTransportationManager->ActivateNavigator(fGhostNavigator); }
146 else
147 {
148 G4Exception("G4WeightCutOffProcess::StartTracking",
149 "ProcParaWorld000",FatalException,
150 "G4WeightCutOffProcess is used for tracking without having a parallel world assigned");
151 }
152//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
153
154// G4cout << "G4ParallelWorldScoringProcess::StartTracking <<<<<<<<<<<<<<<<<< " << G4endl;
155//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
156// Let PathFinder initialize
157//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
158 fPathFinder->PrepareNewTrack(trk->GetPosition(),trk->GetMomentumDirection());
159//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
160
161//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
162// Setup initial touchables for the first step
163//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
164 fOldGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
165 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
166 fNewGhostTouchable = fOldGhostTouchable;
167 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
168
169 // Initialize
170 fGhostSafety = -1.;
171 fOnBoundary = false;
172 }
173}
174
175
179{
180// *condition = Forced;
181// return kInfinity;
182
183// *condition = StronglyForced;
184 *condition = Forced;
185 return DBL_MAX;
186}
187
190 const G4Step &aStep)
191{
192 fParticleChange->Initialize(aTrack);
193
194 if(paraflag) {
195 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
196 //xbug? fOnBoundary = false;
197 CopyStep(aStep);
198
199 if(fOnBoundary)
200 {
201//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
202// Locate the point and get new touchable
203//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
204 //?? fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
205 //?? step.GetPostStepPoint()->GetMomentumDirection());
206 fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
207//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
208 }
209 else
210 {
211 // Do I need this ??????????????????????????????????????????????????????????
212 // fGhostNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
213 // ?????????????????????????????????????????????????????????????????????????
214
215 // fPathFinder->ReLocate(track.GetPosition());
216
217 // reuse the touchable
218 fNewGhostTouchable = fOldGhostTouchable;
219 }
220
221 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
222 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
223
224 }
225
226 if(paraflag) {
227 G4GeometryCell postCell(*(fGhostPostStepPoint->GetPhysicalVolume()),
228 fGhostPostStepPoint->GetTouchable()->GetReplicaNumber());
229
230
231 // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(aStep);
232 // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(fGhostStep);
233 G4double R = fSourceImportance;
234 if (fIStore)
235 {
236 G4double i = fIStore->GetImportance(postCell);
237 if (i>0)
238 {
239 R/=i;
240 }
241 }
242 G4double w = aTrack.GetWeight();
243 if (w<R*fWeightLimit)
244 {
245 G4double ws = fWeightSurvival*R;
246 G4double p = w/(ws);
247 if (G4UniformRand()<p)
248 {
249 fParticleChange->ProposeTrackStatus(fStopAndKill);
250 }
251 else
252 {
253 fParticleChange->ProposeWeight(ws);
254 }
255 }
256 } else {
257
258 G4GeometryCell postCell(*(aStep.GetPostStepPoint()->GetPhysicalVolume()),
260
261 // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(aStep);
262 // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(fGhostStep);
263 G4double R = fSourceImportance;
264 if (fIStore)
265 {
266 G4double i = fIStore->GetImportance(postCell);
267 if (i>0)
268 {
269 R/=i;
270 }
271 }
272 G4double w = aTrack.GetWeight();
273 if (w<R*fWeightLimit)
274 {
275 G4double ws = fWeightSurvival*R;
276 G4double p = w/(ws);
277 if (G4UniformRand()<p)
278 {
279 fParticleChange->ProposeTrackStatus(fStopAndKill);
280 }
281 else
282 {
283 fParticleChange->ProposeWeight(ws);
284 }
285 }
286 }
287
288 return fParticleChange;
289
290}
291
293{
294 return theProcessName;
295}
296
299 const G4Track& track, G4double previousStepSize, G4double currentMinimumStep,
300 G4double& proposedSafety, G4GPILSelection* selection)
301{
302 if(paraflag) {
303 static G4FieldTrack endTrack('0');
304 static ELimited eLimited;
305
306 *selection = NotCandidateForSelection;
307 G4double returnedStep = DBL_MAX;
308
309 if (previousStepSize > 0.)
310 { fGhostSafety -= previousStepSize; }
311 // else
312 // { fGhostSafety = -1.; }
313 if (fGhostSafety < 0.) fGhostSafety = 0.0;
314
315 // ------------------------------------------
316 // Determination of the proposed STEP LENGTH:
317 // ------------------------------------------
318 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
319 {
320 // I have no chance to limit
321 returnedStep = currentMinimumStep;
322 fOnBoundary = false;
323 proposedSafety = fGhostSafety - currentMinimumStep;
324 }
325 else // (currentMinimumStep > fGhostSafety: I may limit the Step)
326 {
327 G4FieldTrackUpdator::Update(&fFieldTrack,&track);
328 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
329 // ComputeStep
330 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
331 returnedStep
332 = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
333 track.GetCurrentStepNumber(),fGhostSafety,eLimited,
334 endTrack,track.GetVolume());
335 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
336 if(eLimited == kDoNot)
337 {
338 // Track is not on the boundary
339 fOnBoundary = false;
340 fGhostSafety = fGhostNavigator->ComputeSafety(endTrack.GetPosition());
341 }
342 else
343 {
344 // Track is on the boundary
345 fOnBoundary = true;
346 proposedSafety = fGhostSafety;
347 }
348 //xbug? proposedSafety = fGhostSafety;
349 if(eLimited == kUnique || eLimited == kSharedOther) {
350 *selection = CandidateForSelection;
351 }else if (eLimited == kSharedTransport) {
352 returnedStep *= (1.0 + 1.0e-9);
353 // Expand to disable its selection in Step Manager comparison
354 }
355
356 }
357
358 // ----------------------------------------------
359 // Returns the fGhostSafety as the proposedSafety
360 // The SteppingManager will take care of keeping
361 // the smallest one.
362 // ----------------------------------------------
363 return returnedStep;
364
365 } else {
366 return DBL_MAX;
367 //not sensible! return -1.0;
368 }
369
370}
371
372
376{
377 return -1.0;
378}
379
382{
383 return 0;
384}
385
388{
389 // Dummy ParticleChange ie: does nothing
390 // Expecting G4Transportation to move the track
392 return pParticleChange;
393
394 // return 0;
395}
396
397void G4WeightCutOffProcess::CopyStep(const G4Step & step)
398{
399 fGhostStep->SetTrack(step.GetTrack());
400 fGhostStep->SetStepLength(step.GetStepLength());
402 fGhostStep->SetControlFlag(step.GetControlFlag());
403
404 *fGhostPreStepPoint = *(step.GetPreStepPoint());
405 *fGhostPostStepPoint = *(step.GetPostStepPoint());
406
407//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
408// Set StepStatus for ghost world
409//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
410 if(fOnBoundary)
411 { fGhostPostStepPoint->SetStepStatus(fGeomBoundary); }
412 else if(fGhostPostStepPoint->GetStepStatus()==fGeomBoundary)
413 { fGhostPostStepPoint->SetStepStatus(fPostStepDoItProc); }
414//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
415}
G4double condition(const G4ErrorSymMatrix &m)
@ FatalException
G4ForceCondition
@ Forced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
ELimited
@ kDoNot
@ kUnique
@ kSharedOther
@ kSharedTransport
@ 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
#define G4UniformRand()
Definition: Randomize.hh:53
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
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
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
static G4TransportationManager * GetTransportationManager()
G4int ActivateNavigator(G4Navigator *aNavigator)
G4Navigator * GetNavigator(const G4String &worldName)
virtual G4double GetImportance(const G4GeometryCell &gCell) const =0
void ProposeTrackStatus(G4TrackStatus status)
void ProposeWeight(G4double finalWeight)
virtual void Initialize(const G4Track &)
const G4String & GetName() const
G4int verboseLevel
Definition: G4VProcess.hh:368
G4String theProcessName
Definition: G4VProcess.hh:335
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 G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
const G4String & GetName() const
G4WeightCutOffProcess(G4double wsurvival, G4double wlimit, G4double isource, G4VIStore *istore, const G4String &aName="WeightCutOffProcess", G4bool para=false)
virtual G4VParticleChange * AtRestDoIt(const G4Track &, const G4Step &)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void SetParallelWorld(G4String parallelWorldName)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define DBL_MAX
Definition: templates.hh:83