Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes | List of all members
GB06BOptrSplitAndKillByImportance Class Reference

#include <Doxymodules_biasing.h>

Inheritance diagram for GB06BOptrSplitAndKillByImportance:
G4VBiasingOperator

Public Member Functions

 GB06BOptrSplitAndKillByImportance (G4String particleToBias, G4String name="SplitAndKillByImportance")
 
virtual ~GB06BOptrSplitAndKillByImportance ()
 
virtual void StartRun ()
 
void SetParallelWorld (G4VPhysicalVolume *parallelWorld)
 
std::map< G4int, G4int > & GetImportanceMap ()
 

Private Member Functions

virtual G4VBiasingOperationProposeOccurenceBiasingOperation (const G4Track *, const G4BiasingProcessInterface *) final
 
virtual G4VBiasingOperationProposeFinalStateBiasingOperation (const G4Track *, const G4BiasingProcessInterface *) final
 
virtual G4VBiasingOperationProposeNonPhysicsBiasingOperation (const G4Track *, const G4BiasingProcessInterface *) final
 

Private Attributes

GB06BOptnSplitAndKillByImportancefSplitAndKillByImportance
 
const G4ParticleDefinitionfParticleToBias
 
G4VPhysicalVolumefParallelWorld
 
G4int fParallelWorldIndex
 
const G4BiasingProcessSharedDatafBiasingSharedData
 
const G4ParallelGeometriesLimiterProcessfBiasingLimiterProcess
 
std::map< G4int, G4int > fImportanceMap
 

Detailed Description

Definition at line 139 of file Doxymodules_biasing.h.

Constructor & Destructor Documentation

◆ GB06BOptrSplitAndKillByImportance()

GB06BOptrSplitAndKillByImportance::GB06BOptrSplitAndKillByImportance ( G4String  particleToBias,
G4String  name = "SplitAndKillByImportance" 
)

Definition at line 46 of file GB06BOptrSplitAndKillByImportance.cc.

49 : G4VBiasingOperator(name),
50 fParallelWorld ( nullptr ),
52 fBiasingLimiterProcess ( nullptr )
53{
54 fParticleToBias = G4ParticleTable::GetParticleTable()->FindParticle(particleName);
55
56 if ( fParticleToBias == 0 )
57 {
58 G4ExceptionDescription ed;
59 ed << "Particle `" << particleName << "' not found !" << G4endl;
60 G4Exception("GB06BOptrSplitAndKillByImportance(...)",
61 "exGB06.01",
62 JustWarning,
63 ed);
64 }
65
67 new GB06BOptnSplitAndKillByImportance("splitterFor_"+particleName);
68
69}
const G4ParallelGeometriesLimiterProcess * fBiasingLimiterProcess
GB06BOptnSplitAndKillByImportance * fSplitAndKillByImportance

◆ ~GB06BOptrSplitAndKillByImportance()

GB06BOptrSplitAndKillByImportance::~GB06BOptrSplitAndKillByImportance ( )
virtual

Definition at line 73 of file GB06BOptrSplitAndKillByImportance.cc.

74{
76}

Member Function Documentation

◆ StartRun()

void GB06BOptrSplitAndKillByImportance::StartRun ( )
virtual

Definition at line 80 of file GB06BOptrSplitAndKillByImportance.cc.

81{
82 // ---------------
83 // -- Setup stage:
84 // ---------------
85 // -- get the particle process manager...
86 const G4ProcessManager* processManager = fParticleToBias->GetProcessManager();
87 // -- ... to obtain the biasing information shared among this particle processes:
88 fBiasingSharedData = G4BiasingProcessInterface::GetSharedData( processManager );
89 // -- Remember the index of the parallel world:
90 fBiasingLimiterProcess = fBiasingSharedData->GetParallelGeometriesLimiterProcess();
92
93 // -- Setup the biasing operation:
94 fSplitAndKillByImportance-> SetBiasingSharedData( fBiasingSharedData );
96 fSplitAndKillByImportance-> SetImportanceMap( &fImportanceMap );
97
98}
const G4BiasingProcessSharedData * fBiasingSharedData

◆ ProposeOccurenceBiasingOperation()

virtual G4VBiasingOperation * GB06BOptrSplitAndKillByImportance::ProposeOccurenceBiasingOperation ( const G4Track ,
const G4BiasingProcessInterface  
)
inlinefinalprivatevirtual

Definition at line 78 of file GB06BOptrSplitAndKillByImportance.hh.

80 { return 0; }

◆ ProposeFinalStateBiasingOperation()

virtual G4VBiasingOperation * GB06BOptrSplitAndKillByImportance::ProposeFinalStateBiasingOperation ( const G4Track ,
const G4BiasingProcessInterface  
)
inlinefinalprivatevirtual

Definition at line 84 of file GB06BOptrSplitAndKillByImportance.hh.

86 { return 0; }

◆ ProposeNonPhysicsBiasingOperation()

G4VBiasingOperation * GB06BOptrSplitAndKillByImportance::ProposeNonPhysicsBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface  
)
finalprivatevirtual

Definition at line 103 of file GB06BOptrSplitAndKillByImportance.cc.

106{
107 // -- Check if current particle type is the one to bias:
108 if ( track->GetDefinition() != fParticleToBias ) return 0;
109
110 // -- if so, request biasing:
112
113}

◆ SetParallelWorld()

void GB06BOptrSplitAndKillByImportance::SetParallelWorld ( G4VPhysicalVolume parallelWorld)
inline

Definition at line 99 of file GB06BOptrSplitAndKillByImportance.hh.

100 { fParallelWorld = parallelWorld; }

◆ GetImportanceMap()

std::map< G4int, G4int > & GB06BOptrSplitAndKillByImportance::GetImportanceMap ( )
inline

Definition at line 103 of file GB06BOptrSplitAndKillByImportance.hh.

103{ return fImportanceMap; }

Member Data Documentation

◆ fSplitAndKillByImportance

GB06BOptnSplitAndKillByImportance* GB06BOptrSplitAndKillByImportance::fSplitAndKillByImportance
private

Definition at line 108 of file GB06BOptrSplitAndKillByImportance.hh.

◆ fParticleToBias

const G4ParticleDefinition* GB06BOptrSplitAndKillByImportance::fParticleToBias
private

Definition at line 109 of file GB06BOptrSplitAndKillByImportance.hh.

◆ fParallelWorld

G4VPhysicalVolume* GB06BOptrSplitAndKillByImportance::fParallelWorld
private

Definition at line 110 of file GB06BOptrSplitAndKillByImportance.hh.

◆ fParallelWorldIndex

G4int GB06BOptrSplitAndKillByImportance::fParallelWorldIndex
private

Definition at line 111 of file GB06BOptrSplitAndKillByImportance.hh.

◆ fBiasingSharedData

const G4BiasingProcessSharedData* GB06BOptrSplitAndKillByImportance::fBiasingSharedData
private

Definition at line 112 of file GB06BOptrSplitAndKillByImportance.hh.

◆ fBiasingLimiterProcess

const G4ParallelGeometriesLimiterProcess* GB06BOptrSplitAndKillByImportance::fBiasingLimiterProcess
private

Definition at line 113 of file GB06BOptrSplitAndKillByImportance.hh.

◆ fImportanceMap

std::map< G4int, G4int > GB06BOptrSplitAndKillByImportance::fImportanceMap
private

Definition at line 114 of file GB06BOptrSplitAndKillByImportance.hh.


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

Applications | User Support | Publications | Collaboration