API reference
TwoBody.Basis
TwoBody.BasisSet
TwoBody.ConstantPotential
TwoBody.ContractedBasis
TwoBody.CoulombPotential
TwoBody.DeltaPotential
TwoBody.ExponentialPotential
TwoBody.FiniteDifferenceMethod
TwoBody.FunctionPotential
TwoBody.GaussianBasis
TwoBody.GaussianPotential
TwoBody.GeometricBasisSet
TwoBody.Hamiltonian
TwoBody.KineticTerm
TwoBody.Laplacian
TwoBody.LinearPotential
TwoBody.NonRelativisticKinetic
TwoBody.Operator
TwoBody.PotentialTerm
TwoBody.PowerLawPotential
TwoBody.PrimitiveBasis
TwoBody.RelativisticCorrection
TwoBody.RelativisticKinetic
TwoBody.RestEnergy
TwoBody.SimpleGaussianBasis
TwoBody.UniformGridPotential
TwoBody.YukawaPotential
TwoBody.element
TwoBody.element
TwoBody.element
TwoBody.element
TwoBody.element
TwoBody.element
TwoBody.element
TwoBody.element
TwoBody.element
TwoBody.element
TwoBody.geometric
TwoBody.matrix
TwoBody.matrix
TwoBody.matrix
TwoBody.matrix
TwoBody.matrix
TwoBody.matrix
TwoBody.optimize
TwoBody.optimize
TwoBody.optimize
TwoBody.solve
TwoBody.solve
TwoBody.solve
TwoBody.solve
TwoBody.solve
TwoBody.Basis
— TypeBasis
is an abstract type.
TwoBody.BasisSet
— TypeBasisSet(basis1, basis2, ...)
\[\{ \phi_1, \phi_2, \phi_3, \cdots \}\]
The basis set is the input for Rayleigh–Ritz method. You can define the basis set like this:
\[\begin{aligned} \phi_1(r) &= \exp(-13.00773 ~r^2), \\ \phi_2(r) &= \exp(-1.962079 ~r^2), \\ \phi_3(r) &= \exp(-0.444529 ~r^2), \\ \phi_4(r) &= \exp(-0.1219492 ~r^2). \end{aligned}\]
BS = BasisSet(
SimpleGaussianBasis(13.00773),
SimpleGaussianBasis(1.962079),
SimpleGaussianBasis(0.444529),
SimpleGaussianBasis(0.1219492),
)
TwoBody.ConstantPotential
— TypeConstantPotential(constant=1)
\[+ c\]
Arguments | Symbol |
---|---|
constant | $c$ |
TwoBody.ContractedBasis
— TypeContractedBasis([c1, c2, ...], [basis1, basis2, ...])
\[\phi' = \sum_i c_i \phi_i\]
TwoBody.CoulombPotential
— TypeCoulombPotential(coefficient=1)
\[+ \frac{a}{r}\]
Arguments | Symbol |
---|---|
coefficient | $a$ |
TwoBody.DeltaPotential
— TypeDeltaPotential(coefficient=1)
\[+ a δ(r)\]
Arguments | Symbol |
---|---|
coefficient | $a$ |
TwoBody.ExponentialPotential
— TypeExponentialPotential(coefficient=1, exponent=1)
\[+ a \exp(- b r)\]
Arguments | Symbol |
---|---|
coefficient | $a$ |
exponent | $b$ |
TwoBody.FiniteDifferenceMethod
— TypeFiniteDifferenceMethod(Δr=0.1, rₘₐₓ=50.0, R=Δr:Δr:rₘₐₓ, l=0, direction=:c, solver=:LinearAlgebra)
Arguments | Symbol |
---|---|
Δr::Real | Radial grid spacing. A uniform grid spacing is used, $r_{i+1} = r_{i} + \Delta r$. |
rₘₐₓ::Real | This value is not directly used in the calculation, but it is used to determine the R . |
R::StepRangeLen | Radial grid |
l::Int | Angular momentum quantum number |
direction::Symbol | :c for central, :f for forward, :b for backward |
solver::Symbol | :LinearAlgebra or :ArnoldiMethod |
TwoBody.FunctionPotential
— TypeFunctionPotential(f)
\[+ f(r)\]
TwoBody.GaussianBasis
— TypeGaussianBasis(a=1, l=0, m=0)
\[\phi_i(r, θ, φ) = N _{il} r^l \exp(-a_i r^2) Y_l^m(θ, φ)\]
TwoBody.GaussianPotential
— TypeGaussianPotential(coefficient=1, exponent=1)
\[+ a \exp(- b r^2)\]
Arguments | Symbol |
---|---|
coefficient | $a$ |
exponent | $b$ |
TwoBody.GeometricBasisSet
— TypeGeometricBasisSet(basistype, r₁, rₙ, n; nₘᵢₙ=1, nₘₐₓ=n)
This is a basis set with exponentials generated by geometric()
.
TwoBody.Hamiltonian
— TypeHamiltonian(operator1, operator2, ...)
\[\hat{H} = \sum_i \hat{o}_i\]
The Hamiltonian is the input for each solver. This is an example for the non-relativistic Hamiltonian of hydrogen atom in atomic units:
\[\hat{H} = - \frac{1}{2} \nabla^2 - \frac{1}{r}\]
H = Hamiltonian(
NonRelativisticKinetic(ℏ =1 , m = 1),
CoulombPotential(coefficient = -1),
)
TwoBody.KineticTerm
— TypeKineticTerm <: Operator
is an abstract type.
TwoBody.Laplacian
— TypeLaplacian(coefficient=1)
\[+ a\nabla^2\]
Arguments | Symbol |
---|---|
coefficient | $a$ |
TwoBody.LinearPotential
— TypeLinearPotential(coefficient=1)
\[+ ar\]
Arguments | Symbol |
---|---|
coefficient | $a$ |
TwoBody.NonRelativisticKinetic
— TypeNonRelativisticKinetic(ℏ=1, m=1)
\[-\frac{\hbar^2}{2m} \nabla^2\]
TwoBody.Operator
— TypeOperator
is an abstract type.
TwoBody.PotentialTerm
— TypePotentialTerm <: Operator
is an abstract type.
TwoBody.PowerLawPotential
— TypePowerLawPotential(coefficient=1, exponent=1)
\[+ ar^n\]
Arguments | Symbol |
---|---|
coefficient | $a$ |
exponent | $n$ |
TwoBody.PrimitiveBasis
— TypePrimitiveBasis <: Basis
is an abstract type.
TwoBody.RelativisticCorrection
— TypeRelativisticCorrection(c=1, m=1, n=2)
The p^{2n} term of the Taylor expansion:
\[\begin{aligned} \sqrt{p^2 c^2 + m^2 c^4} =& m \times c^2 \\ &+ 1 / 2 / m \times p^2 (n=1) \\ &- 1 / 8 / m^3 / c^2 \times p^4 (n=2) \\ &+ 1 / 16 / m^5 / c^4 \times p^6 (n=3) \\ &- 5 / 128 / m^7 / c^6 \times p^8 (n=4) \\ &+ \cdots \end{aligned}\]
Use c = 137.035999177
(from 2022 CODATA) in the atomic units.
TwoBody.RelativisticKinetic
— TypeRelativisticKinetic(c=1, m=1)
\[\sqrt{p^2 c^2 + m^2 c^4} - m c^2\]
Use c = 137.035999177
(from 2022 CODATA) in the atomic units.
TwoBody.RestEnergy
— TypeTwoBody.SimpleGaussianBasis
— TypeSimpleGaussianBasis(a=1)
\[\phi_i(r) = \exp(-a_i r^2)\]
Note: This basis is not normalized. This is only for s-wave.
TwoBody.UniformGridPotential
— TypeUniformGridPotential(R, V)
TwoBody.YukawaPotential
— TypeYukawaPotential(coefficient=1, exponent=1)
\[+ \frac{a}{r} \exp(- b r)\]
Arguments | Symbol |
---|---|
coefficient | $a$ |
exponent | $b$ |
TwoBody.element
— Methodelement(o::ConstantPotential, SGB1::SimpleGaussianBasis, SGB2::SimpleGaussianBasis)
\[\begin{aligned} \langle \phi_{i} | c | \phi_{j} \rangle &= c \langle \phi_{i} | \phi_{j} \rangle \\ &= c \iiint \phi_{i}^*(r) \phi_{j}(r) ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= c \int_0^{2\pi} \mathrm{d}\varphi \int_0^\pi \sin\theta ~\mathrm{d}\theta \int_0^\infty r^{2} \mathrm{e}^{-(\alpha_i + \alpha_j) r^2} ~\mathrm{d}r \\ &= c \times 2\pi \times 2 \times \frac{1!!}{2^{2}} \sqrt{\frac{\pi}{a^{3}}} \\ &= \underline{c \left( \frac{\pi}{\alpha_i + \alpha_j} \right)^{3/2}} \end{aligned}\]
Integral Formula:
\[\int_0^{\infty} r^{2n} \exp \left(-a r^2\right) ~\mathrm{d}r = \frac{(2n-1)!!}{2^{n+1}} \sqrt{\frac{\pi}{a^{2n+1}}}\]
TwoBody.element
— Methodelement(o::CoulombPotential, SGB1::SimpleGaussianBasis, SGB2::SimpleGaussianBasis)
\[\begin{aligned} \langle \phi_{i} | \frac{1}{r} | \phi_{j} \rangle &= \iiint \phi_{i}^*(r) \times \frac{1}{r} \times \phi_{j}(r) ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= \int_0^{2\pi} \mathrm{d}\varphi \int_0^\pi \sin\theta ~\mathrm{d}\theta \int_0^\infty r \mathrm{e}^{-(\alpha_i + \alpha_j) r^2} ~\mathrm{d}r \\ &= 2\pi \times 2 \times \frac{0!}{2 (\alpha_i + \alpha_j)} \\ &= \underline{\frac{2\pi}{\alpha_i + \alpha_j}} \end{aligned}\]
Integral Formula:
\[\int_0^{\infty} r^{2n+1} \exp \left(-a r^2\right) ~\mathrm{d}r = \frac{n!}{2 a^{n+1}}\]
TwoBody.element
— Methodelement(o::GaussianPotential, SGB1::SimpleGaussianBasis, SGB2::SimpleGaussianBasis)
\[\begin{aligned} \langle \phi_{i} | \exp(-br^2) | \phi_{j} \rangle &= \iiint \phi_{i}^*(r) \times \exp(-br^2) \times \phi_{j}(r) ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= \int_0^{2\pi} \mathrm{d}\varphi \int_0^\pi \sin\theta ~\mathrm{d}\theta \int_0^\infty r^2 \mathrm{e}^{-(b+\alpha_i + \alpha_j) r^2} ~\mathrm{d}r \\ &= 2\pi \times 2 \times \frac{1!!}{2^{2}} \sqrt{\frac{\pi}{(b + \alpha_i + \alpha_j)^{2\cdot1+1}}} \\ &= \underline{\left( \frac{\pi}{b + \alpha_i + \alpha_j} \right)^{3/2}} \end{aligned}\]
Integral Formula:
\[\int_0^{\infty} r^{2n} \exp \left(-a r^2\right) ~\mathrm{d}r = \frac{(2n-1)!!}{2^{n+1}} \sqrt{\frac{\pi}{a^{2n+1}}}\]
TwoBody.element
— Methodelement(o::Hamiltonian, SGB1::SimpleGaussianBasis, SGB2::SimpleGaussianBasis)
\[\begin{aligned} H_{ij} &= \langle \phi_{i} | \hat{H} | \phi_{j} \rangle \\ &= \langle \phi_{i} | \sum_k \hat{o}_k | \phi_{j} \rangle \\ &= \sum_k \langle \phi_{i} | \hat{o}_k | \phi_{j} \rangle \\ \end{aligned}\]
TwoBody.element
— Methodelement(o::LinearPotential, SGB1::SimpleGaussianBasis, SGB2::SimpleGaussianBasis)
\[\begin{aligned} \langle \phi_{i} | r | \phi_{j} \rangle &= \iiint \phi_{i}^*(r) \times r \times \phi_{j}(r) ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= \int_0^{2\pi} \mathrm{d}\varphi \int_0^\pi \sin\theta ~\mathrm{d}\theta \int_0^\infty r^3 \mathrm{e}^{-(\alpha_i + \alpha_j) r^2} ~\mathrm{d}r \\ &= 2\pi \times 2 \times \frac{1!}{2 (\alpha_i + \alpha_j)^{2}} \\ &= \underline{\frac{2\pi}{(\alpha_i + \alpha_j)^2}} \end{aligned}\]
Integral Formula:
\[\int_0^{\infty} r^{2n+1} \exp \left(-a r^2\right) ~\mathrm{d}r = \frac{n!}{2 a^{n+1}}\]
TwoBody.element
— Methodelement(o::NonRelativisticKinetic, SGB1::SimpleGaussianBasis, SGB2::SimpleGaussianBasis)
\[\begin{aligned} T_{ij} = \langle \phi_{i} | \hat{T} | \phi_{j} \rangle &= \iiint \mathrm{e}^{-\alpha_i r^2} \left[ -\frac{\hbar^2}{2\mu} \nabla^2 \right] \mathrm{e}^{-\alpha_j r^2} ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= -\frac{\hbar^2}{2\mu} \iiint \mathrm{e}^{-\alpha_i r^2} \left[ \nabla^2 \right] \mathrm{e}^{-\alpha_j r^2} ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= -\frac{\hbar^2}{2\mu} \iiint \mathrm{e}^{-\alpha_i r^2} \left[ -6\alpha_j + 4\alpha_j^2 r^2 \right] \mathrm{e}^{-\alpha_j r^2} ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= -\frac{\hbar^2}{2\mu} \iint \sin\theta ~\mathrm{d}\theta \mathrm{d}\varphi \int \left[ -6\alpha_j + 4\alpha_j^2 r^2 \right] r^2 \mathrm{e}^{-(\alpha_i + \alpha_j) r^2} ~\mathrm{d}r \\ &= -\frac{\hbar^2}{2\mu} \cdot 4\pi \left[ -6\alpha_j \mathrm{GGI}(2, \alpha_i + \alpha_j) +4\alpha_j^2 \mathrm{GGI}(4, \alpha_i + \alpha_j) \right] \\ &= -\frac{\hbar^2}{2\mu} \cdot 4\pi \left[ -6\alpha_j \frac{\Gamma\left( \frac{3}{2} \right)}{2 (\alpha_i + \alpha_j)^{\frac{3}{2}}} +4\alpha_j^2 \frac{\Gamma\left( \frac{5}{2} \right)}{2 (\alpha_i + \alpha_j)^{\frac{5}{2}}} \right] \\ &= -\frac{\hbar^2}{2\mu} \cdot 4\pi \left[ -6\alpha_j \frac{ \sqrt{\pi}/2}{2 (\alpha_i + \alpha_j)^{\frac{3}{2}}} +4\alpha_j^2 \frac{3\sqrt{\pi}/4}{2 (\alpha_i + \alpha_j)^{\frac{5}{2}}} \right] \\ &= -\frac{\hbar^2}{2\mu} \cdot 4\pi \left[ \frac{\alpha_j}{\alpha_i + \alpha_j} - 1 \right] \cdot 6 \alpha_j \cdot \frac{\sqrt{\pi}/2}{2 (\alpha_i + \alpha_j)^{\frac{3}{2}}} \\ &= -\frac{\hbar^2}{2\mu} \cdot 4\pi \left[ - \frac{\alpha_i}{\alpha_i + \alpha_j} \right] \cdot 6 \alpha_j \cdot \frac{\sqrt{\pi}/2}{2 (\alpha_i + \alpha_j)^{\frac{3}{2}}} \\ &= \underline{ \frac{\hbar^2}{2\mu} \cdot 6 \cdot \frac{\alpha_i \alpha_j \pi^{\frac{3}{2}}}{(\alpha_i + \alpha_j)^{\frac{5}{2}}} } \end{aligned}\]
or
\[\begin{aligned} T_{ij} = \langle \phi_{i} | \hat{T} | \phi_{j} \rangle &= \iiint \mathrm{e}^{-\alpha_i r^2} \left[ -\frac{\hbar^2}{2\mu} \nabla^2 \right] \mathrm{e}^{-\alpha_j r^2} ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= -\frac{\hbar^2}{2\mu} \iiint \mathrm{e}^{-\alpha_i r^2} \nabla^2 \mathrm{e}^{-\alpha_j r^2} ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= \frac{\hbar^2}{2\mu} \iiint \left[ \nabla \mathrm{e}^{-\alpha_i r^2} \right] \left[ \nabla \mathrm{e}^{-\alpha_j r^2} \right] ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= \frac{\hbar^2}{2\mu} \iiint \left[ -2 \alpha_i r \mathrm{e}^{-\alpha_i r^2} \right] \left[ -2 \alpha_j r \mathrm{e}^{-\alpha_j r^2} \right] ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= \frac{\hbar^2}{2\mu} \cdot 4 \alpha_i \alpha_j \iiint \left[ r \mathrm{e}^{-\alpha_i r^2} \right] \left[ r \mathrm{e}^{-\alpha_j r^2} \right] ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= \frac{\hbar^2}{2\mu} \cdot 4 \alpha_i \alpha_j \iint \sin\theta ~\mathrm{d}\theta \mathrm{d}\varphi \int r^4 \mathrm{e}^{- (\alpha_i + \alpha_j) r^2} ~\mathrm{d}r \\ &= \frac{\hbar^2}{2\mu} \cdot 4 \alpha_i \alpha_j \cdot 4 \pi \cdot \mathrm{GGI}(4, \alpha_i + \alpha_j) \\ &= \frac{\hbar^2}{2\mu} \cdot 4 \alpha_i \alpha_j \cdot 4 \pi \cdot \frac{\Gamma\left( \frac{5}{2} \right)}{2 (\alpha_i + \alpha_j)^{\frac{5}{2}}} \\ &= \frac{\hbar^2}{2\mu} \cdot 4 \alpha_i \alpha_j \cdot 4 \pi \cdot \frac{3\sqrt{\pi}/4}{2 (\alpha_i + \alpha_j)^{\frac{5}{2}}} \\ &= \underline{ \frac{\hbar^2}{2\mu} \cdot 6 \cdot \frac{\alpha_i \alpha_j \pi^{\frac{3}{2}}}{(\alpha_i + \alpha_j)^{\frac{5}{2}}} } \end{aligned}\]
TwoBody.element
— Methodelement(o::PowerLawPotential, SGB1::SimpleGaussianBasis, SGB2::SimpleGaussianBasis)
\[\begin{aligned} \langle \phi_{i} | r^n | \phi_{j} \rangle &= \iiint \phi_{i}^*(r) \times r^n \times \phi_{j}(r) ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= \int_0^{2\pi} \mathrm{d}\varphi \int_0^\pi \sin\theta ~\mathrm{d}\theta \int_0^\infty r^{n+2} \mathrm{e}^{-(\alpha_i + \alpha_j) r^2} ~\mathrm{d}r \\ &= 2\pi \times 2 \times \frac{\Gamma\left( \frac{n+3}{2} \right)}{2 (\alpha_i + \alpha_j)^{\frac{n+3}{2}}} \\ &= \underline{2\pi\frac{\Gamma\left( \frac{n+3}{2} \right)}{(\alpha_i + \alpha_j)^{\frac{n+3}{2}}}} \end{aligned}\]
Integral Formula:
\[\int_0^{\infty} r^{n} \exp \left(-a r^2\right) ~\mathrm{d}r = \frac{\Gamma\left( \frac{n+1}{2} \right)}{2 a^{\frac{n+1}{2}}}\]
TwoBody.element
— Methodelement(o::RestEnergy, SGB1::SimpleGaussianBasis, SGB2::SimpleGaussianBasis)
\[\begin{aligned} \langle \phi_{i} | mc^2 | \phi_{j} \rangle &= mc^2 \langle \phi_{i} | \phi_{j} \rangle \\ &= mc^2 \iiint \phi_{i}^*(r) \phi_{j}(r) ~r^2 \sin\theta ~\mathrm{d}r \mathrm{d}\theta \mathrm{d}\varphi \\ &= mc^2 \int_0^{2\pi} \mathrm{d}\varphi \int_0^\pi \sin\theta ~\mathrm{d}\theta \int_0^\infty r^{2} \mathrm{e}^{-(\alpha_i + \alpha_j) r^2} ~\mathrm{d}r \\ &= mc^2 \times 2\pi \times 2 \times \frac{1!!}{2^{2}} \sqrt{\frac{\pi}{a^{3}}} \\ &= \underline{mc^2 \left( \frac{\pi}{\alpha_i + \alpha_j} \right)^{3/2}} \end{aligned}\]
Integral Formula:
\[\int_0^{\infty} r^{2n} \exp \left(-a r^2\right) ~\mathrm{d}r = \frac{(2n-1)!!}{2^{n+1}} \sqrt{\frac{\pi}{a^{2n+1}}}\]
TwoBody.element
— Methodelement(SGB1::SimpleGaussianBasis, SGB2::SimpleGaussianBasis)
\[\begin{aligned} S_{ij} = \langle \phi_{i} | \phi_{j} \rangle &= \int \phi_{i}^*(r) \phi_{j}(r) \mathrm{d} \pmb{r} \\ &= \iiint \mathrm{e}^{-\alpha_i r^2} \mathrm{e}^{-\alpha_j r^2} ~r^2 \sin\theta ~ \mathrm{d} r \mathrm{d} \theta \mathrm{d} \varphi \\ &= \int_0^{2\pi} \mathrm{d}\varphi \int_0^\pi \sin\theta ~\mathrm{d}\theta \int_0^\infty r^{2} \mathrm{e}^{-(\alpha_i + \alpha_j) r^2} ~\mathrm{d}r \\ &= 2\pi \times 2 \times \frac{1!!}{2^{2}} \sqrt{\frac{\pi}{a^{3}}} \\ &= \underline{\left( \frac{\pi}{\alpha_i + \alpha_j} \right)^{3/2}} \end{aligned}\]
Integral Formula:
\[\int_0^{\infty} r^{2n} \exp \left(-a r^2\right) ~\mathrm{d}r = \frac{(2n-1)!!}{2^{n+1}} \sqrt{\frac{\pi}{a^{2n+1}}}\]
TwoBody.element
— Methodelement(o::Laplacian, SGB1::SimpleGaussianBasis, SGB2::SimpleGaussianBasis)
\[\begin{aligned} \langle \phi_{i} | \nabla^2 | \phi_{j} \rangle = \underline{ -6 \frac{\alpha_i \alpha_j \pi^{\frac{3}{2}}}{(\alpha_i + \alpha_j)^{\frac{5}{2}}} } \end{aligned}\]
```
TwoBody.geometric
— MethodExponents of Gaussian basis functions are given by geometric progression:
\[\begin{aligned} & v_i = \frac{1}{r_i^2}, \\ & r_i = r_1 a^{i-1}. \end{aligned}\]
This function return array of $\nu_i$:
\[(r_1, r_{n}, n, n_\mathrm{max}) \mapsto (\nu_1, \nu_2, \cdots, \nu_{n-1}, \nu_n, \nu_{n+1}, \cdots, \nu_{n_\mathrm{max}})\]
Usually $n = n_\mathrm{max}$. Set $n<n_\mathrm{max}$ if you want to extend the geometric progression.
Examples:
julia> ν = TwoBody.geometric(0.1, 10.0, 5)
5-element Vector{Float64}:
100.0
10.0
0.9999999999999997
0.09999999999999996
0.009999999999999995
julia> ν = TwoBody.geometric(0.1, 10.0, 5, nₘₐₓ = 10)
10-element Vector{Float64}:
100.0
10.0
0.9999999999999997
0.09999999999999996
0.009999999999999995
0.0009999999999999994
9.999999999999994e-5
9.999999999999992e-6
9.999999999999991e-7
9.999999999999988e-8
TwoBody.matrix
— Methodmatrix(basisset::BasisSet)
This function returns the overlap matrix $\pmb{S}$. The element is written as $S_{ij} = \langle \phi_{i} | \phi_{j} \rangle$.
TwoBody.matrix
— Methodmatrix(hamiltonian::Hamiltonian, basisset::BasisSet)
This function returns the Hamiltonian matrix $\pmb{H}$. The element is written as $H_{ij} = \langle \phi_{i} | \hat{H} | \phi_{j} \rangle$.
TwoBody.matrix
— Methodmatrix(o::Hamiltonian, method::FiniteDifferenceMethod)
The matrix for the Hamiltonian is a sum of matrices for each term,
\[\pmb{H} = \sum_i \pmb{O}_i.\]
TwoBody.matrix
— Methodmatrix(o::NonRelativisticKinetic, method::FiniteDifferenceMethod)
We use the shorthand notation $\psi'(r) = \frac{\mathrm{d}\psi}{\mathrm{d}r}(r)$ and $\psi''(r) = \frac{\mathrm{d}^{2}\psi}{\mathrm{d}r^{2}}(r)$. For the uniform grid spacing ($r_{i+1} = r_{i} + \Delta r$), the finite difference for the first derivative,
\[\frac{\mathrm{d}\psi}{\mathrm{d}r}(r) = \frac{\psi(r+\Delta r) - \psi(r-\Delta r)}{2\Delta r} + O(\Delta r^{2})\]
is written as
\[\left(\begin{array}{ccccc} \psi'(r_1) \\ \psi'(r_2) \\ \psi'(r_3) \\ \psi'(r_4) \\ \vdots \end{array}\right) \simeq \frac{1}{2\Delta r} \left(\begin{array}{ccccc} 0 & 1 & 0 & 0 &\ldots \\ -1 & 0 & 1 & 0 &\ldots \\ 0 & -1 & 0 & 1 &\ldots \\ 0 & 0 & -1 & 0 &\ldots \\ \vdots & \vdots & \vdots & \vdots & \ddots \\ \end{array}\right) \left(\begin{array}{ccccccc} \psi(r_1) \\ \psi(r_2) \\ \psi(r_3) \\ \psi(r_4) \\ \vdots \end{array}\right),\]
and the finite difference for the second derivative,
\[\frac{\mathrm{d}^{2}}{\mathrm{d}r^{2}}(r) = \frac{\psi(r+\Delta r) - 2f(r) + \psi(r-\Delta r)}{\Delta r^{2}} + O(\Delta r^{2}).\]
is written as
\[\left(\begin{array}{ccccc} \psi''(r_1) \\ \psi''(r_2) \\ \psi''(r_3) \\ \psi''(r_4) \\ \vdots \end{array}\right) \simeq \frac{1}{\Delta r^2} \left(\begin{array}{ccccccc} -2 & 1 & 0 & 0 & \ldots \\ 1 & -2 & 1 & 0 & \ldots \\ 0 & 1 & -2 & 1 & \ldots \\ 0 & 0 & 1 & -2 & \ldots \\ \vdots & \vdots & \vdots & \vdots & \ddots \end{array}\right) \left(\begin{array}{ccccccc} \psi(r_1) \\ \psi(r_2) \\ \psi(r_3) \\ \psi(r_4) \\ \vdots \end{array}\right).\]
Similarly, the matrix for the kinetic energy,
\[\hat{T} = -\frac{\hbar^2}{2\mu} \left[ \frac{\partial^2}{\partial r^2} + \frac{2}{r} \frac{\partial}{\partial r} - \frac{l(l+1)}{r^2} \right]\]
is written as
\[\pmb{T} = - \frac{\hbar^2}{2\mu} \left[ \frac{1}{{\Delta r}^2} \left(\begin{array}{ccccccc} -2 & 1 & 0 & \ldots \\ 1 & -2 & 1 & \ldots \\ 0 & 1 & -2 & \ldots \\ \vdots & \vdots & \vdots & \ddots \\ \end{array}\right) + \left(\begin{array}{ccccccc} 2/r_1 & 0 & 0 & \ldots \\ 0 & 2/r_2 & 0 & \ldots \\ 0 & 0 & 2/r_3 & \ldots \\ \vdots & \vdots & \vdots & \ddots \\ \end{array}\right) \frac{1}{\Delta r} \left(\begin{array}{ccccccc} 0 & 1 & 0 & \ldots \\ -1 & 0 & 1 & \ldots \\ 0 & -1 & -1 & \ldots \\ \vdots & \vdots & \vdots & \ddots \\ \end{array}\right) - l(l+1) \left(\begin{array}{ccccccc} {r_1}^2 & 0 & 0 & \ldots \\ 0 & {r_2}^2 & 0 & \ldots \\ 0 & 0 & {r_3}^2 & \ldots \\ \vdots & \vdots & \vdots & \ddots \\ \end{array}\right) \right].\]
TwoBody.matrix
— Methodmatrix(o::RestEnergy, method::FiniteDifferenceMethod)
The matrix for the rest energy $mc^2$ is a diagonal matrix,
\[mc^2 \left(\begin{array}{ccccccc} 1 & 0 & 0 & \ldots \\ 0 & 1 & 0 & \ldots \\ 0 & 0 & 1 & \ldots \\ \vdots & \vdots & \vdots & \ddots \\ \end{array}\right).\]
TwoBody.matrix
— Methodmatrix(o::PotentialTerm, method::FiniteDifferenceMethod)
The matrix for the potential energy $V(r)$ is a diagonal matrix,
\[\pmb{V} = \left(\begin{array}{ccccccc} V(r_1) & 0 & 0 & \ldots \\ 0 & V(r_2) & 0 & \ldots \\ 0 & 0 & V(r_3) & \ldots \\ \vdots & \vdots & \vdots & \ddots \\ \end{array}\right).\]
TwoBody.optimize
— Methodfunction optimize(hamiltonian::Hamiltonian, basisset::BasisSet; perturbation=Hamiltonian(), info=4, progress=true, optimizer=Optim.NelderMead(), options...)
This function minimizes the energy by changing the exponents of the basis functions using Optim.jl.
\[\frac{\partial E}{\partial a_i} = 0\]
TwoBody.optimize
— Methodoptimize(hamiltonian::Hamiltonian, basis::Basis; perturbation=Hamiltonian(), info=4, optimizer=Optim.NelderMead())
This a optimizer for 1-basis calculations. This function returns optimize(hamiltonian, BasisSet(basis); perturbation=perturbation, info=info, progress=progress, optimizer=optimizer, options...)
.
TwoBody.optimize
— Methodoptimize(hamiltonian::Hamiltonian, basisset::GeometricBasisSet; perturbation=Hamiltonian(), info=4, optimizer=Optim.NelderMead())
This function minimizes the energy by optimizing $r_1$ and $r_n$ using Optim.jl.
\[\frac{\partial E}{\partial r_1} = \frac{\partial E}{\partial r_n} = 0\]
TwoBody.solve
— Functionsolve(hamiltonian::Hamiltonian, wavefunction::Function, method::FiniteDifferenceMethod, info=4, nₘₐₓ=4)
TwoBody.solve
— Methodsolve(hamiltonian::Hamiltonian, basisset::BasisSet)
This function returns the eigenvalues $E$ and eigenvectors $\pmb{c}$ for
\[\pmb{H} \pmb{c} = E \pmb{S} \pmb{c}.\]
The Hamiltonian matrix is defined as $H_{ij} = \langle \phi_{i} | \hat{H} | \phi_{j} \rangle$. The overlap matrix is defined as $S_{ij} = \langle \phi_{i} | \phi_{j} \rangle$.
TwoBody.solve
— Methodsolve(hamiltonian::Hamiltonian, basis::Basis; perturbation=Hamiltonian(), info=4)
This a solver for 1-basis calculations. This function returns solve(hamiltonian, BasisSet(basis); perturbation=perturbation, info=info)
.
TwoBody.solve
— Methodsolve(hamiltonian::Hamiltonian, method::FiniteDifferenceMethod; perturbation=Hamiltonian(), info=4, nₘₐₓ=4)
This method solve the eigenvalue problem for the Hamiltonian discretized as a sparse matrix with finite difference approximation,
\[\pmb{H} \pmb{\psi} = E \pmb{\psi}.\]
The eigenvalue $E$ is an approximation of the exact energy and the eigenvector $\pmb{\psi}$ is a vector of the approximated values of the exact wavefunction $\psi(r)$ on points of the grid,
\[\pmb{\psi} = \left(\begin{array}{c} \psi(r_1) \\ \psi(r_2) \\ \psi(r_3) \\ \vdots \\ \end{array}\right).\]
TwoBody.solve
— Methodsolve(hamiltonian::Hamiltonian, basisset::GeometricBasisSet; perturbation=Hamiltonian(), info=4)
This function is a wrapper for solve(hamiltonian::Hamiltonian, basisset::BasisSet, ...)
.