Base.MatrixMethod
Matrix(a, op::QuantumOperator, b[, overlaps])

Generate the matrix corresponding to the quantum operator op, between the configurations (e.g. SlaterDeterminants) a and b, i.e ⟨a|op|b⟩. It is possible to specify non-orthogonalities between single-particle orbitals in overlaps.

Base.MatrixMethod
Matrix(op::QuantumOperator, slater_determinants[, overlaps])

Generate the matrix corresponding to the quantum operator op, between the different slater_determinants. It is possible to specify non-orthogonalities between single-particle orbitals in overlaps.

EnergyExpressions.BitConfigurationsType
BitConfigurations(orbitals, configurations)

Represent collection of configurations as bit vectors, where true values indicate that specific orbitals are occupied.

Example

julia> bcs = BitConfigurations([[:a,:b,:c], [:x,:b,:c], [:a,:y,:c], [:a,:b,:z]])
6-orbital 4-configuration BitConfigurations

1: a b c
2: a -> x
3: b -> y
4: c -> z

julia> h = FieldFreeOneBodyHamiltonian()
ĥ₀

julia> Matrix(bcs, h)
4×4 SparseMatrixCSC{NBodyMatrixElement, Int64} with 10 stored entries:
 (a|a) + (b|b) + (c|c)  (a|x)                  - (b|y)                (c|z)
 (x|a)                  (b|b) + (c|c) + (x|x)  ⋅                      ⋅
 - (y|b)                ⋅                      (a|a) + (c|c) + (y|y)  ⋅
 (z|c)                  ⋅                      ⋅                      (a|a) + (b|b) + (z|z)
EnergyExpressions.ConjugateType
Conjugate(orbital)

Type representing the conjugation of an orbital.

Examples

julia> Conjugate(:a)
:a†
EnergyExpressions.ContractedOperatorType
ContractedOperator(a, o, b)

An NBodyOperator representing the contraction of the operator o over the orbital sets a and b. The lengths of a and b have to equal, and they cannot exceed the dimension of o.

EnergyExpressions.CoulombInteractionType
CoulombInteraction

Two-body Hamiltonian, representing the mutual Coulombic repulsion between two electrons. Is diagonal in spin, i.e. the spin of the orbitals associated with the same coordinate must be the same.

Examples

julia> EnergyExpressions.OrbitalMatrixElement((:a,:b), CoulombInteraction(), (:c,:d))
[a b|c d]

julia> EnergyExpressions.OrbitalMatrixElement((:a,:b), CoulombInteraction(), (:b,:a))
G(a,b)
EnergyExpressions.EnergyExpressionType
EnergyExpression

An energy expression is given by an energy matrix, or interaction matrix, sandwiched between a vector of mixing coefficients: E = c'H*c, where c are the mixing coefficients and H the energy matrix.

EnergyExpressions.MCEquationSystemType
MCEquationSystem(equations, integrals)

Represents a coupled system of integro-differential equations, resulting from the variation of a multi-configurational EnergyExpression, with respect to all constituent orbitals. All integrals that are in common between the equations need only be computed once per iteration, for efficiency.

EnergyExpressions.MCTermType
MCTerm(i, j, coeff, operator, source_orbital, integrals=[])

Represents one term in the multi-configurational expansion. i and j are indices in the mixing-coefficient vector c (which is subject to optimization, and thus has to be referred to), coeff is an additional coefficient, and integrals is a list of indices into the vector of common integrals, the values of which should be multiplied to form the overall coefficient.

EnergyExpressions.NBodyEquationType
NBodyEquation{N,O}(orbital, operator::NBodyOperator[, factor::NBodyTerm])

Equation for an orbital, acted upon by an operator, which may be a single-particle operator, or an N-body operator, contracted over all coordinates but one, and optionally multiplied by an NBodyTerm, corresponding to overlaps/matrix elements between other orbitals.

EnergyExpressions.NBodyMatrixElementMethod
NBodyMatrixElement(a, op, b, overlap)

Generate the matrix element of op, a linear combination of NBodyOperator, between the configurations (e.g. Slater determinants) a and b, according to the Löwdin rules. The matrix overlap contains the mutual overlaps between all single-particle orbitals in the Slater determinants. If the orbitals are all orthogonal, the Löwdin rules collapse to the Slater–Condon rules.

EnergyExpressions.NBodyMatrixElementMethod
NBodyMatrixElement(a, op, b, nzcofactors)

Generate the matrix element of the N-body operator op, between the orbital sets a and b, where nzcofactors list all N-tuples for which the determinantal cofactor of the orbital overlap matrix is non-vanishing.

EnergyExpressions.NBodyMatrixElementMethod
NBodyMatrixElement(a, op, b, overlap)

Generate the matrix element of the N-body operator op, between the Slater determinants a and b, according to the Löwdin rules. The matrix overlap contains the mutual overlaps between all single-particle orbitals in the Slater determinants. If the orbitals are all orthogonal, the Löwdin rules collapse to the Slater–Condon rules.

EnergyExpressions.NBodyTermType
NBodyTerm(factors, coeff)

Structure representing one term in the expansion of a N-body matrix element.

EnergyExpressions.OneBodyHamiltonianType
OneBodyHamiltonian

The one-body Hamiltonian, may include external fields. It is diagonal in spin, i.e. it does not couple orbitals of opposite spin.

EnergyExpressions.OrbitalEquationType
OrbitalEquation(orbital, equation,
                one_body, direct_terms, exchange_terms, source_terms)

Represents the integro-differential equation for orbital, expressed as a linear combination of the different terms, with pointers to the list of common integrals that is stored by the encompassing MCEquationSystem object.

EnergyExpressions.OrbitalMatrixElementType
OrbitalMatrixElement(a,o,b)

Represents the N-body matrix element between the sets of orbitals a and b.

Examples

julia> struct MyTwoBodyOperator <: TwoBodyOperator end

julia> EnergyExpressions.OrbitalMatrixElement((:a,:b), MyTwoBodyOperator(), (:c,:d))
⟨a b|MyTwoBodyOperator()|c d⟩
EnergyExpressions.OrbitalOverlapType
OrbitalOverlap(a,b)

Represents the overlap between the orbitals a and b in a N-body matrix element expansion.

Examples

julia> EnergyExpressions.OrbitalOverlap(:a,:b)
⟨a|b⟩
EnergyExpressions.OrbitalsType
Orbitals(orbitals, overlaps, has_overlap, non_orthogonalities)

Structure storing a common set of orbitals, along with possible overlaps between them, in case of non-orthogonalities. has_overlap is a boolean matrix indicates if a pair of orbitals have overlap, either due to non-orthogonality or if they are the same orbital. non_orthogonalities is a boolean vector that indicates if a specific orbital is non-orthogonal to any other orbital in the set of orbitals. This structure is used internally by BitConfigurations.

EnergyExpressions.SlaterDeterminantType
SlaterDeterminant(orbitals::Vector{O})

Constructs a Slater determinant from a set of spin-orbitals.

Examples

julia> SlaterDeterminant([:a, :b])
a(1)b(2) - a(2)b(1)

julia> SlaterDeterminant([:a, :b, :c])
a(1)b(2)c(3) - a(1)b(3)c(2) - a(2)b(1)c(3) + a(2)b(3)c(1) + a(3)b(1)c(2) - a(3)b(2)c(1)
EnergyExpressions.SlaterDeterminantMethod
SlaterDeterminant(configuration::Configuration{<:SpinOrbital})

Constructs a Slater determinant from the spin-orbitals of the spin-configuration configuration.

Examples

julia> SlaterDeterminant(spin_configurations(c"1s2")[1])
1s₀α(1)1s₀β(2) - 1s₀α(2)1s₀β(1)
Base.:==Method
Base.:(==)(a::NBodyMatrixElement, b::NBodyMatrixElement; kwargs...)

Test if a and b are exactly equal to each other, i.e. their terms all agree exactly, as well as the expansion coefficients. The actual comparison is performed by compare.

Base.adjointMethod
adjoint(slater_determinant)

Construct the adjoint of slater_determinant

Examples

julia> SlaterDeterminant([:a, :b])'
[a(1)b(2) - a(2)b(1)]†
Base.conjMethod
conj(o::AbstractOrbital)

Convenience function to conjugate an AbstractOrbital.

Examples

julia> conj(o"1s")
1s†
Base.conjMethod
conj(o::Conjugate)

Convenience function to unconjugate a conjugated orbital.

Examples

julia> conj(Conjugate(:a))
:a
Base.diffMethod
diff(E::Matrix{NBodyMatrixElement}, o::O)

Vary the matrix of NBodyMatrixElements with respect to the orbital o.

Examples

julia> E = Matrix(OneBodyHamiltonian()+CoulombInteraction(),
                  SlaterDeterminant.([[:a, :b], [:c, :d]]))
2×2 Array{EnergyExpressions.NBodyMatrixElement,2}:
 (a|a) + (b|b) - G(a,b) + F(a,b)  - [a b|d c] + [a b|c d]
 - [c d|b a] + [c d|a b]          (c|c) + (d|d) - G(c,d) + F(c,d)

julia> diff(E, :a)
2×2 SparseArrays.SparseMatrixCSC{LinearCombinationEquation,Int64} with 2 stored entries:
  [1, 1]  =  ⟨a|ĥ + -⟨b|[a|b] + ⟨a|[b|b]
  [2, 1]  =  -⟨d|[c|b] + ⟨c|[d|b]

julia> diff(E, Conjugate(:b))
2×2 SparseArrays.SparseMatrixCSC{LinearCombinationEquation,Int64} with 2 stored entries:
  [1, 1]  =  ĥ|b⟩ + -[a|b]|a⟩ + [a|a]|b⟩
  [1, 2]  =  -[a|d]|c⟩ + [a|c]|d⟩
Base.diffMethod
diff(ome::OrbitalMatrixElement, o::Conjugate{O})

Vary the orbital matrix element ⟨abc...|Ω|xyz...⟩ with respect to ⟨o|.

Base.diffMethod
diff(ome::OrbitalMatrixElement, o::O)

Vary the orbital matrix element ⟨abc...|Ω|xyz...⟩ with respect to |o⟩.

Base.diffMethod
diff(ab::OrbitalOverlap, o::Conjugate{O})

Vary the orbital overlap ⟨a|b⟩ with respect to ⟨o|.

Base.diffMethod
diff(ab::OrbitalOverlap, o::O)

Vary the orbital overlap ⟨a|b⟩ with respect to |o⟩.

Base.diffMethod
diff(energy_expression, orbitals)

Derive the integro-differential equations for all orbitals, from energy_expression. Returns a MCEquationSystem, that gathers information on which integrals are common to all equations, for efficient equation solving.

Base.diffMethod
diff(fun!, energy_expression, orbitals)

Derive the integro-differential equations for all orbitals, from energy_expression; after each orbital equation has been generated fun! is applied to energy_expression with the current orbital as argument, which allows gradual modification of the energy expression. Returns a MCEquationSystem, that gathers information on which integrals are common to all equations, for efficient equation solving.

Base.inMethod
in(corbital::Conjugate, co::ContractedOperator)

Test if corbital is among the left set of orbitals of the ContractedOperator co. Useful to test if co is an integral operator with respect to corbital.

Base.inMethod
in(orbital, co::ContractedOperator)

Test if orbital is among the right set of orbitals of the ContractedOperator co. Useful to test if co is an integral operator with respect to orbital.

Base.isapproxMethod
Base.isapprox(a::NBodyMatrixElement, b::NBodyMatrixElement; kwargs...)

Test if a and b are approximately equal to each other, i.e. their terms all agree exactly, and the expansion coefficients are approximately equal. The actual comparison is performed by compare.

Base.iszeroMethod
iszero(me::EnergyExpressions.OrbitalMatrixElement{1,<:SpinOrbital{<:Orbital},OneBodyHamiltonian,<:SpinOrbital{<:Orbital}})

The matrix element vanishes if the spin-orbitals do not have the same spin.

Base.iszeroMethod
iszero(me::EnergyExpressions.OrbitalMatrixElement{2,<:SpinOrbital{<:Orbital},CoulombInteraction,<:SpinOrbital{<:Orbital}})

The matrix element vanishes if the (non-relativistic) spin-orbitals associated with the same coordinate do not have the same spin.

Base.lengthMethod
length(adjoint_slater_determinant)

Return the number of spin-orbitals in the adjoint Slater determinant.

Base.lengthMethod
length(slater_determinant)

Return the number of spin-orbitals in the Slater determinant.

EnergyExpressions.cofactorMethod
cofactor(k, l, A)

Calculate the cofactor of A, where the rows k and the columns l have been stricken out. The cofactor is calculated recursively, by expanding the minor determinants in cofactors, so this function should only be used in case it is known that the cofactor is non-zero.

EnergyExpressions.compareMethod
compare(a::NBodyMatrixElement, op, b::NBodyMatrixElement; kwargs...)

Compare the NBodyMatrixElements a and b for similarity; all the terms of a need to be present in b, and vice versa, and their expansion coefficients have to agree when compared using op.

This function is mainly designed for testing purposes, i.e. to compare an expression with a reference, generated otherwise. It may not be performant. It may also fail on edge cases.

EnergyExpressions.contractFunction
contract(orbital_matrix_element, i...)

Contract the orbital_matrix_element over all coordinates i....

EnergyExpressions.contractMethod
contract(ome::OrbitalMatrixElement{N}, i...)

Contract ome over all coordinates i.... length(i) cannot be larger than N.

EnergyExpressions.coupled_statesMethod
coupled_states(E[; i₀=1])

Find all states coupled by the energy expression E, starting from the state with index i₀. This can be useful to reduce the necessary basis or to generate invariant sets for split-operator propagation.

EnergyExpressions.detaxisMethod
detaxis(i::CartesianIndex{N})

Generate the axis index vector for the determinant minor, whose rows or columns represented by the CartesianIndex i should be omitted. Implemented via complement.

EnergyExpressions.distinct_permutationsMethod
distinct_permutations(fun::Function, ::NBodyOperator{N}, b)

Generate all distinct permutations p of b (which is expected to be of length N), and call fun(σ, b[p]) where σ=(-1)^p is the sign of the permutation p.

EnergyExpressions.invariant_setsMethod
invariant_sets(E)

Generate a list of all invariant sets, i.e. configurations that are coupled through the matrix elements of E.

Example

julia> E = sparse([1 1 0; 1 1 0; 0 0 1])
3×3 SparseMatrixCSC{Int64, Int64} with 5 stored entries:
 1  1  ⋅
 1  1  ⋅
 ⋅  ⋅  1

julia> invariant_sets(E)
2-element Vector{Vector{Int64}}:
 [1, 2]
 [3]
EnergyExpressions.isdependentMethod
isdependent(nbt::NBodyTerm, o)

Returns true if any of the factors comprising nbt is dependent on the orbital o. Not that the result is dependent on whether o is conjugated or not.

EnergyExpressions.isdependentMethod
isdependent(o::OrbitalMatrixElement, orbital)

Returns true if the OrbitalMatrixElement o depends on orbital.

Examples

julia> isdependent(EnergyExpressions.OrbitalMatrixElement((:a,), OneBodyHamiltonian(), (:b,)), :a)
false

julia> isdependent(EnergyExpressions.OrbitalMatrixElement((:a,), OneBodyHamiltonian(), (:b,)), Conjugate(:a))
true

julia> isdependent(EnergyExpressions.OrbitalMatrixElement((:a,), OneBodyHamiltonian(), (:b,)), :b)
true

julia> isdependent(EnergyExpressions.OrbitalMatrixElement((:a,:b,), CoulombInteraction(), (:c,:d)), :c)
true
EnergyExpressions.isdependentMethod
isdependent(o::OrbitalOverlap, orbital)

Returns true if the OrbitalOverlap o depends on orbital.

Examples

julia> isdependent(OrbitalOverlap(:a,:b), :a)
false

julia> isdependent(OrbitalOverlap(:a,:b), Conjugate(:a))
true

julia> isdependent(OrbitalOverlap(:a,:b), :b)
true
EnergyExpressions.non_zero_cofactorsMethod
non_zero_cofactors(sd, N, i, j)

Find all non-zero cofactors of the orbital overlap matrix between the Slater determinants i & j of sd, generated when striking out N rows & columns. This routine is tailored towards the case when few non-orthogonalities are present, e.g. approximately proportional to the number of orbitals.

Non-orthogonality between spin-orbitals is handled by dividing the them into two subspaces:

  1. The orthogonal spin-orbitals that are common to both Slater determinants (core orbitals),

  2. All non-orthogonal orbitals, and the orbitals which differ between the Slater determinants (i.e. holes of i and particles of j).

The relative phase between the Slater determinants is determined by group 2 alone, by permuting the particles to the positions of the holes, we find this phase. We can then formally permute them together to a diagonal block at lower-right corner of the orbital overlap matrix without incurring a phase change, since we need to permute the same number of rows and columns. We thus get this structure:

                 ╷       ╷
                 │ 1 │   │
 det(Sᵢⱼ) = (-)ᵏ │───┼───│
                 │   │ 2 │
                 ╵       ╵

where k is decided by the permutation necessary to put the particles in the positions of the holes.

Obviously, the determinant of the orbital matrix is now given by det(Sᵢⱼ) = (-)ᵏ*det(2), since we trivially have det(1)==1.

Depending on the rank of 2 (determined by the number of hole–particle pairs and which spin-orbitals are non-orthogonal), we need to strike out at least size(2,1)-rank(2) rows/columns from 2, and at most min(N,size(2,1)), i.e. for each value of n ∈ size(2,1)-rank(2):min(N,size(2,1)), we need to additionally strike out m = N - n rows from 1, but since the determinant of subspace 1 is unity, regardless of how many rows/columns we've stricken out, this is a trivial excercise. Of course, we also require that m ≤ size(1,1).

EnergyExpressions.nonzero_minorsMethod
nonzero_minors(N, overlap) -> (ks,ls)

Find all (distinct) minor determinants of order N of the orbital overlap matrix that do not vanish, i.e. all non-vanishing minors are guaranteed to be present, but not all of the returned minors are guaranteed to be non-zero. Vanishing minors returned arise when the overlap matrix is rank deficient, which is unlikely to happen when computing energy expressions, but must still be guarded against. This is most easily checked by actually calculating the cofactor, which is most likely desired anyway.

EnergyExpressions.numbodiesMethod
numbodies(::OrbitalOverlap)

Returns the number of bodies coupled by the zero-body operator in the orbital overlap, i.e. 0.

EnergyExpressions.numbodiesMethod
numbodies(::OrbitalMatrixElement{N})

Returns the number of bodies coupled by the operator, i.e. N.

EnergyExpressions.numbodiesMethod
numbodies(::NBodyOperator{N})

Returns the number of bodies coupled by the N-body operator, i.e. N.

EnergyExpressions.overlap_matrixMethod
overlap_matrix(a::Cfg, b::Cfg[, overlaps=[]]) where Cfg

Generate the single-particle orbital overlap matrix, between the orbitals in the configurations (e.g. Slater determinants) a and b. All orbitals are assumed to be orthogonal, except for those which are given in overlaps.

Examples

First we define two Slater determinants that have some orbitals in common:

julia> sa = SlaterDeterminant([:i, :j, :l,:k̃])
i(1)j(2)l(3)k̃(4) - i(1)j(2)l(4)k̃(3) - i(1)j(3)l(2)k̃(4) + i(1)j(3)l(4)k̃(2) + …  + i(4)j(1)l(3)k̃(2) + i(4)j(2)l(1)k̃(3) - i(4)j(2)l(3)k̃(1) - i(4)j(3)l(1)k̃(2) + i(4)j(3)l(2)k̃(1)

julia> sb = SlaterDeterminant([:i, :j, :k, :l̃])
i(1)j(2)k(3)l̃(4) - i(1)j(2)k(4)l̃(3) - i(1)j(3)k(2)l̃(4) + i(1)j(3)k(4)l̃(2) + …  + i(4)j(1)k(3)l̃(2) + i(4)j(2)k(1)l̃(3) - i(4)j(2)k(3)l̃(1) - i(4)j(3)k(1)l̃(2) + i(4)j(3)k(2)l̃(1)

The orbital overlap matrix by default is

julia> overlap_matrix(sa, sb)
4×4 SparseArrays.SparseMatrixCSC{EnergyExpressions.NBodyTerm,Int64} with 2 stored entries:
  [1, 1]  =  1
  [2, 2]  =  1

which has only two non-zero entries, since only two of the orbitals are common between the Slater determinants sa and sb.

We can then define that the orbitals and are non-orthogonal:

julia> overlap_matrix(sa, sb, [OrbitalOverlap(:k̃,:l̃)])
4×4 SparseArrays.SparseMatrixCSC{EnergyExpressions.NBodyTerm,Int64} with 3 stored entries:
  [1, 1]  =  1
  [2, 2]  =  1
  [4, 4]  =  ⟨k̃|l̃⟩

We can even specify that the orbital is non-orthogonal to itself (this can be useful when the is a linear combination of orthogonal orbitals):

julia> overlap_matrix(sa, sa, [OrbitalOverlap(:k̃,:k̃)])
4×4 SparseArrays.SparseMatrixCSC{EnergyExpressions.NBodyTerm,Int64} with 4 stored entries:
  [1, 1]  =  1
  [2, 2]  =  1
  [3, 3]  =  1
  [4, 4]  =  ⟨k̃|k̃⟩

Notice that this overlap matrix was calculated between the Slater determinant sa and itself.

EnergyExpressions.pushifmissing!Method
pushifmissing!(vector, element)

Push element to the end of vector, if not already present. Returns the index of element in vector.

EnergyExpressions.transformMethod
transform(f::Function, nbme::NBodyMatrixElement)

Transform integrals of the the N-body matrix element nbme according to the function f, which should accept a single NBodyTermFactor as its argument, and return a NBodyMatrixElement. This is useful for adapting energy expressions to specific symmetries of the system under consideration.

EnergyExpressions.transformMethod
transform(f::Function, nbt::NBodyTerm)

Transform integrals of the the N-body matrix element expansion term nbt according to the function f, which should accept a single NBodyTermFactor as its argument.

LinearAlgebra.detMethod
det(A)

Calculate the determinant of the matrix A whose elements are of the NBodyTerm type, by expanding the determinant along the first column. This is an expensive operation, and should only be done with relatively sparse matrices.

LinearAlgebra.ishermitianMethod
ishermitian(op::QuantumOperator)

By default, all QuantumOperators are Hermitian; this can be overridden for subtypes to explicitly declare an operator non-Hermitian.

EnergyExpressions.@above_diagonal_loopMacro
above_diagonal_loop(N, itersym, imax, args...)

Generate N Cartesian loops for the iteration variables itersym_{1:N}, where itersym_N ∈ 1:imax, itersym_{N-1} ∈ itersym_N+1:imax, etc, i.e. above the hyper-diagonal of the N-dimensional hypercube with the side imax. args... is passed on to Base.Cartesian._nloops. above_diagonal_loop is nestable.

Examples

julia> @above_diagonal_loop 2 i 3 begin
           println("==================================")
           println("i = ", Base.Cartesian.@ntuple 2 i)
           @above_diagonal_loop 2 j 3 begin
               println("j = ", Base.Cartesian.@ntuple 2 j)
           end
       end
==================================
i = (2, 1)
j = (2, 1)
j = (3, 1)
j = (3, 2)
==================================
i = (3, 1)
j = (2, 1)
j = (3, 1)
j = (3, 2)
==================================
i = (3, 2)
j = (2, 1)
j = (3, 1)
j = (3, 2)
EnergyExpressions.@anti_diagonal_loopMacro
anti_diagonal_loop(N, itersym, imax, args...)

Generate N Cartesian loops for the iteration variables itersym_{1:N}, where itersym_N ∈ 1:imax, itersym_{N-1} ∈ 1:imax\itersym_N, etc, i.e. no two iteration variables have the same values simultaneously. args... is passed on to Base.Cartesian._nloops; however, preexpr is already used to skip the diagonal elements. anti_diagonal_loop is nestable.

Examples

julia> @anti_diagonal_loop 3 i 3 begin
           println("-----------------------------")
           t = (Base.Cartesian.@ntuple 3 i)
           println("$t: ", allunique(t))
           @anti_diagonal_loop 2 j 2 begin
               u = (Base.Cartesian.@ntuple 2 j)
               println("$u: ", allunique(u))
           end
       end
-----------------------------
(3, 2, 1): true
(2, 1): true
(1, 2): true
-----------------------------
(2, 3, 1): true
(2, 1): true
(1, 2): true
-----------------------------
(3, 1, 2): true
(2, 1): true
(1, 2): true
-----------------------------
(1, 3, 2): true
(2, 1): true
(1, 2): true
-----------------------------
(2, 1, 3): true
(2, 1): true
(1, 2): true
-----------------------------
(1, 2, 3): true
(2, 1): true
(1, 2): true