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

#include <G4EmStandardPhysicsGS.hh>

Inheritance diagram for G4EmStandardPhysicsGS:

Public Member Functions

 G4EmStandardPhysicsGS (G4int ver=1, const G4String &name="")
 ~G4EmStandardPhysicsGS () override
void ConstructParticle () override
void ConstructProcess () override
Public Member Functions inherited from G4VPhysicsConstructor
 G4VPhysicsConstructor (const G4String &="")
 G4VPhysicsConstructor (const G4String &name, G4int physics_type)
virtual ~G4VPhysicsConstructor ()
void SetPhysicsName (const G4String &="")
const G4StringGetPhysicsName () const
void SetPhysicsType (G4int)
G4int GetPhysicsType () const
G4int GetInstanceID () const
virtual void TerminateWorker ()
void SetVerboseLevel (G4int value)
G4int GetVerboseLevel () const

Additional Inherited Members

Static Public Member Functions inherited from G4VPhysicsConstructor
static const G4VPCManagerGetSubInstanceManager ()
Protected Types inherited from G4VPhysicsConstructor
using PhysicsBuilder_V = G4VPCData::PhysicsBuilders_V
Protected Member Functions inherited from G4VPhysicsConstructor
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
PhysicsBuilder_V GetBuilders () const
void AddBuilder (G4PhysicsBuilderInterface *bld)
Protected Attributes inherited from G4VPhysicsConstructor
G4int verboseLevel = 0
G4String namePhysics = ""
G4int typePhysics = 0
G4ParticleTabletheParticleTable = nullptr
G4int g4vpcInstanceID = 0
Static Protected Attributes inherited from G4VPhysicsConstructor
static G4RUN_DLL G4VPCManager subInstanceManager

Detailed Description

Definition at line 51 of file G4EmStandardPhysicsGS.hh.

Constructor & Destructor Documentation

◆ G4EmStandardPhysicsGS()

G4EmStandardPhysicsGS::G4EmStandardPhysicsGS ( G4int ver = 1,
const G4String & name = "" )
explicit

Definition at line 108 of file G4EmStandardPhysicsGS.cc.

109 : G4VPhysicsConstructor("G4EmStandardGS")
110{
111 SetVerboseLevel(ver);
112 G4EmParameters* param = G4EmParameters::Instance();
113 param->SetDefaults();
114 param->SetVerbose(ver);
115 // use a denser discrete kinetic energy grid for more accurate interpolation
116 param->SetNumberOfBinsPerDecade(16);
117 // set the continuous step limit to: 0.2*Range that goes to Range below 10 um
118 param->SetStepFunction(0.2, 10*CLHEP::um);
119 // set the GS MSC model for e-/e+ to be used below 1.0 GeV with its (Mott,
120 // screening, scattering power) corrections activated and with the accurate
121 // stepping and boundary crossing algorithms (no other options since 11.4)
122 // with a skin of 3 elastic MFP near boundary
123 param->SetMscEnergyLimit(1.0*CLHEP::GeV);
124 param->SetUseMottCorrection(true);
126 param->SetMscSkin(3);
127 param->SetMscRangeFactor(0.08);
128 // activate fluoresence, i.e. emission of characteristic X-ray
129 param->SetFluo(true);
130 // set the energy loss fluctuation type
133}
@ bElectromagnetic
@ fUrbanFluctuation
@ fUseSafetyPlus
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetStepFunction(G4double v1, G4double v2)
void SetFluo(G4bool val)
void SetFluctuationType(G4EmFluctuationType val)
void SetVerbose(G4int val)
void SetMscSkin(G4double val)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetMscEnergyLimit(G4double val)
void SetUseMottCorrection(G4bool val)
void SetMscRangeFactor(G4double val)
G4VPhysicsConstructor(const G4String &="")
void SetVerboseLevel(G4int value)

◆ ~G4EmStandardPhysicsGS()

G4EmStandardPhysicsGS::~G4EmStandardPhysicsGS ( )
override

Definition at line 137 of file G4EmStandardPhysicsGS.cc.

138{}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysicsGS::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 142 of file G4EmStandardPhysicsGS.cc.

143{
144 // minimal set of particles for EM physics
146}
static void ConstructMinimalEmSet()

◆ ConstructProcess()

void G4EmStandardPhysicsGS::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 150 of file G4EmStandardPhysicsGS.cc.

151{
152 if(verboseLevel > 1) {
153 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
154 }
156 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
157
158 // processes used by several particles
159 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
160 G4NuclearStopping* pnuc(nullptr);
161
162 // high energy limit for e+- scattering models and bremsstrahlung
164
165 // Add gamma EM processes
166
167 // gamma
168 G4ParticleDefinition* particle = G4Gamma::Gamma();
169
170 G4PhotoElectricEffect* pe = new G4PhotoElectricEffect();
171 pe->SetEmModel(new G4LivermorePhotoElectricModel());
172
173 G4ComptonScattering* cs = new G4ComptonScattering;
174 G4GammaConversion* gc = new G4GammaConversion;
175 G4RayleighScattering* rs = new G4RayleighScattering;
176
177 if (G4EmParameters::Instance()->GeneralProcessActive()) {
178 G4GammaGeneralProcess* sp = new G4GammaGeneralProcess();
179 sp->AddEmProcess(pe);
180 sp->AddEmProcess(cs);
181 sp->AddEmProcess(gc);
182 sp->AddEmProcess(rs);
184 ph->RegisterProcess(sp, particle);
185 } else {
186 ph->RegisterProcess(pe, particle);
187 ph->RegisterProcess(cs, particle);
188 ph->RegisterProcess(gc, particle);
189 ph->RegisterProcess(rs, particle);
190 }
191
192 // e-
193 particle = G4Electron::Electron();
194
195 // msc: GS[:100 MeV] + WentzelVI[100 MeV:]
196 G4GoudsmitSaundersonMscModel* msc1 = new G4GoudsmitSaundersonMscModel();
197 G4WentzelVIModel* msc2 = new G4WentzelVIModel();
198 msc1->SetHighEnergyLimit(mscEnergyLimit);
199 msc2->SetLowEnergyLimit(mscEnergyLimit);
200 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
201 // (WVI is a mixed model, i.e. needs single scattering)
202 G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
203 G4CoulombScattering* ss = new G4CoulombScattering();
204 ss->SetEmModel(ssm);
205 ss->SetMinKinEnergy(mscEnergyLimit);
206 ssm->SetLowEnergyLimit(mscEnergyLimit);
207 ssm->SetActivationLowEnergyLimit(mscEnergyLimit);
208
209 // ionisation: Penelope[:1.0 GeV] + Moller[1.0 GeV:]
210 G4eIonisation* eioni = new G4eIonisation();
212 G4VEmModel* theIoniMod = new G4PenelopeIonisationModel();
213 theIoniMod->SetHighEnergyLimit(1.0*CLHEP::GeV);
214 eioni->AddEmModel(0, theIoniMod);
215
216 // bremsstrahlung: Seltzer-Berger[:1.0 GeV] + extended Bethe–Heitler[1.0 GeV:]
217 G4eBremsstrahlung* brem = new G4eBremsstrahlung();
218 G4SeltzerBergerModel* br1 = new G4SeltzerBergerModel();
219 G4eBremsstrahlungRelModel* br2 = new G4eBremsstrahlungRelModel();
220 br1->SetAngularDistribution(new G4Generator2BS());
221 br2->SetAngularDistribution(new G4Generator2BS());
222 brem->SetEmModel(br1);
223 brem->SetEmModel(br2);
224 br1->SetHighEnergyLimit(1.0*CLHEP::GeV);
225
226 ph->RegisterProcess(eioni, particle);
227 ph->RegisterProcess(brem, particle);
228 ph->RegisterProcess(ss, particle);
229
230 // e+
231 particle = G4Positron::Positron();
232
233 // msc: GS[:100 MeV] + WentzelVI[100 MeV:]
234 msc1 = new G4GoudsmitSaundersonMscModel();
235 msc2 = new G4WentzelVIModel();
236 msc1->SetHighEnergyLimit(mscEnergyLimit);
237 msc2->SetLowEnergyLimit(mscEnergyLimit);
238 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
239 // (WVI is a mixed model, i.e. needs single scattering)
240 ssm = new G4eCoulombScatteringModel();
241 ss = new G4CoulombScattering();
242 ss->SetEmModel(ssm);
243 ss->SetMinKinEnergy(mscEnergyLimit);
244 ssm->SetLowEnergyLimit(mscEnergyLimit);
245 ssm->SetActivationLowEnergyLimit(mscEnergyLimit);
246
247 // ionisation: Penelope[:1.0 GeV] + Bhabha[1.0 GeV:]
248 eioni = new G4eIonisation();
250 G4VEmModel* pen = new G4PenelopeIonisationModel();
251 pen->SetHighEnergyLimit(1.0*CLHEP::GeV);
252 eioni->AddEmModel(0, pen);
253
254 // bremsstrahlung: Seltzer-Berger[:1.0 GeV] + extended Bethe–Heitler[1.0 GeV:]
255 brem = new G4eBremsstrahlung();
256 br1 = new G4SeltzerBergerModel();
257 br2 = new G4eBremsstrahlungRelModel();
258 br1->SetAngularDistribution(new G4Generator2BS());
259 br2->SetAngularDistribution(new G4Generator2BS());
260 brem->SetEmModel(br1);
261 brem->SetEmModel(br2);
262 br1->SetHighEnergyLimit(1.0*CLHEP::GeV);
263
264 ph->RegisterProcess(eioni, particle);
265 ph->RegisterProcess(brem, particle);
266 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
267 ph->RegisterProcess(ss, particle);
268
269 // generic ion
270 particle = G4GenericIon::GenericIon();
271 G4ionIonisation* ionIoni = new G4ionIonisation();
272 ph->RegisterProcess(hmsc, particle);
273 ph->RegisterProcess(ionIoni, particle);
274
275 // muons, hadrons, ions
277
278 // extra configuration
279 G4EmModelActivator mact(GetPhysicsName());
280}
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition G4Electron.cc:91
static void ConstructCharged(G4hMultipleScattering *hmsc, G4NuclearStopping *nucStopping, G4bool isWVI=true)
static void ConstructElectronMscProcess(G4VMscModel *msc1, G4VMscModel *msc2, G4ParticleDefinition *particle)
static void PrepareEMPhysics()
G4double MscEnergyLimit() const
static G4VEmFluctuationModel * ModelOfFluctuations(G4bool isIon=false)
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4GenericIon * GenericIon()
static G4LossTableManager * Instance()
void SetGammaGeneralProcess(G4VEmProcess *)
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition G4Positron.cc:90
void SetHighEnergyLimit(G4double)
void SetActivationLowEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
void SetAngularDistribution(G4VEmAngularDistribution *)
void SetMinKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=0)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=0)
const G4String & GetPhysicsName() const

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