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

#include <G4SDParticleFilter.hh>

+ Inheritance diagram for G4SDParticleFilter:

Public Member Functions

 G4SDParticleFilter (G4String name)
 
 G4SDParticleFilter (G4String name, const G4String &particleName)
 
 G4SDParticleFilter (G4String name, const std::vector< G4String > &particleNames)
 
 G4SDParticleFilter (G4String name, const std::vector< G4ParticleDefinition * > &particleDef)
 
 ~G4SDParticleFilter () override=default
 
G4bool Accept (const G4Step *) const override
 
void add (const G4String &particleName)
 
void addIon (G4int Z, G4int A)
 
void show ()
 
- Public Member Functions inherited from G4VSDFilter
 G4VSDFilter (G4String name)
 
virtual ~G4VSDFilter ()
 
G4String GetName () const
 

Additional Inherited Members

- Protected Attributes inherited from G4VSDFilter
G4String filterName
 

Detailed Description

Definition at line 52 of file G4SDParticleFilter.hh.

Constructor & Destructor Documentation

◆ G4SDParticleFilter() [1/4]

G4SDParticleFilter::G4SDParticleFilter ( G4String name)

Definition at line 45 of file G4SDParticleFilter.cc.

46 : G4VSDFilter(name)
47{}
G4VSDFilter(G4String name)

◆ G4SDParticleFilter() [2/4]

G4SDParticleFilter::G4SDParticleFilter ( G4String name,
const G4String & particleName )

Definition at line 49 of file G4SDParticleFilter.cc.

51 : G4VSDFilter(name)
52{
55 if(pd == nullptr)
56 {
57 G4String msg = "Particle <";
58 msg += particleName;
59 msg += "> not found.";
60 G4Exception("G4SDParticleFilter::G4SDParticleFilter", "DetPS0101",
61 FatalException, msg);
62 }
63 thePdef.push_back(pd);
64}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()

◆ G4SDParticleFilter() [3/4]

G4SDParticleFilter::G4SDParticleFilter ( G4String name,
const std::vector< G4String > & particleNames )

Definition at line 66 of file G4SDParticleFilter.cc.

68 : G4VSDFilter(name)
69{
70 for(const auto & particleName : particleNames)
71 {
74 if(pd == nullptr)
75 {
76 G4String msg = "Particle <";
77 msg += particleName;
78 msg += "> not found.";
79 G4Exception("G4SDParticleFilter::G4SDParticleFilter", "DetPS0102",
80 FatalException, msg);
81 }
82 thePdef.push_back(pd);
83 }
84}

◆ G4SDParticleFilter() [4/4]

G4SDParticleFilter::G4SDParticleFilter ( G4String name,
const std::vector< G4ParticleDefinition * > & particleDef )

Definition at line 86 of file G4SDParticleFilter.cc.

88 : G4VSDFilter(name)
89 , thePdef(particleDef)
90{
91 for(auto i : particleDef)
92 {
93 if(i == nullptr)
94 G4Exception("G4SDParticleFilter::G4SDParticleFilter", "DetPS0103",
96 "NULL pointer is found in the given particleDef vector.");
97 }
98}

◆ ~G4SDParticleFilter()

G4SDParticleFilter::~G4SDParticleFilter ( )
overridedefault

Member Function Documentation

◆ Accept()

G4bool G4SDParticleFilter::Accept ( const G4Step * aStep) const
overridevirtual

Implements G4VSDFilter.

Definition at line 100 of file G4SDParticleFilter.cc.

101{
102 for(auto i : thePdef)
103 {
104 if(i == aStep->GetTrack()->GetDefinition())
105 return true;
106 }
107
108 // Ions by Z,A
109 for(size_t i = 0; i < theIonZ.size(); i++)
110 {
111 if(theIonZ[i] == aStep->GetTrack()->GetDefinition()->GetAtomicNumber() &&
112 theIonA[i] == aStep->GetTrack()->GetDefinition()->GetAtomicMass())
113 {
114 return true;
115 }
116 }
117
118 return false;
119}
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
G4Track * GetTrack() const
G4ParticleDefinition * GetDefinition() const

Referenced by G4SDParticleWithEnergyFilter::Accept().

◆ add()

void G4SDParticleFilter::add ( const G4String & particleName)

Definition at line 121 of file G4SDParticleFilter.cc.

122{
125 if(pd == nullptr)
126 {
127 G4String msg = "Particle <";
128 msg += particleName;
129 msg += "> not found.";
130 G4Exception("G4SDParticleFilter::add()", "DetPS0104", FatalException, msg);
131 }
132 for(auto & i : thePdef)
133 {
134 if(i == pd)
135 return;
136 }
137 thePdef.push_back(pd);
138}

Referenced by G4SDParticleWithEnergyFilter::add().

◆ addIon()

void G4SDParticleFilter::addIon ( G4int Z,
G4int A )

Definition at line 140 of file G4SDParticleFilter.cc.

141{
142 for(size_t i = 0; i < theIonZ.size(); i++)
143 {
144 if(theIonZ[i] == Z && theIonA[i] == A)
145 {
146 G4cout << "G4SDParticleFilter:: Ion has been already registered."
147 << G4endl;
148 return;
149 }
150 }
151 theIonZ.push_back(Z);
152 theIonA.push_back(A);
153}
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

◆ show()

void G4SDParticleFilter::show ( )

Definition at line 155 of file G4SDParticleFilter.cc.

156{
157 G4cout << "----G4SDParticleFileter particle list------" << G4endl;
158 for(auto & i : thePdef)
159 {
160 G4cout << i->GetParticleName() << G4endl;
161 }
162 for(size_t i = 0; i < theIonZ.size(); i++)
163 {
164 G4cout << " Ion PrtclDef (" << theIonZ[i] << "," << theIonA[i] << ")"
165 << G4endl;
166 }
167 G4cout << "-------------------------------------------" << G4endl;
168}

Referenced by G4SDParticleWithEnergyFilter::show().


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