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

#include <G4EmStandardPhysics_option2.hh>

Inheritance diagram for G4EmStandardPhysics_option2:
G4VPhysicsConstructor

Public Member Functions

 G4EmStandardPhysics_option2 (G4int ver=1)
 
 G4EmStandardPhysics_option2 (G4int ver, const G4String &name)
 
virtual ~G4EmStandardPhysics_option2 ()
 
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 55 of file G4EmStandardPhysics_option2.hh.

Constructor & Destructor Documentation

G4EmStandardPhysics_option2::G4EmStandardPhysics_option2 ( G4int  ver = 1)

Definition at line 120 of file G4EmStandardPhysics_option2.cc.

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

121  : G4VPhysicsConstructor("G4EmStandard_opt2"), verbose(ver)
122 {
125 }
static G4LossTableManager * Instance()
G4VPhysicsConstructor(const G4String &="")
G4EmStandardPhysics_option2::G4EmStandardPhysics_option2 ( G4int  ver,
const G4String name 
)

Definition at line 129 of file G4EmStandardPhysics_option2.cc.

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

130  : G4VPhysicsConstructor("G4EmStandard_opt2"), verbose(ver)
131 {
134 }
static G4LossTableManager * Instance()
G4VPhysicsConstructor(const G4String &="")
G4EmStandardPhysics_option2::~G4EmStandardPhysics_option2 ( )
virtual

Definition at line 138 of file G4EmStandardPhysics_option2.cc.

139 {}

Member Function Documentation

void G4EmStandardPhysics_option2::ConstructParticle ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 143 of file G4EmStandardPhysics_option2.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().

144 {
145 // gamma
146  G4Gamma::Gamma();
147 
148 // leptons
153 
154 // mesons
159 
160 // barions
163 
164 // ions
167  G4He3::He3();
168  G4Alpha::Alpha();
170 }
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_option2::ConstructProcess ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 174 of file G4EmStandardPhysics_option2.cc.

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

175 {
176  if(verbose > 1) {
177  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
178  }
180 
181  // muon & hadron bremsstrahlung and pair production
190 
191  // muon & hadron multiple scattering
193  mumsc->AddEmModel(0, new G4WentzelVIModel());
195  pmsc->AddEmModel(0, new G4WentzelVIModel());
197  pimsc->AddEmModel(0, new G4WentzelVIModel());
199  kmsc->AddEmModel(0, new G4WentzelVIModel());
200  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
201 
202  // high energy limit for e+- scattering models and bremsstrahlung
203  G4double highEnergyLimit = 100*MeV;
204 
205  // Add standard EM Processes
206  aParticleIterator->reset();
207  while( (*aParticleIterator)() ){
208  G4ParticleDefinition* particle = aParticleIterator->value();
209  G4String particleName = particle->GetParticleName();
210 
211  if (particleName == "gamma") {
212 
213  ph->RegisterProcess(new G4PhotoElectricEffect(), particle);
214  ph->RegisterProcess(new G4ComptonScattering(), particle);
215  ph->RegisterProcess(new G4GammaConversion(), particle);
216  ph->RegisterProcess(new G4RayleighScattering(), particle);
217 
218  } else if (particleName == "e-") {
219 
220  G4eIonisation* eioni = new G4eIonisation();
221  eioni->SetStepFunction(0.8, 1.0*mm);
222 
225  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
226  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
227  msc1->SetHighEnergyLimit(highEnergyLimit);
228  msc2->SetLowEnergyLimit(highEnergyLimit);
229  msc->AddEmModel(0, msc1);
230  msc->AddEmModel(0, msc2);
231 
234  ss->SetEmModel(ssm, 1);
235  ss->SetMinKinEnergy(highEnergyLimit);
236  ssm->SetLowEnergyLimit(highEnergyLimit);
237  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
238 
239  G4eBremsstrahlung* brem = new G4eBremsstrahlung();
244  brem->SetEmModel(br1,1);
245  brem->SetEmModel(br2,2);
246  br2->SetLowEnergyLimit(GeV);
247 
248  ph->RegisterProcess(msc, particle);
249  ph->RegisterProcess(eioni, particle);
250  ph->RegisterProcess(brem, particle);
251  ph->RegisterProcess(ss, particle);
252 
253  } else if (particleName == "e+") {
254 
255  G4eIonisation* eioni = new G4eIonisation();
256  eioni->SetStepFunction(0.8, 1.0*mm);
257 
260  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
261  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
262  msc1->SetHighEnergyLimit(highEnergyLimit);
263  msc2->SetLowEnergyLimit(highEnergyLimit);
264  msc->AddEmModel(0, msc1);
265  msc->AddEmModel(0, msc2);
266 
269  ss->SetEmModel(ssm, 1);
270  ss->SetMinKinEnergy(highEnergyLimit);
271  ssm->SetLowEnergyLimit(highEnergyLimit);
272  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
273 
274  ph->RegisterProcess(msc, particle);
275  ph->RegisterProcess(eioni, particle);
276  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
277  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
278  ph->RegisterProcess(ss, particle);
279 
280  } else if (particleName == "mu+" ||
281  particleName == "mu-" ) {
282 
283  ph->RegisterProcess(mumsc, particle);
284  ph->RegisterProcess(new G4MuIonisation(), particle);
285  ph->RegisterProcess(mub, particle);
286  ph->RegisterProcess(mup, particle);
287  ph->RegisterProcess(new G4CoulombScattering(), particle);
288 
289  } else if (particleName == "alpha" ||
290  particleName == "He3") {
291 
292  //ph->RegisterProcess(hmsc, particle);
293  ph->RegisterProcess(new G4hMultipleScattering(), particle);
294  ph->RegisterProcess(new G4ionIonisation(), particle);
295 
296  } else if (particleName == "GenericIon") {
297 
298  G4ionIonisation* ionIoni = new G4ionIonisation();
299  //ionIoni->SetEmModel(new G4IonParametrisedLossModel());
300  //ionIoni->SetStepFunction(0.1, 20*um);
301 
302  ph->RegisterProcess(hmsc, particle);
303  ph->RegisterProcess(ionIoni, particle);
304 
305  } else if (particleName == "pi+" ||
306  particleName == "pi-" ) {
307 
308  ph->RegisterProcess(pimsc, particle);
309  ph->RegisterProcess(new G4hIonisation(), particle);
310  ph->RegisterProcess(pib, particle);
311  ph->RegisterProcess(pip, particle);
312  ph->RegisterProcess(new G4CoulombScattering(), particle);
313 
314  } else if (particleName == "kaon+" ||
315  particleName == "kaon-" ) {
316 
317  ph->RegisterProcess(kmsc, particle);
318  ph->RegisterProcess(new G4hIonisation(), particle);
319  ph->RegisterProcess(kb, particle);
320  ph->RegisterProcess(kp, particle);
321  ph->RegisterProcess(new G4CoulombScattering(), particle);
322 
323  } else if (particleName == "proton" ||
324  particleName == "anti_proton") {
325 
326  ph->RegisterProcess(pmsc, particle);
327  ph->RegisterProcess(new G4hIonisation(), particle);
328  ph->RegisterProcess(pb, particle);
329  ph->RegisterProcess(pp, particle);
330  ph->RegisterProcess(new G4CoulombScattering(), particle);
331 
332  } else if (particleName == "B+" ||
333  particleName == "B-" ||
334  particleName == "D+" ||
335  particleName == "D-" ||
336  particleName == "Ds+" ||
337  particleName == "Ds-" ||
338  particleName == "anti_He3" ||
339  particleName == "anti_alpha" ||
340  particleName == "anti_deuteron" ||
341  particleName == "anti_lambda_c+" ||
342  particleName == "anti_omega-" ||
343  particleName == "anti_sigma_c+" ||
344  particleName == "anti_sigma_c++" ||
345  particleName == "anti_sigma+" ||
346  particleName == "anti_sigma-" ||
347  particleName == "anti_triton" ||
348  particleName == "anti_xi_c+" ||
349  particleName == "anti_xi-" ||
350  particleName == "deuteron" ||
351  particleName == "lambda_c+" ||
352  particleName == "omega-" ||
353  particleName == "sigma_c+" ||
354  particleName == "sigma_c++" ||
355  particleName == "sigma+" ||
356  particleName == "sigma-" ||
357  particleName == "tau+" ||
358  particleName == "tau-" ||
359  particleName == "triton" ||
360  particleName == "xi_c+" ||
361  particleName == "xi-" ) {
362 
363  ph->RegisterProcess(hmsc, particle);
364  ph->RegisterProcess(new G4hIonisation(), particle);
365  }
366  }
367 
368  // Em options
369  //
370  G4EmProcessOptions opt;
371  opt.SetVerbose(verbose);
372 
373  // Scattering options
374  //
375  opt.SetPolarAngleLimit(CLHEP::pi);
376 
377  // Ionization
378  //
379  //opt.SetSubCutoff(true);
380 
381  // Deexcitation
382  //
385 
386 }
static G4LossTableManager * Instance()
void SetStepFunction(G4double v1, G4double v2)
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)
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetEmModel(G4VEmModel *, G4int index=1)
#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 SetStepLimitType(G4MscStepLimitType val)
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: