Geant4 11.4.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 (const G4String &particleName)
void AddChargedParticles ()
void StartTracking (const G4Track *track)
Public Member Functions inherited from G4VBiasingOperator
 G4VBiasingOperator (const G4String &name)
virtual ~G4VBiasingOperator ()=default
virtual void Configure ()
virtual void ConfigureForWorker ()
virtual void StartRun ()
virtual void EndTracking ()
const G4StringGetName () 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 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}
G4VBiasingOperator(const G4String &name)

◆ ~G4ChannelingOptrMultiParticleChangeCrossSection()

virtual G4ChannelingOptrMultiParticleChangeCrossSection::~G4ChannelingOptrMultiParticleChangeCrossSection ( )
inlinevirtual

Definition at line 53 of file G4ChannelingOptrMultiParticleChangeCrossSection.hh.

53{}

Member Function Documentation

◆ AddChargedParticles()

void G4ChannelingOptrMultiParticleChangeCrossSection::AddChargedParticles ( )

Definition at line 70 of file G4ChannelingOptrMultiParticleChangeCrossSection.cc.

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

Referenced by G4ChannelingOptrMultiParticleChangeCrossSection().

◆ AddParticle()

void G4ChannelingOptrMultiParticleChangeCrossSection::AddParticle ( const 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
62 G4ChannelingOptrChangeCrossSection* optr =
63 new G4ChannelingOptrChangeCrossSection(particleName);
64 fParticlesToBias.push_back( particle );
65 fBOptrForParticle[ particle ] = optr;
66}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
#define G4endl
Definition G4ios.hh:67
G4ParticleDefinition * FindParticle(G4int PDGEncoding)

Referenced by AddChargedParticles().

◆ StartTracking()

void G4ChannelingOptrMultiParticleChangeCrossSection::StartTracking ( const G4Track * track)
virtual

Reimplemented from G4VBiasingOperator.

Definition at line 98 of file G4ChannelingOptrMultiParticleChangeCrossSection.cc.

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

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