Geant4  10.02.p03
GB01BOptrMultiParticleChangeCrossSection Class Reference

#include <GB01BOptrMultiParticleChangeCrossSection.hh>

Inheritance diagram for GB01BOptrMultiParticleChangeCrossSection:
Collaboration diagram for GB01BOptrMultiParticleChangeCrossSection:

Public Member Functions

 GB01BOptrMultiParticleChangeCrossSection ()
 
virtual ~GB01BOptrMultiParticleChangeCrossSection ()
 
void AddParticle (G4String particleName)
 
void StartTracking (const G4Track *track)
 
- 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 G4VBiasingOperationProposeOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
virtual G4VBiasingOperationProposeFinalStateBiasingOperation (const G4Track *, const G4BiasingProcessInterface *)
 
virtual G4VBiasingOperationProposeNonPhysicsBiasingOperation (const G4Track *, const G4BiasingProcessInterface *)
 
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 

Private Attributes

std::map< const G4ParticleDefinition *, GB01BOptrChangeCrossSection *> fBOptrForParticle
 
std::vector< const G4ParticleDefinition *> fParticlesToBias
 
GB01BOptrChangeCrossSectionfCurrentOperator
 
G4int fnInteractions
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VBiasingOperator
static const std::vector< G4VBiasingOperator *> & GetBiasingOperators ()
 
static G4VBiasingOperatorGetBiasingOperator (const G4LogicalVolume *)
 
- Protected Member Functions inherited from G4VBiasingOperator
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void ExitBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 

Detailed Description

Definition at line 51 of file GB01BOptrMultiParticleChangeCrossSection.hh.

Constructor & Destructor Documentation

◆ GB01BOptrMultiParticleChangeCrossSection()

GB01BOptrMultiParticleChangeCrossSection::GB01BOptrMultiParticleChangeCrossSection ( )

Definition at line 35 of file GB01BOptrMultiParticleChangeCrossSection.cc.

36  : G4VBiasingOperator("TestManyExponentialTransform")
37 {}
G4VBiasingOperator(G4String name)

◆ ~GB01BOptrMultiParticleChangeCrossSection()

virtual GB01BOptrMultiParticleChangeCrossSection::~GB01BOptrMultiParticleChangeCrossSection ( )
inlinevirtual

Definition at line 54 of file GB01BOptrMultiParticleChangeCrossSection.hh.

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

Member Function Documentation

◆ AddParticle()

void GB01BOptrMultiParticleChangeCrossSection::AddParticle ( G4String  particleName)

Definition at line 39 of file GB01BOptrMultiParticleChangeCrossSection.cc.

40 {
41  const G4ParticleDefinition* particle =
43 
44  if ( particle == 0 )
45  {
47  ed << "Particle `" << particleName << "' not found !" << G4endl;
48  G4Exception("GB01BOptrMultiParticleChangeCrossSection::AddParticle(...)",
49  "exGB01.02",
51  ed);
52  return;
53  }
54 
56  fParticlesToBias.push_back( particle );
57  fBOptrForParticle[ particle ] = optr;
58 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
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
std::vector< const G4ParticleDefinition *> fParticlesToBias
std::map< const G4ParticleDefinition *, GB01BOptrChangeCrossSection *> fBOptrForParticle
Here is the call graph for this function:
Here is the caller graph for this function:

◆ OperationApplied()

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

Reimplemented from G4VBiasingOperator.

Definition at line 94 of file GB01BOptrMultiParticleChangeCrossSection.cc.

100 {
101  // -- count number of biased interactions:
102  fnInteractions++;
103 
104  // -- inform the underneath biasing operator that a biased interaction occured:
106  biasingCase,
107  occurenceOperationApplied,
108  weightForOccurenceInteraction,
109  finalStateOperationApplied,
110  particleChangeProduced );
111 }
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:

◆ ProposeFinalStateBiasingOperation()

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

Implements G4VBiasingOperator.

Definition at line 75 of file GB01BOptrMultiParticleChangeCrossSection.hh.

76  {return 0;}

◆ ProposeNonPhysicsBiasingOperation()

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

Implements G4VBiasingOperator.

Definition at line 78 of file GB01BOptrMultiParticleChangeCrossSection.hh.

79  {return 0;}
Here is the call graph for this function:

◆ ProposeOccurenceBiasingOperation()

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

Implements G4VBiasingOperator.

Definition at line 62 of file GB01BOptrMultiParticleChangeCrossSection.cc.

64 {
65  // -- examples of limitations imposed to apply the biasing:
66  // -- limit application of biasing to primary particles only:
67  if ( track->GetParentID() != 0 ) return 0;
68  // -- limit to at most 5 biased interactions:
69  if ( fnInteractions > 4 ) return 0;
70  // -- and limit to a weight of at least 0.05:
71  if ( track->GetWeight() < 0.05 ) return 0;
72 
73  if ( fCurrentOperator ) return fCurrentOperator->
75  else return 0;
76 }
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 GB01BOptrMultiParticleChangeCrossSection::StartTracking ( const G4Track *  track)
virtual

Reimplemented from G4VBiasingOperator.

Definition at line 79 of file GB01BOptrMultiParticleChangeCrossSection.cc.

80 {
81  // -- fetch the underneath biasing operator, if any, for the current particle type:
82  const G4ParticleDefinition* definition = track->GetParticleDefinition();
83  std::map < const G4ParticleDefinition*, GB01BOptrChangeCrossSection* > :: iterator
84  it = fBOptrForParticle.find( definition );
85  fCurrentOperator = 0;
86  if ( it != fBOptrForParticle.end() ) fCurrentOperator = (*it).second;
87 
88  // -- reset count for number of biased interactions:
89  fnInteractions = 0;
90 }
std::map< const G4ParticleDefinition *, GB01BOptrChangeCrossSection *> fBOptrForParticle
Here is the caller graph for this function:

Member Data Documentation

◆ fBOptrForParticle

std::map< const G4ParticleDefinition*, GB01BOptrChangeCrossSection* > GB01BOptrMultiParticleChangeCrossSection::fBOptrForParticle
private

Definition at line 105 of file GB01BOptrMultiParticleChangeCrossSection.hh.

◆ fCurrentOperator

GB01BOptrChangeCrossSection* GB01BOptrMultiParticleChangeCrossSection::fCurrentOperator
private

Definition at line 107 of file GB01BOptrMultiParticleChangeCrossSection.hh.

◆ fnInteractions

G4int GB01BOptrMultiParticleChangeCrossSection::fnInteractions
private

Definition at line 110 of file GB01BOptrMultiParticleChangeCrossSection.hh.

◆ fParticlesToBias

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

Definition at line 106 of file GB01BOptrMultiParticleChangeCrossSection.hh.


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