Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChannelingOptrMultiParticleChangeCrossSection Class Reference

#include <G4ChannelingOptrMultiParticleChangeCrossSection.hh>

+ Inheritance diagram for G4ChannelingOptrMultiParticleChangeCrossSection:

Public Member Functions

 G4ChannelingOptrMultiParticleChangeCrossSection ()
 
virtual ~G4ChannelingOptrMultiParticleChangeCrossSection ()
 
void AddParticle (G4String particleName)
 
void AddChargedParticles ()
 
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 StartTracking (const G4Track *)
 
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 ()
 

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 G4VBiasingOperationProposeNonPhysicsBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
 
virtual G4VBiasingOperationProposeOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
 
virtual G4VBiasingOperationProposeFinalStateBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0
 
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void ExitBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 

Detailed Description

Constructor & Destructor Documentation

◆ G4ChannelingOptrMultiParticleChangeCrossSection()

G4ChannelingOptrMultiParticleChangeCrossSection::G4ChannelingOptrMultiParticleChangeCrossSection ( )

Definition at line 38 of file G4ChannelingOptrMultiParticleChangeCrossSection.cc.

38 :
39G4VBiasingOperator("ChannelingChangeXS-Many"),
40fCurrentOperator(0),
41fnInteractions(0){
43}

◆ ~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()

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
G4ParticleDefinition * FindParticle(G4int PDGEncoding)

Referenced by AddChargedParticles().

◆ 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 );
101 fCurrentOperator = 0;
102 if ( it != fBOptrForParticle.end() ) fCurrentOperator = (*it).second;
103 fnInteractions = 0;
104}
const G4ParticleDefinition * GetParticleDefinition() const

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