Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4EmStandardPhysics Class Reference

#include <G4EmStandardPhysics.hh>

Inheritance diagram for G4EmStandardPhysics:
G4VPhysicsConstructor

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
 

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)
 
- 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::G4EmStandardPhysics ( G4int  ver = 0)

Definition at line 110 of file G4EmStandardPhysics.cc.

References bElectromagnetic, G4LossTableManager::Instance(), and G4VPhysicsConstructor::SetPhysicsType().

111  : G4VPhysicsConstructor("G4EmStandard"), verbose(ver)
112 {
115 }
static G4LossTableManager * Instance()
G4VPhysicsConstructor(const G4String &="")
G4EmStandardPhysics::G4EmStandardPhysics ( G4int  ver,
const G4String name 
)

Definition at line 119 of file G4EmStandardPhysics.cc.

References bElectromagnetic, G4LossTableManager::Instance(), and G4VPhysicsConstructor::SetPhysicsType().

120  : G4VPhysicsConstructor("G4EmStandard"), verbose(ver)
121 {
124 }
static G4LossTableManager * Instance()
G4VPhysicsConstructor(const G4String &="")
G4EmStandardPhysics::~G4EmStandardPhysics ( )
virtual

Definition at line 128 of file G4EmStandardPhysics.cc.

129 {}

Member Function Documentation

void G4EmStandardPhysics::ConstructParticle ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 133 of file G4EmStandardPhysics.cc.

References G4Alpha::Alpha(), G4AntiProton::AntiProton(), G4Deuteron::Deuteron(), G4Electron::Electron(), G4Gamma::Gamma(), G4GenericIon::GenericIonDefinition(), G4He3::He3(), G4KaonMinus::KaonMinusDefinition(), G4KaonPlus::KaonPlusDefinition(), G4MuonMinus::MuonMinus(), G4MuonPlus::MuonPlus(), G4PionMinus::PionMinusDefinition(), G4PionPlus::PionPlusDefinition(), G4Positron::Positron(), G4Proton::Proton(), and G4Triton::Triton().

134 {
135 // gamma
136  G4Gamma::Gamma();
137 
138 // leptons
143 
144 // mesons
149 
150 // barions
153 
154 // ions
157  G4He3::He3();
158  G4Alpha::Alpha();
160 }
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
void G4EmStandardPhysics::ConstructProcess ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 164 of file G4EmStandardPhysics.cc.

References G4VMultipleScattering::AddEmModel(), aParticleIterator, G4cout, G4endl, G4ParticleDefinition::GetParticleName(), G4PhysicsListHelper::GetPhysicsListHelper(), G4VPhysicsConstructor::GetPhysicsName(), G4LossTableManager::Instance(), python.hepunit::MeV, G4InuclParticleNames::mup, G4InuclParticleNames::pip, G4InuclParticleNames::pp, G4PhysicsListHelper::RegisterProcess(), G4VEmModel::SetActivationLowEnergyLimit(), G4LossTableManager::SetAtomDeexcitation(), G4VEmProcess::SetEmModel(), G4VEmModel::SetHighEnergyLimit(), G4VEmModel::SetLowEnergyLimit(), G4VEmProcess::SetMinKinEnergy(), G4EmProcessOptions::SetPolarAngleLimit(), and G4EmProcessOptions::SetVerbose().

165 {
166  if(verbose > 1) {
167  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
168  }
170 
171  // muon & hadron bremsstrahlung and pair production
180 
181  // muon & hadron multiple scattering
183  mumsc->AddEmModel(0, new G4WentzelVIModel());
185  pimsc->AddEmModel(0, new G4WentzelVIModel());
187  kmsc->AddEmModel(0, new G4WentzelVIModel());
189  pmsc->AddEmModel(0, new G4WentzelVIModel());
190  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
191 
192  // high energy limit for e+- scattering models
193  G4double highEnergyLimit = 100*MeV;
194 
195  // Add standard EM Processes
196  aParticleIterator->reset();
197  while( (*aParticleIterator)() ){
198  G4ParticleDefinition* particle = aParticleIterator->value();
199  G4String particleName = particle->GetParticleName();
200 
201  if (particleName == "gamma") {
202 
203  ph->RegisterProcess(new G4PhotoElectricEffect(), particle);
204  ph->RegisterProcess(new G4ComptonScattering(), particle);
205  ph->RegisterProcess(new G4GammaConversion(), particle);
206 
207  } else if (particleName == "e-") {
208 
210  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
211  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
212  msc1->SetHighEnergyLimit(highEnergyLimit);
213  msc2->SetLowEnergyLimit(highEnergyLimit);
214  msc->AddEmModel(0, msc1);
215  msc->AddEmModel(0, msc2);
216 
219  ss->SetEmModel(ssm, 1);
220  ss->SetMinKinEnergy(highEnergyLimit);
221  ssm->SetLowEnergyLimit(highEnergyLimit);
222  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
223 
224  ph->RegisterProcess(msc, particle);
225  ph->RegisterProcess(new G4eIonisation(), particle);
226  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
227  ph->RegisterProcess(ss, particle);
228 
229  } else if (particleName == "e+") {
230 
232  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
233  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
234  msc1->SetHighEnergyLimit(highEnergyLimit);
235  msc2->SetLowEnergyLimit(highEnergyLimit);
236  msc->AddEmModel(0, msc1);
237  msc->AddEmModel(0, msc2);
238 
241  ss->SetEmModel(ssm, 1);
242  ss->SetMinKinEnergy(highEnergyLimit);
243  ssm->SetLowEnergyLimit(highEnergyLimit);
244  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
245 
246  ph->RegisterProcess(msc, particle);
247  ph->RegisterProcess(new G4eIonisation(), particle);
248  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
249  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
250  ph->RegisterProcess(ss, particle);
251 
252  } else if (particleName == "mu+" ||
253  particleName == "mu-" ) {
254 
255  ph->RegisterProcess(mumsc, particle);
256  ph->RegisterProcess(new G4MuIonisation(), particle);
257  ph->RegisterProcess(mub, particle);
258  ph->RegisterProcess(mup, particle);
259  ph->RegisterProcess(new G4CoulombScattering(), particle);
260 
261  } else if (particleName == "alpha" ||
262  particleName == "He3") {
263 
264  //ph->RegisterProcess(hmsc, particle);
265  ph->RegisterProcess(new G4hMultipleScattering(), particle);
266  ph->RegisterProcess(new G4ionIonisation(), particle);
267 
268  } else if (particleName == "GenericIon") {
269 
270  ph->RegisterProcess(hmsc, particle);
271  ph->RegisterProcess(new G4ionIonisation(), particle);
272 
273  } else if (particleName == "pi+" ||
274  particleName == "pi-" ) {
275 
276  //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
277  ph->RegisterProcess(pimsc, particle);
278  ph->RegisterProcess(new G4hIonisation(), particle);
279  ph->RegisterProcess(pib, particle);
280  ph->RegisterProcess(pip, particle);
281  ph->RegisterProcess(new G4CoulombScattering(), particle);
282 
283  } else if (particleName == "kaon+" ||
284  particleName == "kaon-" ) {
285 
286  //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
287  ph->RegisterProcess(kmsc, particle);
288  ph->RegisterProcess(new G4hIonisation(), particle);
289  ph->RegisterProcess(kb, particle);
290  ph->RegisterProcess(kp, particle);
291  ph->RegisterProcess(new G4CoulombScattering(), particle);
292 
293  } else if (particleName == "proton" ||
294  particleName == "anti_proton") {
295 
296  //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
297  ph->RegisterProcess(pmsc, particle);
298  ph->RegisterProcess(new G4hIonisation(), particle);
299  ph->RegisterProcess(pb, particle);
300  ph->RegisterProcess(pp, particle);
301  ph->RegisterProcess(new G4CoulombScattering(), particle);
302 
303  } else if (particleName == "B+" ||
304  particleName == "B-" ||
305  particleName == "D+" ||
306  particleName == "D-" ||
307  particleName == "Ds+" ||
308  particleName == "Ds-" ||
309  particleName == "anti_He3" ||
310  particleName == "anti_alpha" ||
311  particleName == "anti_deuteron" ||
312  particleName == "anti_lambda_c+" ||
313  particleName == "anti_omega-" ||
314  particleName == "anti_sigma_c+" ||
315  particleName == "anti_sigma_c++" ||
316  particleName == "anti_sigma+" ||
317  particleName == "anti_sigma-" ||
318  particleName == "anti_triton" ||
319  particleName == "anti_xi_c+" ||
320  particleName == "anti_xi-" ||
321  particleName == "deuteron" ||
322  particleName == "lambda_c+" ||
323  particleName == "omega-" ||
324  particleName == "sigma_c+" ||
325  particleName == "sigma_c++" ||
326  particleName == "sigma+" ||
327  particleName == "sigma-" ||
328  particleName == "tau+" ||
329  particleName == "tau-" ||
330  particleName == "triton" ||
331  particleName == "xi_c+" ||
332  particleName == "xi-" ) {
333 
334  ph->RegisterProcess(hmsc, particle);
335  ph->RegisterProcess(new G4hIonisation(), particle);
336  }
337  }
338  G4EmProcessOptions opt;
339  opt.SetVerbose(verbose);
340  // opt.SetApplyCuts(true);
341  opt.SetPolarAngleLimit(CLHEP::pi);
342 
343  // Deexcitation
344  //
347 }
static G4LossTableManager * Instance()
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
void SetEmModel(G4VEmModel *, G4int index=1)
G4GLOB_DLL std::ostream G4cout
#define aParticleIterator
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
const G4String & GetPhysicsName() const
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:704
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
static G4PhysicsListHelper * GetPhysicsListHelper()
#define G4endl
Definition: G4ios.hh:61
void SetMinKinEnergy(G4double e)
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetVerbose(G4int val, const G4String &name="all", G4bool worker=false)
void SetPolarAngleLimit(G4double val)

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