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

#include <Doxymodules_biasing.h>

Inheritance diagram for GB02BOptrMultiParticleForceCollision:
G4VBiasingOperator

Public Member Functions

 GB02BOptrMultiParticleForceCollision ()
 
virtual ~GB02BOptrMultiParticleForceCollision ()
 
void AddParticle (G4String particleName)
 
virtual void StartTracking (const G4Track *track) final
 

Private Member Functions

virtual G4VBiasingOperationProposeNonPhysicsBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess) final
 
virtual G4VBiasingOperationProposeOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess) final
 
virtual G4VBiasingOperationProposeFinalStateBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess) final
 
void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced) final
 
void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced) final
 
void ExitBiasing (const G4Track *, const G4BiasingProcessInterface *) final
 

Private Attributes

std::map< const G4ParticleDefinition *, G4BOptrForceCollision * > fBOptrForParticle
 
std::vector< const G4ParticleDefinition * > fParticlesToBias
 
G4BOptrForceCollisionfCurrentOperator
 

Detailed Description

Definition at line 81 of file Doxymodules_biasing.h.

Constructor & Destructor Documentation

◆ GB02BOptrMultiParticleForceCollision()

GB02BOptrMultiParticleForceCollision::GB02BOptrMultiParticleForceCollision ( )

Definition at line 38 of file GB02BOptrMultiParticleForceCollision.cc.

39: G4VBiasingOperator("TestManyForceCollision")
40{
41}

◆ ~GB02BOptrMultiParticleForceCollision()

virtual GB02BOptrMultiParticleForceCollision::~GB02BOptrMultiParticleForceCollision ( )
inlinevirtual

Definition at line 41 of file GB02BOptrMultiParticleForceCollision.hh.

41{}

Member Function Documentation

◆ AddParticle()

void GB02BOptrMultiParticleForceCollision::AddParticle ( G4String  particleName)

Definition at line 45 of file GB02BOptrMultiParticleForceCollision.cc.

46{
47 const G4ParticleDefinition* particle =
48 G4ParticleTable::GetParticleTable()->FindParticle( particleName );
49
50 if ( particle == 0 )
51 {
52 G4ExceptionDescription ed;
53 ed << "Particle `" << particleName << "' not found !" << G4endl;
54 G4Exception("GB02BOptrMultiParticleForceCollision::AddParticle(...)",
55 "exGB02.01",
56 JustWarning,
57 ed);
58 return;
59 }
60
61 G4BOptrForceCollision* optr = new G4BOptrForceCollision(particleName,
62 "ForceCollisionFor"+particleName);
63 fParticlesToBias.push_back( particle );
64 fBOptrForParticle[ particle ] = optr;
65
66}
std::map< const G4ParticleDefinition *, G4BOptrForceCollision * > fBOptrForParticle
std::vector< const G4ParticleDefinition * > fParticlesToBias

◆ ProposeNonPhysicsBiasingOperation()

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

Definition at line 83 of file GB02BOptrMultiParticleForceCollision.cc.

86{
88 GetProposedNonPhysicsBiasingOperation(track, callingProcess);
89 else return 0;
90}

◆ ProposeOccurenceBiasingOperation()

G4VBiasingOperation * GB02BOptrMultiParticleForceCollision::ProposeOccurenceBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
finalprivatevirtual

Definition at line 71 of file GB02BOptrMultiParticleForceCollision.cc.

74{
76 GetProposedOccurenceBiasingOperation(track, callingProcess);
77 else return 0;
78}

◆ ProposeFinalStateBiasingOperation()

G4VBiasingOperation * GB02BOptrMultiParticleForceCollision::ProposeFinalStateBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
finalprivatevirtual

Definition at line 95 of file GB02BOptrMultiParticleForceCollision.cc.

98{
100 GetProposedFinalStateBiasingOperation(track, callingProcess);
101 else return 0;
102}

◆ OperationApplied() [1/2]

void GB02BOptrMultiParticleForceCollision::OperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation operationApplied,
const G4VParticleChange particleChangeProduced 
)
finalprivate

Definition at line 118 of file GB02BOptrMultiParticleForceCollision.cc.

123{
124 if ( fCurrentOperator ) fCurrentOperator->ReportOperationApplied( callingProcess,
125 biasingCase,
126 operationApplied,
127 particleChangeProduced );
128
129}

◆ OperationApplied() [2/2]

void GB02BOptrMultiParticleForceCollision::OperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation occurenceOperationApplied,
G4double  weightForOccurenceInteraction,
G4VBiasingOperation finalStateOperationApplied,
const G4VParticleChange particleChangeProduced 
)
finalprivate

Definition at line 134 of file GB02BOptrMultiParticleForceCollision.cc.

141{
142 if ( fCurrentOperator ) fCurrentOperator->ReportOperationApplied( callingProcess,
143 biasingCase,
144 occurenceOperationApplied,
145 weightForOccurenceInteraction,
146 finalStateOperationApplied,
147 particleChangeProduced );
148}

◆ ExitBiasing()

void GB02BOptrMultiParticleForceCollision::ExitBiasing ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
finalprivate

Definition at line 153 of file GB02BOptrMultiParticleForceCollision.cc.

156{
157 if ( fCurrentOperator ) fCurrentOperator->ExitingBiasing( track, callingProcess );
158}

◆ StartTracking()

void GB02BOptrMultiParticleForceCollision::StartTracking ( const G4Track track)
finalvirtual

Definition at line 106 of file GB02BOptrMultiParticleForceCollision.cc.

107{
108 const G4ParticleDefinition* definition = track->GetParticleDefinition();
109 std::map < const G4ParticleDefinition*, G4BOptrForceCollision* > :: iterator
110 it = fBOptrForParticle.find( definition );
112 if ( it != fBOptrForParticle.end() ) fCurrentOperator = (*it).second;
113}

Member Data Documentation

◆ fBOptrForParticle

std::map< const G4ParticleDefinition*, G4BOptrForceCollision* > GB02BOptrMultiParticleForceCollision::fBOptrForParticle
private

Definition at line 86 of file GB02BOptrMultiParticleForceCollision.hh.

◆ fParticlesToBias

std::vector< const G4ParticleDefinition* > GB02BOptrMultiParticleForceCollision::fParticlesToBias
private

Definition at line 87 of file GB02BOptrMultiParticleForceCollision.hh.

◆ fCurrentOperator

G4BOptrForceCollision* GB02BOptrMultiParticleForceCollision::fCurrentOperator
private

Definition at line 88 of file GB02BOptrMultiParticleForceCollision.hh.


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

Applications | User Support | Publications | Collaboration