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

#include <G4EmStandardPhysics_option4.hh>

+ Inheritance diagram for G4EmStandardPhysics_option4:

Public Member Functions

 G4EmStandardPhysics_option4 (G4int ver=1, const G4String &name="")
 
 ~G4EmStandardPhysics_option4 () 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 ()
 
virtual void ConstructParticle ()=0
 
virtual void ConstructProcess ()=0
 
void SetPhysicsName (const G4String &="")
 
const G4StringGetPhysicsName () const
 
void SetPhysicsType (G4int)
 
G4int GetPhysicsType () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
G4int GetInstanceID () const
 
virtual void TerminateWorker ()
 

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
 
G4String namePhysics
 
G4int typePhysics
 
G4ParticleTabletheParticleTable
 
G4int g4vpcInstanceID
 
- Static Protected Attributes inherited from G4VPhysicsConstructor
static G4RUN_DLL G4VPCManager subInstanceManager
 

Detailed Description

Definition at line 52 of file G4EmStandardPhysics_option4.hh.

Constructor & Destructor Documentation

◆ G4EmStandardPhysics_option4()

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

Definition at line 104 of file G4EmStandardPhysics_option4.cc.

106 : G4VPhysicsConstructor("G4EmStandard_opt4"), verbose(ver)
107{
109 param->SetDefaults();
110 param->SetVerbose(verbose);
111 param->SetMinEnergy(100*CLHEP::eV);
112 param->SetLowestElectronEnergy(100*CLHEP::eV);
113 param->SetNumberOfBinsPerDecade(20);
115 param->SetStepFunction(0.2, 10*CLHEP::um);
116 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
117 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
118 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
119 param->SetUseMottCorrection(true); // use Mott-correction for e-/e+ msc gs
120 param->SetMscStepLimitType(fUseSafetyPlus); // for e-/e+ msc gs
121 param->SetMscSkin(3); // error-free stepping for e-/e+ msc gs
122 param->SetMscRangeFactor(0.08); // error-free stepping for e-/e+ msc gs
123 param->SetMuHadLateralDisplacement(true);
124 param->SetFluo(true);
125 param->SetUseICRU90Data(true);
126 param->SetMaxNIELEnergy(1*CLHEP::MeV);
128}
@ bElectromagnetic
@ fUseSafetyPlus
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetMuHadLateralDisplacement(G4bool val)
void ActivateAngularGeneratorForIonisation(G4bool val)
void SetStepFunction(G4double v1, G4double v2)
void SetFluo(G4bool val)
void SetStepFunctionMuHad(G4double v1, G4double v2)
void SetVerbose(G4int val)
void SetMscSkin(G4double val)
void SetMaxNIELEnergy(G4double val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetUseICRU90Data(G4bool val)
void SetUseMottCorrection(G4bool val)
void SetMscRangeFactor(G4double val)

◆ ~G4EmStandardPhysics_option4()

G4EmStandardPhysics_option4::~G4EmStandardPhysics_option4 ( )
override

Definition at line 132 of file G4EmStandardPhysics_option4.cc.

133{}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysics_option4::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 137 of file G4EmStandardPhysics_option4.cc.

138{
139 // minimal set of particles for EM physics
141}
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:234

◆ ConstructProcess()

void G4EmStandardPhysics_option4::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 145 of file G4EmStandardPhysics_option4.cc.

146{
147 if(verbose > 1) {
148 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
149 }
151
153
154 // processes used by several particles
156 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
157
158 // nuclear stopping
159 G4double nielEnergyLimit = G4EmParameters::Instance()->MaxNIELEnergy();
161 pnuc->SetMaxKinEnergy(nielEnergyLimit);
162
163 // high energy limit for e+- scattering models and bremsstrahlung
164 G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
165
166 // Add gamma EM Processes
168
169 // Photoelectric
171 G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
172 pe->SetEmModel(theLivermorePEModel);
173
174 // Compton scattering
177 G4VEmModel* theLowEPComptonModel = new G4LowEPComptonModel();
178 theLowEPComptonModel->SetHighEnergyLimit(20*CLHEP::MeV);
179 cs->AddEmModel(0, theLowEPComptonModel);
180
181 // Gamma conversion
183 G4VEmModel* conv = new G4BetheHeitler5DModel();
184 gc->SetEmModel(conv);
185
186 if(G4EmParameters::Instance()->GeneralProcessActive()) {
188 sp->AddEmProcess(pe);
189 sp->AddEmProcess(cs);
190 sp->AddEmProcess(gc);
191 sp->AddEmProcess(new G4RayleighScattering());
193 ph->RegisterProcess(sp, particle);
194 } else {
195 ph->RegisterProcess(pe, particle);
196 ph->RegisterProcess(cs, particle);
197 ph->RegisterProcess(gc, particle);
198 ph->RegisterProcess(new G4RayleighScattering(), particle);
199 }
200
201 // e-
202 particle = G4Electron::Electron();
203
204 // multiple scattering
206 // e-/e+ msc gs with Mott-correction
207 // (Mott-correction is set through G4EmParameters)
210 msc1->SetHighEnergyLimit(highEnergyLimit);
211 msc2->SetLowEnergyLimit(highEnergyLimit);
212 msc->SetEmModel(msc1);
213 msc->SetEmModel(msc2);
214
217 ss->SetEmModel(ssm);
218 ss->SetMinKinEnergy(highEnergyLimit);
219 ssm->SetLowEnergyLimit(highEnergyLimit);
220 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
221
222 // ionisation
223 G4eIonisation* eioni = new G4eIonisation();
224 G4VEmModel* theIoniLiv = new G4LivermoreIonisationModel();
225 theIoniLiv->SetHighEnergyLimit(0.1*CLHEP::MeV);
226 eioni->AddEmModel(0, theIoniLiv, new G4UniversalFluctuation() );
227
228 // bremsstrahlung
234 brem->SetEmModel(br1);
235 brem->SetEmModel(br2);
236 br1->SetHighEnergyLimit(CLHEP::GeV);
237
238 // register processes
239 ph->RegisterProcess(msc, particle);
240 ph->RegisterProcess(eioni, particle);
241 ph->RegisterProcess(brem, particle);
242 ph->RegisterProcess(ee, particle);
243 ph->RegisterProcess(ss, particle);
244
245 // e+
246 particle = G4Positron::Positron();
247
248 // multiple scattering
249 msc = new G4eMultipleScattering();
250 // e-/e+ msc gs with Mott-correction
251 // (Mott-correction is set through G4EmParameters)
252 msc1 = new G4GoudsmitSaundersonMscModel();
253 msc2 = new G4WentzelVIModel();
254 msc1->SetHighEnergyLimit(highEnergyLimit);
255 msc2->SetLowEnergyLimit(highEnergyLimit);
256 msc->SetEmModel(msc1);
257 msc->SetEmModel(msc2);
258
259 ssm = new G4eCoulombScatteringModel();
260 ss = new G4CoulombScattering();
261 ss->SetEmModel(ssm);
262 ss->SetMinKinEnergy(highEnergyLimit);
263 ssm->SetLowEnergyLimit(highEnergyLimit);
264 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
265
266 // ionisation
267 eioni = new G4eIonisation();
269 pen->SetHighEnergyLimit(0.1*CLHEP::MeV);
270 eioni->AddEmModel(0, pen, new G4UniversalFluctuation());
271
272 // bremsstrahlung
273 brem = new G4eBremsstrahlung();
274 br1 = new G4SeltzerBergerModel();
275 br2 = new G4eBremsstrahlungRelModel();
278 brem->SetEmModel(br1);
279 brem->SetEmModel(br2);
280 br1->SetHighEnergyLimit(CLHEP::GeV);
281
282 // register processes
283 ph->RegisterProcess(msc, particle);
284 ph->RegisterProcess(eioni, particle);
285 ph->RegisterProcess(brem, particle);
286 ph->RegisterProcess(ee, particle);
287 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
288 ph->RegisterProcess(ss, particle);
289
290 // generic ion
291 particle = G4GenericIon::GenericIon();
292 G4ionIonisation* ionIoni = new G4ionIonisation();
294 ph->RegisterProcess(hmsc, particle);
295 ph->RegisterProcess(ionIoni, particle);
296 ph->RegisterProcess(pnuc, particle);
297
298 // muons, hadrons, ions
300
301 // extra configuration
303}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition: G4Electron.cc:93
static void ConstructCharged(G4hMultipleScattering *hmsc, G4NuclearStopping *nucStopping, G4bool isWVI=true)
Definition: G4EmBuilder.cc:170
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:259
G4double MaxNIELEnergy() const
G4double MscEnergyLimit() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
static G4LossTableManager * Instance()
void SetGammaGeneralProcess(G4VEmProcess *)
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition: G4Positron.cc:93
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:778
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
void SetMinKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=0)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
void SetEmModel(G4VMscModel *, size_t index=0)
const G4String & GetPhysicsName() const

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