Geant4 9.6.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 53 of file G4SDParticleFilter.hh.

Constructor & Destructor Documentation

◆ G4SDParticleFilter() [1/4]

G4SDParticleFilter::G4SDParticleFilter ( G4String  name)

Definition at line 46 of file G4SDParticleFilter.cc.

47 :G4VSDFilter(name)
48{
49 thePdef.clear();
50 theIonZ.clear();
51 theIonA.clear();
52}

◆ G4SDParticleFilter() [2/4]

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

Definition at line 54 of file G4SDParticleFilter.cc.

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

◆ G4SDParticleFilter() [3/4]

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

Definition at line 73 of file G4SDParticleFilter.cc.

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

◆ G4SDParticleFilter() [4/4]

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

Definition at line 94 of file G4SDParticleFilter.cc.

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

◆ ~G4SDParticleFilter()

G4SDParticleFilter::~G4SDParticleFilter ( )
virtual

Definition at line 107 of file G4SDParticleFilter.cc.

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

Member Function Documentation

◆ Accept()

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

Implements G4VSDFilter.

Definition at line 114 of file G4SDParticleFilter.cc.

115{
116
117 for ( size_t i = 0; i < thePdef.size(); i++){
118 if ( thePdef[i] == aStep->GetTrack()->GetDefinition() ) return TRUE;
119 }
120
121 // Ions by Z,A
122 for ( size_t i = 0; i < theIonZ.size(); i++){
123 if ( theIonZ[i] == aStep->GetTrack()->GetDefinition()->GetAtomicNumber()
124 && theIonA[i] == aStep->GetTrack()->GetDefinition()->GetAtomicMass() ){
125 return TRUE;
126 }
127 }
128
129 return FALSE;
130}
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
G4Track * GetTrack() const
G4ParticleDefinition * GetDefinition() const
#define TRUE
Definition: globals.hh:55
#define FALSE
Definition: globals.hh:52

Referenced by G4SDParticleWithEnergyFilter::Accept().

◆ add()

void G4SDParticleFilter::add ( const G4String particleName)

Definition at line 132 of file G4SDParticleFilter.cc.

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

Referenced by G4SDParticleWithEnergyFilter::add().

◆ addIon()

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

Definition at line 150 of file G4SDParticleFilter.cc.

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

◆ show()

void G4SDParticleFilter::show ( )

Definition at line 161 of file G4SDParticleFilter.cc.

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

Referenced by G4SDParticleWithEnergyFilter::show().


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