#include <G4VEmModel.hh>
Inherited by G4BetheBlochModel, G4BetheHeitlerModel, G4BoldyshevTripletModel, G4BraggIonModel, G4BraggModel, G4DNABornExcitationModel1, G4DNABornExcitationModel2, G4DNABornIonisationModel1, G4DNABornIonisationModel2, G4DNAChampionElasticModel, G4DNADingfelderChargeDecreaseModel, G4DNADingfelderChargeIncreaseModel, G4DNAEmfietzoglouExcitationModel, G4DNAEmfietzoglouIonisationModel, G4DNAIonElasticModel, G4DNAMeltonAttachmentModel, G4DNAMillerGreenExcitationModel, G4DNARuddIonisationExtendedModel, G4DNARuddIonisationModel, G4DNASancheExcitationModel, G4DNAScreenedRutherfordElasticModel, G4DNATransformElectronModel, G4DNAUeharaScreenedRutherfordElasticModel, G4eBremParametrizedModel, G4eBremsstrahlungRelModel, G4eCoulombScatteringModel, G4eeToHadronsModel, G4eeToHadronsMultiModel, G4eeToTwoGammaModel, G4EmMultiModel, G4eSingleCoulombScatteringModel, G4hCoulombScatteringModel, G4ICRU49NuclearStoppingModel, G4ICRU73QOModel, G4IonCoulombScatteringModel, G4IonParametrisedLossModel, G4KleinNishinaCompton, G4KleinNishinaModel, G4LivermoreComptonModel, G4LivermoreComptonModifiedModel, G4LivermoreGammaConversionModel, G4LivermoreGammaConversionModelRC, G4LivermoreIonisationModel, G4LivermoreNuclearGammaConversionModel, G4LivermorePhotoElectricModel, G4LivermorePolarizedComptonModel, G4LivermorePolarizedGammaConversionModel, G4LivermorePolarizedPhotoElectricGDModel, G4LivermorePolarizedPhotoElectricModel, G4LivermorePolarizedRayleighModel, G4LivermoreRayleighModel, G4LowEPComptonModel, G4LowEPPolarizedComptonModel, G4MicroElecElasticModel, G4MicroElecInelasticModel, G4MollerBhabhaModel, G4mplIonisationModel, G4mplIonisationWithDeltaModel, G4MuBetheBlochModel, G4MuBremsstrahlungModel, G4MuElecElasticModel, G4MuElecInelasticModel, G4MuPairProductionModel, G4PAIModel, G4PAIPhotModel, G4PairProductionRelModel, G4PEEffectFluoModel, G4PenelopeAnnihilationModel, G4PenelopeBremsstrahlungModel, G4PenelopeComptonModel, G4PenelopeGammaConversionModel, G4PenelopeIonisationModel, G4PenelopePhotoElectricModel, G4PenelopeRayleighModel, G4TDNAOneStepThermalizationModel< MODEL >, G4VDNAPTBModel, G4VLEPTSModel, G4VMscModel, and G4XrayRayleighModel.
|
| G4VEmModel (const G4String &nam) |
|
virtual | ~G4VEmModel () |
|
virtual void | Initialise (const G4ParticleDefinition *, const G4DataVector &)=0 |
|
virtual void | SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0 |
|
virtual void | InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) |
|
virtual void | InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *) |
|
virtual void | InitialiseForElement (const G4ParticleDefinition *, G4int Z) |
|
virtual G4double | ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX) |
|
virtual G4double | CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) |
|
virtual G4double | GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy) |
|
virtual G4double | ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) |
|
virtual G4double | ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) |
|
virtual G4double | ChargeSquareRatio (const G4Track &) |
|
virtual G4double | GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) |
|
virtual G4double | GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) |
|
virtual void | StartTracking (G4Track *) |
|
virtual void | CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length) |
|
virtual G4double | Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy) |
|
virtual G4double | MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0) |
|
virtual G4double | MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *) |
|
virtual void | SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy) |
|
virtual void | DefineForRegion (const G4Region *) |
|
virtual void | ModelDescription (std::ostream &outFile) const |
|
void | InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &) |
|
std::vector
< G4EmElementSelector * > * | GetElementSelectors () |
|
void | SetElementSelectors (std::vector< G4EmElementSelector * > *) |
|
virtual G4double | ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX) |
|
G4double | CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) |
|
G4double | ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) |
|
G4double | ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) |
|
G4int | SelectIsotopeNumber (const G4Element *) |
|
const G4Element * | SelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) |
|
const G4Element * | SelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) |
|
G4int | SelectRandomAtomNumber (const G4Material *) |
|
void | SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr) |
|
void | SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal) |
|
G4ElementData * | GetElementData () |
|
G4PhysicsTable * | GetCrossSectionTable () |
|
G4VEmFluctuationModel * | GetModelOfFluctuations () |
|
G4VEmAngularDistribution * | GetAngularDistribution () |
|
void | SetAngularDistribution (G4VEmAngularDistribution *) |
|
G4double | HighEnergyLimit () const |
|
G4double | LowEnergyLimit () const |
|
G4double | HighEnergyActivationLimit () const |
|
G4double | LowEnergyActivationLimit () const |
|
G4double | PolarAngleLimit () const |
|
G4double | SecondaryThreshold () const |
|
G4bool | LPMFlag () const |
|
G4bool | DeexcitationFlag () const |
|
G4bool | ForceBuildTableFlag () const |
|
G4bool | UseAngularGeneratorFlag () const |
|
void | SetAngularGeneratorFlag (G4bool) |
|
void | SetHighEnergyLimit (G4double) |
|
void | SetLowEnergyLimit (G4double) |
|
void | SetActivationHighEnergyLimit (G4double) |
|
void | SetActivationLowEnergyLimit (G4double) |
|
G4bool | IsActive (G4double kinEnergy) |
|
void | SetPolarAngleLimit (G4double) |
|
void | SetSecondaryThreshold (G4double) |
|
void | SetLPMFlag (G4bool val) |
|
void | SetDeexcitationFlag (G4bool val) |
|
void | SetForceBuildTable (G4bool val) |
|
void | SetMasterThread (G4bool val) |
|
G4bool | IsMaster () const |
|
G4double | MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle) |
|
const G4String & | GetName () const |
|
void | SetCurrentCouple (const G4MaterialCutsCouple *) |
|
const G4Element * | GetCurrentElement () const |
|
const G4Isotope * | GetCurrentIsotope () const |
|
G4bool | IsLocked () const |
|
void | SetLocked (G4bool) |
|
Definition at line 104 of file G4VEmModel.hh.
G4VEmModel::G4VEmModel |
( |
const G4String & |
nam | ) |
|
|
explicit |
Definition at line 68 of file G4VEmModel.cc.
69 flucModel(
nullptr),anglModel(
nullptr),
name(nam), lowLimit(0.1*
CLHEP::keV),
72 theLPMflag(
false),flagDeexcitation(
false),flagForceBuildTable(
false),
75 fCurrentElement(
nullptr),fCurrentIsotope(
nullptr),nsec(5)
79 elmSelectors =
nullptr;
80 localElmSelectors =
true;
82 useAngularGenerator =
false;
static G4LossTableManager * Instance()
static constexpr double keV
const std::vector< G4int > * theDensityIdx
void Register(G4VEnergyLossProcess *p)
const std::vector< G4double > * theDensityFactor
G4VParticleChange * pParticleChange
static constexpr double TeV
G4PhysicsTable * xSectionTable
G4ElementData * fElementData
static constexpr double pi
G4VEmModel::~G4VEmModel |
( |
| ) |
|
|
virtual |
Definition at line 92 of file G4VEmModel.cc.
94 if(localElmSelectors) {
96 for(
G4int i=0; i<nSelectors; ++i) {
97 delete (*elmSelectors)[i];
void DeRegister(G4VEnergyLossProcess *p)
G4PhysicsTable * xSectionTable
G4ElementData * fElementData
Reimplemented in G4BraggIonGasModel, and G4BetheBlochIonGasModel.
Definition at line 345 of file G4VEmModel.cc.
G4double GetKineticEnergy() const
const G4ParticleDefinition * GetParticleDefinition() const
G4Material * GetMaterial() const
virtual G4double GetChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Reimplemented in G4IonParametrisedLossModel, G4LowEPPolarizedComptonModel, G4LowEPComptonModel, G4BetheBlochModel, G4eCoulombScatteringModel, G4MuBetheBlochModel, G4MuBremsstrahlungModel, G4BraggIonModel, G4BraggModel, G4ICRU73QOModel, G4UrbanAdjointMscModel, G4MuPairProductionModel, G4UrbanMscModel, G4PenelopeComptonModel, G4eBremParametrizedModel, G4WentzelVIModel, G4eeToHadronsMultiModel, G4PenelopeIonisationModel, G4IonCoulombScatteringModel, G4WentzelVIRelModel, G4mplIonisationWithDeltaModel, G4eBremsstrahlungRelModel, G4hCoulombScatteringModel, G4MollerBhabhaModel, G4eSingleCoulombScatteringModel, G4EmMultiModel, G4eeToTwoGammaModel, G4eeToHadronsModel, G4PenelopeBremsstrahlungModel, G4PairProductionRelModel, G4BetheHeitlerModel, G4KleinNishinaCompton, G4KleinNishinaModel, G4PenelopePhotoElectricModel, G4LivermoreIonisationModel, G4PenelopeRayleighModel, G4PEEffectFluoModel, G4LivermorePolarizedComptonModel, G4PenelopeAnnihilationModel, G4PenelopeGammaConversionModel, G4PolarizedComptonModel, G4LivermorePhotoElectricModel, G4LivermorePolarizedPhotoElectricModel, G4LivermorePolarizedRayleighModel, G4LivermoreComptonModel, G4LivermoreGammaConversionModel, G4LivermoreComptonModifiedModel, G4LivermorePolarizedGammaConversionModel, G4LivermoreNuclearGammaConversionModel, G4BoldyshevTripletModel, G4LivermoreGammaConversionModelRC, G4LivermorePolarizedPhotoElectricGDModel, G4LivermoreRayleighModel, and G4XrayRayleighModel.
Definition at line 321 of file G4VEmModel.cc.
Definition at line 527 of file G4VEmModel.hh.
535 cutEnergy,maxEnergy);
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void SetCurrentElement(const G4Element *)
Reimplemented in G4EmMultiModel.
Definition at line 489 of file G4VEmModel.hh.
virtual G4double ComputeDEDXPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
void SetCurrentCouple(const G4MaterialCutsCouple *)
const G4Material * GetMaterial() const
Reimplemented in G4IonParametrisedLossModel, G4BetheBlochModel, G4MuBetheBlochModel, G4PenelopeIonisationModel, G4BraggIonModel, G4BraggModel, G4ICRU73QOModel, G4PenelopeBremsstrahlungModel, G4MuBremsstrahlungModel, G4MuPairProductionModel, G4MollerBhabhaModel, G4PAIModel, G4LivermoreIonisationModel, G4PAIPhotModel, G4eBremParametrizedModel, G4eBremsstrahlungRelModel, G4ICRU49NuclearStoppingModel, G4mplIonisationModel, G4mplIonisationWithDeltaModel, G4BetheBlochNoDeltaModel, G4BraggNoDeltaModel, and G4ICRU73NoDeltaModel.
Definition at line 248 of file G4VEmModel.cc.
Definition at line 514 of file G4VEmModel.hh.
521 return cross > 0.0 ? 1./cross :
DBL_MAX;
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
static const G4double emax
Definition at line 500 of file G4VEmModel.hh.
508 cutEnergy,maxEnergy);
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void SetCurrentCouple(const G4MaterialCutsCouple *)
const G4Material * GetMaterial() const
Reimplemented in G4IonParametrisedLossModel, G4TDNAOneStepThermalizationModel< MODEL >, G4BetheBlochModel, G4MuBetheBlochModel, G4BraggIonModel, G4BraggModel, G4ICRU73QOModel, G4PenelopeIonisationModel, G4PAIModel, G4MollerBhabhaModel, G4eeToTwoGammaModel, G4PenelopeBremsstrahlungModel, G4PAIPhotModel, G4eeToHadronsMultiModel, G4PenelopeComptonModel, G4PEEffectFluoModel, G4MicroElecInelasticModel, G4MuElecInelasticModel, G4eeToHadronsModel, G4DNATransformElectronModel, G4BetheBlochNoDeltaModel, G4BraggNoDeltaModel, G4ICRU73NoDeltaModel, G4DNAEmfietzoglouIonisationModel, G4VDNAPTBModel, G4DNAEmfietzoglouExcitationModel, G4MicroElecElasticModel, G4MuElecElasticModel, G4DNAIonElasticModel, G4LivermorePhotoElectricModel, G4LivermorePolarizedPhotoElectricModel, G4DNABornIonisationModel1, G4DNABornIonisationModel2, G4DNARuddIonisationExtendedModel, G4DNARuddIonisationModel, G4DNABornExcitationModel2, G4DNAMillerGreenExcitationModel, G4DNASancheExcitationModel, G4DNABornExcitationModel1, G4DNAChampionElasticModel, G4DNAMeltonAttachmentModel, G4DNADingfelderChargeIncreaseModel, G4DNADingfelderChargeDecreaseModel, G4DNAScreenedRutherfordElasticModel, G4LEPTSAttachmentModel, G4LEPTSExcitationModel, G4LEPTSDissociationModel, G4LEPTSElasticModel, G4LEPTSIonisationModel, G4LEPTSPositroniumModel, G4LEPTSRotExcitationModel, G4LEPTSVibExcitationModel, and G4DNAUeharaScreenedRutherfordElasticModel.
Definition at line 257 of file G4VEmModel.cc.
266 const G4double* theAtomNumDensityVector =
273 for (
G4int i=0; i<nelm; ++i) {
274 cross += theAtomNumDensityVector[i]*
std::vector< G4Element * > G4ElementVector
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
const G4ElementVector * GetElementVector() const
const G4double * GetVecNbOfAtomsPerVolume() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
static const G4double emax
size_t GetNumberOfElements() const
G4bool G4VEmModel::DeexcitationFlag |
( |
| ) |
const |
|
inline |
G4bool G4VEmModel::ForceBuildTableFlag |
( |
| ) |
const |
|
inline |
const G4Element * G4VEmModel::GetCurrentElement |
( |
| ) |
const |
|
inline |
const G4Isotope * G4VEmModel::GetCurrentIsotope |
( |
| ) |
const |
|
inline |
const G4String & G4VEmModel::GetName |
( |
| ) |
const |
|
inline |
Definition at line 132 of file G4VEmModel.cc.
G4VParticleChange * pParticleChange
Definition at line 118 of file G4VEmModel.cc.
G4VParticleChange * pParticleChange
G4double G4VEmModel::HighEnergyActivationLimit |
( |
| ) |
const |
|
inline |
G4double G4VEmModel::HighEnergyLimit |
( |
| ) |
const |
|
inline |
Implemented in G4GoudsmitSaundersonMscModel, G4TDNAOneStepThermalizationModel< MODEL >, G4IonParametrisedLossModel, G4LowEPPolarizedComptonModel, G4LowEPComptonModel, G4eCoulombScatteringModel, G4UrbanAdjointMscModel, G4MuBremsstrahlungModel, G4MuPairProductionModel, G4BetheBlochModel, G4BraggModel, G4ICRU73QOModel, G4UrbanMscModel, G4IonCoulombScatteringModel, G4MuBetheBlochModel, G4PAIModel, G4BraggIonModel, G4WentzelVIRelModel, G4PenelopeIonisationModel, G4hCoulombScatteringModel, G4WentzelVIModel, G4eSingleCoulombScatteringModel, G4eeToHadronsMultiModel, G4MicroElecInelasticModel, G4MuElecInelasticModel, G4PAIPhotModel, G4PenelopeBremsstrahlungModel, G4PenelopeComptonModel, G4PolarizedAnnihilationModel, G4MollerBhabhaModel, G4PolarizedPEEffectModel, G4eBremsstrahlungRelModel, G4ICRU49NuclearStoppingModel, G4LivermoreBremsstrahlungModel, G4eBremParametrizedModel, G4PairProductionRelModel, G4SeltzerBergerModel, G4BetheHeitlerModel, G4eeToTwoGammaModel, G4EmMultiModel, G4eeToHadronsModel, G4PenelopePhotoElectricModel, G4mplIonisationModel, G4mplIonisationWithDeltaModel, G4LivermoreIonisationModel, G4ePolarizedBremsstrahlungModel, G4PolarizedGammaConversionModel, G4KleinNishinaCompton, G4KleinNishinaModel, G4PEEffectFluoModel, G4DNATransformElectronModel, G4PenelopeGammaConversionModel, G4PenelopeRayleighModel, G4DNAEmfietzoglouIonisationModel, G4PenelopeAnnihilationModel, G4LivermorePolarizedComptonModel, G4DummyModel, G4DNAEmfietzoglouExcitationModel, G4MicroElecElasticModel, G4MuElecElasticModel, G4LivermorePolarizedRayleighModel, G4LivermorePhotoElectricModel, G4LivermorePolarizedPhotoElectricModel, G4DNABornIonisationModel1, G4DNABornIonisationModel2, G4DNAIonElasticModel, G4LivermoreComptonModifiedModel, G4LivermoreComptonModel, G4DNARuddIonisationExtendedModel, G4DNARuddIonisationModel, G4DNAMillerGreenExcitationModel, G4DNABornExcitationModel2, G4DNASancheExcitationModel, G4LivermorePolarizedGammaConversionModel, G4LivermorePolarizedPhotoElectricGDModel, G4DNAMeltonAttachmentModel, G4XrayRayleighModel, G4DNABornExcitationModel1, G4DNAChampionElasticModel, G4DNADingfelderChargeIncreaseModel, G4LivermoreGammaConversionModel, G4DNADingfelderChargeDecreaseModel, G4LivermoreGammaConversionModelRC, G4LivermoreNuclearGammaConversionModel, G4LivermoreRayleighModel, G4DNAScreenedRutherfordElasticModel, G4VDNAPTBModel, G4BoldyshevTripletModel, G4DNAUeharaScreenedRutherfordElasticModel, G4LEPTSAttachmentModel, G4LEPTSDissociationModel, G4LEPTSElasticModel, G4LEPTSExcitationModel, G4LEPTSIonisationModel, G4LEPTSPositroniumModel, G4LEPTSRotExcitationModel, and G4LEPTSVibExcitationModel.
Definition at line 146 of file G4VEmModel.cc.
159 if(highLimit <= lowLimit) {
return; }
169 elmSelectors =
new std::vector<G4EmElementSelector*>;
171 if(numOfCouples > nSelectors) {
172 for(
G4int i=nSelectors; i<numOfCouples; ++i) {
173 elmSelectors->push_back(
nullptr);
175 nSelectors = numOfCouples;
179 for(
G4int i=0; i<numOfCouples; ++i) {
182 if(cuts[i] ==
DBL_MAX) {
continue; }
189 if((*elmSelectors)[i]) {
190 if(material == ((*elmSelectors)[i])->GetMaterial()) { create =
false; }
191 else {
delete (*elmSelectors)[i]; }
204 ((*elmSelectors)[i])->
Initialise(part, cuts[i]);
G4int NumberOfBinsPerDecade() const
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
size_t GetTableSize() const
static const G4double emax
G4double G4Log(G4double x)
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static G4EmParameters * Instance()
const G4Material * GetMaterial() const
Reimplemented in G4LowEPPolarizedComptonModel, G4LowEPComptonModel, G4LivermorePhotoElectricModel, G4LivermorePolarizedPhotoElectricModel, G4LivermoreBremsstrahlungModel, G4SeltzerBergerModel, G4LivermorePolarizedPhotoElectricGDModel, G4LivermorePolarizedComptonModel, G4LivermorePolarizedRayleighModel, G4LivermoreComptonModel, G4LivermoreGammaConversionModel, G4LivermorePolarizedGammaConversionModel, G4LivermoreGammaConversionModelRC, G4LivermoreNuclearGammaConversionModel, G4BoldyshevTripletModel, and G4LivermoreRayleighModel.
Definition at line 243 of file G4VEmModel.cc.
Definition at line 224 of file G4VEmModel.cc.
230 for(
G4int i=0; i<
n; ++i) {
231 G4int Z = ((*theElementVector)[i])->GetZasInt();
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
size_t GetNumberOfElements() const
Reimplemented in G4LowEPPolarizedComptonModel, G4LowEPComptonModel, G4eCoulombScatteringModel, G4MuBremsstrahlungModel, G4MuPairProductionModel, G4PAIModel, G4hCoulombScatteringModel, G4WentzelVIModel, G4eSingleCoulombScatteringModel, G4PenelopeIonisationModel, G4PAIPhotModel, G4PenelopeBremsstrahlungModel, G4PenelopeComptonModel, G4eBremsstrahlungRelModel, G4PairProductionRelModel, G4BetheHeitlerModel, G4eBremParametrizedModel, G4KleinNishinaCompton, G4KleinNishinaModel, G4PenelopePhotoElectricModel, G4PenelopeGammaConversionModel, G4PenelopeRayleighModel, G4PenelopeAnnihilationModel, G4LivermorePolarizedComptonModel, G4LivermorePolarizedRayleighModel, G4LivermoreComptonModel, G4LivermoreGammaConversionModel, G4LivermorePolarizedGammaConversionModel, G4LivermoreGammaConversionModelRC, G4LivermoreNuclearGammaConversionModel, G4BoldyshevTripletModel, and G4LivermoreRayleighModel.
Definition at line 218 of file G4VEmModel.cc.
Definition at line 752 of file G4VEmModel.hh.
754 return (kinEnergy >= eMinActive && kinEnergy <= eMaxActive);
G4bool G4VEmModel::IsLocked |
( |
| ) |
const |
|
inline |
G4bool G4VEmModel::IsMaster |
( |
| ) |
const |
|
inline |
G4double G4VEmModel::LowEnergyActivationLimit |
( |
| ) |
const |
|
inline |
G4double G4VEmModel::LowEnergyLimit |
( |
| ) |
const |
|
inline |
G4bool G4VEmModel::LPMFlag |
( |
| ) |
const |
|
inline |
Definition at line 481 of file G4VEmModel.hh.
484 dynPart->GetKineticEnergy());
virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition *, G4double kineticEnergy)
void G4VEmModel::ModelDescription |
( |
std::ostream & |
outFile | ) |
const |
|
virtual |
Definition at line 440 of file G4VEmModel.cc.
442 outFile <<
"The description for this model has not been written yet.\n";
G4double G4VEmModel::PolarAngleLimit |
( |
| ) |
const |
|
inline |
Implemented in G4IonParametrisedLossModel, G4BetheBlochModel, G4TDNAOneStepThermalizationModel< MODEL >, G4LowEPPolarizedComptonModel, G4MuBetheBlochModel, G4BraggIonModel, G4BraggModel, G4ICRU73QOModel, G4LowEPComptonModel, G4MuBremsstrahlungModel, G4eCoulombScatteringModel, G4MuPairProductionModel, G4MollerBhabhaModel, G4PenelopeIonisationModel, G4PAIModel, G4PenelopeBremsstrahlungModel, G4PenelopeComptonModel, G4eeToTwoGammaModel, G4eeToHadronsModel, G4IonCoulombScatteringModel, G4PAIPhotModel, G4eeToHadronsMultiModel, G4eBremParametrizedModel, G4hCoulombScatteringModel, G4mplIonisationWithDeltaModel, G4eSingleCoulombScatteringModel, G4eBremsstrahlungRelModel, G4EmMultiModel, G4VMscModel, G4PolarizedAnnihilationModel, G4PairProductionRelModel, G4BetheHeitlerModel, G4PEEffectFluoModel, G4MicroElecInelasticModel, G4MuElecInelasticModel, G4KleinNishinaCompton, G4KleinNishinaModel, G4PenelopePhotoElectricModel, G4ICRU49NuclearStoppingModel, G4PolarizedComptonModel, G4LivermoreIonisationModel, G4PenelopeRayleighModel, G4LivermorePolarizedComptonModel, G4PenelopeGammaConversionModel, G4mplIonisationModel, G4PenelopeAnnihilationModel, G4PolarizedPEEffectModel, G4DNATransformElectronModel, G4LivermoreBremsstrahlungModel, G4LivermorePhotoElectricModel, G4LivermorePolarizedPhotoElectricModel, G4LivermorePolarizedRayleighModel, G4PolarizedMollerBhabhaModel, G4SeltzerBergerModel, G4DNAEmfietzoglouIonisationModel, G4VDNAPTBModel, G4LivermoreComptonModel, G4DNAEmfietzoglouExcitationModel, G4ePolarizedBremsstrahlungModel, G4PolarizedGammaConversionModel, G4LivermoreGammaConversionModel, G4MicroElecElasticModel, G4MuElecElasticModel, G4LivermorePolarizedGammaConversionModel, G4DNABornExcitationModel2, G4DNAIonElasticModel, G4DNAMillerGreenExcitationModel, G4LivermoreComptonModifiedModel, G4HeatedKleinNishinaCompton, G4DNABornIonisationModel1, G4DNABornIonisationModel2, G4LivermoreNuclearGammaConversionModel, G4DNABornExcitationModel1, G4BoldyshevTripletModel, G4LivermoreGammaConversionModelRC, G4DummyModel, G4DNARuddIonisationExtendedModel, G4DNARuddIonisationModel, G4LivermorePolarizedPhotoElectricGDModel, G4LivermoreRayleighModel, G4DNASancheExcitationModel, G4XrayRayleighModel, G4DNAChampionElasticModel, G4DNAMeltonAttachmentModel, G4DNADingfelderChargeIncreaseModel, G4DNADingfelderChargeDecreaseModel, G4DNAScreenedRutherfordElasticModel, G4DNAUeharaScreenedRutherfordElasticModel, G4LEPTSAttachmentModel, G4LEPTSExcitationModel, G4LEPTSDissociationModel, G4LEPTSElasticModel, G4LEPTSIonisationModel, G4LEPTSPositroniumModel, G4LEPTSRotExcitationModel, and G4LEPTSVibExcitationModel.
G4double G4VEmModel::SecondaryThreshold |
( |
| ) |
const |
|
inline |
Definition at line 584 of file G4VEmModel.hh.
589 fCurrentIsotope =
nullptr;
595 for(; idx<ni; ++idx) {
597 if (x <= 0.0) {
break; }
599 if(idx >= ni) { idx = ni - 1; }
602 N = fCurrentIsotope->
GetN();
size_t GetNumberOfIsotopes() const
G4double * GetRelativeAbundanceVector() const
const G4Isotope * GetIsotope(G4int iso) const
void SetCurrentElement(const G4Element *)
Definition at line 541 of file G4VEmModel.hh.
547 fCurrentCouple = couple;
553 cutEnergy,maxEnergy);
555 fCurrentIsotope =
nullptr;
556 return fCurrentElement;
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
const G4Material * GetMaterial() const
Definition at line 297 of file G4VEmModel.cc.
305 fCurrentElement = (*theElementVector)[
n];
309 for(
G4int i=0; i<
n; ++i) {
311 fCurrentElement = (*theElementVector)[i];
316 return fCurrentElement;
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
size_t GetNumberOfElements() const
Definition at line 561 of file G4VEmModel.hh.
567 G4int Z = (*elmv)[0]->GetZasInt();
571 for(
size_t i=0; i<
nn; ++i) {
574 Z = (*elmv)[i]->GetZasInt();
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
const G4double * GetVecNbOfAtomsPerVolume() const
G4double GetTotNbOfAtomsPerVolume() const
size_t GetNumberOfElements() const
void G4VEmModel::SetActivationHighEnergyLimit |
( |
G4double |
val | ) |
|
|
inline |
void G4VEmModel::SetActivationLowEnergyLimit |
( |
G4double |
val | ) |
|
|
inline |
void G4VEmModel::SetAngularGeneratorFlag |
( |
G4bool |
val | ) |
|
|
inline |
Definition at line 426 of file G4VEmModel.cc.
435 localTable = isLocal;
G4PhysicsTable * xSectionTable
Definition at line 458 of file G4VEmModel.hh.
460 fCurrentElement = elm;
461 fCurrentIsotope =
nullptr;
void G4VEmModel::SetDeexcitationFlag |
( |
G4bool |
val | ) |
|
|
inline |
Definition at line 809 of file G4VEmModel.hh.
812 if(elmSelectors) { nSelectors = elmSelectors->size(); }
813 localElmSelectors =
false;
Definition at line 418 of file G4VEmModel.cc.
G4VParticleChange * pParticleChange
Definition at line 759 of file G4VEmModel.hh.
761 if(!isLocked) { polarAngleLimit = val; }
G4bool G4VEmModel::UseAngularGeneratorFlag |
( |
| ) |
const |
|
inline |
Definition at line 377 of file G4VEmModel.cc.
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void SetCurrentCouple(const G4MaterialCutsCouple *)
const G4Material * GetMaterial() const
size_t G4VEmModel::idxTable |
|
protected |
const std::vector<G4double>* G4VEmModel::theDensityFactor |
|
protected |
const std::vector<G4int>* G4VEmModel::theDensityIdx |
|
protected |
The documentation for this class was generated from the following files:
- geant4.10.03.p01/source/processes/electromagnetic/utils/include/G4VEmModel.hh
- geant4.10.03.p01/source/processes/electromagnetic/utils/src/G4VEmModel.cc