Angular momentum

JenaAtomicCalculator.AngularMomentum.ChengIMethod

`AngularMomentum.ChengI

  • (kapa::Int64, ma::AngularM64, kapb::Int64, mb::AngularM64, L::AngularJ64, M::AngularM64)` ... evaluates the angular I (kappa m, kappa' m', LM) integral as defined by Cheng, NATO summerschool (198x), Eq. (A4.5), and including the full magnetic (orientational) dependence. A value::Float64 is returned.
source
JenaAtomicCalculator.AngularMomentum.ChengIMethod
  • (kapa::Int64, kapb::Int64, L::AngularJ64)` ... evaluates the same angular I (kappa m, kappa' m', LM) integral but without the magnetic (orientational) dependence. A value::Float64 is returned.
source
JenaAtomicCalculator.AngularMomentum.ClebschGordanMethod

AngularMomentum.ClebschGordan(ja, ma, jb, mb, Jab, Mab) ... calculates the Clebsch-Gordan coefficient <ja, ma, jb, mb; Jab, Mab> for given quantum numbers by a proper call to a Wigner 3-j symbol. A value::Float64 is returned.

source
JenaAtomicCalculator.AngularMomentum.ClebschGordan_oldMethod

AngularMomentum.ClebschGordan_old(ja::AngularJ64, ma::AngularM64, jb::AngularJ64, mb::AngularM64, Jab::AngularJ64, Mab::AngularM64) ... calculates the Clebsch-Gordan coefficient <ja, ma, jb, mb; Jab, Mab> for given quantum numbers by a proper call to a Wigner 3-j symbol; a value::Float64 is returned.

source
JenaAtomicCalculator.AngularMomentum.JohnsonIMethod

AngularMomentum.JohnsonI(kapa::Int64, ma::AngularM64, kapb::Int64, mb::AngularM64, L::AngularJ64, M::AngularM64) ... evaluates the angular CL (kappa m, kappa' m', L M) integral as defined in his book. A value::Float64 is returned.

source
JenaAtomicCalculator.AngularMomentum.Wigner_3jMethod

AngularMomentum.Wigner_3j(a, b, c, m_a, m_b, m_c) ... calculates the value of a Wigner 3-j symbol for given quantum numbers as displayed in many texts on the theory of angular momentum (see R. D. Cowan, The Theory of Atomic Structure and Spectra; University of California Press, 1981, p. 142); it calls the corresponding function from the GNU Scientific Library. A value::Float64 is returned.

source
JenaAtomicCalculator.AngularMomentum.Wigner_3j_oldMethod

AngularMomentum.Wigner_3j_old(a::AngularJ64, b::AngularJ64, c::AngularJ64, m_a::AngularM64, m_b::AngularM64, m_c::AngularM64) ... calculates the value of a Wigner 3-j symbol for given quantum numbers by its algebraic formulae as displayed in many texts on the theory of angular momentum (see R. D. Cowan, The Theory of Atomic Structure and Spectra; University of California Press, 1981, p. 142); a value::Float64 is returned.

source
JenaAtomicCalculator.AngularMomentum.Wigner_6jMethod

AngularMomentum.Wigner_6j(a, b, c, d, e, f) ... calculates the value of a Wigner 6-j symbol for given quantum numbers as displayed in many texts on the theory of angular momentum (see R. D. Cowan, The Theory of Atomic Structure and Spectra; University of California Press, 1981, p. 142); it calls the corresponding function from the GNU Scientific Library. A value::Float64 is returned.

source
JenaAtomicCalculator.AngularMomentum.Wigner_6j_oldMethod

AngularMomentum.Wigner_6j_old(a::AngularJ64, b::AngularJ64, c::AngularJ64, d::AngularJ64, e::AngularJ64, f::AngularJ64) ... calculates the value of a Wigner 6-j symbol for given quantum numbers by its algebraic formulae as displayed in many texts on the theory of angular momentum (see R. D. Cowan, The Theory of Atomic Structure and Spectra; University of California Press, 1981, p. 142); a value::Float64 is returned.

source
JenaAtomicCalculator.AngularMomentum.Wigner_9jMethod

AngularMomentum.Wigner_9j(a, b, c, d, e, f, g, h, i) ... calculates the value of a Wigner 3-j symbol for given quantum numbers as displayed in many texts on the theory of angular momentum (see R. D. Cowan, The Theory of Atomic Structure and Spectra; University of California Press, 1981, p. 142); it calls the corresponding function from the GNU Scientific Library. A value::Float64 is returned.

source
JenaAtomicCalculator.AngularMomentum.Wigner_9j_oldMethod

AngularMomentum.Wigner_9j_old(a::AngularJ64, b::AngularJ64, c::AngularJ64, d::AngularJ64, e::AngularJ64, f::AngularJ64, g::AngularJ64, h::AngularJ64, i::AngularJ64) ... calculates the value of a Wigner 6-j symbol for given quantum numbers by its algebraic formulae as displayed in many texts on the theory of angular momentum (see R. D. Cowan, The Theory of Atomic Structure and Spectra; University of California Press, 1981, p. 142); a value::Float64 is returned.

source
JenaAtomicCalculator.AngularMomentum.Wigner_DFunctionMethod

AngularMomentum.Wigner_DFunction(j, p, q, alpha::Float64, beta::Float64, gamma::Float64) ... calculates the value of a Wigner D^j_pq (alpha, beta, gamma) for given quantum numbers and (Euler) angles (alpha, beta, gamma). It makes use of the small Wigner d(beta) matrix as the key part. A value::Float64 is returned.

source
JenaAtomicCalculator.AngularMomentum.Wigner_dmatrixMethod

AngularMomentum.Wigner_dmatrix(j, mp, m, beta::Float64) ... calculates the value of the small Wigner d^j_m',m (beta) for given quantum numbers and the angle beta. Wigner's formula is applied to evaluate the small Wigner matrix; a value::Float64 is returned.

source
JenaAtomicCalculator.AngularMomentum.allowedDoubleKappaCouplingSequenceMethod

AngularMomentum.allowedDoubleKappaCouplingSequence(syma::LevelSymmetry, symb::LevelSymmetry, maxKappa::Int64) ... to determine all allowed coupling sequences that fulfill syma + (kappa1, symx, kappa2) –> symb == syma + kappa1 –> symx + kappa2 –> symb, and where kappa1|, |kappa2| <= maxKappa. A list::Array{Tuple{Int64,LevelSymmetry,Int64},1}) of (kappa1, symx, kappa2) is returned.

source
JenaAtomicCalculator.AngularMomentum.allowedDoubleKappaSymmetriesMethod

AngularMomentum.allowedDoubleKappaSymmetries(syma::LevelSymmetry, kappa1::Int64, kappa2::Int64, symb::LevelSymmetry) ... to determine all allowed level symmetries symx that can be coupled to the sequence syma + kappa1 –> {symx} + kappa2 –> symb. A list::Array{LevelSymmetry,1} of symx is returned.

source
JenaAtomicCalculator.AngularMomentum.allowedDoubleKappasMethod

AngularMomentum.allowedDoubleKappas(syma::LevelSymmetry, symb::LevelSymmetry, maxKappa::Int64) ... to determine all allowed pairs of kappa that fulfill syma + (kappa1, kappa2) –> symb, and where |kappa1|, |kappa2| <= maxKappa. A list::Array{Tuple{Int64,Int64},1}) of (kappa1, kappa2) is returned.

source
JenaAtomicCalculator.AngularMomentum.allowedKappaSymmetriesMethod

AngularMomentum.allowedKappaSymmetries(syma::LevelSymmetry, symb::LevelSymmetry) ... to determine all allowed single-electron symmetries/partial waves kappa (l,j) that can be coupled to the given level symmetries. A list::Array{Int64,1} of kappa-values is returned.

source
JenaAtomicCalculator.AngularMomentum.allowedTotalSymmetriesMethod

AngularMomentum.allowedTotalSymmetries(symf::LevelSymmetry, mp2::EmMultipole, mp1::EmMultipole, symi::LevelSymmetry) ... to determine all allowed total symmetries J^P that can be constructed by coupling a multipole wave mp1 to the initial symmetry symi, and which can be further coupled with mp2 to the final symmetry symf. A list::Array{LevelSymmetry,1} of total symmetries is returned.

source
JenaAtomicCalculator.AngularMomentum.allowedTotalSymmetriesMethod

AngularMomentum.allowedTotalSymmetries(syma::LevelSymmetry, kappa::Int64) ... to determine all allowed total symmetries J^P that can be constructed by coupling a partial wave kappa (l,j) to the given level symmetry syma. A list::Array{LevelSymmetry,1} of total symmetries is returned.

source
JenaAtomicCalculator.AngularMomentum.isTriangleMethod

AngularMomentum.isTriangle(ja::AngularJ64, jb::AngularJ64, jc::AngularJ64) ... evaluates to true if Delta(ja,jb,jc) = 1, ie. if the angular momenta ja, jb and jc can couple to each other, and false otherwise.

source
JenaAtomicCalculator.AngularMomentum.parityEmMultipolePiMethod

AngularMomentum.parityEmMultipolePi(pa::Parity, multipole::EmMultipole, pb::Parity) ... evaluates to true if the given multipole fullfills the parity selection rule pi(a, multipole, b) = 1, and false otherwise. This includes a proper test for both, electric and magnetic multipoles, based on multipole.electric.

source
JenaAtomicCalculator.AngularMomentum.phaseFactorMethod

AngularMomentum.phaseFactor(list::Array{Any,1}) ... checks and calculates the phase factor (-1)^(ja + mb -jc ...) that occur frequently in angular momentum theory; a value +1. or -1. is returned. Use phaseFactor([ja::Union{AngularJ64,AngularM64), -1, mb::Union{AngularJ64,AngularM64), ..., 1, jc::Union{AngularJ64,AngularM64)]) to specify the phase.

source
JenaAtomicCalculator.AngularMomentum.sphericalYlmMethod

AngularMomentum.sphericalYlm(l::Int64, m::Int64, theta::Float64, phi::Float64) ... calculates the spherical harmonics for low l-values explicitly. A value::Complex{Float64} is returned. Ylm = sqrt( (2l+1) / (twotwo*pi) ) * spherical_Clm(l,m,theta,phi).

source
JenaAtomicCalculator.AngularMomentum.triangularDeltaMethod

AngularMomentum.triangularDelta(ia2::Int64, ib2::Int64, ic2::Int64) ... calculates the tringular Delta(ja,jb,jc). The arguments in this integer function are i2a = 2*ja+1, ... The result is 0 if the triangular condition failes and 1 otherwise.

source

Atomic interaction strength

JenaAtomicCalculator.InteractionStrength.XLCoefficientType

struct InteractionStrength.XLCoefficient ... defines a type for coefficients of the two-electron (Breit) interaction

+ kind      ::Char       ... Kind of integral, either 'S' or 'T'
+ nu        ::Int64      ... Rank of the integral.
+ a         ::Orbital    ... Orbitals a, b, c, d.
+ b         ::Orbital
+ c         ::Orbital
+ d         ::Orbital
+ coeff     ::Float64    ... corresponding coefficient.
source
JenaAtomicCalculator.InteractionStrength.MabEmissionJohnsonyMethod

InteractionStrength.MabEmissionJohnsony(mp::EmMultipole, gauge::EmGauge, omega::Float64, a::Orbital, b::Orbital, grid::Radial.Grid) ... to compute the (single-electron reduced matrix element) interaction strength <a || O^(Mp,emission) || b> for the interaction with the Mp multipole component of the radiation field and the transition frequency omega, and within the given gauge. A value::Float64 is returned. This procedure has been re-worked due to the book by Johnson (2007).

source
JenaAtomicCalculator.InteractionStrength.MabEmissionJohnsony_WuMethod

InteractionStrength.MabEmissionJohnsony_Wu(mp::EmMultipole, gauge::EmGauge, omega::Float64, a::Orbital, b::Orbital, grid::Radial.Grid) ... to compute the (single-electron reduced matrix element) interaction strength <a || O^(Mp,emission) || b> for the interaction with the Mp multipole component of the radiation field and the transition frequency omega, and within the given gauge. The caluclation is performed by using the function AngularMomentum.CLreducedme_sms. A value::Float64 is returned. This procedure has been re-worked due to the book by Johnson (2007).

source
JenaAtomicCalculator.InteractionStrength.MbaAbsorptionChengMethod

InteractionStrength.MbaAbsorptionCheng(mp::EmMultipole, gauge::EmGauge, omega::Float64, b::Orbital, a::Orbital, grid::Radial.Grid) ... to compute the (single-electron reduced matrix element) interaction strength <b || O^(Mp, absorption) || a> for the interaction with the Mp multipole component of the radiation field and the transition frequency omega, and within the given gauge. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.MbaEmissionAndreyOldMethod

InteractionStrength.MbaEmissionAndrey(mp::EmMultipole, gauge::EmGauge, omega::Float64, b::Orbital, a::Orbital, grid::Radial.Grid) ... to compute the (single-electron reduced matrix element) interaction strength <b || O^(Mp, emission) || a> for the interaction with the Mp multipole component of the radiation field and the transition frequency omega, and within the given gauge. A value::Float64 is returned. At present, only the magnetic matrix elements are implemented. This procedure has been worked out with Andrey but is currently not in use.

source
JenaAtomicCalculator.InteractionStrength.MbaEmissionChengMethod

InteractionStrength.MbaEmissionCheng(mp::EmMultipole, gauge::EmGauge, omega::Float64, b::Orbital, a::Orbital, grid::Radial.Grid) ... to compute the (single-electron reduced matrix element) interaction strength <b || O^(Mp,emission) || a> for the interaction with the Mp multipole component of the radiation field and the transition frequency omega, and within the given gauge. A value::Float64 is returned. This procedure has been first worked out with Andrey; in this case, however, the phases are not under good control, and this gives rise to wrong amplitudes and rates. The procedure is currently not in use.

source
JenaAtomicCalculator.InteractionStrength.MbaEmissionJohnsonxMethod

InteractionStrength.MbaEmissionJohnsonx(mp::EmMultipole, gauge::EmGauge, omega::Float64, b::Orbital, a::Orbital, grid::Radial.Grid) ... to compute the (single-electron reduced matrix element) interaction strength <b || O^(Mp,emission) || a> for the interaction with the Mp multipole component of the radiation field and the transition frequency omega, and within the given gauge. A value::Float64 is returned. This procedure has been adapted from Jiri's work but modified for the Coulomb gauge which was apparently wrong following some former implementation with RATIP.

source
JenaAtomicCalculator.InteractionStrength.MbaEmissionMigdalekMethod

InteractionStrength.MbaEmissionMigdalek(cp::CorePolarization, a::Orbital, b::Orbital, grid::Radial.Grid) ... to compute the (single-electron reduced matrix element) interaction strength <b || O^(E1, emission with core-polarization) || a> in length gauge. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.XL_BreitMethod

InteractionStrength.XL_Breit(L::Int64, a::Orbital, b::Orbital, c::Orbital, d::Orbital, grid::Radial.Grid, eeint::Union{BreitInteraction, CoulombBreit, CoulombGaunt}; keep::Bool=false) ... computes the the effective Breit interaction strengths X^LBreit (abcd) or Gaunt interaction strengths X^LGaunt (abcd) for given rank L and orbital functions a, b, c and d at the given grid. For keep=true, the procedure looks up the (global) directory GBLStorageXL_Coulomb and returns the corresponding value without re-calculation of the interaction strength; it also 'stores' the calculated value if not yet included. For keep=false, the interaction strength is always computed on-fly. A value::Float64 is returned. At present, only the zero-frequency Breit or Gaunt interaction is taken into account.

source
JenaAtomicCalculator.InteractionStrength.XL_BreitDampedMethod

InteractionStrength.XL_BreitDamped(tau::Float64, L::Int64, a::Orbital, b::Orbital, c::Orbital, d::Orbital, grid::Radial.Grid) ... computes the the effective Breit interaction strengths X^L_Breit (abcd) for given rank L and orbital functions a, b, c and d at the given grid. A value::Float64 is returned. At present, only the zero-frequency Breit interaction is taken into account.

source
JenaAtomicCalculator.InteractionStrength.XL_Breit_coefficientsMethod

InteractionStrength.XL_Breit_coefficients(L::Int64, a::Orbital, b::Orbital, c::Orbital, d::Orbital; onlyGaunt::Bool=false) ... evaluates the combinations and pre-coefficients for the zero-frequency Breit interaction X^L_Breit (omega=0.; abcd) for given rank L and orbital functions a, b, c and d. A list of coefficients xcList::Array{XLCoefficient,1} is returned.

source
JenaAtomicCalculator.InteractionStrength.XL_Breit_densitiesMethod

InteractionStrength.XL_Breit_densities(xcList::Array{XLCoefficient,1}, factor::Float64, grid::Radial.Grid) ... computes the the effective Breit interaction strengths X^L,0_Breit (abcd) for given rank L and a list of orbital functions a, b, c, d and angular coefficients at the given grid. A value::Float64 is returned. At present, only the zero-frequency Breit interaction is taken into account.

source
JenaAtomicCalculator.InteractionStrength.XL_CoulombMethod

InteractionStrength.XL_Coulomb(L::Int64, a::Orbital, b::Orbital, c::Orbital, d::Orbital, grid::Radial.Grid; keep::Bool=false) ... computes the the effective Coulomb interaction strengths X^LCoulomb (abcd) for given rank L and orbital functions a, b, c and d at the given grid. For keep=true, the procedure looks up the (global) directory GBLStorageXLCoulomb and returns the corresponding value without re-calculation of the interaction strength; it also 'stores' the calculated value if not yet included. For keep=false, the interaction strength is always computed on-fly. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.XL_CoulombMethod

InteractionStrength.XL_Coulomb(L::Int64, a::Subshell, b::Orbital, c::Orbital, d::Subshell, primitives::BsplinesN.Primitives) ... computes the (exchange) Coulomb interaction strengths X^L_Coulomb (.bc.) for given rank L and orbital functions as well as the given primitives. A (nsL+nsS) x (nsL+nsS) matrixV::Array{Float64,2} is returned.

source
JenaAtomicCalculator.InteractionStrength.XL_CoulombMethod

InteractionStrength.XL_Coulomb(L::Int64, a::Subshell, b::Orbital, c::Subshell, d::Orbital, primitives::BsplinesN.Primitives) ... computes the (direct) Coulomb interaction strengths X^L_Coulomb (.b.d) for given rank L and orbital functions as well as the given primitives. A (nsL+nsS) x (nsL+nsS) matrixV::Array{Float64,2} is returned.

source
JenaAtomicCalculator.InteractionStrength.XL_CoulombDampedMethod

InteractionStrength.XL_CoulombDamped(tau::Float64, L::Int64, a::Orbital, b::Orbital, c::Orbital, d::Orbital, grid::Radial.Grid) ... computes the the effective Coulomb interaction strengths X^L_Coulomb (abcd) for given rank L and orbital functions a, b, c and d at the given grid. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.XL_Coulomb_DHMethod

InteractionStrength.XL_Coulomb_DH(L::Int64, a::Orbital, b::Orbital, c::Orbital, d::Orbital, grid::Radial.Grid, lambda::Float64) ... computes the the effective Coulomb-Debye-Hückel interaction strengths X^LCoulombDH (abcd) for given rank L and orbital functions a, b, c and d at the given grid and for the given screening parameter lambda. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.XL_Coulomb_WOMethod

InteractionStrength.XL_Coulomb_WO(L::Int64, a::Orbital, b::Orbital, c::Orbital, d::Orbital, grid::Radial.Grid) ... computes the the effective Coulomb interaction strengths X^L_Coulomb (abcd) for given rank L and orbital functions a, b, c and d at the given grid but without optimization. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.XL_plasma_ionSphereMethod

InteractionStrength.XL_plasma_ionSphere(L::Int64, a::Orbital, b::Orbital, c::Orbital, d::Orbital, lambda::Float64) ... computes the effective interaction strengths X^L_ion-sphere (abcd) for given rank L and orbital functions a, b, c and d and for the plasma parameter lambda. A value::Float64 is returned. Not yet implemented !

source
JenaAtomicCalculator.InteractionStrength.X_smsAMethod

InteractionStrength.X_smsA(a::Orbital, b::Orbital, c::Orbital, d::Orbital, nm::Nuclear.Model, grid::Radial.Grid) ... computes the the effective interaction strengths X^1_sms,A (abcd) for fixed rank 1 and orbital functions a, b, c and d at the given grid. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.X_smsBMethod

InteractionStrength.X_smsB(a::Orbital, b::Orbital, c::Orbital, d::Orbital, nm::Nuclear.Model, grid::Radial.Grid) ... computes the the effective interaction strengths X^1_sms,B (abcd) for fixed rank 1 and orbital functions a, b, c and d at the given grid. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.X_smsCMethod

InteractionStrength.X_smsC(a::Orbital, b::Orbital, c::Orbital, d::Orbital, nm::Nuclear.Model, grid::Radial.Grid) ... computes the the effective interaction strengths X^1_sms,C (abcd) for fixed rank 1 and orbital functions a, b, c and d at the given grid. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.bosonShiftMethod

InteractionStrength.bosonShift(a::Orbital, b::Orbital, potential::Array{Float64,1}, grid::Radial.Grid) ... computes the <a|| h^(boson-field) ||b> reduced matrix element of the boson-field shift Hamiltonian for orbital functions a, b. This boson-field shift Hamiltonian just refers to the effective potential of the given isotope due to the (assumed) boson mass. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.dipoleMethod

InteractionStrength.dipole(a::Orbital, b::Orbital, grid::Radial.Grid) ... computes the <a|| d ||b> reduced matrix element of the dipole operator for orbital functions a, b. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.eMultipoleMethod

InteractionStrength.eMultipole(k::Int64, a::Orbital, b::Orbital, grid::Radial.Grid) ... computes the <a|| t^(Ek) ||b> reduced matrix element of the dipole operator for orbital functions a, b. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.fieldShiftMethod

InteractionStrength.fieldShift(a::Orbital, b::Orbital, deltaPotential::Array{Float64,1}, grid::Radial.Grid) ... computes the <a|| h^(field-shift) ||b> reduced matrix element of the field-shift Hamiltonian for orbital functions a, b. This field-shift Hamiltonian just refers to the difference of the nuclear potential deltaPotential of two isotopes, and which is already divided by the difference of the mean-square radii. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.hamiltonian_nmsMethod

InteractionStrength.hamiltonian_nms(a::Orbital, b::Orbital, nm::Nuclear.Model, grid::Radial.Grid) ... computes the <a|| h_nms ||b> reduced matrix element of the normal-mass-shift Hamiltonian for orbital functions a, b. A value::Float64 is returned. For details, see Naze et al., CPC 184 (2013) 2187, Eq. (37).

source
JenaAtomicCalculator.InteractionStrength.hfs_tE1Method

InteractionStrength.hfs_tE1(a::Orbital, b::Orbital, grid::Radial.Grid) ... computes the <a|| t^(E1) ||b> reduced matrix element for the HFS coupling to the electric-quadrupole moment of the nucleus for orbital functions a, b. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.hfs_tE2Method

InteractionStrength.hfs_tE2(a::Orbital, b::Orbital, grid::Radial.Grid) ... computes the <a|| t^(2) ||b> reduced matrix element for the HFS coupling to the electric-quadrupole moment of the nucleus for orbital functions a, b. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.hfs_tE3Method

InteractionStrength.hfs_tE3(a::Orbital, b::Orbital, grid::Radial.Grid) ... computes the <a|| t^(E3) ||b> reduced matrix element for the HFS coupling to the electric-quadrupole moment of the nucleus for orbital functions a, b. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.hfs_tM1Method

InteractionStrength.hfs_tM1(a::Orbital, b::Orbital, grid::Radial.Grid) ... computes the <a|| t^(1) ||b> reduced matrix element for the HFS coupling to the magnetic-dipole moment of the nucleus for orbital functions a, b. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.hfs_tM2Method

InteractionStrength.hfs_tM2(a::Orbital, b::Orbital, grid::Radial.Grid) ... computes the <a|| t^(M2) ||b> reduced matrix element for the HFS coupling to the magnetic-dipole moment of the nucleus for orbital functions a, b. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.hfs_tM3Method

InteractionStrength.hfs_tM3(a::Orbital, b::Orbital, grid::Radial.Grid) ... computes the <a|| t^(M3) ||b> reduced matrix element for the HFS coupling to the magnetic-dipole moment of the nucleus for orbital functions a, b. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.matrixL_CoulombMethod

InteractionStrength.matrixL_Coulomb(L::Int64, a::Orbital, b::Orbital, c::Orbital, d::Orbital, primitives::BsplinesN.Primitives) ... computes the partly-contracted (effective) Coulomb interaction matrices M^L_Coulomb (abcd) for given rank L and orbital functions a, b, c and d at the given grid. The matrix M^L is defined for the primitives and contracted over the two orbitals b, d (for a=c) or b, c (for a=d). An error message is issued if a != c && a != d. A matrix::Array{Float64,2} is returned.

source
JenaAtomicCalculator.InteractionStrength.multipoleTransitionMethod

InteractionStrength.multipoleTransition(mp::EmMultipole, gauge::EmGauge, omega::Float64, b::Orbital, a::Orbital, grid::Radial.Grid) ... to compute the (single-electron reduced matrix element) multipole-transition interaction strength <b || T^(Mp, absorption) || a> due to Johnson (2007) for the interaction with the Mp multipole component of the radiation field and the transition frequency omega, and within the given gauge. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.schiffMomentMethod

InteractionStrength.schiffMoment(a::Orbital, b::Orbital, nm::Nuclear.Model, grid::Radial.Grid) ... computes the <a|| h^(Schiff-moment) ||b> reduced matrix element of the Schiff-moment Hamiltonian for orbital functions a, b and for the nuclear density as given by the nuclear model. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.weakChargeMethod

InteractionStrength.weakCharge(a::Orbital, b::Orbital, nm::Nuclear.Model, grid::Radial.Grid) ... computes the <a|| h^(weak-charge) ||b> reduced matrix element of the weak-charge Hamiltonian for orbital functions a, b and for the nuclear density as given by the nuclear model. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.zeeman_Delta_n1Method

InteractionStrength.zeeman_Delta_n1(a::Orbital, b::Orbital, grid::Radial.Grid) ... computes the <a|| Delta n^(1) ||b> reduced matrix element for the Zeeman-Schwinger contribution to the coupling to an external magnetic field for orbital functions a, b. A value::Float64 is returned.

source
JenaAtomicCalculator.InteractionStrength.zeeman_n1Method

InteractionStrength.zeeman_n1(a::Orbital, b::Orbital, grid::Radial.Grid) ... computes the <a|| n^(1) ||b> reduced matrix element for the Zeeman coupling to an external magnetic field for orbital functions a, b. A value::Float64 is returned.

source

Basics: data types

JenaAtomicCalculator.Basics.AbstractAngularMomentumType

abstract type Basics.AbstractAngularMomentum ... defines an abstract and a number of concrete types for dealing with angular momentum variables.

+ AngularJ64        ... to deal with (total) angular momenta J >= 0
+ AngularM64        ... to deal with magnetic quantum numbers -J <= M <= J
+ LevelSymmetry     ... to deal with an overall J^P symmetry of a level.
source
JenaAtomicCalculator.Basics.AbstractContinuumSolutionsType

abstract type Basics.AbstractContinuumSolutions ... defines an abstract and a number of singleton types for solving the continuum orbitals in a given potential.

+ ContBessel              ... generate a pure Bessel function for the large component together with kinetic balance.
+ ContSine                ... generate a pure Sine function for the large component together with kinetic balance.
+ NonrelativisticCoulomb  ... generate a non-relativistic Coulomb function for the large component together with kinetic balance.
+ AsymptoticCoulomb       ... generate a pure (asymptotic) Coulombic function for both components.
+ BsplineGalerkin         ... generate a continuum orbital with the Galerkin method.dealing with warnings that are made during a run or REPL session.
source
JenaAtomicCalculator.Basics.AbstractEeInteractionType

abstract type Basics.AbstractEeInteraction ... defines an abstract and a number of singleton types for specifying the electron-electron interaction.

+ struct DiagonalCoulomb                   ... to represent the Coulomb part of the e-e interaction for just diagonal ME.        
+ struct CoulombInteraction                ... to represent the Coulomb part of the electron-electron interaction.        
+ struct CoulombGaunt                      ... to represent the Coulomb part of the electron-electron interaction.        
+ struct BreitInteraction(factor::Float64) 
    ... to represent the (frequency-dependent) Breit part of the electron-electron interaction with factor*omega,
        including the zero-Breit approximation (if factor=0.).
+ struct CoulombBreit(factor::Float64)     ... to represent the Coulomb+Breit part of the electron-electron interaction.
source
JenaAtomicCalculator.Basics.AbstractEmFieldType

abstract type Basics.AbstractEmField ... defines an abstract type to distinguish between different static and time-dependent electric and magnetic fields; this is useful for atomic compass simulations and for the computation of Stark shifts.

+ NoEmField               ... No electric and magnetic field is defined.
+ StaticField             ... A static field that is characterized by its (real) amplitude into a given direction.
+ TimeHarmonicField       ... Define a time-harmonic field that is characterized by its (real) amplitude into a 
                              given direction and a (harmonic) frequency.
source
JenaAtomicCalculator.Basics.AbstractExcitationSchemeType

abstract type Basics.AbstractExcitationScheme ... defines an abstract and a number of singleton types to distinguish between different schemes for generating configuration lists as they frequently occur in Green function and cascade computations.

+ struct NoExcitationScheme        
... dummy scheme for (unsupported) initialization of this abstract tpye.

+ struct DeExciteSingleElectron        
... generates all excitations and de-excitations of a single electron from a given list of bound
    electron configurations. The number of electrons of the generated configurations is the same as 
    for the given bound configurations.

+ struct DeExciteTwoElectrons       
... generates all excitations and de-excitations of one or two electrons from a given list of bound
    electron configurations. The number of electrons of the generated configurations is the same as 
    for the given bound configurations.
    
+ struct AddSingleElectron             
... generates configurations by just adding a single electrons to a given list of bound
    electron configurations. The number of electrons of the generated configurations is N+1.
    
+ struct ExciteByCapture             
... generates all excitations and de-excitations of one or more electron from a given list of bound
    electron configurations, together with an capture of an additional electron. 
    The number of electrons of the generated configurations is N+1.
source
JenaAtomicCalculator.Basics.AbstractFieldValueType

abstract type Basics.AbstractFieldValue ... to specify and deal with the function values of different physical fields, such as f(x), f(x,y,z), f(rho,phi), f(r,theta,phi), and for both, scalar and vector fields. Often, such field values are the outcome of some computation which can be used for integration, display, etc.

+ struct Cartesian2DFieldValue{Type}     ... to specify a field value of type T in terms of x, y.
+ struct Cartesian3DFieldValue{Type}     ... to specify a field value of type T in terms of x, y, z.
+ struct PolarFieldValue{Type}           ... to specify a field value of type T in terms of rho, phi.
+ struct SphericalFieldValue{Type}       ... to specify a field value of type T in terms of r, theta, phi.
source
JenaAtomicCalculator.Basics.AbstractFieldVectorType

abstract type Basics.AbstractFieldVector ... to specify and deal with different physical field vectors, such as (Ax, Ay), (Ax, Ay, Az), (A1, A0, A_-1). Often, such field values are used to charaterize electric or magnetic fields (vector potentials).

+ struct Cartesian2DFieldVector{Type}    ... to specify a field vector of type T in terms of (Ax, Ay).
+ struct Cartesian3DFieldVector{Type}    ... to specify a field vector of type T in terms of (Ax, Ay, Az).
+ struct SphericalFieldVector{Type}      ... to specify a field vector of type T in terms of (A_1, A_0, A_-1).
source
JenaAtomicCalculator.Basics.AbstractLevelPopulationType

abstract type Basics.AbstractLevelPopulation ... defines an abstract and a number of singleton types to distinguish between different level population (models).

+ struct BoltzmannLevelPopulation     ... to represent a Boltzmann level population.       
+ struct SahaLevelPopulation          ... to represent a Saha level population.
source
JenaAtomicCalculator.Basics.AbstractLineShiftSettingsType

abstract type Basics.AbstractLineShiftSettings ... defines an abstract and a number of singleton types for the the (allowed) plasma models.

+ Basics.NoLineShiftSettings         ... No line-shift settings are defined.
+ AutoIonization.PlasmaSettings      ... Settings for Auger-line computations.
+ PhotoIonization.PlasmaSettings     ... Settings for photoionization-line computations.
source
JenaAtomicCalculator.Basics.AbstractMeshType

abstract type Basics.AbstractMesh ... to specify different mesh types in terms of their basic parameters; these meshes can be used, for example, for integration or for defining the representation of observables.

+ struct Cartesian2DMesh    ... to specify a 2D Cartesian mesh in terms of x and y.
+ struct GLegenreMesh       ... to specify a Gauss-Legendre mesh in terms of [a,b] and number of zeros.
+ struct LinearMesh         ... to specify a linear mesh in terms of [a,b] and number of points.
+ struct PolarMesh          ... to specify a 2D mesh for rho and phi.
+ struct SphercialMesh      ... to specify a 3D mesh in terms of r, theta, phi.
source
JenaAtomicCalculator.Basics.AbstractPlasmaModelType

abstract type Basics.AbstractPlasmaModel ... defines an abstract and a number of singleton types for the the (allowed) plasma models.

+ NoPlasmaModel                 ... No plasma model defined.
+ DebyeHueckelModel             ... Debye-Hueckel plasma model.
+ IonSphereModel                ... Ion-sphere (not yet supported).
+ StewartPyattModel             ... Stewart-Pyatt model (not yet supported).
+ WithoutAutoionizationModel    ... Just excludes all autoionizing levels; no original plasma model.
source
JenaAtomicCalculator.Basics.AbstractPolarizationType

abstract type Basics.AbstractPolarization ... defines an abstract type to comprise various polarizations of light and electron beams.

+ LinearPolarization        ... to specify a linearly-polarized pulse/beam.
+ LeftCircular              ... to specify a left-circularly polarized pulse/beam.
+ RightCircular             ... to specify a right-circularly polarized pulse/beam.
+ LeftElliptical            ... to specify an elliptically polarized pulse/beam.
+ RightElliptical           ... to specify an elliptically polarized pulse/beam.
+ NonePolarization          ... to specify an upolarized pulse/beam.
+ DensityMatrixPolarization ... to specify the polarization of a pulse/beam by its (2x2) 
                                density matrix (not yet).
source
JenaAtomicCalculator.Basics.AbstractPotentialType

abstract type Basics.AbstractPotential ... defines an abstract and a number of singleton types to distinguish between different (electronic) atomic potentials.

+ struct DFSpotential     ... to represent a Dirac-Fock-Slater potential.        
+ struct CoreHartree      ... to represent a core-Hartree potential.        
+ struct KohnSham         ... to represent a Kohn-Sham potential.        
+ struct HartreeSlater    ... to represent a Hartree-Slater potential.
source
JenaAtomicCalculator.Basics.AbstractProcessType

abstract type Basics.AbstractProcess ... defines an abstract and a number of singleton types to distinguish different atomic processes.

+ struct Auger            ... Auger transitions, i.e. single autoionization or the emission of a single free electron into the continuum.
+ struct AugerInPlasma    ... Auger transitions but calculated for a specified plasma model.
+ struct Compton          ... Rayleigh-Compton scattering cross sections.
+ struct Coulex           ... Coulomb-excitation of target or projeticle electrons by fast, heavy ions.
+ struct Coulion          ... Coulomb-ionization of target or projeticle electrons by fast, heavy ions.
+ struct Dierec           ... di-electronic recombination, i.e. the dielectronic capture of a free electron and the subsequent emission of a photon.
+ struct DoubleAuger      ... Double Auger rates.
+ struct ImpactExcAuto    ... di-electronic recombination, i.e. the dielectronic capture of a free electron and the subsequent emission of a photon.
+ struct MultiPhotonDE    ... multi-photon excitation and decay rates, including 2-photon, etc. processes.
+ struct MultiPI          ... multi-photon (single-electron) ionization.
+ struct MultiPDI         ... multi-photon (single-electron) double ionization.
+ struct Photo            ... Photoionization processes, i.e. the emission of a single free electron into the continuum due to an external light field.
+ struct PhotoDouble      ... Photo-double ionization rates.
+ struct PhotoExc         ... Photoexcitation rates.
+ struct PhotoExcFluor    ... photoexcitation fluorescence rates and cross sections.
+ struct PhotoExcAuto     ... photoexcitation autoionization cross sections and collision strengths.
+ struct PhotoInPlasma    ... Photoionization processes but calculated for a specified plasma model.
+ struct PhotoIonFluor    ... photoionization fluorescence rates and cross sections.
+ struct PhotoIonAuto     ... photoionization autoionization cross sections and collision strengths.
+ struct Radiative        ... Radiative (multipole) transitions between bound-state levels of the same charge state.
+ struct Rec              ... radiative electron capture, i.e. the capture of a free electron with the simultaneous emission of a photon.
+ struct Eimex            ... electron-impact excitation cross sections and collision strengths.
+ struct RAuger           ... Radiative Auger rates.
source
JenaAtomicCalculator.Basics.AbstractQuantizationAxisType

abstract type Basics.AbstractQuantizationAxis ... defines an abstract type to distinguish between different choices of the quantization axis in light-atom interactions.

+ DefaultQuantizationAxis      ... Use the (default) z-axis for quantization of atomic levels.
+ StaticQuantizationAxis       ... Define a static quantization axis in terms of a unit vector.
+ HarmonicQuantizationAxis     ... Define a time-harmonic quanization axis in terms of a unit vector and a frequency.
source
JenaAtomicCalculator.Basics.AbstractScFieldType

abstract type Basics.AbstractScField ... defines an abstract and a number of singleton types to distinguish between different self-consistent fields

+ struct ALField          ... to represent an average-level field.        
+ struct EOLField         ... to represent an (extended) optimized-level field.        
+ struct DFSField         ... to represent an mean Dirac-Fock-Slater field.        
+ struct DFSwCPField      ... to represent an mean Dirac-Fock-Slater with core-polarization field.        
+ struct HSField          ... to represent an mean Hartree-Slater field.        
+ struct EHField          ... to represent an mean extended-Hartree field.        
+ struct KSField          ... to represent an mean Kohn-Sham field.        
+ struct HartreeField     ... to represent an mean Hartree field.        
+ struct CHField          ... to represent an mean core-Hartree field.        
+ struct NuclearField     ... to represent a pure nuclear (potential) field.
source
JenaAtomicCalculator.Basics.AbstractSelectionType

abstract type Basics.AbstractSelection ... defines an abstract and a number of concrete types to distinguish between level- and line-selectors

+ struct LevelSelection   ... to specify a list of levels by means of their (level) indices or level symmetries.         
+ struct LineSelection    ... to specify a list of lines by means of their (level) indices or level symmetries.      
+ struct PathwaySelection ... to specify a list of lines by means of their (level) indices or level symmetries.      
+ struct ShellSelection   ... to specify a list of lines by means of their (level) indices or level symmetries.
source
JenaAtomicCalculator.Basics.AbstractSpectrumKindType

abstract type Basics.AbstractSpectrumKind ... defines an abstract and a number of singleton types for specifying different kinds of spectra that need to be displayed

+ struct BarIntensities                         
    ... to display just intensity-bars at given x-positions.  
+ struct DiscreteLines                         
    ... to display lines (for guiding the eyes) but which are only defined at discrete x-points.
+ struct DiscretePoints                         
    ... to display discrete points that are defined at discrete x-points.
+ struct LorentzianIntensitiesConstantWidths    
    ... to display the superposition of Lorentzians with given position and intensity but constant widths.        
+ struct LorentzianIntensities  
    ... to display the superposition of Lorentzians with given position, intensity and individual widths.
source
JenaAtomicCalculator.Basics.AbstractWarningType

abstract type Basics.AbstractWarning ... defines an abstract and a number of singleton types for dealing with warnings that are made during a run or REPL session. Cf. Defaults.warn().

+ AddWarning        ... add a Warning to a warningList.
+ PrintWarnings     ... print all warnings into a jac-warn.report file.
+ ResetWarnings     ... reset (empty) the warningList, usually at the beginning of a new run.to distinguish between different warnings
source
JenaAtomicCalculator.Basics.AngularM64Type

struct Basics.AngularM64 <: AbstractAngularMomentum ... defines a type for angular momentum projections m.

+ num  ::Int64              ... numerator
+ den  ::Int64              ... denominator, must be 1 or 2
source
JenaAtomicCalculator.Basics.Cartesian2DFieldValueType

struct Basics.Cartesian2DFieldValue{Type} <: Basics.AbstractFieldValue ... to specify a scalar field value of type Type in terms of x, y.

+ x       ::Float64       ... x-coordinate.
+ y       ::Float64       ... y-coordinate.
+ val     ::Type          ... field value of type Type.
source
JenaAtomicCalculator.Basics.Cartesian2DMeshType

struct Basics.Cartesian2DMesh <: Basics.AbstractMesh ... to specify a 2D Cartesian mesh in terms of x and y.

+ xMesh       ::Basics.AbstractMesh       ... mesh for the x-coordinate.
+ yMesh       ::Basics.AbstractMesh       ... mesh for the y-coordinate.
source
JenaAtomicCalculator.Basics.Cartesian3DFieldValueType

struct Basics.Cartesian3DFieldValue{Type} <: Basics.AbstractFieldValue ... to specify a scalar field value of type Type in terms of x, y, z.

+ x       ::Float64       ... x-coordinate.
+ y       ::Float64       ... y-coordinate.
+ z       ::Float64       ... z-coordinate.
+ val     ::Type          ... field value of type Type.
source
JenaAtomicCalculator.Basics.Cartesian3DFieldVectorType

struct Basics.Cartesian3DFieldVector{Type} <: Basics.AbstractFieldVector ... to specify a scalar field vector of type T in terms of (Ax, Ay, Az).

+ x       ::Type       ... x-component.
+ y       ::Type       ... y-component.
+ z       ::Type       ... z-component.
source
JenaAtomicCalculator.Basics.ContinuumNormalizationType

abstract type Basics.AbstractContinuumNormalization ... defines an abstract and a number of singleton types for dealing with the normalization of continuum orbitals.

+ PureSineNorm       ... normalize with regard to an (asymtotic) pure sine funtion, sin(kr).
+ CoulombSineNorm    ... normalize with regard to an (asymtotic) Coulombic-sine funtion, sin(kr + ...).
+ OngRussekNorm      ... normalize by following Ong & Russek (1973).
+ AlokNorm           ... normalize following Salvats Radial code
source
JenaAtomicCalculator.Basics.DFSFieldType

struct Basics.DFSField <: AbstractScField ... defines a type to describe a mean Dirac-Fock-Slater field.

+ strength           ::Float64   ... Strength factor of the DFS potential (default: strength=1.0)
source
JenaAtomicCalculator.Basics.DFSwCPFieldType

struct Basics.DFSwCPField <: AbstractScField ... defines a type to describe a mean Dirac-Fock-Slater field with core-polarization.

+ corePolarization  ::CorePolarization   ... Parametrization of the core-polarization potential/contribution.
source
JenaAtomicCalculator.Basics.DebyeHueckelModelType

struct Basics.DebyeHueckelModel <: AbstractPlasmaModel ... to specify (the parameters of) a Debye-Hückel potential.

+ debyeLength  ::Float64               ... the Debye length D [a_o].
+ radius       ::Float64               ... the Debye radius R_D [a_o].
source
JenaAtomicCalculator.Basics.GLegenreMeshType

struct Basics.GLegenreMesh <: Basics.AbstractMesh ... to specify a Gauss-Legendre mesh in terms of [a,b] and number of zeros.

+ a           ::Float64       ... mesh as defined in the interval [a,b].
+ b           ::Float64
+ NoZeros     ::Int64         ... Number of GL zeros of the mesh.
source
JenaAtomicCalculator.Basics.HarmonicQuantizationAxisType

struct Basics.HarmonicQuantizationAxis <: AbstractQuantizationAxis ... to specify a time-harmonic quantization axis in terms of a unit vector for its direction and a frequency omega. !!! It need to be explained how omega is related to the components of nVector; perhaps, some further further specification is required to make this axis unique.

+ nVector  ::Basics.Cartesian3DFieldValue{Float64} ... 3D unit vector that specifies the quantization axis.
+ omega    ::Float64                               ... Frequency of the time-harmonic motion.
source
JenaAtomicCalculator.Basics.IonSphereModelType

struct Basics.IonSphereModel <: AbstractPlasmaModel ... to specify (the parameters of) a ion-sphere potential.

+ radius          ::Float64               ... the ion-sphere radius R [a_o].
+ electronDensity ::Float64               ... electron density n_e (T).
source
JenaAtomicCalculator.Basics.LevelSelectionMethod

Basics.LevelSelection(active::Bool; indices::Array{Int64,1}=Int64[], symmetries::Array{LevelSymmetry,1}=LevelSymmetry[]) ... constructor for specifying the details of a LevelSelection.

source
JenaAtomicCalculator.Basics.LevelSelectionType

struct Basics.LevelSelection < Basics.AbstractSelection ... defines a struct to specify a list of levels by means of their (level) indices or level symmetries.

+ active       ::Bool                     ... true, if some selection has been made.
+ indices      ::Array{Int64,1}           ... List of selected indices.
+ symmetries   ::Array{LevelSymmetry,1}   ... List of selected symmetries
source
JenaAtomicCalculator.Basics.LevelSymmetryType

struct Basics.LevelSymmetry < AbstractAngularMomentum ... defines a struct for defining the overall J^P symmetry of a level.

+ J          ::AngularJ64  ... total angular momentum of a level
+ parity     ::Parity      ... total parity of the level
source
JenaAtomicCalculator.Basics.LineSelectionMethod

Basics.LineSelection(active::Bool; indexPairs::Array{Tuple{Int64,Int64},1}=Tuple{Int64,Int64}[], symmetryPairs::Array{Tuple{LevelSymmetry,LevelSymmetry},1}=Tuple{LevelSymmetry,LevelSymmetry}[]) ... constructor for specifying the details of a LineSelection.

source
JenaAtomicCalculator.Basics.LineSelectionType

struct Basics.LineSelection < Basics.AbstractSelection ... defines a struct to specify a list of level pair by means of their (level) indices or level symmetries.

+ active        ::Bool                                          ... true, if some selection has been made.
+ indexPairs    ::Array{Tuple{Int64,Int64},1}                   ... List of selected index pairs.
+ symmetryPairs ::Array{Tuple{LevelSymmetry,LevelSymmetry},1}   ... List of selected symmetry pairs.
source
JenaAtomicCalculator.Basics.LinearMeshType

struct Basics.LinearMesh <: Basics.AbstractMesh ... to specify a linear mesh in terms of [a,b] and number of mesh points.

+ a           ::Float64       ... mesh as defined in the interval [a,b].
+ b           ::Float64
+ NoPoints    ::Int64         ... Number of mesh points, including a, b.
source
JenaAtomicCalculator.Basics.PathwaySelectionMethod

Basics.PathwaySelection(active::Bool; indexTriples::Array{Tuple{Int64,Int64,Int64},1}=Tuple{Int64,Int64,Int64}[], symmetryTriples::Array{Tuple{LevelSymmetry,LevelSymmetry,LevelSymmetry},1}=Tuple{LevelSymmetry,LevelSymmetry,LevelSymmetry}[]) ... constructor for specifying the details of a PathwaySelection.

source
JenaAtomicCalculator.Basics.PathwaySelectionType

struct Basics.PathwaySelection < Basics.AbstractSelection ... defines a struct to specify a list of level triple (pathways) by means of their (level) indices or level symmetries.

+ active          ::Bool                                          ... true, if some selection has been made.
+ indexTriples    ::Array{Tuple{Int64,Int64,Int64},1}             ... List of selected index triples.
+ symmetryTriples ::Array{Tuple{LevelSymmetry,LevelSymmetry,LevelSymmetry},1}  ... List of selected symmetry triples.
source
JenaAtomicCalculator.Basics.PolarFieldValueType

struct Basics.PolarFieldValue{Type} <: Basics.AbstractFieldValue ... to specify a field value of type Type in terms of rho, phi.

+ rho     ::Float64       ... rho-coordinate.
+ phi     ::Float64       ... phi-coordinate.
+ val     ::Type          ... field value of type Type.
source
JenaAtomicCalculator.Basics.PolarMeshType

struct Basics.PolarMesh <: Basics.AbstractMesh ... to specify a 2D polar mesh in terms of rho and phi.

+ rhoMesh     ::Basics.AbstractMesh       ... mesh for the rho-coordinate.
+ phiMesh     ::Basics.AbstractMesh       ... mesh for the phi-coordinate.
source
JenaAtomicCalculator.Basics.ShellSelectionType

struct Basics.ShellSelection < Basics.AbstractSelection ... defines a struct to specify a list of shells by means of different constructors.

+ active       ::Bool                     ... true, if some selection has been made.
+ shells       ::Array{Shell,1}           ... List of explicitly selected shells.
+ lSymmetries  ::Array{Int64,1}           ... List of selected l-symmetries
source
JenaAtomicCalculator.Basics.SphericalFieldValueType

struct Basics.SphericalFieldValue{Type} <: Basics.AbstractFieldValue ... to specify a field value of type Type in terms of r, theta, phi.

+ r       ::Float64       ... r-coordinate.
+ theta   ::Float64       ... theta-coordinate.
+ phi     ::Float64       ... phi-coordinate.
+ val     ::Type          ... field value of type Type.
source
JenaAtomicCalculator.Basics.SphericalMeshType

struct Basics.SphericalMesh <: Basics.AbstractMesh ... to specify a 3D polar mesh in terms of r, theta and phi.

+ rMesh       ::Basics.AbstractMesh         ... mesh for the r-coordinate.
+ thetaMesh   ::Basics.AbstractMesh         ... mesh for the theta-coordinate.
+ phiMesh     ::Basics.AbstractMesh         ... mesh for the phi-coordinate.
source
JenaAtomicCalculator.Basics.StaticFieldType

struct Basics.StaticField <: AbstractEmField ... to specify a static – electric or magnetic – field that is characterized by its (real) amplitude into a given direction.

+ amplitude    ::Basics.Cartesian3DFieldVector{Float64}  ... 3D vector that represents the amplitude A_o of the field.
source
JenaAtomicCalculator.Basics.StaticQuantizationAxisType

struct Basics.StaticQuantizationAxis <: AbstractQuantizationAxis ... to specify a static quantization axis in terms of a unit vector for its direction.

+ nVector      ::Basics.Cartesian3DFieldValue{Float64}
    ... 3D unit vector that specifies the quantization axis.
source
JenaAtomicCalculator.Basics.StewartPyattModelType

struct Basics.StewartPyattModel <: AbstractPlasmaModel ... to specify (the parameters of) a Stewart-Pyatt plasma model.

+ radius          ::Float64               ... the Stewart-Pyatt radius R [a_o].
+ electronDensity ::Float64               ... electron density n_e (T).
+ lambda          ::Float64               ... lambda^(ST) (T, ne, ni)
source
JenaAtomicCalculator.Basics.TimeHarmonicFieldType

struct Basics.TimeHarmonicField <: AbstractEmField ... to specify a time-harmonic – electric or magnetic – field that is characterized by its (real) amplitude into a given direction and a (harmonic) frequency.

+ amplitude    ::Basics.Cartesian3DFieldVector{Float64}  ... 3D vector that represents the amplitude A_o of the field.
+ omega        ::Float64                                ... Frequency of the time-harmonic field.
source
JenaAtomicCalculator.Basics.oplusMethod

Basics.oplus(ja::AngularJ64, jb::AngularJ64) ... adds the angular momenta ja oplus jb and returns a list::Array{AngularJ64,1} of j-valuescin the interval |ja - jb| <= j <= ja + jb.

source
JenaAtomicCalculator.Basics.projectionsMethod

Basics.projections(ja::AngularJ64) ... returns all allowed projections of a given angular momenta ja as a list::Array{AngularM64,1} of m-values, i.e. -ja, -ja+1, ..., j.

source
JenaAtomicCalculator.Basics.CorePolarizationType

struct Basics.CorePolarization ... defines a type to describe the influence of core-polarization upon the potential or interaction operator.

+ doApply          ::Bool                ... True if core-polarization is to be applied, and false otherwise.
+ coreAlpha        ::Float64             ... Polarizibility of the ionic core.
+ coreRadius       ::Float64             ... Radius of the ionic core.
+ coreShells       ::Array{Shell,1}      ... Shells of the ionic core
source
JenaAtomicCalculator.Basics.EmGaugeType

@enum Basics.EmGauge ... defines a enumeration for the (allowed) gauges to deal with the radiation field; this differs from the UseGauge which can be selected for explicit computations.

    + NoGauge      ... No gauge define (yet).
    + Coulomb      ... Coulomb gauge (= velocity gauge in the non-relativistic limit; Grant, 1974).
    + Babushkin    ... Babushkin gauge (= length gauge in the non-relativistic limit; Grant, 1974). 
    + Magnetic     ... Magnetic gauge (transverse velocity gauge).
    + Velocity     ... Transverse velocity gauge (Johnson, 2007).
    + Length       ... Length gauge (Johnson, 2007).
source
JenaAtomicCalculator.Basics.EmMultipoleType

struct Basics.EmMultipole ... defines a struct for the em multipoles, E1, M1, E2, ...

+ L           ::Int64    ... multipolarity
+ electric    ::Bool     ... True if the multipole is an electric one, and false otherwise.
source
JenaAtomicCalculator.Basics.EmPropertyType

struct EmProperty ... defines a type to maintain two gauge forms of a computed result that depends on the radiation field.

+ Coulomb         ::Float64    ... Value for the Coulomb gauge of the radiation field.
+ Babushkin       ::Float64    ... Value for the Coulomb gauge of the radiation field.
source
JenaAtomicCalculator.Basics.IsotopicFractionType

struct Basics.IsotopicFraction ... defines a type to deal with different – individual or averaged – isotopic fractions in a computation.

+ Z        ::Float64    ... nuclear charge
+ A        ::Float64    ... (mean) mass number of the isotope/isotopic mixture.
+ x        ::Float64    ... fraction 0 <= x <=1 of the given isotope/isotopic mixture.
source
JenaAtomicCalculator.Basics.ShellType

struct Basics.Shell ... defines a enumeration for the allowed values of a non-relativistic shell.

+ n     ::Int64  ... principal quantum number 
+ l     ::Int64  ... orbital angular quantum number
source
JenaAtomicCalculator.Basics.SubshellType

struct Basics.Subshell ... defines a type for the allowed values of a relativistic subhell.

+ n        ::Int64  ... principal quantum number 
+ kappa    ::Int64  ... relativistic angular quantum number
source
JenaAtomicCalculator.Basics.SubshellStateRType

struct Basics.SubshellStateR ... defines a struct for the relativistic antisymmetric subshell states in the seniority scheme; this struct can be accessed only internally and, therefore, the only the standard constructor is supported.

+ subshell     ::Subshell        ... to refer to a particular subshell. 
+ occ          ::Int64           ... occupation of the subshell
+ nu           ::Int64           ... seniority of the subshell state, often called nu
+ Jsub2        ::Int64           ... 2*Jsub
source
JenaAtomicCalculator.Basics.EigenType

struct Basics.Eigen ... defines a simple struct to communicate eigenvalues and eigenvectors if different diagonalization procedures are used.

+ values   ::Array{Float64,1}            ... List of eigenvalues.
+ vectors  ::Array{Vector{Float64},1}    ... List of eigenvectors.
source
JenaAtomicCalculator.Basics.EmStokesType

struct Basics.EmStokes ... defines a type to maintain the (computed) given Stokes parameter for the polarization of emitted photons or electrons.

+ P1              ::EmProperty    ... Stokes P1 parameter in Coulomb and Babushkin gauge
+ P2              ::EmProperty    ... Stokes P2 parameter.
+ P3              ::EmProperty    ... Stokes P3 parameter.
source
JenaAtomicCalculator.Basics.ExpStokesType

struct Basics.ExpStokes ... defines a type to maintain the (experimentally) given Stokes parameter for the polarization of incoming photons or electrons.

+ P1              ::Float64    ... Stokes P1 parameter.
+ P2              ::Float64    ... Stokes P2 parameter.
+ P3              ::Float64    ... Stokes P3 parameter.
source
JenaAtomicCalculator.Basics.LevelKeyType

struct Basics.LevelKey ... defines a struct for identifying a level by its symmetry, energy, etc.

+ sym          ::LevelSymmetry    ... Level symmetry.
+ index        ::Int64            ... Index of this level in its original multiplet, or 0.
+ energy       ::Float64          ... Total energy.
+ relativeOcc  ::Float64          ... Relative occupation of this level, if part of some multiplet/cascade.
source
JenaAtomicCalculator.Basics.LineKeyType

struct Basics.LineKey ... defines a struct for identifying a line by the keys of the initial and final level, its transition energy, decay strength, etc.

+ iLevelKey    ::LevelKey         ... Key of the initial level.
+ fLevelKey    ::LevelKey         ... Key of the final level.
+ energy       ::Float64          ... Total energy.
+ strength     ::Float64          ... Strength of the given line.
source
JenaAtomicCalculator.Basics.ScalarPropertyType

struct ScalarProperty{T} ... defines a type to maintain a scalar function dependence f(x) but for different types input and function values.

+ arg         ::Float64    ... Argument of function f(x).
+ value       ::T          ... Function value f(x).
source
JenaAtomicCalculator.Basics.SolidAngleType

struct Basics.SolidAngle ... defines a type for a solid angle Omega = (theta, phi).

+ theta  ::Float64       ... polar angle theta with regard to the z-axis; 0 <= theta <= pi.
+ phi    ::Float64       ... azimuthal angle phi (with regard to the x-z plane in mathematical positive direction); 0 <= phi <= 2 pi.
source
JenaAtomicCalculator.Basics.TensorCompType

struct Basics.TensorComp ... defines a type for a component of the statistical tensor as associated with an atomic or ionic level.

+ k               ::Int64          ... rank of the tensor component
+ q               ::Int64          ... q-component
+ comp            ::ComplexF64     ... value of the tensor component
source
JenaAtomicCalculator.Basics.isEqualMethod

Basics.isEqual(a::Eigen, b::Eigen, accuracy::Float64) ... determines whether a == b with the requested accuracy; a value::Bool is returned. The test is made element-wise on abs(a.number - b.number) > accuracy.

source
JenaAtomicCalculator.Basics.tensorCompMethod

Basics.tensorComp(k::Int64, q::Int64, tcomps::Array{TensorComp,1}; withZeros::Bool = false) ... returns the value of the statistical tensor component rho_kq if it is contained in tcomps.

source

Basics: functions

Basics: compute & generate

B-Spline basis

JenaAtomicCalculator.BsplinesN.BsplineType

struct BsplinesN.Bspline ... defines a type for a (single) B-spline that is defined on a given radial grid from r[lower:upper]. Note that only the non-zero values are specified for the B-spline function and its derivative.

+ lower        ::Int64               ... lower radial index (on the radial grid.r) from where the functions is nonzero.
+ upper        ::Int64               ... upper radial index up to which the functions is nonzero.
+ bs           ::Array{Float64,1}    ... radial B-spline functions as defined on the predefined grid.r[lower:upper]
+ bp           ::Array{Float64,1}    ... derivative of bs on the predefined grid grid.r[lower:upper]
source
JenaAtomicCalculator.BsplinesN.PrimitivesType

struct BsplinesN.Primitives ... defines a type for a set of primitive functions which typically belongs to a well-defined grid.

+ grid         ::Radial.Grid         ... radial grid on which the states are represented.
+ bsplinesL    ::Array{Bspline,1}    ... set of B-splines for the large components on the given radial grid.
+ bsplinesS    ::Array{Bspline,1}    ... set of B-splines for the small components on the given radial grid.
source
JenaAtomicCalculator.BsplinesN.computeNondiagonalDMethod

BsplinesN.computeNondiagonalD(pm::Int64, kappa::Int64, bspline1::BsplinesN.Bspline, bspline2::BsplinesN.Bspline, grid::Radial.Grid) ... computes the (radial and non-diagonal) D_kappa^+/- integral two the bsplines, all defined on grid <bspline1| +/- d/dr + kappa/r | bspline2>. – pm = +1/-1 provides the phase for taking the derivative.

source
JenaAtomicCalculator.BsplinesN.computeOverlapMethod

BsplinesN.computeOverlap(bspline1::BsplinesN.Bspline, bspline2::BsplinesN.Bspline, grid::Radial.Grid) ... computes the (radial) overlap integral <bspline1|bsplines> for two bpslines as defined on grid.

source
JenaAtomicCalculator.BsplinesN.computeVlocalMethod

BsplinesN.computeVlocal(bspline1::BsplinesN.Bspline, bspline2::BsplinesN.Bspline, pot::Radial.Potential, grid::Radial.Grid) ... computes the (radial) integral <bspline1| V_pot |bsplines> for two bpslines and the given radial potential as defined on grid.

source
JenaAtomicCalculator.BsplinesN.extractBsplineCoefficientsMethod

BsplinesN.extractBsplineCoefficients(sh::Subshell, wc::Basics.Eigen, grid::Radial.Grid) ... Here, it is assumed that the matrix wc contains the (column) eigenvectors as associated with a single-electron Dirac Hamiltonian matrix for symmetry kappa in Subshell(n, kappa). The procedure then extracts the (full) vector of B-spline coefficients for the radial orbital of subshell sh by applying the standard rules of atomic physics for the principal quantum number n. A vector::Array{Float64,1} is returned, whose length is nsL+nsS in the original basis of B-spline functions/primitives.

source
JenaAtomicCalculator.BsplinesN.generateGalerkinMatrixMethod

BsplinesN.generateGalerkinMatrix(sh::Subshell, energy::Float64, pot::Radial.Potential, primitives::BsplinesN.Primitives) ... generates the Galerkin-A matrix for the given potential and B-spline primitives; a matrix::Array{Float64,2} is returned.

source
JenaAtomicCalculator.BsplinesN.generateOrbitalFromPrimitivesMethod

BsplinesN.generateOrbitalFromPrimitives(sh::Subshell, energy::Float64, mtp::Int64, ev::Array{Float64,1}, primitives::BsplinesN.Primitives) ... generates the large and small components of a (relativistic) orbital for the subshell sh from the given primitives and the eigenvector ev. A (non-normalized) orbital::Orbital is returned.

source
JenaAtomicCalculator.BsplinesN.generateOrbitalFromPrimitivesMethod

BsplinesN.generateOrbitalFromPrimitives(sh::Subshell, wc::Basics.Eigen, primitives::BsplinesN.Primitives) ... generates the large and small components for the subshell sh from the primitives and their eigenvalues & eigenvectors. A (normalized) orbital::Orbital is returned.

source
JenaAtomicCalculator.BsplinesN.generateOrbitalsMethod

BsplinesN.generateOrbitals(subshells::Array{Subshell,1}, pot::Radial.Potential, nm::Nuclear.Model, primitives::BsplinesN.Primitives; printout::Bool=true) ... generates all single-electron orbitals from subshell list for the radial potential pot. A set of orbitals::Dict{Subshell, Orbital} is returned.

source
JenaAtomicCalculator.BsplinesN.generateOrbitalsHydrogenicMethod

BsplinesN.generateOrbitalsHydrogenic(subshells::Array{Subshell,1}, nm::Nuclear.Model, primitives::BsplinesN.Primitives; printout::Bool=true) ... generates all single-electron orbitals from subshell list for the nuclear potential as specified by nm. A set of orbitals::Dict{Subshell, Orbital} is returned.

source
JenaAtomicCalculator.BsplinesN.generatePrimitivesMethod

BsplinesN.generatePrimitives(grid::Radial.Grid) ... generates the breaks, knots and the B-spline primitives of order k, both for the large and small components. The function applies the given grid parameters; no primitive is defined beyond grid[n_max]. The definition of the primitives follow the work of Zatsarinny and Froese Fischer, CPC 202 (2016) 287. –- A (set of) primitives::BsplinesN.Primitives is returned.

source
JenaAtomicCalculator.BsplinesN.generateTTpMatrix!Method

BsplinesN.generateTTpMatrix!(TTp::String, kappa::Int64, primitives::BsplinesN.Primitives, storage::Dict{String,Array{Float64,2}}) ... returns the TTp block of the (single-electron) Dirac Hamiltonian matrix for an electron with symmetry kappa without any potential. The following TTp strings are allowed: ["LL-overlap", "SS-overlap", "LS-Dkappa^-", "LS-Dkappa^+"].

    Two modes are distinguished owing to the values that are available in the storage (Dict).
        * The TTp matrix block from the storage is returned, if an entry is known; it is assumed that this matrix
          block belong to the given set of primitives.
        * The TTp matrix is computed and set to the storage otherwise; from the TTp string, the key string
          key = string(kappa) * ":" * TTp is generated an applied in the storage dictionary.
          
    All B-splines are supposed to be defined for the same (radial) grid; a  matrix::Array{Float64,2}  is returned which 
    is quadratic for 'LL-overlap' and 'SS-overlap' and whose dimension depends on the number of B-splines for the large 
    and small component, otherwise.
source
JenaAtomicCalculator.BsplinesN.setupLocalMatrixMethod

`BsplinesN.setupLocalMatrix(kappa::Int64, primitives::BsplinesN.Primitives, pot::Radial.Potential, storage::Dict{String,Array{Float64,2}}) ...set-up the local parts of the generalized eigenvalue problem for the symmetry block kappa and the given (local) potential pot. The B-spline (basis) functions are defined by primitivesL for the large component and primitivesS for the small one, respectively.

source

Electron continuum

JenaAtomicCalculator.Continuum.SettingsType

struct Continuum.Settings ... defines a type for the parameters for computing continuum orbitals.

+ includeExchange         ::Bool    ... True, if the exchange is to be included and false otherwise.
+ mtp                     ::Int64   ... No of grid points for which the continuum orbital(s) are to be computed.
source
JenaAtomicCalculator.Continuum.generateOrbitalAsymptoticCoulombMethod

Continuum.generateOrbitalAsymptoticCoulomb(energy::Float64, sh::Subshell, pot::Radial.Potential, settings::Continuum.Settings) ... to generate a simple relativistic and non-normalized continuum orbital but by assuming an asymptotic sin/cos behaviour for both components at all mesh points. For this behaviour, both components and their derivatives can be calculated analytically. The effective charge Zbar is determined for the given potential and the non-Coulombic phase shifts is set to zero. An non-normalized cOrbital::Orbital is returned.

    Warning: At present, only the real part eta.re is taken into account; check for consistency !!
source
JenaAtomicCalculator.Continuum.generateOrbitalBesselMethod

Continuum.generateOrbitalBessel(energy::Float64, sh::Subshell, grid::Radial.Grid, settings::Continuum.Settings) ... to generate a simple, non-relativistic and non-normalized Bessel wave for the large component of the continuum orbital and a small component due to the kinetic-balance condition. A (non-normalized) orbital::Orbital is returned. While Pprime is obtained from the analytical expression of the large component, Qprime is set to zero here.

source
JenaAtomicCalculator.Continuum.generateOrbitalForLevelMethod

Continuum.generateOrbitalForLevel(energy::Float64, sh::Subshell, level::Level, nm::Nuclear.Model, grid::Radial.Grid, basis::Basis, settings::Continuum.Settings) ... to generate a continuum orbital for the (continuum) subshell sh, the energy and the effective charge within the given potential. The continuum orbital is generated orthogonal with regard to all subshells of the same symmetry in the basis. All further specifications about this generations are made by proper settings. A tupel of a (continuum) (orbital::Orbital, phase::Float64) is returned.

source
JenaAtomicCalculator.Continuum.generateOrbitalGalerkinMethod

Continuum.generateOrbitalGalerkin(energy::Float64, sh::Subshell, pot::Radial.Potential, settings::Continuum.Settings) ... to generate a non-normalized continuum orbital within the given local potential by using the Galerkin method and a given B-spline basis. A (non-normalized) orbital::Orbital is returned.

source
JenaAtomicCalculator.Continuum.generateOrbitalLocalPotentialMethod

Continuum.generateOrbitalLocalPotential(energy::Float64, sh::Subshell, pot::Radial.Potential, settings::Continuum.Settings) ... to generate a continuum orbital for the (continuum) subshell sh, the energy and the given potential. The effective charge for the normalization of the continuum orbital is derived from the potential. All further specifications about this generations are made by proper settings; however, the function termintates if the settings.includeExchange = true. A tupel of a (continuum) (orbital::Orbital, phase::Float64, normFactor::Float64) is returned.

source
JenaAtomicCalculator.Continuum.generateOrbitalNonrelativisticCoulombMethod

Continuum.generateOrbitalNonrelativisticCoulomb(energy::Float64, sh::Subshell, Zeff::Float64, grid::Radial.Grid, settings::Continuum.Settings) ... to generate a simple, non-relativistic and non-normalized Coulomb wave for the (continuum) subshell sh, the small component is simply obtained by the kinetic-balance condition. A orbital::Orbital is returned. ***** This function does not yet work because there is no Julia implementation of the hypergeometric function with complex arguments available.

source
JenaAtomicCalculator.Continuum.generateOrbitalNonrelativisticCoulombMethod
  • (energy::Float64, sh::Subshell, pot::Radial.Potential, settings::Continuum.Settings) ... to generate a simple, non-relativistic and non-normalized Coulomb wave for the (continuum) subshell sh, the small component is simply obtained by the kinetic-balance condition. A orbital::Orbital is returned.
source
JenaAtomicCalculator.Continuum.generateOrbitalPureSineMethod

Continuum.generateOrbitalPureSine(energy::Float64, sh::Subshell, grid::Radial.Grid, settings::Continuum.Settings) ... to generate a simple, non-relativistic and non-normalized Bessel wave for the large component of the continuum orbital and a small component due to the kinetic-balance condition. A (non-normalized) orbital::Orbital is returned. While Pprime is obtained from the analytical expression of the large component, Qprime is set to zero here.

source
JenaAtomicCalculator.Continuum.gridConsistencyMethod

Continuum.gridConsistency(maxEnergy::Float64, grid::Radial.Grid) ... to check the consistency of the given grid with the maximum energy of the required continuum electrons; an error message is issued if the grid.hp = 0. or 15 * grid.hp < wavelength(maxEnergy) or if the grid has less than 600 grid points. The function also returns the recommended grid point where the normalization and phase is to be determined. This number is currently set to nrContinuum = grid.NoPoints - 200 ... to correct for the wrong 'phase behaviour' at large r-values.

source
JenaAtomicCalculator.Continuum.normalizeOrbitalAlokMethod

Continuum.normalizeOrbitalAlok(cOrbital::Orbital, pot::Radial.Potential, settings::Continuum.Settings) ... to normalize the given continuum orbital with regard to a (asymptotic) wave function as per Salvat Code. The orbitasl are normalized to unit amplitude. An ( orbital::Orbital, (δ + Δ)::Float64, N::Float64 ) is returned.

source
JenaAtomicCalculator.Continuum.normalizeOrbitalCoulombSineMethod

Continuum.normalizeOrbitalCoulombSine(cOrbital::Orbital, pot::Radial.Potential, settings::Continuum.Settings) ... to normalize the given continuum orbital with regard to a (asymptotic) pure sine-function. An (on-energy-scale-normalized) orbital::Orbital is returned.

source
JenaAtomicCalculator.Continuum.normalizeOrbitalOngRussekMethod

Continuum.normalizeOrbitalOngRussek(cOrbital::Orbital, pot::Radial.Potential, settings::Continuum.Settings) ... to normalize the given continuum orbital with regard to a (asymptotic) pure sine-function. An (on-energy-scale-normalized) orbital::Orbital is returned.

source
JenaAtomicCalculator.Continuum.normalizeOrbitalPureSineMethod

Continuum.normalizeOrbitalPureSine(cOrbital::Orbital, grid::Radial.Grid, settings::Continuum.Settings) ... to normalize the given continuum orbital with regard to a (asymptotic) pure sine-function. An (on-energy-scale-normalized) orbital::Orbital and its relative phase phi w.r.t. sin(kr + phi) is returned.

source
JenaAtomicCalculator.Continuum.twoFzeroMethod

function twoFzero(CA::ComplexF64, CB::ComplexF64, CZ::ComplexF64) ... Calculates the Hypergeometric function 2F0(CA,CB;1/CZ) hypergeometric asymptotic series. Taken from Radial package by Salvat et al. A ComplexF64 value is returned.

source

Hydrogenic ions

JenaAtomicCalculator.HydrogenicIon.energyMethod

HydrogenicIon.energy(sh::Shell, Z::Float64) ... to compute the (non-relativistic) energy for the given shell and for a point-like nuclear charge Z; the energy is printed in the current energy units to screen but is returned in Hartree. A energy::Float64 is returned.

source
JenaAtomicCalculator.HydrogenicIon.energyMethod

HydrogenicIon.energy(sh::Subshell, Z::Float64) ... to computes the Dirac energy for the hydrogenic subshell sh and for point-like nucleus with nuclear charge Z; a energy::Float64 in atomic units and without the rest energy of the electron is returned. That is the binding energy of a 1s_1/2 electron for Z=1 is -0.50000665.

source
JenaAtomicCalculator.HydrogenicIon.orbitalMethod

HydrogenicIon.orbital(sh::Subshell, nm::Nuclear.Model, grid::Radial.Grid) ... to compute a relativstic hydrogenic Dirac orbital for the given nuclear model by using an explicit diagonalization of the Dirac Hamiltonian in a B-spline basis; an orbital::Radial.Orbital is returned.

source
JenaAtomicCalculator.HydrogenicIon.radialOrbitalMethod

HydrogenicIon.radialOrbital(sh::Subshell, Z::Float64, r::Float64) ... to compute a relativistic hydrogenic Dirac orbital for the given subshell and nuclear charge Z; a value::Float64 is returned; contributed by C Naumann (2022). Implementation of the analytical solution of the Dirac equation for a point-like nucleus following the book by Johnson: Atomic Structure Theory.

source
JenaAtomicCalculator.HydrogenicIon.radialOrbitalMethod

HydrogenicIon.radialOrbital(sh::Subshell, Z::Float64, grid::Radial.Grid) ... to compute a relativistic hydrogenic Dirac orbital on the given grid; an Array with the large and small component is returned; contributed by C Naumann (2022).

source
JenaAtomicCalculator.HydrogenicIon.radialOrbitalMethod

HydrogenicIon.radialOrbital(sh::Subshell, nm::Nuclear.Model, grid::Radial.Grid) ... to compute a relativstic hydrogenic Dirac orbital for the given nuclear model by using an explicit diagonalization of the Dirac Hamiltonian in a B-spline basis; an orbital::Radial.Orbital is returned; contributed by C Naumann (2022).

source
JenaAtomicCalculator.HydrogenicIon.radialOrbital_old2022Method

HydrogenicIon.radialOrbital_old2022(sh::Subshell, Z::Float64, grid::Radial.Grid) ... to compute a relativstic hydrogenic Dirac orbital on the given grid by applying the kinetic-balance to a corresponding non-relavistic orbital; an orbital::Radial.Orbital is returned.

source
JenaAtomicCalculator.HydrogenicIon.radialOrbital_old2022Method

HydrogenicIon.radialOrbital_old2022(sh::Subshell, nm::Nuclear.Model, grid::Radial.Grid) ... to compute a relativstic hydrogenic Dirac orbital for the given nuclear model by using an explicit diagonalization of the Dirac Hamiltonian in a B-spline basis; an orbital::Radial.Orbital is returned.

source
JenaAtomicCalculator.HydrogenicIon.rkExpectationMethod

HydrogenicIon.rkExpectation(srk::String, sh::Shell, Z::Float64) ... to compute the (non-relativistic) r^k expectation value for the shell sh of an ion with charge Z; a value::Float64 [in a_o^k] is returned. The string can takes values srk = ["r", "r^2", "1/r"]

source

Spin angular coefficients

JenaAtomicCalculator.SpinAngular.AbstractAngularTypeType

abstract type SpinAngular.AbstractAngularType ... defines an abstract type and a number of data types to work with one- and two-particle operators of given rank, see also:

+ struct OneParticleOperator    ... to represent a single-particle operator with well-defined spherical tensor rank.
+ struct TwoParticleOperator    ... to represent a two-particle operator with well-defined spherical tensor rank.
source
JenaAtomicCalculator.SpinAngular.Coefficient1pType

struct SpinAngular.Coefficient1p ... a struct for defining a single spin-angular coefficient for a reduced one-particle matrix element <a || o^(L) || b>.

+ nu       ::Int64             ... Rank (L or nu) of the single-particle interaction strength.
+ a        ::Subshell          ... Left-hand subshell (orbital).
+ b        ::Subshell          ... Right-hand subshell (orbital).
+ T        ::Float64           ... (Value of) spin-angular coefficient.
source
JenaAtomicCalculator.SpinAngular.Coefficient2pType

struct SpinAngular.Coefficient2p ... a struct for defining a single spin-angular coefficient for a reduced two-particle matrix element <ab || o^(L) || cd>, such as the Coulomb interaction strength.

+ nu       ::Int64             ... Rank (L or nu) of the single-particle interaction strength.
+ a        ::Subshell          ... Left-hand subshell (orbital).
+ b        ::Subshell          ... Left-hand subshell (orbital).
+ c        ::Subshell          ... Right-hand subshell (orbital).
+ d        ::Subshell          ... Right-hand subshell (orbital).
+ V        ::Float64           ... (Value of) spin-angular coefficient.
source
JenaAtomicCalculator.SpinAngular.OneParticleOperatorType

struct SpinAngular.OneParticleOperator <: AbstractAngularType ... a struct for defining the spherial tensor structure of a single-particle operator.

+ rank            ::Int64             ... Rank of the operator.
+ parity          ::Basics.Parity     ... Parity of the operator (if needed ??)
+ sameOrbitalSet  ::Bool              ... True if orbitals for both CSF are taken from the same orbital set (if needed ??)
source
JenaAtomicCalculator.SpinAngular.QspaceTermType

struct SpinAngular.QspaceTerm ... a struct for defining a subshell term/state |j (nu) alpha Q J> == |j (nu) Q J Nr> for a subshell with well-defined j.

+ j        ::AngularJ64   ... subshell j
+ Q        ::AngularJ64   ... quasi-spin
+ J        ::AngularJ64   ... total J of subshell term
+ Nr       ::Int64        ... Additional quantum number Nr = 0,1,2.
+ min_odd  ::Int64        ... the minimal limits of the subshell terms for odd number operators in second quantization
+ max_odd  ::Int64        ... the maximal limits of the subshell terms for odd number operators in second quantization
+ min_even ::Int64        ... the minimal limits of the subshell terms for even number operators in second quantization
+ max_even ::Int64        ... the maximal limits of the subshell terms for even number operators in second quantization
source
JenaAtomicCalculator.SpinAngular.SchemeEta_aType

struct SpinAngular.SchemeEta ... to define various singleton structs in order to distinguish between different irreducible tensors and matrix elements.

+ Eta_a         ... a^(qj)_mq
+ Eta_W         ... W^(k_12) = [a x a]
+ Eta_aW        ... W^(k_12, k_2) = [a1 x [a2 x a3]]
+ Eta_Wa        ... W^(k_12, k_2) = [[a1 x a2] x a3]
+ Eta_WW        ... W^(kk0) = [[a1 x a2]^k x [a3 x a4]^k]^0
source
JenaAtomicCalculator.SpinAngular.TwoParticleOperatorType

struct SpinAngular.TwoParticleOperator ... a struct for defining the spherial tensor structure of a two-particle operator.

+ rank            ::Int64             ... Rank of the operator (if needed ??).
+ parity          ::Basics.Parity     ... Parity of the operator (if needed ??)
+ sameOrbitalSet  ::Bool              ... True if orbitals for both CSF are taken from the same orbital set (if needed ??)
source