Geant4 9.6.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)
 
 G4EmStandardPhysics_option4 (G4int ver, const G4String &name)
 
virtual ~G4EmStandardPhysics_option4 ()
 
virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
- 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
 

Additional Inherited Members

- Protected Member Functions inherited from G4VPhysicsConstructor
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 
- Protected Attributes inherited from G4VPhysicsConstructor
G4int verboseLevel
 
G4String namePhysics
 
G4int typePhysics
 
G4ParticleTabletheParticleTable
 
G4ParticleTable::G4PTblDicIteratortheParticleIterator
 
G4PhysicsListHelperthePLHelper
 

Detailed Description

Definition at line 53 of file G4EmStandardPhysics_option4.hh.

Constructor & Destructor Documentation

◆ G4EmStandardPhysics_option4() [1/2]

G4EmStandardPhysics_option4::G4EmStandardPhysics_option4 ( G4int  ver = 1)

Definition at line 122 of file G4EmStandardPhysics_option4.cc.

123 : G4VPhysicsConstructor("G4EmStandard_opt4"), verbose(ver)
124{
127}
@ bElectromagnetic
static G4LossTableManager * Instance()

◆ G4EmStandardPhysics_option4() [2/2]

G4EmStandardPhysics_option4::G4EmStandardPhysics_option4 ( G4int  ver,
const G4String name 
)

Definition at line 131 of file G4EmStandardPhysics_option4.cc.

132 : G4VPhysicsConstructor("G4EmStandard_opt3"), verbose(ver)
133{
136}

◆ ~G4EmStandardPhysics_option4()

G4EmStandardPhysics_option4::~G4EmStandardPhysics_option4 ( )
virtual

Definition at line 140 of file G4EmStandardPhysics_option4.cc.

141{}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysics_option4::ConstructParticle ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 145 of file G4EmStandardPhysics_option4.cc.

146{
147// gamma
149
150// leptons
155
156// mesons
161
162// barions
165
166// ions
169 G4He3::He3();
172}
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:87
static G4He3 * He3()
Definition: G4He3.cc:94
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Triton * Triton()
Definition: G4Triton.cc:95

◆ ConstructProcess()

void G4EmStandardPhysics_option4::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 176 of file G4EmStandardPhysics_option4.cc.

177{
179
180 // muon & hadron bremsstrahlung and pair production
189
190 // muon & hadron multiple scattering
192 mumsc->AddEmModel(0, new G4WentzelVIModel());
194 pimsc->AddEmModel(0, new G4WentzelVIModel());
196 kmsc->AddEmModel(0, new G4WentzelVIModel());
198 pmsc->AddEmModel(0, new G4WentzelVIModel());
199 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
200
201 // high energy limit for e+- scattering models
202 G4double highEnergyLimit = 100*MeV;
203
204 // nuclear stopping
205 G4NuclearStopping* ionnuc = new G4NuclearStopping();
207
208 // Add standard EM Processes
210 while( (*theParticleIterator)() ){
212 G4String particleName = particle->GetParticleName();
213 if(verbose > 1)
214 G4cout << "### " << GetPhysicsName() << " instantiates for "
215 << particleName << G4endl;
216
217 if (particleName == "gamma") {
218
219 // Compton scattering
221 cs->SetEmModel(new G4KleinNishinaModel(),1);
222 G4VEmModel* theLowEPComptonModel = new G4LowEPComptonModel();
223 theLowEPComptonModel->SetHighEnergyLimit(20*MeV);
224 cs->AddEmModel(0, theLowEPComptonModel);
225 ph->RegisterProcess(cs, particle);
226
227 // Photoelectric
229 G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
230 theLivermorePEModel->SetHighEnergyLimit(10*GeV);
231 pe->SetEmModel(theLivermorePEModel,1);
232 ph->RegisterProcess(pe, particle);
233
234 // Gamma conversion
236 G4VEmModel* thePenelopeGCModel = new G4PenelopeGammaConversionModel();
237 thePenelopeGCModel->SetHighEnergyLimit(1*GeV);
238 gc->SetEmModel(thePenelopeGCModel,1);
239 ph->RegisterProcess(gc, particle);
240
241 // Rayleigh scattering
242 ph->RegisterProcess(new G4RayleighScattering(), particle);
243
244 } else if (particleName == "e-") {
245
246 // multiple scattering
251 msc1->SetHighEnergyLimit(highEnergyLimit);
252 msc2->SetLowEnergyLimit(highEnergyLimit);
253 msc->AddEmModel(0, msc1);
254 msc->AddEmModel(0, msc2);
255
258 ss->SetEmModel(ssm, 1);
259 ss->SetMinKinEnergy(highEnergyLimit);
260 ssm->SetLowEnergyLimit(highEnergyLimit);
261 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
262
263 // ionisation
264 G4eIonisation* eIoni = new G4eIonisation();
265 eIoni->SetStepFunction(0.2, 100*um);
266 G4VEmModel* theIoniPenelope = new G4PenelopeIonisationModel();
267 theIoniPenelope->SetHighEnergyLimit(0.1*MeV);
268 eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
269
270 // bremsstrahlung
272
273 ph->RegisterProcess(msc, particle);
274 ph->RegisterProcess(eIoni, particle);
275 ph->RegisterProcess(eBrem, particle);
276 ph->RegisterProcess(ss, particle);
277
278 } else if (particleName == "e+") {
279
280 // multiple scattering
285 msc1->SetHighEnergyLimit(highEnergyLimit);
286 msc2->SetLowEnergyLimit(highEnergyLimit);
287 msc->AddEmModel(0, msc1);
288 msc->AddEmModel(0, msc2);
289
292 ss->SetEmModel(ssm, 1);
293 ss->SetMinKinEnergy(highEnergyLimit);
294 ssm->SetLowEnergyLimit(highEnergyLimit);
295 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
296
297 // ionisation
298 G4eIonisation* eIoni = new G4eIonisation();
299 eIoni->SetStepFunction(0.2, 100*um);
300 G4VEmModel* theIoniPenelope = new G4PenelopeIonisationModel();
301 theIoniPenelope->SetHighEnergyLimit(0.1*MeV);
302 eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
303
304 // bremsstrahlung
306
307 ph->RegisterProcess(msc, particle);
308 ph->RegisterProcess(eIoni, particle);
309 ph->RegisterProcess(eBrem, particle);
310 ph->RegisterProcess(ss, particle);
311
312 // annihilation at rest and in flight
313 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
314
315 } else if (particleName == "mu+" ||
316 particleName == "mu-" ) {
317
318 G4MuIonisation* muIoni = new G4MuIonisation();
319 muIoni->SetStepFunction(0.2, 50*um);
320
321 ph->RegisterProcess(mumsc, particle);
322 ph->RegisterProcess(muIoni, particle);
323 ph->RegisterProcess(mub, particle);
324 ph->RegisterProcess(mup, particle);
325 ph->RegisterProcess(new G4CoulombScattering(), particle);
326
327 } else if (particleName == "alpha" ||
328 particleName == "He3") {
329
331 G4ionIonisation* ionIoni = new G4ionIonisation();
332 ionIoni->SetStepFunction(0.1, 10*um);
333
334 ph->RegisterProcess(msc, particle);
335 ph->RegisterProcess(ionIoni, particle);
336 ph->RegisterProcess(ionnuc, particle);
337
338 } else if (particleName == "GenericIon") {
339
340 G4ionIonisation* ionIoni = new G4ionIonisation();
342 ionIoni->SetStepFunction(0.1, 1*um);
343
344 ph->RegisterProcess(hmsc, particle);
345 ph->RegisterProcess(ionIoni, particle);
346 ph->RegisterProcess(ionnuc, particle);
347
348 } else if (particleName == "pi+" ||
349 particleName == "pi-" ) {
350
351 //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
352 G4hIonisation* hIoni = new G4hIonisation();
353 hIoni->SetStepFunction(0.2, 50*um);
354
355 ph->RegisterProcess(pimsc, particle);
356 ph->RegisterProcess(hIoni, particle);
357 ph->RegisterProcess(pib, particle);
358 ph->RegisterProcess(pip, particle);
359
360 } else if (particleName == "kaon+" ||
361 particleName == "kaon-" ) {
362
363 //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
364 G4hIonisation* hIoni = new G4hIonisation();
365 hIoni->SetStepFunction(0.2, 50*um);
366
367 ph->RegisterProcess(kmsc, particle);
368 ph->RegisterProcess(hIoni, particle);
369 ph->RegisterProcess(kb, particle);
370 ph->RegisterProcess(kp, particle);
371
372 } else if (particleName == "proton" ||
373 particleName == "anti_proton") {
374
375 //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
376 G4hIonisation* hIoni = new G4hIonisation();
377 hIoni->SetStepFunction(0.2, 50*um);
378
379 ph->RegisterProcess(pmsc, particle);
380 ph->RegisterProcess(hIoni, particle);
381 ph->RegisterProcess(pb, particle);
382 ph->RegisterProcess(pp, particle);
383 ph->RegisterProcess(pnuc, particle);
384
385 } else if (particleName == "B+" ||
386 particleName == "B-" ||
387 particleName == "D+" ||
388 particleName == "D-" ||
389 particleName == "Ds+" ||
390 particleName == "Ds-" ||
391 particleName == "anti_He3" ||
392 particleName == "anti_alpha" ||
393 particleName == "anti_deuteron" ||
394 particleName == "anti_lambda_c+" ||
395 particleName == "anti_omega-" ||
396 particleName == "anti_sigma_c+" ||
397 particleName == "anti_sigma_c++" ||
398 particleName == "anti_sigma+" ||
399 particleName == "anti_sigma-" ||
400 particleName == "anti_triton" ||
401 particleName == "anti_xi_c+" ||
402 particleName == "anti_xi-" ||
403 particleName == "deuteron" ||
404 particleName == "lambda_c+" ||
405 particleName == "omega-" ||
406 particleName == "sigma_c+" ||
407 particleName == "sigma_c++" ||
408 particleName == "sigma+" ||
409 particleName == "sigma-" ||
410 particleName == "tau+" ||
411 particleName == "tau-" ||
412 particleName == "triton" ||
413 particleName == "xi_c+" ||
414 particleName == "xi-" ) {
415
416 ph->RegisterProcess(hmsc, particle);
417 ph->RegisterProcess(new G4hIonisation(), particle);
418 ph->RegisterProcess(pnuc, particle);
419 }
420 }
421
422 // Em options
423 //
425 opt.SetVerbose(verbose);
426
427 // Multiple Coulomb scattering
428 //
429 opt.SetPolarAngleLimit(CLHEP::pi);
430
431 // Physics tables
432 //
433 opt.SetMinEnergy(10*eV);
434 opt.SetMaxEnergy(10*TeV);
435 opt.SetDEDXBinning(240);
436 opt.SetLambdaBinning(240);
437
438 // Nuclear stopping
439 pnuc->SetMaxKinEnergy(MeV);
440
441 // Ionization
442 //
443 //opt.SetSubCutoff(true);
444
445 // Deexcitation
448 de->SetFluo(true);
449}
@ fUseDistanceToBoundary
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void SetMaxEnergy(G4double val)
void SetDEDXBinning(G4int val)
void SetLambdaBinning(G4int val)
void SetPolarAngleLimit(G4double val)
void SetVerbose(G4int val, const G4String &name="all")
void SetMinEnergy(G4double val)
void SetAtomDeexcitation(G4VAtomDeexcitation *)
const G4String & GetParticleName() const
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:606
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void SetMinKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=1)
void SetMaxKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=1)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
void SetStepFunction(G4double v1, G4double v2)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
void SetStepLimitType(G4MscStepLimitType val)
const G4String & GetPhysicsName() const
G4ParticleTable::G4PTblDicIterator * theParticleIterator

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