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

#include <PhysListEmStandardNR.hh>

Inheritance diagram for PhysListEmStandardNR:
G4VPhysicsConstructor

Public Member Functions

 PhysListEmStandardNR (const G4String &name="standardNR")
 
 ~PhysListEmStandardNR ()
 
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 42 of file PhysListEmStandardNR.hh.

Constructor & Destructor Documentation

PhysListEmStandardNR::PhysListEmStandardNR ( const G4String name = "standardNR")
PhysListEmStandardNR::~PhysListEmStandardNR ( )

Definition at line 92 of file PhysListEmStandardNR.cc.

93 {}

Member Function Documentation

virtual void PhysListEmStandardNR::ConstructParticle ( void  )
inlinevirtual

Implements G4VPhysicsConstructor.

Definition at line 50 of file PhysListEmStandardNR.hh.

50 {};
void PhysListEmStandardNR::ConstructProcess ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 97 of file PhysListEmStandardNR.cc.

References aParticleIterator, python.hepunit::eV, G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGCharge(), G4PhysicsListHelper::GetPhysicsListHelper(), python.hepunit::GeV, G4LossTableManager::Instance(), G4ParticleDefinition::IsShortLived(), python.hepunit::MeV, G4InuclParticleNames::mup, G4PhysicsListHelper::RegisterProcess(), G4VEmModel::SetActivationLowEnergyLimit(), G4LossTableManager::SetAtomDeexcitation(), G4EmProcessOptions::SetDEDXBinning(), G4VMultipleScattering::SetEmModel(), G4VEmProcess::SetEmModel(), G4VEnergyLossProcess::SetEmModel(), G4VAtomDeexcitation::SetFluo(), G4VEmModel::SetHighEnergyLimit(), G4EmProcessOptions::SetLambdaBinning(), G4EmProcessOptions::SetMaxEnergy(), G4ScreenedNuclearRecoil::SetMaxEnergyForScattering(), G4EmProcessOptions::SetMinEnergy(), G4EmProcessOptions::SetPolarAngleLimit(), G4VEnergyLossProcess::SetStepFunction(), and python.hepunit::TeV.

98 {
100 
101  // muon & hadron bremsstrahlung and pair production
104 
106  G4double energyLimit = 100.*MeV;
107  nucr->SetMaxEnergyForScattering(energyLimit);
109  csm->SetActivationLowEnergyLimit(energyLimit);
110 
111  aParticleIterator->reset();
112  while( (*aParticleIterator)() ){
113  G4ParticleDefinition* particle = aParticleIterator->value();
114  G4String particleName = particle->GetParticleName();
115 
116  if (particleName == "gamma") {
117 
118  // Compton scattering
120  cs->SetEmModel(new G4KleinNishinaModel(),1);
121  ph->RegisterProcess(cs, particle);
122 
123  // Photoelectric
125  G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
126  theLivermorePEModel->SetHighEnergyLimit(10*GeV);
127  pe->SetEmModel(theLivermorePEModel,1);
128  ph->RegisterProcess(pe, particle);
129 
130  // Gamma conversion
132  G4VEmModel* thePenelopeGCModel = new G4PenelopeGammaConversionModel();
133  thePenelopeGCModel->SetHighEnergyLimit(1*GeV);
134  gc->SetEmModel(thePenelopeGCModel,1);
135  ph->RegisterProcess(gc, particle);
136 
137  // Rayleigh scattering
138  ph->RegisterProcess(new G4RayleighScattering(), particle);
139 
140  } else if (particleName == "e-") {
141 
142  // ionisation
143  G4eIonisation* eIoni = new G4eIonisation();
144  eIoni->SetStepFunction(0.2, 100*um);
145 
146  // bremsstrahlung
147  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
148 
149  ph->RegisterProcess(new G4eMultipleScattering(), particle);
150  ph->RegisterProcess(eIoni, particle);
151  ph->RegisterProcess(eBrem, particle);
152 
153  } else if (particleName == "e+") {
154  // ionisation
155  G4eIonisation* eIoni = new G4eIonisation();
156  eIoni->SetStepFunction(0.2, 100*um);
157 
158  // bremsstrahlung
159  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
160 
161  ph->RegisterProcess(new G4eMultipleScattering(), particle);
162  ph->RegisterProcess(eIoni, particle);
163  ph->RegisterProcess(eBrem, particle);
164 
165  // annihilation at rest and in flight
166  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
167 
168  } else if (particleName == "mu+" ||
169  particleName == "mu-" ) {
170 
171  G4MuIonisation* muIoni = new G4MuIonisation();
172  muIoni->SetStepFunction(0.2, 50*um);
173 
174  ph->RegisterProcess(muIoni, particle);
175  ph->RegisterProcess(mub, particle);
176  ph->RegisterProcess(mup, particle);
177  ph->RegisterProcess(new G4CoulombScattering(), particle);
178 
179  } else if (particleName == "alpha" || particleName == "He3") {
180 
183  model->SetActivationLowEnergyLimit(energyLimit);
184  msc->SetEmModel(model, 1);
185  ph->RegisterProcess(msc, particle);
186 
187  G4ionIonisation* ionIoni = new G4ionIonisation();
188  ionIoni->SetStepFunction(0.1, 10*um);
189  ph->RegisterProcess(ionIoni, particle);
190 
191  ph->RegisterProcess(nucr, particle);
192 
193  } else if (particleName == "GenericIon" ) {
194 
196  G4UrbanMscModel* model = new G4UrbanMscModel();
197  model->SetActivationLowEnergyLimit(energyLimit);
198  msc->SetEmModel(model, 1);
199  ph->RegisterProcess(msc, particle);
200 
201  G4ionIonisation* ionIoni = new G4ionIonisation();
202  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
203  ionIoni->SetStepFunction(0.1, 1*um);
204  ph->RegisterProcess(ionIoni, particle);
205 
206  ph->RegisterProcess(nucr, particle);
207 
208  } else if (particleName == "proton" ||
209  particleName == "deuteron" ||
210  particleName == "triton") {
211 
213  G4UrbanMscModel* model = new G4UrbanMscModel();
214  model->SetActivationLowEnergyLimit(energyLimit);
215  msc->SetEmModel(model, 1);
216  ph->RegisterProcess(msc, particle);
217 
218  G4hIonisation* hIoni = new G4hIonisation();
219  hIoni->SetStepFunction(0.05, 1*um);
220  ph->RegisterProcess(hIoni, particle);
221 
222  ph->RegisterProcess(nucr, particle);
223 
224  } else if ((!particle->IsShortLived()) &&
225  (particle->GetPDGCharge() != 0.0) &&
226  (particle->GetParticleName() != "chargedgeantino")) {
227  //all others charged particles except geantino
228 
229  ph->RegisterProcess(new G4hMultipleScattering(), particle);
230  ph->RegisterProcess(new G4hIonisation(), particle);
231  }
232  }
233 
234  // Em options
235  //
236  // Main options and setting parameters are shown here.
237  // Several of them have default values.
238  //
239  G4EmProcessOptions emOptions;
240 
241  //physics tables
242  //
243  emOptions.SetMinEnergy(10*eV);
244  emOptions.SetMaxEnergy(10*TeV);
245  emOptions.SetDEDXBinning(12*20);
246  emOptions.SetLambdaBinning(12*20);
247 
248  // scattering
249  emOptions.SetPolarAngleLimit(0.0);
250 
251  // Deexcitation
254  de->SetFluo(true);
255 }
static G4LossTableManager * Instance()
void SetMinEnergy(G4double val)
void SetStepFunction(G4double v1, G4double v2)
void SetEmModel(G4VMscModel *, G4int index=1)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
void SetDEDXBinning(G4int val)
void SetEmModel(G4VEmModel *, G4int index=1)
void SetLambdaBinning(G4int val)
#define aParticleIterator
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
const XML_Char XML_Content * model
void SetMaxEnergyForScattering(G4double energy)
set the upper energy beyond which this process has no cross section
void SetMaxEnergy(G4double val)
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:704
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetEmModel(G4VEmModel *, G4int index=1)
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
A process which handles screened Coulomb collisions between nuclei.
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetPolarAngleLimit(G4double val)

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