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

#include <G4EmBuilder.hh>

Static Public Member Functions

static void ConstructBasicEmPhysics (G4hMultipleScattering *hmsc, const std::vector< G4int > &listHadrons)
 
static void ConstructLightHadrons (G4ParticleDefinition *part1, G4ParticleDefinition *part2, G4bool isHEP, G4bool isProton, G4bool isWVI)
 
static void ConstructIonEmPhysics (G4hMultipleScattering *hmsc, G4NuclearStopping *nucStopping)
 
static void ConstructCharged (G4hMultipleScattering *hmsc, G4NuclearStopping *nucStopping, G4bool isWVI=true)
 
static void ConstructMinimalEmSet ()
 
static void PrepareEMPhysics ()
 

Detailed Description

Definition at line 44 of file G4EmBuilder.hh.

Member Function Documentation

◆ ConstructBasicEmPhysics()

void G4EmBuilder::ConstructBasicEmPhysics ( G4hMultipleScattering hmsc,
const std::vector< G4int > &  listHadrons 
)
static

Definition at line 90 of file G4EmBuilder.cc.

92{
95
96 for( auto & pdg : partList ) {
97 auto part = table->FindParticle( pdg );
98 if ( part == nullptr || part->GetPDGCharge() == 0.0 ) { continue; }
99 ph->RegisterProcess(hmsc, part);
100 ph->RegisterProcess(new G4hIonisation(), part);
101 }
102}
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()

Referenced by ConstructCharged().

◆ ConstructCharged()

void G4EmBuilder::ConstructCharged ( G4hMultipleScattering hmsc,
G4NuclearStopping nucStopping,
G4bool  isWVI = true 
)
static

Definition at line 170 of file G4EmBuilder.cc.

173{
177 G4bool isHEP = ( param->MaxKinEnergy() > hpar->EnergyThresholdForHeavyHadrons() );
178
179 // muon bremsstrahlung and pair production
180 G4MuBremsstrahlung* mub = ( isHEP ) ? new G4MuBremsstrahlung() : nullptr;
181 G4MuPairProduction* mup = ( isHEP ) ? new G4MuPairProduction() : nullptr;
182
183 // muon multiple and single scattering
185 if(isWVI) { mumsc->SetEmModel(new G4WentzelVIModel()); }
186 G4CoulombScattering* muss = ( isWVI ) ? new G4CoulombScattering() : nullptr;
187
188 // Add standard EM Processes
189 // mu+-
191 ph->RegisterProcess(mumsc, part);
192 ph->RegisterProcess(new G4MuIonisation(), part);
193 if( isHEP ) {
194 ph->RegisterProcess(mub, part);
195 ph->RegisterProcess(mup, part);
196 }
197 if( isWVI ) { ph->RegisterProcess(muss, part); }
198
199 part = G4MuonMinus::MuonMinus();
200 ph->RegisterProcess(mumsc, part);
201 ph->RegisterProcess(new G4MuIonisation(), part);
202 if( isHEP ) {
203 ph->RegisterProcess(mub, part);
204 ph->RegisterProcess(mup, part);
205 }
206 if( isWVI ) { ph->RegisterProcess(muss, part); }
207
208 // pi+-
210
211 // K+-
213
214 // p, pbar
216 if( nucStopping != nullptr ) {
217 ph->RegisterProcess(nucStopping, G4Proton::Proton());
218 }
219
220 // ions
221 ConstructIonEmPhysics(hmsc, nucStopping);
222
223 // hyperons and anti particles
224 if( isHEP ) {
226
227 // b- and c- charged particles
228 if( hpar->EnableBCParticles() ) {
230 }
231 }
232}
bool G4bool
Definition: G4Types.hh:86
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:92
static void ConstructLightHadrons(G4ParticleDefinition *part1, G4ParticleDefinition *part2, G4bool isHEP, G4bool isProton, G4bool isWVI)
Definition: G4EmBuilder.cc:132
static void ConstructIonEmPhysics(G4hMultipleScattering *hmsc, G4NuclearStopping *nucStopping)
Definition: G4EmBuilder.cc:104
static void ConstructBasicEmPhysics(G4hMultipleScattering *hmsc, const std::vector< G4int > &listHadrons)
Definition: G4EmBuilder.cc:90
static G4EmParameters * Instance()
G4double MaxKinEnergy() const
static const std::vector< G4int > & GetBCChargedHadrons()
static const std::vector< G4int > & GetHeavyChargedParticles()
static G4HadronicParameters * Instance()
G4bool EnableBCParticles() const
G4double EnergyThresholdForHeavyHadrons() const
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:112
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:112
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:98
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:97
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:97
static G4Proton * Proton()
Definition: G4Proton.cc:92
void SetEmModel(G4VMscModel *, size_t index=0)

Referenced by G4EmLivermorePhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysicsGS::ConstructProcess(), and G4EmStandardPhysicsWVI::ConstructProcess().

◆ ConstructIonEmPhysics()

void G4EmBuilder::ConstructIonEmPhysics ( G4hMultipleScattering hmsc,
G4NuclearStopping nucStopping 
)
static

Definition at line 104 of file G4EmBuilder.cc.

106{
108
110 ph->RegisterProcess(hmsc, part);
111 ph->RegisterProcess(new G4hIonisation(), part);
112
113 part = G4Triton::Triton();
114 ph->RegisterProcess(hmsc, part);
115 ph->RegisterProcess(new G4hIonisation(), part);
116
117 part = G4He3::He3();
118 ph->RegisterProcess(new G4hMultipleScattering(), part);
119 ph->RegisterProcess(new G4ionIonisation(), part);
120 if( nucStopping != nullptr ) {
121 ph->RegisterProcess(nucStopping, part);
122 }
123
124 part = G4Alpha::Alpha();
125 ph->RegisterProcess(new G4hMultipleScattering(), part);
126 ph->RegisterProcess(new G4ionIonisation(), part);
127 if( nucStopping != nullptr ) {
128 ph->RegisterProcess(nucStopping, part);
129 }
130}
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
static G4He3 * He3()
Definition: G4He3.cc:93
static G4Triton * Triton()
Definition: G4Triton.cc:94

Referenced by ConstructCharged().

◆ ConstructLightHadrons()

void G4EmBuilder::ConstructLightHadrons ( G4ParticleDefinition part1,
G4ParticleDefinition part2,
G4bool  isHEP,
G4bool  isProton,
G4bool  isWVI 
)
static

Definition at line 132 of file G4EmBuilder.cc.

136{
138
139 G4hBremsstrahlung* brem = ( isHEP ) ? new G4hBremsstrahlung() : nullptr;
140 G4hPairProduction* pair = ( isHEP ) ? new G4hPairProduction() : nullptr;
141
143 if(isWVI) { msc->SetEmModel(new G4WentzelVIModel()); }
144 G4CoulombScattering* ss = ( isWVI ) ? new G4CoulombScattering() : nullptr;
145
146 ph->RegisterProcess(msc, part1);
147 ph->RegisterProcess(new G4hIonisation(), part1);
148 if( isHEP ) {
149 ph->RegisterProcess(brem, part1);
150 ph->RegisterProcess(pair, part1);
151 }
152 if( isWVI ) { ph->RegisterProcess(ss, part1); }
153
154 if( isProton ) {
155 msc = new G4hMultipleScattering();
156 if(isWVI) {
157 msc->SetEmModel(new G4WentzelVIModel());
158 ss = new G4CoulombScattering();
159 }
160 }
161 ph->RegisterProcess(msc, part2);
162 ph->RegisterProcess(new G4hIonisation(), part2);
163 if( isHEP ) {
164 ph->RegisterProcess(brem, part2);
165 ph->RegisterProcess(pair, part2);
166 }
167 if( isWVI ) { ph->RegisterProcess(ss, part2); }
168}

Referenced by ConstructCharged().

◆ ConstructMinimalEmSet()

◆ PrepareEMPhysics()


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