Geant4  10.02.p03
G4EmStandardPhysics Class Reference

#include <G4EmStandardPhysics.hh>

Inheritance diagram for G4EmStandardPhysics:
Collaboration diagram for G4EmStandardPhysics:

Public Member Functions

 G4EmStandardPhysics (G4int ver=0)
 
 G4EmStandardPhysics (G4int ver, const G4String &name)
 
virtual ~G4EmStandardPhysics ()
 
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
 

Private Attributes

G4int verbose
 

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 52 of file G4EmStandardPhysics.hh.

Constructor & Destructor Documentation

◆ G4EmStandardPhysics() [1/2]

G4EmStandardPhysics::G4EmStandardPhysics ( G4int  ver = 0)

Definition at line 111 of file G4EmStandardPhysics.cc.

112  : G4VPhysicsConstructor("G4EmStandard"), verbose(ver)
113 {
115  param->SetDefaults();
116  param->SetVerbose(verbose);
118 }
void SetVerbose(G4int val)
static G4EmParameters * Instance()
G4VPhysicsConstructor(const G4String &="")
Here is the call graph for this function:

◆ G4EmStandardPhysics() [2/2]

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

Definition at line 122 of file G4EmStandardPhysics.cc.

123  : G4VPhysicsConstructor("G4EmStandard"), verbose(ver)
124 {
126  param->SetDefaults();
127  param->SetVerbose(verbose);
129 }
void SetVerbose(G4int val)
static G4EmParameters * Instance()
G4VPhysicsConstructor(const G4String &="")
Here is the call graph for this function:

◆ ~G4EmStandardPhysics()

G4EmStandardPhysics::~G4EmStandardPhysics ( )
virtual

Definition at line 133 of file G4EmStandardPhysics.cc.

134 {}

Member Function Documentation

◆ ConstructParticle()

void G4EmStandardPhysics::ConstructParticle ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 138 of file G4EmStandardPhysics.cc.

139 {
140  // gamma
141  G4Gamma::Gamma();
142 
143  // leptons
148 
149  // mesons
154 
155  // barions
158 
159  // ions
162  G4He3::He3();
163  G4Alpha::Alpha();
165 
166  // dna
167  G4EmModelActivator mact;
168  mact.ConstructParticle();
169 }
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:

◆ ConstructProcess()

void G4EmStandardPhysics::ConstructProcess ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 173 of file G4EmStandardPhysics.cc.

174 {
175  if(verbose > 1) {
176  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
177  }
179 
180  // muon & hadron bremsstrahlung and pair production
189 
190  // muon & hadron multiple scattering
192  mumsc->AddEmModel(0, new G4WentzelVIModel());
194 
196  pimsc->AddEmModel(0, new G4WentzelVIModel());
198 
200  kmsc->AddEmModel(0, new G4WentzelVIModel());
202 
204  pmsc->AddEmModel(0, new G4WentzelVIModel());
206 
207  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
208 
209  // high energy limit for e+- scattering models
210  G4double highEnergyLimit = 100*MeV;
211 
212  // Add standard EM Processes
213  auto myParticleIterator=GetParticleIterator();
214  myParticleIterator->reset();
215  while( (*myParticleIterator)() ){
216  G4ParticleDefinition* particle = myParticleIterator->value();
217  G4String particleName = particle->GetParticleName();
218 
219  if (particleName == "gamma") {
220 
221  ph->RegisterProcess(new G4PhotoElectricEffect(), particle);
222  ph->RegisterProcess(new G4ComptonScattering(), particle);
223  ph->RegisterProcess(new G4GammaConversion(), particle);
224 
225  } else if (particleName == "e-") {
226 
228  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
229  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
230  msc1->SetNewDisplacementFlag(false);
231  msc1->SetHighEnergyLimit(highEnergyLimit);
232  msc2->SetLowEnergyLimit(highEnergyLimit);
233  msc->AddEmModel(0, msc1);
234  msc->AddEmModel(0, msc2);
235 
238  ss->SetEmModel(ssm, 1);
239  ss->SetMinKinEnergy(highEnergyLimit);
240  ssm->SetLowEnergyLimit(highEnergyLimit);
241  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
242 
243  ph->RegisterProcess(msc, particle);
244  ph->RegisterProcess(new G4eIonisation(), particle);
245  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
246  ph->RegisterProcess(ss, particle);
247 
248  } else if (particleName == "e+") {
249 
251  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
252  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
253  msc1->SetNewDisplacementFlag(false);
254  msc1->SetHighEnergyLimit(highEnergyLimit);
255  msc2->SetLowEnergyLimit(highEnergyLimit);
256  msc->AddEmModel(0, msc1);
257  msc->AddEmModel(0, msc2);
258 
261  ss->SetEmModel(ssm, 1);
262  ss->SetMinKinEnergy(highEnergyLimit);
263  ssm->SetLowEnergyLimit(highEnergyLimit);
264  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
265 
266  ph->RegisterProcess(msc, particle);
267  ph->RegisterProcess(new G4eIonisation(), particle);
268  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
269  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
270  ph->RegisterProcess(ss, particle);
271 
272  } else if (particleName == "mu+" ||
273  particleName == "mu-" ) {
274 
275  ph->RegisterProcess(mumsc, particle);
276  ph->RegisterProcess(new G4MuIonisation(), particle);
277  ph->RegisterProcess(mub, particle);
278  ph->RegisterProcess(mup, particle);
279  ph->RegisterProcess(muss, particle);
280 
281  } else if (particleName == "alpha" ||
282  particleName == "He3") {
283 
284  //ph->RegisterProcess(hmsc, particle);
285  ph->RegisterProcess(new G4hMultipleScattering(), particle);
286  ph->RegisterProcess(new G4ionIonisation(), particle);
287 
288  } else if (particleName == "GenericIon") {
289 
290  ph->RegisterProcess(hmsc, particle);
291  ph->RegisterProcess(new G4ionIonisation(), particle);
292 
293  } else if (particleName == "pi+" ||
294  particleName == "pi-" ) {
295 
296  //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
297  ph->RegisterProcess(pimsc, particle);
298  ph->RegisterProcess(new G4hIonisation(), particle);
299  ph->RegisterProcess(pib, particle);
300  ph->RegisterProcess(pip, particle);
301  ph->RegisterProcess(piss, particle);
302 
303  } else if (particleName == "kaon+" ||
304  particleName == "kaon-" ) {
305 
306  //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
307  ph->RegisterProcess(kmsc, particle);
308  ph->RegisterProcess(new G4hIonisation(), particle);
309  ph->RegisterProcess(kb, particle);
310  ph->RegisterProcess(kp, particle);
311  ph->RegisterProcess(kss, particle);
312 
313  } else if (particleName == "proton" ||
314  particleName == "anti_proton") {
315 
316  //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
317  ph->RegisterProcess(pmsc, particle);
318  ph->RegisterProcess(new G4hIonisation(), particle);
319  ph->RegisterProcess(pb, particle);
320  ph->RegisterProcess(pp, particle);
321  ph->RegisterProcess(pss, particle);
322 
323  } else if (particleName == "B+" ||
324  particleName == "B-" ||
325  particleName == "D+" ||
326  particleName == "D-" ||
327  particleName == "Ds+" ||
328  particleName == "Ds-" ||
329  particleName == "anti_He3" ||
330  particleName == "anti_alpha" ||
331  particleName == "anti_deuteron" ||
332  particleName == "anti_lambda_c+" ||
333  particleName == "anti_omega-" ||
334  particleName == "anti_sigma_c+" ||
335  particleName == "anti_sigma_c++" ||
336  particleName == "anti_sigma+" ||
337  particleName == "anti_sigma-" ||
338  particleName == "anti_triton" ||
339  particleName == "anti_xi_c+" ||
340  particleName == "anti_xi-" ||
341  particleName == "deuteron" ||
342  particleName == "lambda_c+" ||
343  particleName == "omega-" ||
344  particleName == "sigma_c+" ||
345  particleName == "sigma_c++" ||
346  particleName == "sigma+" ||
347  particleName == "sigma-" ||
348  particleName == "tau+" ||
349  particleName == "tau-" ||
350  particleName == "triton" ||
351  particleName == "xi_c+" ||
352  particleName == "xi-" ) {
353 
354  ph->RegisterProcess(hmsc, particle);
355  ph->RegisterProcess(new G4hIonisation(), particle);
356  }
357  }
358 
359  // Deexcitation
360  //
363 
364  G4EmModelActivator mact;
365  mact.ConstructProcess();
366 }
static const double MeV
Definition: G4SIunits.hh:211
static G4LossTableManager * Instance()
void SetNewDisplacementFlag(G4bool)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
void SetEmModel(G4VEmModel *, G4int index=1)
const G4String & GetParticleName() const
const G4String & GetPhysicsName() const
G4GLOB_DLL std::ostream G4cout
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
static G4PhysicsListHelper * GetPhysicsListHelper()
#define G4endl
Definition: G4ios.hh:61
void SetMinKinEnergy(G4double e)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
void SetAtomDeexcitation(G4VAtomDeexcitation *)
Here is the call graph for this function:

Member Data Documentation

◆ verbose

G4int G4EmStandardPhysics::verbose
private

Definition at line 66 of file G4EmStandardPhysics.hh.


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