Geant4  10.02.p03
GB02BOptrMultiParticleForceCollision Class Reference

#include <GB02BOptrMultiParticleForceCollision.hh>

Inheritance diagram for GB02BOptrMultiParticleForceCollision:
Collaboration diagram for GB02BOptrMultiParticleForceCollision:

Public Member Functions

 GB02BOptrMultiParticleForceCollision ()
 
virtual ~GB02BOptrMultiParticleForceCollision ()
 
void AddParticle (G4String particleName)
 
virtual void StartTracking (const G4Track *track) final
 
- Public Member Functions inherited from G4VBiasingOperator
 G4VBiasingOperator (G4String name)
 
virtual ~G4VBiasingOperator ()
 
virtual void Configure ()
 
virtual void ConfigureForWorker ()
 
virtual void StartRun ()
 
virtual void EndTracking ()
 
const G4String GetName () const
 
void AttachTo (const G4LogicalVolume *)
 
G4BiasingAppliedCase GetPreviousBiasingAppliedCase () const
 
G4VBiasingOperationGetProposedOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
G4VBiasingOperationGetProposedFinalStateBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
G4VBiasingOperationGetProposedNonPhysicsBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
void ExitingBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
void ReportOperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
void ReportOperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 
const G4VBiasingOperationGetPreviousNonPhysicsAppliedOperation ()
 

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
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VBiasingOperator
static const std::vector< G4VBiasingOperator *> & GetBiasingOperators ()
 
static G4VBiasingOperatorGetBiasingOperator (const G4LogicalVolume *)
 

Detailed Description

Definition at line 35 of file GB02BOptrMultiParticleForceCollision.hh.

Constructor & Destructor Documentation

◆ GB02BOptrMultiParticleForceCollision()

GB02BOptrMultiParticleForceCollision::GB02BOptrMultiParticleForceCollision ( )

Definition at line 33 of file GB02BOptrMultiParticleForceCollision.cc.

34 : G4VBiasingOperator("TestManyForceCollision")
35 {
36 }
G4VBiasingOperator(G4String name)

◆ ~GB02BOptrMultiParticleForceCollision()

virtual GB02BOptrMultiParticleForceCollision::~GB02BOptrMultiParticleForceCollision ( )
inlinevirtual

Definition at line 38 of file GB02BOptrMultiParticleForceCollision.hh.

38 {}
Here is the call graph for this function:

Member Function Documentation

◆ AddParticle()

void GB02BOptrMultiParticleForceCollision::AddParticle ( G4String  particleName)

Definition at line 38 of file GB02BOptrMultiParticleForceCollision.cc.

39 {
40  const G4ParticleDefinition* particle =
42 
43  if ( particle == 0 )
44  {
46  ed << "Particle `" << particleName << "' not found !" << G4endl;
47  G4Exception("GB02BOptrMultiParticleForceCollision::AddParticle(...)",
48  "exGB02.01",
50  ed);
51  return;
52  }
53 
54  G4BOptrForceCollision* optr = new G4BOptrForceCollision(particleName,
55  "ForceCollisionFor"+particleName);
56  fParticlesToBias.push_back( particle );
57  fBOptrForParticle[ particle ] = optr;
58 
59 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
std::vector< const G4ParticleDefinition *> fParticlesToBias
std::map< const G4ParticleDefinition *, G4BOptrForceCollision *> fBOptrForParticle
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4ParticleTable * GetParticleTable()
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ExitBiasing()

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

Reimplemented from G4VBiasingOperator.

Definition at line 134 of file GB02BOptrMultiParticleForceCollision.cc.

136 {
137  if ( fCurrentOperator ) fCurrentOperator->ExitingBiasing( track, callingProcess );
138 }
void ExitingBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ OperationApplied() [1/2]

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

Reimplemented from G4VBiasingOperator.

Definition at line 103 of file GB02BOptrMultiParticleForceCollision.cc.

107 {
109  biasingCase,
110  operationApplied,
111  particleChangeProduced );
112 
113 }
void ReportOperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ OperationApplied() [2/2]

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

Reimplemented from G4VBiasingOperator.

Definition at line 117 of file GB02BOptrMultiParticleForceCollision.cc.

123 {
125  biasingCase,
126  occurenceOperationApplied,
127  weightForOccurenceInteraction,
128  finalStateOperationApplied,
129  particleChangeProduced );
130 }
void ReportOperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
Here is the call graph for this function:

◆ ProposeFinalStateBiasingOperation()

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

Implements G4VBiasingOperator.

Definition at line 83 of file GB02BOptrMultiParticleForceCollision.cc.

85 {
86  if ( fCurrentOperator ) return fCurrentOperator->
88  else return 0;
89 }
G4VBiasingOperation * GetProposedFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ProposeNonPhysicsBiasingOperation()

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

Implements G4VBiasingOperator.

Definition at line 73 of file GB02BOptrMultiParticleForceCollision.cc.

75 {
76  if ( fCurrentOperator ) return fCurrentOperator->
78  else return 0;
79 }
G4VBiasingOperation * GetProposedNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ProposeOccurenceBiasingOperation()

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

Implements G4VBiasingOperator.

Definition at line 63 of file GB02BOptrMultiParticleForceCollision.cc.

65 {
66  if ( fCurrentOperator ) return fCurrentOperator->
68  else return 0;
69 }
G4VBiasingOperation * GetProposedOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ StartTracking()

void GB02BOptrMultiParticleForceCollision::StartTracking ( const G4Track *  track)
finalvirtual

Reimplemented from G4VBiasingOperator.

Definition at line 92 of file GB02BOptrMultiParticleForceCollision.cc.

93 {
94  const G4ParticleDefinition* definition = track->GetParticleDefinition();
95  std::map < const G4ParticleDefinition*, G4BOptrForceCollision* > :: iterator
96  it = fBOptrForParticle.find( definition );
97  fCurrentOperator = 0;
98  if ( it != fBOptrForParticle.end() ) fCurrentOperator = (*it).second;
99 }
std::map< const G4ParticleDefinition *, G4BOptrForceCollision *> fBOptrForParticle
Here is the caller graph for this function:

Member Data Documentation

◆ fBOptrForParticle

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

Definition at line 83 of file GB02BOptrMultiParticleForceCollision.hh.

◆ fCurrentOperator

G4BOptrForceCollision* GB02BOptrMultiParticleForceCollision::fCurrentOperator
private

Definition at line 85 of file GB02BOptrMultiParticleForceCollision.hh.

◆ fParticlesToBias

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

Definition at line 84 of file GB02BOptrMultiParticleForceCollision.hh.


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