Geant4-11
Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
G4ChannelingOptrMultiParticleChangeCrossSection Class Reference

#include <G4ChannelingOptrMultiParticleChangeCrossSection.hh>

Inheritance diagram for G4ChannelingOptrMultiParticleChangeCrossSection:
G4VBiasingOperator

Public Member Functions

void AddChargedParticles ()
 
void AddParticle (G4String particleName)
 
void AttachTo (const G4LogicalVolume *)
 
virtual void Configure ()
 
virtual void ConfigureForWorker ()
 
virtual void EndTracking ()
 
void ExitingBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
 G4ChannelingOptrMultiParticleChangeCrossSection ()
 
const G4String GetName () const
 
G4BiasingAppliedCase GetPreviousBiasingAppliedCase () const
 
const G4VBiasingOperationGetPreviousNonPhysicsAppliedOperation ()
 
G4VBiasingOperationGetProposedFinalStateBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
G4VBiasingOperationGetProposedNonPhysicsBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
G4VBiasingOperationGetProposedOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
void ReportOperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 
void ReportOperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void StartRun ()
 
void StartTracking (const G4Track *track)
 
virtual ~G4ChannelingOptrMultiParticleChangeCrossSection ()
 

Static Public Member Functions

static G4VBiasingOperatorGetBiasingOperator (const G4LogicalVolume *)
 
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators ()
 

Protected Member Functions

virtual void ExitBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 

Private Member Functions

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

Private Attributes

std::map< const G4ParticleDefinition *, G4ChannelingOptrChangeCrossSection * > fBOptrForParticle
 
G4ChannelingOptrChangeCrossSectionfCurrentOperator
 
std::map< const G4LogicalVolume *, G4intfDepthInTree
 
G4VBiasingOperationfFinalStateBiasingOperation
 
const G4String fName
 
G4int fnInteractions
 
G4VBiasingOperationfNonPhysicsBiasingOperation
 
G4VBiasingOperationfOccurenceBiasingOperation
 
std::vector< const G4ParticleDefinition * > fParticlesToBias
 
const G4VBiasingOperationfPreviousAppliedFinalStateBiasingOperation
 
const G4VBiasingOperationfPreviousAppliedNonPhysicsBiasingOperation
 
const G4VBiasingOperationfPreviousAppliedOccurenceBiasingOperation
 
G4BiasingAppliedCase fPreviousBiasingAppliedCase
 
const G4VBiasingOperationfPreviousProposedFinalStateBiasingOperation
 
const G4VBiasingOperationfPreviousProposedNonPhysicsBiasingOperation
 
const G4VBiasingOperationfPreviousProposedOccurenceBiasingOperation
 
std::vector< const G4LogicalVolume * > fRootVolumes
 

Static Private Attributes

static G4MapCache< const G4LogicalVolume *, G4VBiasingOperator * > fLogicalToSetupMap
 
static G4VectorCache< G4VBiasingOperator * > fOperators
 
static G4Cache< G4BiasingOperatorStateNotifier * > fStateNotifier
 

Detailed Description

Definition at line 50 of file G4ChannelingOptrMultiParticleChangeCrossSection.hh.

Constructor & Destructor Documentation

◆ G4ChannelingOptrMultiParticleChangeCrossSection()

G4ChannelingOptrMultiParticleChangeCrossSection::G4ChannelingOptrMultiParticleChangeCrossSection ( )

◆ ~G4ChannelingOptrMultiParticleChangeCrossSection()

virtual G4ChannelingOptrMultiParticleChangeCrossSection::~G4ChannelingOptrMultiParticleChangeCrossSection ( )
inlinevirtual

Definition at line 53 of file G4ChannelingOptrMultiParticleChangeCrossSection.hh.

53{}

Member Function Documentation

◆ AddChargedParticles()

void G4ChannelingOptrMultiParticleChangeCrossSection::AddChargedParticles ( )

Definition at line 69 of file G4ChannelingOptrMultiParticleChangeCrossSection.cc.

69 {
70 G4ParticleTable::G4PTblDicIterator* aParticleIterator =
71 (G4ParticleTable::GetParticleTable())->GetIterator();
72
73 aParticleIterator->reset();
74
75 while( (*aParticleIterator)() ){
76 G4ParticleDefinition* particle = aParticleIterator->value();
77 if (particle->GetPDGCharge() !=0) {
78 AddParticle(particle->GetParticleName());
79 }
80 }
81}
G4double GetPDGCharge() const
const G4String & GetParticleName() const
void reset(G4bool ifSkipIon=true)
static G4ParticleTable * GetParticleTable()

References AddParticle(), G4ParticleDefinition::GetParticleName(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGCharge(), G4ParticleTableIterator< K, V >::reset(), and G4ParticleTableIterator< K, V >::value().

Referenced by G4ChannelingOptrMultiParticleChangeCrossSection().

◆ AddParticle()

void G4ChannelingOptrMultiParticleChangeCrossSection::AddParticle ( G4String  particleName)

Definition at line 47 of file G4ChannelingOptrMultiParticleChangeCrossSection.cc.

47 {
48 const G4ParticleDefinition* particle =
50
51 if ( particle == 0 )
52 {
54 ed << "Particle `" << particleName << "' not found !" << G4endl;
55 G4Exception("G4ChannelingOptrMultiParticleChangeCrossSection::AddParticle(...)",
56 "G4Channeling",
58 ed);
59 return;
60 }
61
63 fParticlesToBias.push_back( particle );
64 fBOptrForParticle[ particle ] = optr;
65}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4endl
Definition: G4ios.hh:57
std::map< const G4ParticleDefinition *, G4ChannelingOptrChangeCrossSection * > fBOptrForParticle
G4ParticleDefinition * FindParticle(G4int PDGEncoding)

References fBOptrForParticle, G4ParticleTable::FindParticle(), fParticlesToBias, G4endl, G4Exception(), G4ParticleTable::GetParticleTable(), and JustWarning.

Referenced by AddChargedParticles().

◆ AttachTo()

void G4VBiasingOperator::AttachTo ( const G4LogicalVolume logical)
inherited

Definition at line 59 of file G4VBiasingOperator.cc.

60{
62 it = fLogicalToSetupMap.Find(logical);
63 if ( it == fLogicalToSetupMap.End() ) fLogicalToSetupMap[logical] = this;
64 else if ( (*it).second != this )
65 {
67 ed << "Biasing operator `" << GetName()
68 << "' can not be attached to Logical volume `"
69 << logical->GetName() << "' which is already used by another operator !" << G4endl;
70 G4Exception("G4VBiasingOperator::AttachTo(...)",
71 "BIAS.MNG.01",
73 ed);
74 }
75}
const G4String & GetName() const
iterator Find(const key_type &k)
Definition: G4Cache.hh:459
iterator End()
Definition: G4Cache.hh:453
static G4MapCache< const G4LogicalVolume *, G4VBiasingOperator * > fLogicalToSetupMap
const G4String GetName() const

References G4MapCache< KEYTYPE, VALTYPE >::End(), G4MapCache< KEYTYPE, VALTYPE >::Find(), G4VBiasingOperator::fLogicalToSetupMap, G4endl, G4Exception(), G4LogicalVolume::GetName(), G4VBiasingOperator::GetName(), and JustWarning.

◆ Configure()

virtual void G4VBiasingOperator::Configure ( )
inlinevirtualinherited

Reimplemented in G4BOptrForceCollision.

Definition at line 271 of file G4VBiasingOperator.hh.

271{}

◆ ConfigureForWorker()

virtual void G4VBiasingOperator::ConfigureForWorker ( )
inlinevirtualinherited

Reimplemented in G4BOptrForceCollision.

Definition at line 274 of file G4VBiasingOperator.hh.

274{}

◆ EndTracking()

virtual void G4VBiasingOperator::EndTracking ( )
inlinevirtualinherited

Reimplemented in G4BOptrForceCollision.

Definition at line 279 of file G4VBiasingOperator.hh.

279{}

◆ ExitBiasing()

void G4VBiasingOperator::ExitBiasing ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
protectedvirtualinherited

Reimplemented in G4BOptrForceCollision.

Definition at line 173 of file G4VBiasingOperator.cc.

174{}

Referenced by G4VBiasingOperator::ExitingBiasing().

◆ ExitingBiasing()

void G4VBiasingOperator::ExitingBiasing ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
inherited

Definition at line 152 of file G4VBiasingOperator.cc.

153{
154 ExitBiasing( track, callingProcess );
155
156 // -- reset all data members:
167}
virtual void ExitBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
const G4VBiasingOperation * fPreviousProposedNonPhysicsBiasingOperation
G4VBiasingOperation * fNonPhysicsBiasingOperation
G4BiasingAppliedCase fPreviousBiasingAppliedCase
G4VBiasingOperation * fOccurenceBiasingOperation
G4VBiasingOperation * fFinalStateBiasingOperation
const G4VBiasingOperation * fPreviousProposedOccurenceBiasingOperation
const G4VBiasingOperation * fPreviousAppliedOccurenceBiasingOperation
const G4VBiasingOperation * fPreviousAppliedNonPhysicsBiasingOperation
const G4VBiasingOperation * fPreviousProposedFinalStateBiasingOperation
const G4VBiasingOperation * fPreviousAppliedFinalStateBiasingOperation

References BAC_None, G4VBiasingOperator::ExitBiasing(), G4VBiasingOperator::fFinalStateBiasingOperation, G4VBiasingOperator::fNonPhysicsBiasingOperation, G4VBiasingOperator::fOccurenceBiasingOperation, G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation, G4VBiasingOperator::fPreviousAppliedNonPhysicsBiasingOperation, G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation, G4VBiasingOperator::fPreviousBiasingAppliedCase, G4VBiasingOperator::fPreviousProposedFinalStateBiasingOperation, G4VBiasingOperator::fPreviousProposedNonPhysicsBiasingOperation, and G4VBiasingOperator::fPreviousProposedOccurenceBiasingOperation.

◆ GetBiasingOperator()

G4VBiasingOperator * G4VBiasingOperator::GetBiasingOperator ( const G4LogicalVolume logical)
staticinherited

◆ GetBiasingOperators()

static const std::vector< G4VBiasingOperator * > & G4VBiasingOperator::GetBiasingOperators ( )
inlinestaticinherited

◆ GetName()

const G4String G4VBiasingOperator::GetName ( ) const
inlineinherited

◆ GetPreviousBiasingAppliedCase()

G4BiasingAppliedCase G4VBiasingOperator::GetPreviousBiasingAppliedCase ( ) const
inlineinherited

◆ GetPreviousNonPhysicsAppliedOperation()

const G4VBiasingOperation * G4VBiasingOperator::GetPreviousNonPhysicsAppliedOperation ( )
inlineinherited

◆ GetProposedFinalStateBiasingOperation()

G4VBiasingOperation * G4VBiasingOperator::GetProposedFinalStateBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
inherited

Definition at line 92 of file G4VBiasingOperator.cc.

93{
96}
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

References G4VBiasingOperator::fFinalStateBiasingOperation, and G4VBiasingOperator::ProposeFinalStateBiasingOperation().

◆ GetProposedNonPhysicsBiasingOperation()

G4VBiasingOperation * G4VBiasingOperator::GetProposedNonPhysicsBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
inherited

Definition at line 98 of file G4VBiasingOperator.cc.

99{
102}
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

References G4VBiasingOperator::fNonPhysicsBiasingOperation, and G4VBiasingOperator::ProposeNonPhysicsBiasingOperation().

◆ GetProposedOccurenceBiasingOperation()

G4VBiasingOperation * G4VBiasingOperator::GetProposedOccurenceBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
inherited

Definition at line 86 of file G4VBiasingOperator.cc.

87{
90}
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

References G4VBiasingOperator::fOccurenceBiasingOperation, and G4VBiasingOperator::ProposeOccurenceBiasingOperation().

Referenced by ProposeOccurenceBiasingOperation().

◆ OperationApplied() [1/3]

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

Reimplemented from G4VBiasingOperator.

Definition at line 238 of file G4VBiasingOperator.cc.

182{
183}

◆ OperationApplied() [2/3]

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

Reimplemented from G4VBiasingOperator.

Definition at line 109 of file G4ChannelingOptrMultiParticleChangeCrossSection.cc.

115 {
118 biasingCase,
119 occurenceOperationApplied,
120 weightForOccurenceInteraction,
121 finalStateOperationApplied,
122 particleChangeProduced );
123}
void ReportOperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)

References fCurrentOperator, fnInteractions, and G4VBiasingOperator::ReportOperationApplied().

◆ OperationApplied() [3/3]

void G4VBiasingOperator::OperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation operationApplied,
const G4VParticleChange particleChangeProduced 
)
privatevirtual

Reimplemented from G4VBiasingOperator.

Definition at line 232 of file G4VBiasingOperator.cc.

177{
178}

◆ ProposeFinalStateBiasingOperation()

virtual G4VBiasingOperation * G4ChannelingOptrMultiParticleChangeCrossSection::ProposeFinalStateBiasingOperation ( const G4Track ,
const G4BiasingProcessInterface  
)
inlineprivatevirtual

Implements G4VBiasingOperator.

Definition at line 75 of file G4ChannelingOptrMultiParticleChangeCrossSection.hh.

76 {return 0;}

◆ ProposeNonPhysicsBiasingOperation()

virtual G4VBiasingOperation * G4ChannelingOptrMultiParticleChangeCrossSection::ProposeNonPhysicsBiasingOperation ( const G4Track ,
const G4BiasingProcessInterface  
)
inlineprivatevirtual

Implements G4VBiasingOperator.

Definition at line 78 of file G4ChannelingOptrMultiParticleChangeCrossSection.hh.

79 {return 0;}

◆ ProposeOccurenceBiasingOperation()

G4VBiasingOperation * G4ChannelingOptrMultiParticleChangeCrossSection::ProposeOccurenceBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
privatevirtual

Implements G4VBiasingOperator.

Definition at line 86 of file G4ChannelingOptrMultiParticleChangeCrossSection.cc.

88 {
89
91 GetProposedOccurenceBiasingOperation(track, callingProcess);
92 else return 0;
93}
G4VBiasingOperation * GetProposedOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)

References fCurrentOperator, and G4VBiasingOperator::GetProposedOccurenceBiasingOperation().

◆ ReportOperationApplied() [1/2]

void G4VBiasingOperator::ReportOperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation occurenceOperationApplied,
G4double  weightForOccurenceInteraction,
G4VBiasingOperation finalStateOperationApplied,
const G4VParticleChange particleChangeProduced 
)
inherited

Definition at line 138 of file G4VBiasingOperator.cc.

144{
145 fPreviousBiasingAppliedCase = biasingCase;
146 fPreviousAppliedOccurenceBiasingOperation = occurenceOperationApplied;
147 fPreviousAppliedFinalStateBiasingOperation = finalStateOperationApplied;
148 OperationApplied( callingProcess, biasingCase, occurenceOperationApplied, weightForOccurenceInteraction, finalStateOperationApplied, particleChangeProduced );
149}
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)

References G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation, G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation, G4VBiasingOperator::fPreviousBiasingAppliedCase, and G4VBiasingOperator::OperationApplied().

◆ ReportOperationApplied() [2/2]

void G4VBiasingOperator::ReportOperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation operationApplied,
const G4VParticleChange particleChangeProduced 
)
inherited

Definition at line 104 of file G4VBiasingOperator.cc.

108{
109 fPreviousBiasingAppliedCase = biasingCase;
113 switch ( biasingCase )
114 {
115 case BAC_None:
116 break;
117 case BAC_NonPhysics:
119 break;
120 case BAC_FinalState:
122 break;
123 case BAC_Occurence:
124 G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
125 "BIAS.MNG.02",
127 "Internal logic error, please report !");
128 break;
129 default:
130 G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
131 "BIAS.MNG.03",
133 "Internal logic error, please report !");
134 }
135 OperationApplied( callingProcess, biasingCase, operationApplied, particleChangeProduced );
136}
@ BAC_Occurence
@ BAC_FinalState
@ BAC_NonPhysics

References BAC_FinalState, BAC_None, BAC_NonPhysics, BAC_Occurence, G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation, G4VBiasingOperator::fPreviousAppliedNonPhysicsBiasingOperation, G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation, G4VBiasingOperator::fPreviousBiasingAppliedCase, G4Exception(), JustWarning, and G4VBiasingOperator::OperationApplied().

Referenced by OperationApplied().

◆ StartRun()

virtual void G4VBiasingOperator::StartRun ( )
inlinevirtualinherited

Reimplemented in G4ChannelingOptrChangeCrossSection, and G4BOptrForceCollision.

Definition at line 276 of file G4VBiasingOperator.hh.

276{}

◆ StartTracking()

void G4ChannelingOptrMultiParticleChangeCrossSection::StartTracking ( const G4Track track)
virtual

Reimplemented from G4VBiasingOperator.

Definition at line 97 of file G4ChannelingOptrMultiParticleChangeCrossSection.cc.

97 {
98 const G4ParticleDefinition* definition = track->GetParticleDefinition();
99 std::map < const G4ParticleDefinition*, G4ChannelingOptrChangeCrossSection* > :: iterator
100 it = fBOptrForParticle.find( definition );
102 if ( it != fBOptrForParticle.end() ) fCurrentOperator = (*it).second;
103 fnInteractions = 0;
104}
const G4ParticleDefinition * GetParticleDefinition() const

References fBOptrForParticle, fCurrentOperator, fnInteractions, and G4Track::GetParticleDefinition().

Field Documentation

◆ fBOptrForParticle

std::map< const G4ParticleDefinition*, G4ChannelingOptrChangeCrossSection* > G4ChannelingOptrMultiParticleChangeCrossSection::fBOptrForParticle
private

◆ fCurrentOperator

G4ChannelingOptrChangeCrossSection* G4ChannelingOptrMultiParticleChangeCrossSection::fCurrentOperator
private

◆ fDepthInTree

std::map< const G4LogicalVolume*, G4int > G4VBiasingOperator::fDepthInTree
privateinherited

Definition at line 333 of file G4VBiasingOperator.hh.

◆ fFinalStateBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fFinalStateBiasingOperation
privateinherited

◆ fLogicalToSetupMap

G4MapCache< const G4LogicalVolume *, G4VBiasingOperator * > G4VBiasingOperator::fLogicalToSetupMap
staticprivateinherited

◆ fName

const G4String G4VBiasingOperator::fName
privateinherited

Definition at line 319 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::GetName().

◆ fnInteractions

G4int G4ChannelingOptrMultiParticleChangeCrossSection::fnInteractions
private

◆ fNonPhysicsBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fNonPhysicsBiasingOperation
privateinherited

◆ fOccurenceBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fOccurenceBiasingOperation
privateinherited

◆ fOperators

G4VectorCache< G4VBiasingOperator * > G4VBiasingOperator::fOperators
staticprivateinherited

◆ fParticlesToBias

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

Definition at line 105 of file G4ChannelingOptrMultiParticleChangeCrossSection.hh.

Referenced by AddParticle().

◆ fPreviousAppliedFinalStateBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation
privateinherited

◆ fPreviousAppliedNonPhysicsBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedNonPhysicsBiasingOperation
privateinherited

◆ fPreviousAppliedOccurenceBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation
privateinherited

◆ fPreviousBiasingAppliedCase

G4BiasingAppliedCase G4VBiasingOperator::fPreviousBiasingAppliedCase
privateinherited

◆ fPreviousProposedFinalStateBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedFinalStateBiasingOperation
privateinherited

Definition at line 342 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::ExitingBiasing().

◆ fPreviousProposedNonPhysicsBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedNonPhysicsBiasingOperation
privateinherited

Definition at line 343 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::ExitingBiasing().

◆ fPreviousProposedOccurenceBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedOccurenceBiasingOperation
privateinherited

Definition at line 341 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::ExitingBiasing().

◆ fRootVolumes

std::vector< const G4LogicalVolume* > G4VBiasingOperator::fRootVolumes
privateinherited

Definition at line 332 of file G4VBiasingOperator.hh.

◆ fStateNotifier

G4Cache< G4BiasingOperatorStateNotifier * > G4VBiasingOperator::fStateNotifier
staticprivateinherited

Definition at line 328 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::G4VBiasingOperator().


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