Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EmPenelopePhysics Class Reference

#include <G4EmPenelopePhysics.hh>

Inheritance diagram for G4EmPenelopePhysics:
Collaboration diagram for G4EmPenelopePhysics:

Public Member Functions

 G4EmPenelopePhysics (G4int ver=1, const G4String &name="")
 
virtual ~G4EmPenelopePhysics ()
 
virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
- 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
G4int GetInstanceID () const
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VPhysicsConstructor
static const G4VPCManagerGetSubInstanceManager ()
 
- Protected Member Functions inherited from G4VPhysicsConstructor
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
 
- 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 36 of file G4EmPenelopePhysics.hh.

Constructor & Destructor Documentation

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

Definition at line 133 of file G4EmPenelopePhysics.cc.

134  : G4VPhysicsConstructor("G4EmPenelope"), verbose(ver)
135 {
137  param->SetDefaults();
138  param->SetVerbose(verbose);
139  param->SetMinEnergy(100*eV);
140  param->SetMaxEnergy(10*TeV);
141  param->SetLowestElectronEnergy(100*eV);
142  param->SetNumberOfBinsPerDecade(20);
143  param->SetMscRangeFactor(0.02);
145  param->SetMuHadLateralDisplacement(true);
146  param->SetFluo(true);
147  param->SetPIXEElectronCrossSectionModel("Penelope");
149 }
void SetVerbose(G4int val)
void SetLowestElectronEnergy(G4double val)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetPIXEElectronCrossSectionModel(const G4String &)
void SetMaxEnergy(G4double val)
static constexpr double TeV
Definition: G4SIunits.hh:218
void SetMscRangeFactor(G4double val)
void SetNumberOfBinsPerDecade(G4int val)
static constexpr double eV
Definition: G4SIunits.hh:215
void SetMuHadLateralDisplacement(G4bool val)
void SetMinEnergy(G4double val)
static G4EmParameters * Instance()
G4VPhysicsConstructor(const G4String &="")
void SetFluo(G4bool val)

Here is the call graph for this function:

G4EmPenelopePhysics::~G4EmPenelopePhysics ( )
virtual

Definition at line 153 of file G4EmPenelopePhysics.cc.

154 {}

Member Function Documentation

void G4EmPenelopePhysics::ConstructParticle ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 158 of file G4EmPenelopePhysics.cc.

159 {
160  // gamma
161  G4Gamma::Gamma();
162 
163  // leptons
168 
169  // mesons
174 
175  // baryons
178 
179  // ions
182  G4He3::He3();
183  G4Alpha::Alpha();
185 }
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4He3 * He3()
Definition: G4He3.cc:94

Here is the call graph for this function:

void G4EmPenelopePhysics::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 189 of file G4EmPenelopePhysics.cc.

190 {
191  if(verbose > 1) {
192  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
193  }
195 
196  // muon & hadron bremsstrahlung and pair production
205 
206  // muon & hadron multiple scattering
208  mumsc->AddEmModel(0, new G4WentzelVIModel());
209  //G4MuMultipleScattering* pimsc = new G4MuMultipleScattering();
210  //pimsc->AddEmModel(0, new G4WentzelVIModel());
211  //G4MuMultipleScattering* kmsc = new G4MuMultipleScattering();
212  //kmsc->AddEmModel(0, new G4WentzelVIModel());
213  //G4MuMultipleScattering* pmsc = new G4MuMultipleScattering();
214  //pmsc->AddEmModel(0, new G4WentzelVIModel());
215  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
216 
217  // high energy limit for e+- scattering models
218  G4double highEnergyLimit = 100*MeV;
219 
220  // nuclear stopping
221  G4NuclearStopping* pnuc = new G4NuclearStopping();
222 
223  // Add Penelope EM Processes
224  auto myParticleIterator=GetParticleIterator();
225  myParticleIterator->reset();
226 
227  while( (*myParticleIterator)() ){
228 
229  G4ParticleDefinition* particle = myParticleIterator->value();
230  G4String particleName = particle->GetParticleName();
231 
232  //Applicability range for Penelope models
233  //for higher energies, the Standard models are used
234  G4double PenelopeHighEnergyLimit = 1.0*GeV;
235 
236  if (particleName == "gamma") {
237 
238  //Photo-electric effect
239  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
240  G4PenelopePhotoElectricModel* thePEPenelopeModel = new
242  thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
243  thePhotoElectricEffect->SetEmModel(thePEPenelopeModel, 1);
244  ph->RegisterProcess(thePhotoElectricEffect, particle);
245 
246  //Compton scattering
247  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
248  G4PenelopeComptonModel* theComptonPenelopeModel =
250  theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
251  theComptonScattering->SetEmModel(theComptonPenelopeModel, 1);
252  ph->RegisterProcess(theComptonScattering, particle);
253 
254  //Gamma conversion
255  G4GammaConversion* theGammaConversion = new G4GammaConversion();
256  G4PenelopeGammaConversionModel* theGCPenelopeModel =
258  theGammaConversion->SetEmModel(theGCPenelopeModel,1);
259  ph->RegisterProcess(theGammaConversion, particle);
260 
261  //Rayleigh scattering
262  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
263  G4PenelopeRayleighModel* theRayleighPenelopeModel =
265  //theRayleighPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
266  theRayleigh->SetEmModel(theRayleighPenelopeModel, 1);
267  ph->RegisterProcess(theRayleigh, particle);
268 
269  } else if (particleName == "e-") {
270 
271  // multiple scattering
273  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
274  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
275  msc1->SetHighEnergyLimit(highEnergyLimit);
276  msc2->SetLowEnergyLimit(highEnergyLimit);
277  msc->AddEmModel(0, msc1);
278  msc->AddEmModel(0, msc2);
279 
282  ss->SetEmModel(ssm, 1);
283  ss->SetMinKinEnergy(highEnergyLimit);
284  ssm->SetLowEnergyLimit(highEnergyLimit);
285  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
286 
287  //Ionisation
288  G4eIonisation* eIoni = new G4eIonisation();
289  G4PenelopeIonisationModel* theIoniPenelope =
291  theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
292  eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
293  eIoni->SetStepFunction(0.2, 100*um); //
294 
295  //Bremsstrahlung
296  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
297  G4PenelopeBremsstrahlungModel* theBremPenelope = new
299  theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
300  eBrem->AddEmModel(0,theBremPenelope);
301 
302  // register processes
303  ph->RegisterProcess(msc, particle);
304  ph->RegisterProcess(eIoni, particle);
305  ph->RegisterProcess(eBrem, particle);
306  ph->RegisterProcess(ss, particle);
307 
308  } else if (particleName == "e+") {
309 
310  // multiple scattering
312  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
313  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
314  msc1->SetHighEnergyLimit(highEnergyLimit);
315  msc2->SetLowEnergyLimit(highEnergyLimit);
316  msc->AddEmModel(0, msc1);
317  msc->AddEmModel(0, msc2);
318 
321  ss->SetEmModel(ssm, 1);
322  ss->SetMinKinEnergy(highEnergyLimit);
323  ssm->SetLowEnergyLimit(highEnergyLimit);
324  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
325 
326  //Ionisation
327  G4eIonisation* eIoni = new G4eIonisation();
328  G4PenelopeIonisationModel* theIoniPenelope =
330  theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
331  eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
332  eIoni->SetStepFunction(0.2, 100*um); //
333 
334  //Bremsstrahlung
335  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
336  G4PenelopeBremsstrahlungModel* theBremPenelope = new
338  theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
339  eBrem->AddEmModel(0,theBremPenelope);
340 
341  //Annihilation
343  G4PenelopeAnnihilationModel* theAnnPenelope = new
345  theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
346  eAnni->AddEmModel(0,theAnnPenelope);
347 
348  // register processes
349  ph->RegisterProcess(msc, particle);
350  ph->RegisterProcess(eIoni, particle);
351  ph->RegisterProcess(eBrem, particle);
352  ph->RegisterProcess(eAnni, particle);
353  ph->RegisterProcess(ss, particle);
354 
355  } else if (particleName == "mu+" ||
356  particleName == "mu-" ) {
357 
358  G4MuIonisation* muIoni = new G4MuIonisation();
359  muIoni->SetStepFunction(0.2, 50*um);
360 
361  ph->RegisterProcess(mumsc, particle);
362  ph->RegisterProcess(muIoni, particle);
363  ph->RegisterProcess(mub, particle);
364  ph->RegisterProcess(mup, particle);
365  ph->RegisterProcess(new G4CoulombScattering(), particle);
366 
367  } else if (particleName == "alpha" ||
368  particleName == "He3" ) {
369 
371  G4ionIonisation* ionIoni = new G4ionIonisation();
372  ionIoni->SetStepFunction(0.1, 10*um);
373 
374  ph->RegisterProcess(msc, particle);
375  ph->RegisterProcess(ionIoni, particle);
376  ph->RegisterProcess(pnuc, particle);
377 
378  } else if (particleName == "GenericIon") {
379 
380  G4ionIonisation* ionIoni = new G4ionIonisation();
381  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
382  ionIoni->SetStepFunction(0.1, 1*um);
383 
384  ph->RegisterProcess(hmsc, particle);
385  ph->RegisterProcess(ionIoni, particle);
386  ph->RegisterProcess(pnuc, particle);
387 
388  } else if (particleName == "pi+" ||
389  particleName == "pi-" ) {
390 
392  G4hIonisation* hIoni = new G4hIonisation();
393  hIoni->SetStepFunction(0.2, 50*um);
394 
395  ph->RegisterProcess(pimsc, particle);
396  ph->RegisterProcess(hIoni, particle);
397  ph->RegisterProcess(pib, particle);
398  ph->RegisterProcess(pip, particle);
399 
400  } else if (particleName == "kaon+" ||
401  particleName == "kaon-" ) {
402 
404  G4hIonisation* hIoni = new G4hIonisation();
405  hIoni->SetStepFunction(0.2, 50*um);
406 
407  ph->RegisterProcess(kmsc, particle);
408  ph->RegisterProcess(hIoni, particle);
409  ph->RegisterProcess(kb, particle);
410  ph->RegisterProcess(kp, particle);
411 
412  } else if (particleName == "proton" ||
413  particleName == "anti_proton") {
414 
416  G4hIonisation* hIoni = new G4hIonisation();
417  hIoni->SetStepFunction(0.2, 50*um);
418 
419  ph->RegisterProcess(pmsc, particle);
420  ph->RegisterProcess(hIoni, particle);
421  ph->RegisterProcess(pb, particle);
422  ph->RegisterProcess(pp, particle);
423  ph->RegisterProcess(pnuc, particle);
424 
425  } else if (particleName == "B+" ||
426  particleName == "B-" ||
427  particleName == "D+" ||
428  particleName == "D-" ||
429  particleName == "Ds+" ||
430  particleName == "Ds-" ||
431  particleName == "anti_He3" ||
432  particleName == "anti_alpha" ||
433  particleName == "anti_deuteron" ||
434  particleName == "anti_lambda_c+" ||
435  particleName == "anti_omega-" ||
436  particleName == "anti_sigma_c+" ||
437  particleName == "anti_sigma_c++" ||
438  particleName == "anti_sigma+" ||
439  particleName == "anti_sigma-" ||
440  particleName == "anti_triton" ||
441  particleName == "anti_xi_c+" ||
442  particleName == "anti_xi-" ||
443  particleName == "deuteron" ||
444  particleName == "lambda_c+" ||
445  particleName == "omega-" ||
446  particleName == "sigma_c+" ||
447  particleName == "sigma_c++" ||
448  particleName == "sigma+" ||
449  particleName == "sigma-" ||
450  particleName == "tau+" ||
451  particleName == "tau-" ||
452  particleName == "triton" ||
453  particleName == "xi_c+" ||
454  particleName == "xi-" ) {
455 
456  ph->RegisterProcess(hmsc, particle);
457  ph->RegisterProcess(new G4hIonisation(), particle);
458  }
459  }
460 
461  // Nuclear stopping
462  pnuc->SetMaxKinEnergy(MeV);
463 
464  // Deexcitation
465  //
466  G4VAtomDeexcitation* deexcitation = new G4UAtomicDeexcitation();
468 
470 }
static G4LossTableManager * Instance()
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
void SetEmModel(G4VEmModel *, G4int index=1)
G4GLOB_DLL std::ostream G4cout
static constexpr double um
Definition: G4SIunits.hh:113
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
const G4String & GetPhysicsName() const
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:745
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetMaxKinEnergy(G4double e)
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetEmModel(G4VEmModel *, G4int index=1)
static constexpr double GeV
Definition: G4SIunits.hh:217
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=nullptr)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
void SetMinKinEnergy(G4double e)
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
void SetAtomDeexcitation(G4VAtomDeexcitation *)

Here is the call graph for this function:


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