Geant4-11
Public Member Functions | Static Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | Static Protected Attributes
G4EmLivermorePolarizedPhysics Class Reference

#include <G4EmLivermorePolarizedPhysics.hh>

Inheritance diagram for G4EmLivermorePolarizedPhysics:
G4EmLivermorePhysics G4VPhysicsConstructor

Public Member Functions

void ConstructParticle () override
 
void ConstructProcess () override
 
 G4EmLivermorePolarizedPhysics (G4int ver=1, const G4String &name="")
 
G4int GetInstanceID () const
 
const G4StringGetPhysicsName () const
 
G4int GetPhysicsType () const
 
G4int GetVerboseLevel () const
 
void SetPhysicsName (const G4String &="")
 
void SetPhysicsType (G4int)
 
void SetVerboseLevel (G4int value)
 
virtual void TerminateWorker ()
 
 ~G4EmLivermorePolarizedPhysics () override
 

Static Public Member Functions

static const G4VPCManagerGetSubInstanceManager ()
 

Protected Types

using PhysicsBuilder_V = G4VPCData::PhysicsBuilders_V
 

Protected Member Functions

void AddBuilder (G4PhysicsBuilderInterface *bld)
 
PhysicsBuilder_V GetBuilders () const
 
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
 
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 

Protected Attributes

G4int g4vpcInstanceID = 0
 
G4String namePhysics = ""
 
G4ParticleTabletheParticleTable = nullptr
 
G4int typePhysics = 0
 
G4int verboseLevel = 0
 

Static Protected Attributes

static G4RUN_DLL G4VPCManager subInstanceManager
 

Detailed Description

Definition at line 34 of file G4EmLivermorePolarizedPhysics.hh.

Member Typedef Documentation

◆ PhysicsBuilder_V

Definition at line 149 of file G4VPhysicsConstructor.hh.

Constructor & Destructor Documentation

◆ G4EmLivermorePolarizedPhysics()

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

Definition at line 37 of file G4EmLivermorePolarizedPhysics.cc.

39 : G4EmLivermorePhysics(ver, "G4EmLivermorePolarized")
40{
42 param->SetEnablePolarisation(true);
43}
G4EmLivermorePhysics(G4int ver=1, const G4String &name="G4EmLivermore")
void SetEnablePolarisation(G4bool val)
static G4EmParameters * Instance()

References G4EmParameters::Instance(), and G4EmParameters::SetEnablePolarisation().

◆ ~G4EmLivermorePolarizedPhysics()

G4EmLivermorePolarizedPhysics::~G4EmLivermorePolarizedPhysics ( )
override

Definition at line 47 of file G4EmLivermorePolarizedPhysics.cc.

48{}

Member Function Documentation

◆ AddBuilder()

void G4VPhysicsConstructor::AddBuilder ( G4PhysicsBuilderInterface bld)
protectedinherited

Definition at line 99 of file G4VPhysicsConstructor.cc.

100{
101 (subInstanceManager.offset[g4vpcInstanceID])._builders->push_back(bld);
102}
static G4RUN_DLL G4VPCManager subInstanceManager
G4RUN_DLL G4ThreadLocalStatic T * offset

References G4VPhysicsConstructor::g4vpcInstanceID, G4VUPLSplitter< T >::offset, and G4VPhysicsConstructor::subInstanceManager.

Referenced by G4HadronPhysicsFTFP_BERT::Kaon(), G4HadronPhysicsFTF_BIC::Kaon(), G4HadronPhysicsINCLXX::Kaon(), G4HadronPhysicsFTFP_BERT::Neutron(), G4HadronPhysicsQGSP_BERT::Neutron(), G4HadronPhysicsQGSP_BIC::Neutron(), G4HadronPhysicsFTF_BIC::Neutron(), G4HadronPhysicsFTFP_BERT_HP::Neutron(), G4HadronPhysicsINCLXX::Neutron(), G4HadronPhysicsQGS_BIC::Neutron(), G4HadronPhysicsQGSP_BERT_HP::Neutron(), G4HadronPhysicsQGSP_BIC_HP::Neutron(), G4HadronPhysicsShielding::Neutron(), G4HadronPhysicsFTFP_BERT::Pion(), G4HadronPhysicsQGSP_BERT::Pion(), G4HadronPhysicsQGSP_BIC::Pion(), G4HadronPhysicsFTF_BIC::Pion(), G4HadronPhysicsINCLXX::Pion(), G4HadronPhysicsQGS_BIC::Pion(), G4HadronPhysicsFTFP_BERT::Proton(), G4HadronPhysicsQGSP_BERT::Proton(), G4HadronPhysicsQGSP_BIC::Proton(), G4HadronPhysicsFTF_BIC::Proton(), G4HadronPhysicsINCLXX::Proton(), G4HadronPhysicsNuBeam::Proton(), G4HadronPhysicsQGS_BIC::Proton(), and G4HadronPhysicsQGSP_BIC_AllHP::Proton().

◆ ConstructParticle()

void G4EmLivermorePhysics::ConstructParticle ( )
overridevirtualinherited

Implements G4VPhysicsConstructor.

Definition at line 143 of file G4EmLivermorePhysics.cc.

144{
145 // minimal set of particles for EM physics
147}
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:347

References G4EmBuilder::ConstructMinimalEmSet().

◆ ConstructProcess()

void G4EmLivermorePhysics::ConstructProcess ( )
overridevirtualinherited

Implements G4VPhysicsConstructor.

Definition at line 151 of file G4EmLivermorePhysics.cc.

152{
153 if(verboseLevel > 1) {
154 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
155 }
159
160 // processes used by several particles
161 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
162
163 // high energy limit for e+- scattering models
164 G4double highEnergyLimit= param->MscEnergyLimit();
165 G4double livEnergyLimit = 1*CLHEP::GeV;
166
167 // nuclear stopping
168 G4double nielEnergyLimit = param->MaxNIELEnergy();
169 G4NuclearStopping* pnuc = nullptr;
170 if(nielEnergyLimit > 0.0) {
171 pnuc = new G4NuclearStopping();
172 pnuc->SetMaxKinEnergy(nielEnergyLimit);
173 }
174
175 // Add Livermore EM Processes
176
177 // Add gamma EM Processes
179 G4bool polar = param->EnablePolarisation();
180
181 // photoelectric effect - Livermore model only
184 pe->SetEmModel(peModel);
185 if(polar) {
187 }
188
189 // Compton scattering - Livermore model only
192 G4VEmModel* cModel = nullptr;
193 if(polar) {
195 cModel->SetHighEnergyLimit(20*CLHEP::MeV);
196 } else {
197 cModel = new G4LivermoreComptonModel();
198 cModel->SetHighEnergyLimit(livEnergyLimit);
199 }
200 cs->AddEmModel(0, cModel);
201
202 // gamma conversion
205 gc->SetEmModel(convLiv);
206
207 // Rayleigh scattering
209 if(polar) {
211 }
212
213 ph->RegisterProcess(pe, particle);
214 ph->RegisterProcess(cs, particle);
215 ph->RegisterProcess(gc, particle);
216 ph->RegisterProcess(rl, particle);
217
218 // e-
219 particle = G4Electron::Electron();
220
221 // multiple and single scattering
225 msc1->SetHighEnergyLimit(highEnergyLimit);
226 msc2->SetLowEnergyLimit(highEnergyLimit);
227 msc->SetEmModel(msc1);
228 msc->SetEmModel(msc2);
229
232 ss->SetEmModel(ssm);
233 ss->SetMinKinEnergy(highEnergyLimit);
234 ssm->SetLowEnergyLimit(highEnergyLimit);
235 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
236
237 // Ionisation - Livermore should be used only for low energies
238 G4eIonisation* eioni = new G4eIonisation();
239 G4VEmModel* theIoniLiv = new G4LivermoreIonisationModel();
240 theIoniLiv->SetHighEnergyLimit(0.1*CLHEP::MeV);
241 eioni->AddEmModel(0, theIoniLiv, new G4UniversalFluctuation() );
242
243 // Bremsstrahlung
249 brem->SetEmModel(br1);
250 brem->SetEmModel(br2);
252
254
255 // register processes
256 ph->RegisterProcess(msc, particle);
257 ph->RegisterProcess(eioni, particle);
258 ph->RegisterProcess(brem, particle);
259 ph->RegisterProcess(ee, particle);
260 ph->RegisterProcess(ss, particle);
261
262 // e+
263 particle = G4Positron::Positron();
264
265 // multiple and single scattering
266 msc = new G4eMultipleScattering;
267 msc1 = new G4GoudsmitSaundersonMscModel();
268 msc2 = new G4WentzelVIModel();
269 msc1->SetHighEnergyLimit(highEnergyLimit);
270 msc2->SetLowEnergyLimit(highEnergyLimit);
271 msc->SetEmModel(msc1);
272 msc->SetEmModel(msc2);
273
274 ssm = new G4eCoulombScatteringModel();
275 ss = new G4CoulombScattering();
276 ss->SetEmModel(ssm);
277 ss->SetMinKinEnergy(highEnergyLimit);
278 ssm->SetLowEnergyLimit(highEnergyLimit);
279 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
280
281 // ionisation
282 eioni = new G4eIonisation();
283
284 // Bremsstrahlung from standard
285 brem = new G4eBremsstrahlung();
286 br1 = new G4SeltzerBergerModel();
287 br2 = new G4eBremsstrahlungRelModel();
290 brem->SetEmModel(br1);
291 brem->SetEmModel(br2);
293
294 // register processes
295 ph->RegisterProcess(msc, particle);
296 ph->RegisterProcess(eioni, particle);
297 ph->RegisterProcess(brem, particle);
298 ph->RegisterProcess(ee, particle);
299 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
300 ph->RegisterProcess(ss, particle);
301
302 // generic ion
303 particle = G4GenericIon::GenericIon();
304 G4ionIonisation* ionIoni = new G4ionIonisation();
306 ph->RegisterProcess(hmsc, particle);
307 ph->RegisterProcess(ionIoni, particle);
308 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
309
310 // muons, hadrons, ions
312
313 // extra configuration
315
316}
static constexpr double GeV
Definition: G4SIunits.hh:203
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition: G4Electron.cc:93
static void ConstructCharged(G4hMultipleScattering *hmsc, G4NuclearStopping *nucStopping, G4bool isWVI=true)
Definition: G4EmBuilder.cc:219
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:375
G4bool EnablePolarisation() const
G4double MaxNIELEnergy() const
G4double MscEnergyLimit() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition: G4Positron.cc:93
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:788
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:628
void SetMinKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetEmModel(G4VMscModel *, G4int idx=0)
const G4String & GetPhysicsName() const
static constexpr double GeV
static constexpr double MeV

References G4VEmProcess::AddEmModel(), G4VEnergyLossProcess::AddEmModel(), G4EmBuilder::ConstructCharged(), G4Electron::Electron(), G4EmParameters::EnablePolarisation(), G4cout, G4endl, G4Gamma::Gamma(), G4GenericIon::GenericIon(), G4PhysicsListHelper::GetPhysicsListHelper(), G4VPhysicsConstructor::GetPhysicsName(), CLHEP::GeV, GeV, G4EmParameters::Instance(), G4EmParameters::MaxNIELEnergy(), CLHEP::MeV, G4EmParameters::MscEnergyLimit(), G4Positron::Positron(), G4EmBuilder::PrepareEMPhysics(), G4PhysicsListHelper::RegisterProcess(), G4VEmModel::SetActivationLowEnergyLimit(), G4VEmModel::SetAngularDistribution(), G4VEmProcess::SetEmModel(), G4VEnergyLossProcess::SetEmModel(), G4VMultipleScattering::SetEmModel(), G4VEmModel::SetHighEnergyLimit(), G4VEmModel::SetLowEnergyLimit(), G4VEmProcess::SetMaxKinEnergy(), G4VEmProcess::SetMinKinEnergy(), and G4VPhysicsConstructor::verboseLevel.

◆ GetBuilders()

G4VPhysicsConstructor::PhysicsBuilder_V G4VPhysicsConstructor::GetBuilders ( ) const
protectedinherited

Definition at line 86 of file G4VPhysicsConstructor.cc.

87{
88 const auto& tls = *((subInstanceManager.offset[g4vpcInstanceID])._builders);
89 PhysicsBuilder_V copy(tls.size());
90 G4int i = 0;
91 for(const auto& el : tls)
92 {
93 copy[i++] = el;
94 }
95 return copy;
96}
int G4int
Definition: G4Types.hh:85
G4VPCData::PhysicsBuilders_V PhysicsBuilder_V
void copy(G4double dst[], const G4double src[], size_t size=G4FieldTrack::ncompSVEC)
Definition: G4FieldUtils.cc:98

References field_utils::copy(), G4VPhysicsConstructor::g4vpcInstanceID, G4VUPLSplitter< T >::offset, and G4VPhysicsConstructor::subInstanceManager.

◆ GetInstanceID()

G4int G4VPhysicsConstructor::GetInstanceID ( ) const
inlineinherited

◆ GetParticleIterator()

G4ParticleTable::G4PTblDicIterator * G4VPhysicsConstructor::GetParticleIterator ( ) const
protectedinherited

◆ GetPhysicsName()

const G4String & G4VPhysicsConstructor::GetPhysicsName ( ) const
inlineinherited

Definition at line 191 of file G4VPhysicsConstructor.hh.

192{
193 return namePhysics;
194}

References G4VPhysicsConstructor::namePhysics.

Referenced by G4EmDNAPhysics_option1::ConstructProcess(), G4EmDNAPhysics_option2::ConstructProcess(), G4EmDNAPhysics_option3::ConstructProcess(), G4EmDNAPhysics_option4::ConstructProcess(), G4EmDNAPhysics_option5::ConstructProcess(), G4EmDNAPhysics_option6::ConstructProcess(), G4EmDNAPhysics_option7::ConstructProcess(), G4EmDNAPhysics_option8::ConstructProcess(), G4EmDNAPhysics_stationary_option2::ConstructProcess(), G4EmDNAPhysics_stationary_option4::ConstructProcess(), G4EmDNAPhysics_stationary_option6::ConstructProcess(), G4EmDNAPhysics::ConstructProcess(), G4EmDNAPhysics_stationary::ConstructProcess(), G4EmLivermorePhysics::ConstructProcess(), G4EmLowEPPhysics::ConstructProcess(), G4EmPenelopePhysics::ConstructProcess(), G4EmStandardPhysics::ConstructProcess(), G4EmStandardPhysics_option1::ConstructProcess(), G4EmStandardPhysics_option2::ConstructProcess(), G4EmStandardPhysics_option3::ConstructProcess(), G4EmStandardPhysics_option4::ConstructProcess(), G4EmStandardPhysicsGS::ConstructProcess(), G4EmStandardPhysicsSS::ConstructProcess(), G4EmStandardPhysicsWVI::ConstructProcess(), G4ThermalNeutrons::ConstructProcess(), G4HadronPhysicsFTFP_BERT::DumpBanner(), G4HadronPhysicsQGSP_BERT::DumpBanner(), export_G4VPhysicsConstructor(), G4HadronDElasticPhysics::G4HadronDElasticPhysics(), G4HadronElasticPhysics::G4HadronElasticPhysics(), G4HadronElasticPhysicsHP::G4HadronElasticPhysicsHP(), G4HadronElasticPhysicsLEND::G4HadronElasticPhysicsLEND(), G4HadronElasticPhysicsPHP::G4HadronElasticPhysicsPHP(), G4HadronElasticPhysicsXS::G4HadronElasticPhysicsXS(), G4HadronHElasticPhysics::G4HadronHElasticPhysics(), G4IonElasticPhysics::G4IonElasticPhysics(), G4VModularPhysicsList::RegisterPhysics(), and G4VModularPhysicsList::ReplacePhysics().

◆ GetPhysicsType()

G4int G4VPhysicsConstructor::GetPhysicsType ( ) const
inlineinherited

◆ GetSubInstanceManager()

const G4VPCManager & G4VPhysicsConstructor::GetSubInstanceManager ( )
inlinestaticinherited

◆ GetVerboseLevel()

G4int G4VPhysicsConstructor::GetVerboseLevel ( ) const
inlineinherited

◆ RegisterProcess()

G4bool G4VPhysicsConstructor::RegisterProcess ( G4VProcess process,
G4ParticleDefinition particle 
)
inlineprotectedinherited

◆ SetPhysicsName()

void G4VPhysicsConstructor::SetPhysicsName ( const G4String name = "")
inlineinherited

Definition at line 186 of file G4VPhysicsConstructor.hh.

187{
189}
const char * name(G4int ptype)

References G4InuclParticleNames::name(), and G4VPhysicsConstructor::namePhysics.

Referenced by export_G4VPhysicsConstructor().

◆ SetPhysicsType()

void G4VPhysicsConstructor::SetPhysicsType ( G4int  val)
inlineinherited

Definition at line 196 of file G4VPhysicsConstructor.hh.

197{
198 if(val > 0) { typePhysics = val; }
199}

References G4VPhysicsConstructor::typePhysics.

Referenced by G4DecayPhysics::G4DecayPhysics(), G4EmDNAPhysics::G4EmDNAPhysics(), G4EmDNAPhysics_option1::G4EmDNAPhysics_option1(), G4EmDNAPhysics_option2::G4EmDNAPhysics_option2(), G4EmDNAPhysics_option3::G4EmDNAPhysics_option3(), G4EmDNAPhysics_option4::G4EmDNAPhysics_option4(), G4EmDNAPhysics_option5::G4EmDNAPhysics_option5(), G4EmDNAPhysics_option6::G4EmDNAPhysics_option6(), G4EmDNAPhysics_option7::G4EmDNAPhysics_option7(), G4EmDNAPhysics_option8::G4EmDNAPhysics_option8(), G4EmDNAPhysics_stationary_option2::G4EmDNAPhysics_stationary_option2(), G4EmDNAPhysics_stationary_option4::G4EmDNAPhysics_stationary_option4(), G4EmDNAPhysics_stationary_option6::G4EmDNAPhysics_stationary_option6(), G4EmExtraPhysics::G4EmExtraPhysics(), G4EmLivermorePhysics::G4EmLivermorePhysics(), G4EmLowEPPhysics::G4EmLowEPPhysics(), G4EmPenelopePhysics::G4EmPenelopePhysics(), G4EmStandardPhysics::G4EmStandardPhysics(), G4EmStandardPhysics_option1::G4EmStandardPhysics_option1(), G4EmStandardPhysics_option2::G4EmStandardPhysics_option2(), G4EmStandardPhysics_option3::G4EmStandardPhysics_option3(), G4EmStandardPhysics_option4::G4EmStandardPhysics_option4(), G4EmStandardPhysicsGS::G4EmStandardPhysicsGS(), G4EmStandardPhysicsSS::G4EmStandardPhysicsSS(), G4EmStandardPhysicsWVI::G4EmStandardPhysicsWVI(), G4HadronElasticPhysics::G4HadronElasticPhysics(), G4HadronInelasticQBBC::G4HadronInelasticQBBC(), G4HadronPhysicsFTFP_BERT::G4HadronPhysicsFTFP_BERT(), G4HadronPhysicsQGSP_BERT::G4HadronPhysicsQGSP_BERT(), G4HadronPhysicsQGSP_BIC::G4HadronPhysicsQGSP_BIC(), G4IonINCLXXPhysics::G4IonINCLXXPhysics(), G4IonPhysics::G4IonPhysics(), G4IonPhysicsPHP::G4IonPhysicsPHP(), G4IonQMDPhysics::G4IonQMDPhysics(), G4NeutronTrackingCut::G4NeutronTrackingCut(), G4StepLimiterPhysics::G4StepLimiterPhysics(), G4StoppingPhysics::G4StoppingPhysics(), and G4StoppingPhysicsFritiofWithBinaryCascade::G4StoppingPhysicsFritiofWithBinaryCascade().

◆ SetVerboseLevel()

void G4VPhysicsConstructor::SetVerboseLevel ( G4int  value)
inlineinherited

◆ TerminateWorker()

void G4VPhysicsConstructor::TerminateWorker ( )
virtualinherited

Definition at line 105 of file G4VPhysicsConstructor.cc.

106{
107 if(subInstanceManager.offset[g4vpcInstanceID]._builders != nullptr)
108 {
109 std::for_each(subInstanceManager.offset[g4vpcInstanceID]._builders->begin(),
110 subInstanceManager.offset[g4vpcInstanceID]._builders->end(),
111 [](PhysicsBuilder_V::value_type bld) { delete bld; });
112 subInstanceManager.offset[g4vpcInstanceID]._builders->clear();
113 }
114}

References G4VPhysicsConstructor::g4vpcInstanceID, G4VUPLSplitter< T >::offset, and G4VPhysicsConstructor::subInstanceManager.

Referenced by G4VPhysicsConstructor::~G4VPhysicsConstructor().

Field Documentation

◆ g4vpcInstanceID

G4int G4VPhysicsConstructor::g4vpcInstanceID = 0
protectedinherited

◆ namePhysics

G4String G4VPhysicsConstructor::namePhysics = ""
protectedinherited

◆ subInstanceManager

G4VPCManager G4VPhysicsConstructor::subInstanceManager
staticprotectedinherited

◆ theParticleTable

G4ParticleTable* G4VPhysicsConstructor::theParticleTable = nullptr
protectedinherited

◆ typePhysics

G4int G4VPhysicsConstructor::typePhysics = 0
protectedinherited

◆ verboseLevel

G4int G4VPhysicsConstructor::verboseLevel = 0
protectedinherited

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