Geant4  10.02.p03
PhysListEmStandardSS Class Reference

#include <PhysListEmStandardSS.hh>

Inheritance diagram for PhysListEmStandardSS:
Collaboration diagram for PhysListEmStandardSS:

Public Member Functions

 PhysListEmStandardSS (const G4String &name="standardSS")
 
 ~PhysListEmStandardSS ()
 
virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
 PhysListEmStandardSS (const G4String &name="standardSS")
 
 ~PhysListEmStandardSS ()
 
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)
 
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

Constructor & Destructor Documentation

◆ PhysListEmStandardSS() [1/2]

PhysListEmStandardSS::PhysListEmStandardSS ( const G4String name = "standardSS")

Definition at line 80 of file electromagnetic/TestEm7/src/PhysListEmStandardSS.cc.

81  : G4VPhysicsConstructor(name)
82 {
85 }
static G4LossTableManager * Instance()
G4VPhysicsConstructor(const G4String &="")
Here is the call graph for this function:

◆ ~PhysListEmStandardSS() [1/2]

PhysListEmStandardSS::~PhysListEmStandardSS ( )

Definition at line 89 of file electromagnetic/TestEm7/src/PhysListEmStandardSS.cc.

90 {}

◆ PhysListEmStandardSS() [2/2]

PhysListEmStandardSS::PhysListEmStandardSS ( const G4String name = "standardSS")

◆ ~PhysListEmStandardSS() [2/2]

PhysListEmStandardSS::~PhysListEmStandardSS ( )

Member Function Documentation

◆ ConstructParticle() [1/2]

virtual void PhysListEmStandardSS::ConstructParticle ( void  )
inlinevirtual

Implements G4VPhysicsConstructor.

Definition at line 50 of file electromagnetic/TestEm7/include/PhysListEmStandardSS.hh.

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

◆ ConstructParticle() [2/2]

virtual void PhysListEmStandardSS::ConstructParticle ( void  )
inlinevirtual

Implements G4VPhysicsConstructor.

Definition at line 50 of file medical/electronScattering/include/PhysListEmStandardSS.hh.

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

◆ ConstructProcess() [1/2]

virtual void PhysListEmStandardSS::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

◆ ConstructProcess() [2/2]

void PhysListEmStandardSS::ConstructProcess ( void  )
virtual

Implements G4VPhysicsConstructor.

Definition at line 94 of file electromagnetic/TestEm7/src/PhysListEmStandardSS.cc.

95 {
97 
98  // muon & hadron bremsstrahlung and pair production
101 
103  particleIterator->reset();
104  while( (*particleIterator)() ){
105  G4ParticleDefinition* particle = particleIterator->value();
106  G4String particleName = particle->GetParticleName();
107 
108  if (particleName == "gamma") {
109 
110  // Compton scattering
112  cs->SetEmModel(new G4KleinNishinaModel(),1);
113  G4VEmModel* theLowEPComptonModel = new G4LowEPComptonModel();
114  theLowEPComptonModel->SetHighEnergyLimit(20*MeV);
115  cs->AddEmModel(0, theLowEPComptonModel);
116  ph->RegisterProcess(cs, particle);
117 
118  // Photoelectric
120  G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
121  theLivermorePEModel->SetHighEnergyLimit(10*GeV);
122  pe->SetEmModel(theLivermorePEModel,1);
123  ph->RegisterProcess(pe, particle);
124 
125  // Gamma conversion
127  G4VEmModel* thePenelopeGCModel = new G4PenelopeGammaConversionModel();
128  thePenelopeGCModel->SetHighEnergyLimit(1*GeV);
129  gc->SetEmModel(thePenelopeGCModel,1);
130  ph->RegisterProcess(gc, particle);
131 
132  // Rayleigh scattering
133  ph->RegisterProcess(new G4RayleighScattering(), particle);
134 
135  } else if (particleName == "e-") {
136 
137  // ionisation
138  G4eIonisation* eIoni = new G4eIonisation();
139  eIoni->SetStepFunction(0.2, 100*um);
140  G4VEmModel* theIoniPenelope = new G4PenelopeIonisationModel();
141  theIoniPenelope->SetHighEnergyLimit(0.1*MeV);
142  eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
143 
144  // bremsstrahlung
145  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
146 
147  ph->RegisterProcess(eIoni, particle);
148  ph->RegisterProcess(eBrem, particle);
149  ph->RegisterProcess(new G4CoulombScattering(), particle);
150 
151  } else if (particleName == "e+") {
152  // ionisation
153  G4eIonisation* eIoni = new G4eIonisation();
154  eIoni->SetStepFunction(0.2, 100*um);
155  G4VEmModel* theIoniPenelope = new G4PenelopeIonisationModel();
156  theIoniPenelope->SetHighEnergyLimit(0.1*MeV);
157  eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
158 
159  // bremsstrahlung
160  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
161 
162  ph->RegisterProcess(eIoni, particle);
163  ph->RegisterProcess(eBrem, particle);
164  ph->RegisterProcess(new G4CoulombScattering(), particle);
165 
166  // annihilation at rest and in flight
167  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
168 
169  } else if (particleName == "mu+" ||
170  particleName == "mu-" ) {
171 
172  G4MuIonisation* muIoni = new G4MuIonisation();
173  muIoni->SetStepFunction(0.2, 50*um);
174 
175  ph->RegisterProcess(muIoni, particle);
176  ph->RegisterProcess(mub, particle);
177  ph->RegisterProcess(mup, particle);
178  ph->RegisterProcess(new G4CoulombScattering(), particle);
179 
180  } else if (particleName == "alpha" || particleName == "He3") {
181 
182  G4ionIonisation* ionIoni = new G4ionIonisation();
183  ionIoni->SetStepFunction(0.1, 10*um);
184  ph->RegisterProcess(ionIoni, particle);
185 
187  cs->SetBuildTableFlag(false);
188  ph->RegisterProcess(cs, particle);
189 
190  } else if (particleName == "GenericIon" ) {
191 
192  G4ionIonisation* ionIoni = new G4ionIonisation();
193  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
194  ionIoni->SetStepFunction(0.1, 1*um);
195  ph->RegisterProcess(ionIoni, particle);
196 
198  cs->SetBuildTableFlag(false);
199  ph->RegisterProcess(cs, particle);
200 
201  } else if ((!particle->IsShortLived()) &&
202  (particle->GetPDGCharge() != 0.0) &&
203  (particle->GetParticleName() != "chargedgeantino")) {
204  //all others charged particles except geantino
205 
206  ph->RegisterProcess(new G4hIonisation(), particle);
208  cs->SetBuildTableFlag(false);
209  ph->RegisterProcess(cs, particle);
210 
211  }
212  }
213 
214  // Em options
215  //
216  // Main options and setting parameters are shown here.
217  // Several of them have default values.
218  //
219  G4EmProcessOptions emOptions;
220 
221  //physics tables
222  //
223  emOptions.SetMinEnergy(10*eV);
224  emOptions.SetMaxEnergy(10*TeV);
225  emOptions.SetDEDXBinning(12*20);
226  emOptions.SetLambdaBinning(12*20);
227 
228  // scattering
229  emOptions.SetPolarAngleLimit(0.0);
230 
231  // Deexcitation
234  de->SetFluo(true);
235 }
static const double MeV
Definition: G4SIunits.hh:211
static G4LossTableManager * Instance()
void SetBuildTableFlag(G4bool val)
void SetMinEnergy(G4double val)
void SetStepFunction(G4double v1, G4double v2)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
void SetDEDXBinning(G4int val)
void SetEmModel(G4VEmModel *, G4int index=1)
const G4String & GetParticleName() const
void SetLambdaBinning(G4int val)
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static const double GeV
Definition: G4SIunits.hh:214
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
void SetMaxEnergy(G4double val)
static const double eV
Definition: G4SIunits.hh:212
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
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
static const double TeV
Definition: G4SIunits.hh:215
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
void SetAtomDeexcitation(G4VAtomDeexcitation *)
G4double GetPDGCharge() const
void SetPolarAngleLimit(G4double val)
Here is the call graph for this function:
Here is the caller graph for this function:

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