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

#include <G4WeightCutOffProcess.hh>

+ Inheritance diagram for G4WeightCutOffProcess:

Public Member Functions

 G4WeightCutOffProcess (G4double wsurvival, G4double wlimit, G4double isource, G4VIStore *istore, const G4String &aName="WeightCutOffProcess", G4bool para=false)
 
virtual ~G4WeightCutOffProcess ()
 
 G4WeightCutOffProcess (const G4WeightCutOffProcess &)=delete
 
G4WeightCutOffProcessoperator= (const G4WeightCutOffProcess &)=delete
 
void SetParallelWorld (const G4String &parallelWorldName)
 
void SetParallelWorld (G4VPhysicalVolume *parallelWorld)
 
void StartTracking (G4Track *)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 StartTracking (G4Track *)
 
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 &)
 

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
 

Detailed Description

Definition at line 54 of file G4WeightCutOffProcess.hh.

Constructor & Destructor Documentation

◆ G4WeightCutOffProcess() [1/2]

G4WeightCutOffProcess::G4WeightCutOffProcess ( G4double  wsurvival,
G4double  wlimit,
G4double  isource,
G4VIStore istore,
const G4String aName = "WeightCutOffProcess",
G4bool  para = false 
)

Definition at line 47 of file G4WeightCutOffProcess.cc.

53 : G4VProcess(aName),
54 fParticleChange(new G4ParticleChange),
55 fWeightSurvival(wsurvival),
56 fWeightLimit(wlimit),
57 fSourceImportance(isource),
58 fIStore(istore),
59 fParaflag(para)
60{
61 if (!fParticleChange)
62 {
63 G4Exception("G4WeightCutOffProcess::G4WeightCutOffProcess()",
64 "FatalError", FatalException,
65 "Failed to allocate G4ParticleChange !");
66 }
67
68 G4VProcess::pParticleChange = fParticleChange;
69
70 fGhostStep = new G4Step();
71 fGhostPreStepPoint = fGhostStep->GetPreStepPoint();
72 fGhostPostStepPoint = fGhostStep->GetPostStepPoint();
73
75 fPathFinder = G4PathFinder::GetInstance();
76
77 if (verboseLevel>0)
78 {
79 G4cout << GetProcessName() << " is created " << G4endl;
80 }
81}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4PathFinder * GetInstance()
Definition: G4PathFinder.cc:52
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
static G4TransportationManager * GetTransportationManager()
G4int verboseLevel
Definition: G4VProcess.hh:360
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

◆ ~G4WeightCutOffProcess()

G4WeightCutOffProcess::~G4WeightCutOffProcess ( )
virtual

Definition at line 83 of file G4WeightCutOffProcess.cc.

84{
85 delete fParticleChange;
86 // delete fGhostStep;
87}

◆ G4WeightCutOffProcess() [2/2]

G4WeightCutOffProcess::G4WeightCutOffProcess ( const G4WeightCutOffProcess )
delete

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4WeightCutOffProcess::AlongStepDoIt ( const G4Track track,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 373 of file G4WeightCutOffProcess.cc.

374{
375 // Dummy ParticleChange ie: does nothing
376 // Expecting G4Transportation to move the track
378 return pParticleChange;
379}
virtual void Initialize(const G4Track &)

◆ AlongStepGetPhysicalInteractionLength()

G4double G4WeightCutOffProcess::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
)
virtual

Implements G4VProcess.

Definition at line 285 of file G4WeightCutOffProcess.cc.

289{
290 if(fParaflag) {
291
292 *selection = NotCandidateForSelection;
293 G4double returnedStep = DBL_MAX;
294
295 if (previousStepSize > 0.)
296 { fGhostSafety -= previousStepSize; }
297 // else
298 // { fGhostSafety = -1.; }
299 if (fGhostSafety < 0.) fGhostSafety = 0.0;
300
301 // ------------------------------------------
302 // Determination of the proposed STEP LENGTH:
303 // ------------------------------------------
304 if (currentMinimumStep <= fGhostSafety && currentMinimumStep > 0.)
305 {
306 // I have no chance to limit
307 returnedStep = currentMinimumStep;
308 fOnBoundary = false;
309 proposedSafety = fGhostSafety - currentMinimumStep;
310 }
311 else // (currentMinimumStep > fGhostSafety: I may limit the Step)
312 {
313 G4FieldTrackUpdator::Update(&fFieldTrack,&track);
314 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
315 // ComputeStep
316 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
317 returnedStep
318 = fPathFinder->ComputeStep(fFieldTrack,currentMinimumStep,fNavigatorID,
319 track.GetCurrentStepNumber(),fGhostSafety,feLimited,
320 fEndTrack,track.GetVolume());
321 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
322 if(feLimited == kDoNot)
323 {
324 // Track is not on the boundary
325 fOnBoundary = false;
326 fGhostSafety = fGhostNavigator->ComputeSafety(fEndTrack.GetPosition());
327 }
328 else
329 {
330 // Track is on the boundary
331 fOnBoundary = true;
332 proposedSafety = fGhostSafety;
333 }
334 //xbug? proposedSafety = fGhostSafety;
335 if(feLimited == kUnique || feLimited == kSharedOther) {
336 *selection = CandidateForSelection;
337 }else if (feLimited == kSharedTransport) {
338 returnedStep *= (1.0 + 1.0e-9);
339 // Expand to disable its selection in Step Manager comparison
340 }
341
342 }
343
344 // ----------------------------------------------
345 // Returns the fGhostSafety as the proposedSafety
346 // The SteppingManager will take care of keeping
347 // the smallest one.
348 // ----------------------------------------------
349 return returnedStep;
350
351 } else {
352 return DBL_MAX;
353 //not sensible! return -1.0;
354 }
355
356}
@ CandidateForSelection
@ NotCandidateForSelection
@ kDoNot
@ kUnique
@ kSharedOther
@ kSharedTransport
double G4double
Definition: G4Types.hh:83
static void Update(G4FieldTrack *, const G4Track *)
G4ThreeVector GetPosition() const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
G4double ComputeStep(const G4FieldTrack &pFieldTrack, G4double pCurrentProposedStepLength, G4int navigatorId, G4int stepNo, G4double &pNewSafety, ELimited &limitedStep, G4FieldTrack &EndState, G4VPhysicalVolume *currentVolume)
G4VPhysicalVolume * GetVolume() const
G4int GetCurrentStepNumber() const
#define DBL_MAX
Definition: templates.hh:62

◆ AtRestDoIt()

G4VParticleChange * G4WeightCutOffProcess::AtRestDoIt ( const G4Track ,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 367 of file G4WeightCutOffProcess.cc.

368{
369 return nullptr;
370}

◆ AtRestGetPhysicalInteractionLength()

G4double G4WeightCutOffProcess::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
)
virtual

Implements G4VProcess.

Definition at line 359 of file G4WeightCutOffProcess.cc.

362{
363 return -1.0;
364}

◆ GetName()

const G4String & G4WeightCutOffProcess::GetName ( ) const

Definition at line 280 of file G4WeightCutOffProcess.cc.

281{
282 return theProcessName;
283}
G4String theProcessName
Definition: G4VProcess.hh:345

◆ operator=()

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

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 177 of file G4WeightCutOffProcess.cc.

179{
180 fParticleChange->Initialize(aTrack);
181
182 if(fParaflag) {
183 fOldGhostTouchable = fGhostPostStepPoint->GetTouchableHandle();
184 //xbug? fOnBoundary = false;
185 CopyStep(aStep);
186
187 if(fOnBoundary)
188 {
189//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
190// Locate the point and get new touchable
191//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
192 //?? fPathFinder->Locate(step.GetPostStepPoint()->GetPosition(),
193 //?? step.GetPostStepPoint()->GetMomentumDirection());
194 fNewGhostTouchable = fPathFinder->CreateTouchableHandle(fNavigatorID);
195//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
196 }
197 else
198 {
199 // Do I need this ??????????????????????????????????????????????????????????
200 // fGhostNavigator->LocateGlobalPointWithinVolume(track.GetPosition());
201 // ?????????????????????????????????????????????????????????????????????????
202
203 // fPathFinder->ReLocate(track.GetPosition());
204
205 // reuse the touchable
206 fNewGhostTouchable = fOldGhostTouchable;
207 }
208
209 fGhostPreStepPoint->SetTouchableHandle(fOldGhostTouchable);
210 fGhostPostStepPoint->SetTouchableHandle(fNewGhostTouchable);
211
212 }
213
214 if(fParaflag) {
215 G4GeometryCell postCell(*(fGhostPostStepPoint->GetPhysicalVolume()),
216 fGhostPostStepPoint->GetTouchable()->GetReplicaNumber());
217
218
219 // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(aStep);
220 // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(fGhostStep);
221 G4double R = fSourceImportance;
222 if (fIStore != nullptr)
223 {
224 G4double i = fIStore->GetImportance(postCell);
225 if (i>0)
226 {
227 R/=i;
228 }
229 }
230 G4double w = aTrack.GetWeight();
231 if (w<R*fWeightLimit)
232 {
233 G4double ws = fWeightSurvival*R;
234 G4double p = w/(ws);
235 if (G4UniformRand()<p)
236 {
237 fParticleChange->ProposeTrackStatus(fStopAndKill);
238 }
239 else
240 {
241 fParticleChange->ProposeWeight(ws);
242 }
243 }
244 } else {
245
246 G4GeometryCell postCell(*(aStep.GetPostStepPoint()->GetPhysicalVolume()),
248
249 // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(aStep);
250 // G4GeometryCell postCell = fGCellFinder.GetPostGeometryCell(fGhostStep);
251 G4double R = fSourceImportance;
252 if (fIStore != nullptr)
253 {
254 G4double i = fIStore->GetImportance(postCell);
255 if (i>0)
256 {
257 R/=i;
258 }
259 }
260 G4double w = aTrack.GetWeight();
261 if (w<R*fWeightLimit)
262 {
263 G4double ws = fWeightSurvival*R;
264 G4double p = w/(ws);
265 if (G4UniformRand()<p)
266 {
267 fParticleChange->ProposeTrackStatus(fStopAndKill);
268 }
269 else
270 {
271 fParticleChange->ProposeWeight(ws);
272 }
273 }
274 }
275
276 return fParticleChange;
277
278}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
void Initialize(const G4Track &) override
G4TouchableHandle CreateTouchableHandle(G4int navId) const
const G4VTouchable * GetTouchable() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
const G4TouchableHandle & GetTouchableHandle() const
G4VPhysicalVolume * GetPhysicalVolume() const
G4double GetWeight() const
virtual G4double GetImportance(const G4GeometryCell &gCell) const =0
void ProposeTrackStatus(G4TrackStatus status)
void ProposeWeight(G4double finalWeight)
virtual G4int GetReplicaNumber(G4int depth=0) const
Definition: G4VTouchable.cc:50

◆ PostStepGetPhysicalInteractionLength()

G4double G4WeightCutOffProcess::PostStepGetPhysicalInteractionLength ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 164 of file G4WeightCutOffProcess.cc.

167{
168// *condition = Forced;
169// return kInfinity;
170
171// *condition = StronglyForced;
172 *condition = Forced;
173 return DBL_MAX;
174}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced

◆ SetParallelWorld() [1/2]

void G4WeightCutOffProcess::SetParallelWorld ( const G4String parallelWorldName)

Definition at line 95 of file G4WeightCutOffProcess.cc.

97{
98//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
99// Get pointers of the parallel world and its navigator
100//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
101 fGhostWorldName = parallelWorldName;
102 fGhostWorld = fTransportationManager->GetParallelWorld(fGhostWorldName);
103 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
104//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
105}
G4VPhysicalVolume * GetParallelWorld(const G4String &worldName)
G4Navigator * GetNavigator(const G4String &worldName)

Referenced by G4WeightCutOffConfigurator::Configure().

◆ SetParallelWorld() [2/2]

void G4WeightCutOffProcess::SetParallelWorld ( G4VPhysicalVolume parallelWorld)

Definition at line 107 of file G4WeightCutOffProcess.cc.

109{
110//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
111// Get pointer of navigator
112//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
113 fGhostWorldName = parallelWorld->GetName();
114 fGhostWorld = parallelWorld;
115 fGhostNavigator = fTransportationManager->GetNavigator(fGhostWorld);
116//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
117}
const G4String & GetName() const

◆ StartTracking()

void G4WeightCutOffProcess::StartTracking ( G4Track trk)
virtual

Reimplemented from G4VProcess.

Definition at line 124 of file G4WeightCutOffProcess.cc.

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

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