Geant4  10.02.p03
PhysListEmStandardMP Class Reference

#include <PhysListEmStandardMP.hh>

Inheritance diagram for PhysListEmStandardMP:
Collaboration diagram for PhysListEmStandardMP:

Public Member Functions

 PhysListEmStandardMP (G4int iNumAngles)
 
 ~PhysListEmStandardMP ()
 
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 fNumAngles
 

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 43 of file PhysListEmStandardMP.hh.

Constructor & Destructor Documentation

◆ PhysListEmStandardMP()

PhysListEmStandardMP::PhysListEmStandardMP ( G4int  iNumAngles)

Definition at line 76 of file PhysListEmStandardMP.cc.

77  : G4VPhysicsConstructor("G4EmStandardMP"), fNumAngles(iNumAngles)
78 {
80  param->SetDefaults();
81  param->SetFluo(true);
84 }
void SetMscStepLimitType(G4MscStepLimitType val)
static G4EmParameters * Instance()
G4VPhysicsConstructor(const G4String &="")
void SetFluo(G4bool val)
Here is the call graph for this function:

◆ ~PhysListEmStandardMP()

PhysListEmStandardMP::~PhysListEmStandardMP ( )

Definition at line 88 of file PhysListEmStandardMP.cc.

89 {}

Member Function Documentation

◆ ConstructParticle()

virtual void PhysListEmStandardMP::ConstructParticle ( void  )
inlinevirtual

Implements G4VPhysicsConstructor.

Definition at line 51 of file PhysListEmStandardMP.hh.

51 {};
Here is the call graph for this function:

◆ ConstructProcess()

void PhysListEmStandardMP::ConstructProcess ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 93 of file PhysListEmStandardMP.cc.

94 {
96 
97  // Add standard EM Processes
98  //
100  particleIterator->reset();
101  while( (*particleIterator)() ){
102  G4ParticleDefinition* particle = particleIterator->value();
103  G4String particleName = particle->GetParticleName();
104 
105  if (particleName == "gamma") {
106 
107  ph->RegisterProcess(new G4PhotoElectricEffect, particle);
109  cs->SetEmModel(new G4KleinNishinaModel());
110  ph->RegisterProcess(cs, particle);
111  ph->RegisterProcess(new G4GammaConversion, particle);
113 
114  } else if (particleName == "e-") {
115 
117  particle);
118 
119  // ph->RegisterProcess(new G4eMultipleScattering(), particle);
120  //
121  G4eIonisation* eIoni = new G4eIonisation();
122  eIoni->SetStepFunction(0.1, 100*um);
123  ph->RegisterProcess(eIoni, particle);
124  //
125  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
126 
127  } else if (particleName == "e+") {
128 
129  ph->RegisterProcess(new G4eMultipleScattering(), particle);
130  //
131  G4eIonisation* eIoni = new G4eIonisation();
132  eIoni->SetStepFunction(0.1, 100*um);
133  ph->RegisterProcess(eIoni, particle);
134  //
135  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
136  //
137  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
138 
139  } else if (particleName == "mu+" ||
140  particleName == "mu-" ) {
141 
142  ph->RegisterProcess(new G4MuMultipleScattering(), particle);
143  G4MuIonisation* muIoni = new G4MuIonisation();
144  muIoni->SetStepFunction(0.1, 50*um);
145  ph->RegisterProcess(muIoni, particle);
146  ph->RegisterProcess(new G4MuBremsstrahlung(), particle);
147  ph->RegisterProcess(new G4MuPairProduction(), particle);
148 
149  } else if( particleName == "proton" ||
150  particleName == "pi-" ||
151  particleName == "pi+" ) {
152 
153  ph->RegisterProcess(new G4hMultipleScattering(), particle);
154  G4hIonisation* hIoni = new G4hIonisation();
155  hIoni->SetStepFunction(0.1, 20*um);
156  ph->RegisterProcess(hIoni, particle);
157  ph->RegisterProcess(new G4hBremsstrahlung(), particle);
158  ph->RegisterProcess(new G4hPairProduction(), particle);
159 
160  } else if( particleName == "alpha" ||
161  particleName == "He3" ) {
162 
163  ph->RegisterProcess(new G4hMultipleScattering(), particle);
164  G4ionIonisation* ionIoni = new G4ionIonisation();
165  ionIoni->SetStepFunction(0.1, 1*um);
166  ph->RegisterProcess(ionIoni, particle);
167  ph->RegisterProcess(new G4NuclearStopping(), particle);
168 
169  } else if( particleName == "GenericIon" ) {
170 
171  ph->RegisterProcess(new G4hMultipleScattering(), particle);
172  G4ionIonisation* ionIoni = new G4ionIonisation();
173  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
174  ionIoni->SetStepFunction(0.1, 1*um);
175  ph->RegisterProcess(ionIoni, particle);
176  ph->RegisterProcess(new G4NuclearStopping(), particle);
177 
178  } else if ((!particle->IsShortLived()) &&
179  (particle->GetPDGCharge() != 0.0) &&
180  (particle->GetParticleName() != "chargedgeantino")) {
181 
182  //all others charged particles except geantino
183  ph->RegisterProcess(new G4hMultipleScattering(), particle);
184  ph->RegisterProcess(new G4hIonisation(), particle);
185  }
186  }
187 
188  // Deexcitation
189  //
192 
193  G4EmModelActivator mact;
194  mact.ConstructProcess();
195 }
static G4LossTableManager * Instance()
void SetStepFunction(G4double v1, G4double v2)
void SetEmModel(G4VEmModel *, G4int index=1)
const G4String & GetParticleName() const
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetEmModel(G4VEmModel *, G4int index=1)
static const double um
Definition: G4SIunits.hh:112
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:65
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
void SetAtomDeexcitation(G4VAtomDeexcitation *)
G4double GetPDGCharge() const
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fNumAngles

G4int PhysListEmStandardMP::fNumAngles
private

Definition at line 59 of file PhysListEmStandardMP.hh.


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