AtomicStructure.ManyElectronWavefunctionType
ManyElectronWavefunction

A many-electron wave function configuration can either be given as a CSF (summed over spins) or a configuration of spin-orbitals (where all quantum numbers are specified).

AtomicStructure.AtomType
Atom(radial_orbitals, orbitals, configurations, mix_coeffs, potential)

An atom constitutes a set of single-electron orbitals with associated radial_orbitals, configurations which are ManyElectronWavefunction:s, comprising of anti-symmetrized combinations of such orbitals. The expansion coefficients mix_coeffs determine the linear combination of the configurations for multi-configurational atoms.

The potential can be used to model either the nucleus by itself (a point charge or a nucleus of finite extent) or the core orbitals (i.e. a pseudo-potential).

AtomicStructure.AtomMethod
Atom(init, R::AbstractQuasiMatrix, configurations, potential, ::Type{C})

Create an Atom on the space spanned by R, from the list of electronic configurations, with a nucleus modelled by potential, and initialize the orbitals according to init. C determines the eltype of the mixing coefficients.

AtomicStructure.AtomMethod
Atom(undef, ::Type{T}, R::AbstractQuasiMatrix, configurations, potential, ::Type{C}[, mix_coeffs])

Create an Atom on the space spanned by R, from the list of electronic configurations, with a nucleus modelled by potential, and leave the orbitals uninitialized. T determines the eltype of the radial orbitals and C the mixing coefficients, which by default, are initialized to [1,0,0,...].

AtomicStructure.AtomMethod
Atom(other_atom::Atom, configurations)

Create a new atom using the same basis and nuclear potential as other_atom, but with a different set of configurations. The orbitals of other_atom are copied over as starting guess.

AtomicStructure.AtomMethod
Atom(R::AbstractQuasiMatrix, configurations, potential[, ::Type{C}=eltype(R)])

Create an Atom on the space spanned by R, from the list of electronic configurations, with a nucleus modelled by potential, and initialize the orbitals to their hydrogenic values.

AtomicStructure.AtomMethod
Atom(init, ::Type{T}, R::AbstractQuasiMatrix, configurations, potential, ::Type{C})

Create an Atom on the space spanned by R, from the list of electronic configurations, with a nucleus modelled by potential, and initialize the orbitals according to init. T determines the eltype of the radial orbitals and C the mixing coefficients.

AtomicStructure.AtomicEquationsType
AtomicEquations(atom, equations, integrals)

Structure representing the (e.g. Hartree–Fock) equations for atom, along with all integrals that are shared between the equations.

AtomicStructure.AtomicOrbitalEquationType
AtomicOrbitalEquation(atom, equation, orbital, ϕ, hamiltonian)

Governs the evolution of an atomic orbital belonging to an atom. equation is the symbolic expression, from which hamiltonian is constructed. ϕ is the QuasiVector representing the radial orbital.

AtomicStructure.DirectPotentialType
DirectPotential

Special case of HFPotential for the direct interaction, in which case the potential formed from two orbitals can be precomputed before acting on a third orbital.

AtomicStructure.ExchangePotentialType
ExchangePotential

Special case of HFPotential for the exchange interaction, in which case the potential is formed from the orbital acted upon, along with another orbital, and then applied to a third orbital. Thus this potential cannot be precomputed, but must be recomputed every time the operator is applied. This makes this potential expensive to handle and the number of times it is applied should be minimized, if possible.

AtomicStructure.HFPotentialType
HFPotential(k, a, b, av, bv, V̂, poisson)

Represents the k:th multipole exansion of the Hartree–Fock potential formed by orbitals a and b (av and bv being views of their corresponding radial orbitals). is the resultant one-body potential formed, which can act on a third orbital and poisson computes the potential by solving Poisson's problem.

AtomicStructure.KineticEnergyHamiltonianType
KineticEnergyHamiltonian

The kinetic energy part of the one-body Hamiltonian, ubcluding the centrifugal potential. It is diagonal in spin, i.e. it does not couple orbitals of opposite spin.

AtomicStructure.ObservableType
Observable

Represents a physical quantity that can be observed, which is calculated as the matrix element of an operator between two configurations. All physical observables are real.

AtomicStructure.ObservableMethod
Observable(operator, atom, overlaps, integrals, integral_map,
           [; selector=default_selector(atom), double_counted=false])

Construct an observable corresponding the operator acting on atom; if double_counted, only return those terms that would be double-counted, otherwise return the normal observable equations. overlaps is a list of non-orthogonal, integrals a list of common integrals, and integral_map is a mapping from symbolic integrals to OrbitalIntegrals. selector selects those orbitals not modelled by the potential, e.g. all orbitals in case of a nuclear potential, fewer in case of a pseudo-potential. double_counted can be used to only return those terms of the sum that would be double-counted if simply summing over orbital contributions (applicable to two-body operators only). Half the result of the double-counted term is then subtracted the sum over orbital contributions; this makes the expression slightly more symmetric in orbital space than if the sum is derived avoiding double-counting.

AtomicStructure.OperatorMatrixElementType
OperatorMatrixElement(a, b, Â, coeff, value)

Represents the matrix element coeff*⟨a|Â|b⟩, for the operator and the orbitals a and b, along with the current value of the integral. Typically, is a radial part of an operator and coeff is the associated angular coefficient; coeff can be of any type convertible to a scalar.

AtomicStructure.OrbitalHamiltonianType
OrbitalHamiltonian(R, terms, mix_coeffs, projector, orbital)

The Hamiltonian for orbital is constructed from a radial basis R, a set of OrbitalHamiltonianTerm terms that describe the various interactions between orbitals, mix_coeffs which are the mixing coefficents for the multi-configurational expansion. The projector ensures orthogonality between orbital pairs which have Lagrange multipliers associated with them, by projecting out components of other orbitals every time the OrbitalHamiltonian action on orbital is computed.

AtomicStructure.OrbitalHamiltonianTermType
OrbitalHamiltonianTerm(i, j, coeff, A, integrals)

Represents a term in the orbital Hamiltonian arising from a variation of the energy expressions between configurations i and j in the multi-configurational expansion. coeff is the numeric coefficient, A is the operator acting on the orbital, and integrals is a vector of OrbitalIntegrals arising from the presence of non-orthogonal orbitals and whose values should be multiplied to form the overall coefficient.

AtomicStructure.OrbitalIntegralType
OrbitalIntegral{N}

Abstract type for integrals of rank N of orbitals, whose values need to be recomputed every time the orbitals are updated. Rank 0 corresponds to a scalar value, rank 1 to a diagonal matrix, etc.

AtomicStructure.OrbitalOverlapIntegralType
OrbitalOverlapIntegral(a, b, av, bv, value)

Represents the orbital overlap integral ⟨a|b⟩, for orbitals a and b, along with views of their radial orbitals av and bv and the current value of the integral.

AtomicStructure.ProjectorType
Projector(ϕs, orbitals, S)

Represents the projector on the subspace spanned by the radial orbitals ϕs (corresponding to orbitals).

AtomicStructure.RadialOrbitalType
RadialOrbital

A radial orbital is represented using basis coupled with a vector of expansion coefficients with respect to that basis; the basis implemented as an AbstractQuasimatrix.

AtomicStructure.RadialOrbitalsType
RadialOrbitals

A collection of radial orbitals is instead represented using a matrix of such expansion coefficients, where each matrix column corresponds to a single radial orbital.

AtomicStructure.SCF.KrylovWrapperMethod
SCF.KrylovWrapper(hamiltonian::OrbitalHamiltonian)

Construct a KrylovWrapper such that hamiltonian, that acts on function spaces, can be used in a Krylov solver, which works with linear algebra vector spaces.

AtomicStructure.ShiftTermType
ShiftTerm(λ)

The point of ShiftTerm is to implement an overall energy shift of the Hamiltonian.

AtomicStructure.SourceTermType
SourceTerm(operator, source_orbital, ov)

The point of SourceTerm is to implement inhomogeneous terms that contribute to the equation for an orbital, and whose input is some other source_orbital. This kind of term appears in multi-configurational problems.

AtomicStructure.SCF.energyFunction
SCF.energy(hfeq::AtomicOrbitalEquation[, which=:total])

Compute the orbital energy for the orbital governed by hfeq. Optionally select which contribution is computed (:total, :onebody, :direct, or :exchange).

AtomicStructure.SCF.energy_matrix!Method
energy_matrix!(H, hfeqs::AtomicEquations[, which=:energy])

Compute the energy matrix by computing the energy observable and storing it in H. Requires that hfeqs has the :energy and :kinetic_energy Observables registered (this is the default).

AtomicStructure.SCF.energy_matrix!Method
energy_matrix!(H, hamiltonian, ϕ)

Compute the contribution of hamiltonian to the Hamiltonian matrix H by repeatedly acting on the associated radial orbital ϕ with the different multi-configurational OrbitalHamiltonianTerms of hamiltonian.

AtomicStructure.SCF.orbitalsMethod
SCF.orbitals(atom)

Returns a view of the radial orbital coefficients (NB, it does not return the MulQuasiMatrix, but the actual underlying expansion coefficients, since SCF operates on them in the self-consistent iteration).

AtomicStructure.SCF.update!Method
update!(equations::AtomicEquations[, atom::Atom])

Recompute all integrals using the current values for the radial orbitals (optionally specifying which atom from which the orbitals are taken).

AtomicStructure.SCF.update!Method
SCF.update!(p::DirectPotential)

Update the direct potential p by solving the Poisson problem with the current values of the orbitals forming the mutual density.

AtomicStructure.SCF.update!Method
SCF.update!(p::DirectPotential, atom::Atom)

Update the direct potential p by solving the Poisson problem with the current values of the orbitals of atom forming the mutual density.

AtomicStructure.coefficientMethod
coefficient(term::OrbitalHamiltonianTerm, c::Vector)

Return the multiplicative coefficient pertaining to term, including the conj(c_i)*c_j mixing coefficients, due to the configuration-interaction.

AtomicStructure.coefficientMethod
coefficient(term::OrbitalHamiltonianTerm)

Return the multiplicative coefficient pertaining to term, excluding the conj(c_i)*c_j mixing coefficients, due to the configuration-interaction.

AtomicStructure.diagonalize_one_bodyMethod
diagonalize_one_body(H, nev; method=:arnoldi_shift_invert, tol=1e-10, σ=-1)

Diagonalize the one-body Hamiltonian H and find the nev lowest eigenpairs, using the specified diagonalization method; valid choices are

  • :arnoldi which performs the standard Krylov iteration looking for the eigenvalues with smallest real values,

  • :arnoldi_shift_invert which performs the Krylov iteration but with the shifted and inverted matrix (H - I*σ)⁻¹ looking for the eigenvalues with largest real values,

  • :eigen which uses Julia's built-in eigensolver.

tol sets the Krylov tolerance.

AtomicStructure.find_symmetriesMethod
find_symmetries(orbitals)

Group all orbitals according to their symmetries, e.g. ℓ for Orbitals. This is used to determine which off-diagonal Lagrange multipliers are necessary to maintain orthogonality.

AtomicStructure.hydrogenic!Method
hydrogenic!(atom[; find_lowest=false, find_lowest_ℓmax=Inf, kwargs...])

Initialize the radial orbitals of atom to their unscreened hydrogenic values. This is done via simple diagonalization of the one-body Hamiltonian for each angular symmetry. If find_lowest is true, only the orbital(s) with the lowest energy is kept (out of those with ℓ≤find_lowest_ℓmax). The kwargs are passed on to diagonalize_one_body and can be used to influence how the diagonalization is performed.

AtomicStructure.mⱼMethod
mⱼ(o::SpinOrbital)

Return the total mⱼ that o couples to in LS coupling, i.e. mℓ+mₛ.

AtomicStructure.observe!Method
observe!(A::M, atom::Atom, o::Observable)

Compute the observable o between all configurations and store the results as matrix elements of A. The right-hand side vectors are taken from the equations, whereas the lefth-hand side vectors are taken from atom.

AtomicStructure.observe!Method
observe!(A::M, o::Observable, a::Atom, b::Atom)

Compute the observable o between all configurations, store the results as matrix elements of A, and finally contract with respect to the mixing coefficients of atoms a & b and return the result as a scalar, which corresponds to a transition between the states of the respective atoms. o will be update!d with respect to b.

AtomicStructure.observe!Method
observe!(A::M, o::Observable, atom::Atom)

Compute the observable o between all configurations, store the results as matrix elements of A, and finally contract with respect to the mixing coefficients of atom and return the result as a scalar. o is not update!d with respect to atom.

AtomicStructure.observe!Method
observe!(A::M, o::Observable)

Compute the observable o between all configurations and store the results as matrix elements of A.

AtomicStructure.one_body_hamiltonianMethod
one_body_hamiltonian(::Type{Tuple}, atom, orb)

Return the kinetic and one-body potential energy operators (as a tuple) for the orbital orb of atom.

AtomicStructure.outsidecoremodelMethod
outsidecoremodel(configuration::Configuration, potential::P)

Return the part of the electronic configuration that is not part of the the configuration modelled by the potential. For a point charge, this is the same as the configuration itself, but for pseudo-potentials, typically only the outer shells remain.

AtomicStructure.projectout!Method
projectout!(y, projector)

Project out all components of y parallel to the radial orbitals projector.ϕs.

AtomicStructure.pushterms!Method
pushterms!(terms, operator, equation_terms,
           integrals, integral_map, symbolic_integrals)

For each term in equation_terms, push a term, located at CI coordinates i,j, of the overall orbital Hamiltonian to terms, constructed from operator and a product of orbital integrals, multiplied by an overall factor given by expression and multipole expansions. integrals contain common OrbitalIntegrals and integral_map maps from symbolic_integrals to integrals.

AtomicStructure.screened_hydrogenic!Method
screened_hydrogenic!(atom[; kwargs...])

Initialize the radial orbitals of atom to their screened hydrogenic values. This is done via simple diagonalization of the one-body Hamiltonian for each orbital with screening computed from all the other orbitals of the first configuration of atom. The kwargs are passed on to diagonalize_one_body and can be used to influence how the diagonalization is performed.

AtomicStructure.screeningMethod
screening(i, j)

Compute the amount of screening of orbital i due to orbital j, according to the formula

\[\sigma_{ij} = \left\{ 1 + \left[ \frac{3n_j^2 - \ell_j(\ell_j+1)}{3n_i^2 - \ell_i(\ell_i+1)} \right]^2 \right\}^{-3/2}\]

taken from Eq. (10) of

Bessis, N., & Bessis, G. (1981). Analytic Atomic Shielding
Parameters. The Journal of Chemical Physics, 74(6),
3628–3630. http://dx.doi.org/10.1063/1.441475
AtomicStructure.screeningMethod
screening(i, c)

Compute the screening of orbital i due to all other orbitals of configuration c:

\[\sigma_i = \sum_j (w_j-\delta_{ij})\sigma_{ij}\]

where $w_j$ is the occupancy of orbital j.

Examples

julia> AtomicStructure.screening(o"1s", c"1s2 2s2")
0.3820869935387247

The 1s orbital is only slightly screened by the other 1s electron and the 2 2s electrons, whereas

julia> AtomicStructure.screening(o"2s", c"1s2 2s2")
2.179703979102134

shows that the 2s electron is screened by both the 1s electrons and a little bit of the the other 2s electron.

Base.Broadcast.materialize!Method
materialize!(ma::MulAdd{<:Any, <:Any, <:Any, T, <:DirectPotential, Source, Dest})

Materialize the lazy multiplication–addition of the type y ← α*V̂*x + β*y where is a DirectPotential (with a precomputed direct potential computed via SCF.update!) and x and y are RadialOrbitals.

Base.Broadcast.materialize!Method
materialize!(ma::MulAdd{<:Any, <:Any, <:Any, T, <:ExchangePotential, Source, Dest})

Materialize the lazy multiplication–addition of the type y ← α*V̂*x + β*y where is a ExchangePotential (by solving the Poisson problem with x as one of the constituent source orbitals in the mutual density) and x and y are RadialOrbitals.

Base.copyto!Method
copyto!(dst::Atom, src::Atom)

Copy the radial orbitals of src to the corresponding ones in dst. Orbitals of src missing from dst are skipped. It is assumed that the underlying bases are compatible, i.e. that the radial coordinate of dst encompasses that of src and grid spacing &c. agree.

Base.copyto!Method
copyto!(dest::AbstractMatix, hamiltonian::OrbitalHamiltonian)

Materialize the orbital hamiltonian into matrix form and store it in dest, using the current values of all other orbitals. This is only possible if the orbital hamiltonian does not contain any ExchangePotentials or SourceTerms (which are not diagonal in orbital space), since the former is non-local (and thus not representable as a matrix) and the latter is not a linear operator (but an affine one).

Typical usage is to compute an easily factorizable matrix that can be used for preconditioning the solution of the full equation.

Base.diffMethod
diff(atom; H=atomic_hamiltonian(atom), overlaps=[], selector=outsidecoremodel, verbosity=0)

Differentiate the energy expression of the Hamiltonian H associated with the atom's configurations(s) with respect to the atomic orbitals to derive the Hartree–Fock equations for the orbitals.

By default, the Hamiltonian H=FieldFreeOneBodyHamiltonian()+CoulombInteraction().

Non-orthogonality between orbitals can be specified by providing OrbitalOverlaps between these pairs. Only those electrons not modelled by atom.potential of each configuration are considered for generating the energy expression, this can be changed by choosing another value for selector.

Base.getindexMethod
getindex(atom, js)

Returns a copy of all radial orbitals with index ∈ js.

Base.getindexMethod
getindex(atom, j)

Returns a copy of the j:th radial orbital.

Base.getindexMethod
getindex(atom, orbs)

Returns a copy of the radial orbitals corresponding to orbs.

Base.getindexMethod
getindex(atom, orb)

Returns a copy of the radial orbital corresponding to orb.

Base.hashMethod
hash(atom::Atom, h::UInt)

This computes a hash for atom, mixed with the hash seed h, considering only the fixed quantities, i.e. which radial grid is used, which orbitals are present (but not their radial coefficients), and which potential describes the nucleus. The usage is then not to identify a particular state of the atom, but rather the basis for its Hilbert space. This is useful if one would like to store the radial and mixing coefficients to a uniquely named file, after e.g. a lengthy Hartree–Fock solution; one would then set up the Atom structure as usual, but instead of running the SCF procedure again, compute the hash and load the data from the corresponding file.

Base.viewMethod
view(atom, j)

Returns a view of all radial orbitals with index ∈ js.

Base.viewMethod
view(atom, j)

Returns a view of the j:th radial orbital.

Base.viewMethod
view(atom, orb)

Returns a view of the radial orbital corresponding to orb.

LinearAlgebra.mul!Function
mul!(y, h::OrbitalHamiltonian, x)

Materialize the action of the OrbitalHamiltonian on the linear algebra vector x and store the result in y, by wrapping them both with the QuasiMatrix necessary to transform x and y to the function space of the Hamiltonian.

LinearAlgebra.normMethod
norm(atom)

This calculates the amplitude norm of the atom, i.e. $√N$ where $N$ is the number electrons. Each configuration is weighted by its mixing coefficient.