Geant4 10.7.0
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)
 
virtual ~G4SDParticleFilter ()
 
virtual G4bool Accept (const G4Step *) const
 
void add (const G4String &particleName)
 
void addIon (G4int Z, G4int A)
 
void show ()
 
- Public Member Functions inherited from G4VSDFilter
 G4VSDFilter (G4String name)
 
virtual ~G4VSDFilter ()
 
virtual G4bool Accept (const G4Step *) const =0
 
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{
48 thePdef.clear();
49 theIonZ.clear();
50 theIonA.clear();
51}

◆ G4SDParticleFilter() [2/4]

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

Definition at line 53 of file G4SDParticleFilter.cc.

55 :G4VSDFilter(name)
56{
57 thePdef.clear();
59 if(!pd)
60 {
61 G4String msg = "Particle <";
62 msg += particleName;
63 msg += "> not found.";
64 G4Exception("G4SDParticleFilter::G4SDParticleFilter",
65 "DetPS0101",FatalException,msg);
66 }
67 thePdef.push_back(pd);
68 theIonZ.clear();
69 theIonA.clear();
70}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()

◆ G4SDParticleFilter() [3/4]

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

Definition at line 72 of file G4SDParticleFilter.cc.

74 :G4VSDFilter(name)
75{
76 thePdef.clear();
77 for ( size_t i = 0; i < particleNames.size(); i++){
79 if(!pd)
80 {
81 G4String msg = "Particle <";
82 msg += particleNames[i];
83 msg += "> not found.";
84 G4Exception("G4SDParticleFilter::G4SDParticleFilter",
85 "DetPS0102",FatalException,msg);
86 }
87 thePdef.push_back(pd);
88 theIonZ.clear();
89 theIonA.clear();
90 }
91}

◆ G4SDParticleFilter() [4/4]

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

Definition at line 93 of file G4SDParticleFilter.cc.

95 :G4VSDFilter(name), thePdef(particleDef)
96{
97 for ( size_t i = 0; i < particleDef.size(); i++){
98 if(!particleDef[i]) G4Exception("G4SDParticleFilter::G4SDParticleFilter",
99 "DetPS0103",FatalException,
100 "NULL pointer is found in the given particleDef vector.");
101 }
102 theIonZ.clear();
103 theIonA.clear();
104}

◆ ~G4SDParticleFilter()

G4SDParticleFilter::~G4SDParticleFilter ( )
virtual

Definition at line 106 of file G4SDParticleFilter.cc.

107{
108 thePdef.clear();
109 theIonZ.clear();
110 theIonA.clear();
111 }

Member Function Documentation

◆ Accept()

G4bool G4SDParticleFilter::Accept ( const G4Step aStep) const
virtual

Implements G4VSDFilter.

Definition at line 113 of file G4SDParticleFilter.cc.

114{
115
116 for ( size_t i = 0; i < thePdef.size(); i++){
117 if ( thePdef[i] == aStep->GetTrack()->GetDefinition() ) return TRUE;
118 }
119
120 // Ions by Z,A
121 for ( size_t i = 0; i < theIonZ.size(); i++){
122 if ( theIonZ[i] == aStep->GetTrack()->GetDefinition()->GetAtomicNumber()
123 && theIonA[i] == aStep->GetTrack()->GetDefinition()->GetAtomicMass() ){
124 return TRUE;
125 }
126 }
127
128 return FALSE;
129}
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
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 131 of file G4SDParticleFilter.cc.

132{
135 if(!pd)
136 {
137 G4String msg = "Particle <";
138 msg += particleName;
139 msg += "> not found.";
140 G4Exception("G4SDParticleFilter::add()",
141 "DetPS0104",FatalException,msg);
142 }
143 for ( size_t i = 0; i < thePdef.size(); i++){
144 if ( thePdef[i] == pd ) return;
145 }
146 thePdef.push_back(pd);
147}

Referenced by G4SDParticleWithEnergyFilter::add().

◆ addIon()

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

Definition at line 149 of file G4SDParticleFilter.cc.

149 {
150 for ( size_t i = 0; i < theIonZ.size(); i++){
151 if ( theIonZ[i] == Z && theIonA[i] == A ){
152 G4cout << "G4SDParticleFilter:: Ion has been already registered."<<G4endl;
153 return;
154 }
155 }
156 theIonZ.push_back(Z);
157 theIonA.push_back(A);
158}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ show()

void G4SDParticleFilter::show ( )

Definition at line 160 of file G4SDParticleFilter.cc.

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

Referenced by G4SDParticleWithEnergyFilter::show().


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