Core Types
Quiqbox.DI
— TypeDI{F<:Function} <: StructFunction{F}
A "dressed-up" itself
that carries the information of a function (of type F
). For an instance di=Quiqbox.DI(someFunction)
where someFunction isa Function
, di(anyArgument) === anyArgument
and di() === someFunction
.
Quiqbox.ParamBox
— TypeParamBox{T, V, F<:Function} <: DifferentiableParameter{T, ParamBox}
Parameter container that can enable differentiation.
≡≡≡ Field(s) ≡≡≡
data::Array{Pair{Array{T, 0}, Symbol}, 0}
: The container of the input variable data in the form of a Pair
of its value container and symbol) stored in the ParamBox
. The value of the input variable can be accessed by syntax []
; to modify it, for example for a pb::ParamBox{T}
, use the syntax pb[] = newVal
where newVal
is the new value that is or can be converted into type T
.
map::Union{F,
DI
{F}}
: The mapping of the value of the input variable (i.e. the input value) within the same domain (.map(::T)->T
). The result (i.e., the value of the output variable, or the "output value") can be accessed by syntax ()
.
canDiff::Array{Bool, 0}
: Indicator of whether the output variable is "marked" as differentiable with respect to the input variable in the differentiation process. In other words, it determines whether the output variable represented by the ParamBox
is treated as a dependent variable or an independent variable.
index::Array{Union{Int, Nothing}, 0}
: Additional index assigned to the ParamBox
.
≡≡≡ Initialization Method(s) ≡≡≡
ParamBox(inVar::Union{T, Array{T, 0}}, outSym::Symbol=:undef,
inSym::Symbol=Symbol(IVsymSuffix, outSym);
index::Union{Int, Nothing}=nothing, canDiff::Bool=false) where {T} ->
ParamBox{T, outSym, typeof(itself)}
ParamBox(inVar::Union{T, Array{T, 0}}, outSym::Symbol, mapFunction::Function,
inSym::Symbol=Symbol(IVsymSuffix, outSym);
index::Union{Int, Nothing}=nothing, canDiff::Bool=true) where {T} ->
ParamBox{T, outSym}
=== Positional argument(s) ===
inVar::Union{T, Array{T, 0}}
: The value or the container of the input variable to be stored. If the latter is the type of data
, then it will directly used to construct .data[]
with without any copy.
outSym::Symbol
: The symbol of the output variable represented by the constructed ParamBox
. It's equal to the type parameter V
of the constructed ParamBox
.
inSym::Symbol
: The symbol of the input variable held by the constructed ParamBox
.
mapFunction::Function
: The mapping (mapFunction(::T)->T
) of the input variable, which will be assigned to the field .map
. When mapFunction
is not provided, .map
is set to itself
that maps the input variable to an identical-valued output variable.
=== Keyword argument(s) ===
index::Union{Int, Nothing}
: The index of the constructed ParamBox
. It's should be left with its default value unless the user plans to utilize the index of a ParamBox
for specific application other than differentiation.
canDiff::Bool
: Determine whether the output variable is marked as "differentiable" with respect to the input variable.
≡≡≡ Example(s) ≡≡≡
julia> ParamBox(1.0)
ParamBox{Float64, :undef, …}{0}[∂][undef]⟦=⟧[1.0]
julia> ParamBox(1.0, :a)
ParamBox{Float64, :a, …}{0}[∂][a]⟦=⟧[1.0]
julia> ParamBox(1.0, :a, abs)
ParamBox{Float64, :a, …}{1}[𝛛][x_a]⟦→⟧[1.0]
NOTE 1: The markers "[∂][IV]
" in the printed info indicate the differentiability and the name (the symbol with an assigned index if applied) respectively of the independent variable tied to the ParamBox
. When the ParamBox
is marked as non-differentiable, "[∂]
" is grey and IV
corresponds to the name of the output variable; when the ParamBox
is marked as differentiable, "[∂]
" becomes a green "[𝛛]
", and IV
corresponds to the name of the stored input variable.
NOTE 2: The output variable of a ParamBox
is normally used to differentiate a parameter functional (e.g., the Hartree–Fock energy). However, the derivative with respect to the stored input variable can also be computed to when the ParamBox
is marked as differentiable.
Quiqbox.BasisFunc
— TypeBasisFunc{T, D, 𝑙, GN, PT} <: FloatingGTBasisFuncs{T, D, 𝑙, GN, PT, 1}
A (floating) Gaussian-type basis function with its center assigned to a defined coordinate.
≡≡≡ Field(s) ≡≡≡
center::SpatialPoint{T, D, PT}
: The D
-dimensional center.
gauss::NTuple{GN, GaussFunc{T, <:Any}}
: Gaussian functions within the basis function.
l::Tuple{LTuple{D, 𝑙}}
: Cartesian representation of the angular momentum. E.g., LTuple{3, 1}((1, 0, 0))
(X¹Y⁰Z⁰) would correspond to an specific angular momentum configuration where the sum of all the components is 𝑙=1
.
normalizeGTO::Bool
: Whether each GaussFunc
inside will be normalized in calculations.
param::NTuple{D+GN*2, ParamBox}
: All the tunable parameters stored in the basis function.
≡≡≡ Initialization Method(s) ≡≡≡
BasisFunc(cen::SpatialPoint{T, D, PT},
gs::Tuple{AbstractGaussFunc{T}, Vararg{AbstractGaussFunc{T}, GN-1}},
l::Union{Tuple{LTuple{D, 𝑙}}, LTuple{D, 𝑙}}, normalizeGTO::Bool) where
{T, D, PT, 𝑙, GN} ->
BasisFunc{T, D, 𝑙, GN, PT}
BasisFunc(cen::SpatialPoint{T, D, PT}, gs::AbstractGaussFunc{T},
l::Union{Tuple{LTuple{D, 𝑙}}, LTuple{D, 𝑙}}, normalizeGTO::Bool) where
{T, D, PT, 𝑙, GN} ->
BasisFunc{T, D, 𝑙, 1, PT}
Quiqbox.BasisFuncMix
— TypeBasisFuncMix{T, D, BN, BFT<:BasisFunc{T, D}} <: CompositeGTBasisFuncs{T, D, BN, 1}
Sum of multiple FloatingGTBasisFuncs{<:Any, <:Any, <:Any, <:Any, <:Any, 1}
without any reformulation, treated as one basis function in the calculations.
≡≡≡ Field(s) ≡≡≡
BasisFunc::NTuple{BN, BFT}
: Basis functions used to sum up.
param::Tuple{Vararg{ParamBox}}
: Contained parameters.
≡≡≡ Initialization Method(s) ≡≡≡
BasisFuncMix(bfs::Union{Tuple{Vararg{T}}, AbstractArray{T}}) where
{T<:FloatingGTBasisFuncs{<:Any, <:Any, <:Any, <:Any, <:Any, 1}} ->
BasisFuncMix
Quiqbox.BasisFuncs
— TypeBasisFuncs{T, D, 𝑙, GN, PT, ON} <: FloatingGTBasisFuncs{T, D, 𝑙, GN, PT, ON}
A collection of basis functions with identical parameters except having different orientations within a specified subshell (i.e. same total orbital angular momentum). It has the same fields as BasisFunc
. Specifically, for l
, its size ON
can be no less than 1 and no larger than the size of the corresponding subshell.
Quiqbox.GTBasis
— TypeGTBasis{T, D, BN, BFT<:GTBasisFuncs{T, D, 1}} <: BasisSetData{T, D, BFT}
The container of basis set information.
≡≡≡ Field(s) ≡≡≡
basis::Vector{BFT}
: Stored basis set with its size equal to BN
by default (assuming no deleting or appending to basis set has been made).
S::Matrix{T}
: Overlap matrix.
Te::Matrix{T}
: Kinetic energy part of the electronic core Hamiltonian.
eeI::Array{T, 4}
: Electron-electron interaction.
≡≡≡ Initialization Method(s) ≡≡≡
GTBasis(basis::Union{Tuple{Vararg{GTBasisFuncs{T, D}}},
AbstractVector{<:GTBasisFuncs{T, D}}}) where {T, D} ->
GTBasis{T, D}
Construct a GTBasis
given a basis set.
Quiqbox.GaussFunc
— TypeGaussFunc{T, F1, F2} <: AbstractGaussFunc{T}
A contracted primitive Gaussian-type function.
≡≡≡ Field(s) ≡≡≡
xpn::ParamBox{T, :α, F1}
:The exponent coefficient.
con::ParamBox{T, :d, F2}
: The contraction coefficient.
param::Tuple{ParamBox{T, α}, ParamBox{T, d}}
: The parameter containers inside a GaussFunc
.
≡≡≡ Initialization Method(s) ≡≡≡
GaussFunc(e::Union{T, Array{T, 0}, ParamBox{T}}, d::Union{T, ParamBox{T}}) where
{T<:AbstractFloat} ->
GaussFunc{T}
Quiqbox.SpatialPoint
— TypeSpatialPoint{T, D, PT} <: AbstractSpatialPoint{T, D}
A D
-dimensional spatial point.
≡≡≡ Field(s) ≡≡≡
param::PT
: A Tuple
of ParamBox
s as the components of the spatial coordinate.
marker::Symbol
: A marker that indicates the purpose or meaning of the constructed SpatialPoint
. The default marker is set to :point
.
≡≡≡ Initialization Method(s) ≡≡≡
SpatialPoint(pbs::Union{Tuple{ParamBox{T, :X, Fx}} where Fx, Tuple{ParamBox{T, :X, Fx}, ParamBox{T, :Y, Fy}} where {Fx, Fy}, Tuple{ParamBox{T, :X, Fx}, ParamBox{T, :Y, Fy}, ParamBox{T, :Z, Fz}} where {Fx, Fy, Fz}} where T) -> SpatialPoint
Quiqbox.HFconfig
— TypeHFconfig{T1, HFT, F, T2, L, MS} <: ConfigBox{T1, HFconfig, HFT}
The container of Hartree–Fock method configuration.
≡≡≡ Field(s) ≡≡≡
HF::Val{HFT}
: Hartree–Fock method type. Available values of HFT
are :RHF, :UHF.
C0::InitialC{T1, HFT, F}
: Initial guess of the orbital coefficient matrix(s) C of the canonical orbitals. When C0
is as an argument of HFconfig
's constructor, it can be set to sym::Symbol
where available values of sym
are :GWH, :Hcore, :SAD
; it can also be a Tuple
of prepared orbital coefficient matrix(s) for the corresponding Hartree–Fock method type.
SCF::SCFconfig{T2, L, MS}
: SCF iteration configuration. For more information please refer to SCFconfig
.
maxStep::Int
: Maximum iteration steps allowed regardless if the iteration converges.
earlyStop::Bool
: Whether automatically terminate (or skip) a convergence method early when its performance becomes unstable or poor.
saveTrace::NTuple{4, Bool}
: Determine whether saving (by pushing) the intermediate information from all the iterations steps to the field .temp
of the output HFfinalVars
of runHF
. The types of relevant information are:
Sequence | Information | Corresponding field in HFtempVars |
---|---|---|
1 | orbital coefficient matrix(s) | .Cs |
2 | density matrix(s) | .Ds , .shared.Dtots |
3 | Fock matrix(s) | .Fs |
4 | unconverged Hartree–Fock energy(s) | .Es , .shared.Etots |
≡≡≡ Initialization Method(s) ≡≡≡
HFconfig(;HF::Union{Symbol, Val}=:RHF,
C0::Union{Tuple{AbstractMatrix}, NTuple{2, AbstractMatrix},
Symbol, Val}=:SAD,
SCF::SCFconfig=SCFconfig{Float64, 2, Tuple{Val{:ADIIS}, Val{:DIIS}}}(method, interval=(0.001, 1.0e-10), methodConfig, secondaryConvRatio, oscillateThreshold),
maxStep::Int=150,
earlyStop::Bool=true,
saveTrace::NTuple{4, Bool}=(false, false, false, true)) ->
HFconfig
≡≡≡ Example(s) ≡≡≡
julia> HFconfig();
julia> HFconfig(HF=:UHF);
Quiqbox.HFfinalVars
— TypeHFfinalVars{T, D, HFT, NN, BN, HFTS} <: HartreeFockFinalValue{T, HFT}
The container of the final values after a Hartree–Fock SCF procedure.
≡≡≡ Field(s) ≡≡≡
Ehf::T
: Hartree–Fock energy of the electronic Hamiltonian.
Enn::T
: The nuclear repulsion energy.
Ns::NTuple{HFTS, Int}
: The number(s) of electrons with same spin configurations(s). For restricted closed-shell Hartree–Fock (RHF), the single element in .Ns
represents both spin-up electrons and spin-down electrons.
nuc::NTuple{NN, String}
: The nuclei in the studied system.
nucCoords::NTuple{NN, NTuple{D, T}}
: The coordinates of corresponding nuclei.
C::NTuple{HFTS, Matrix{T}}
: Orbital coefficient matrix(s) for one spin configuration.
D::NTuple{HFTS, Matrix{T}}
: Density matrix(s) for one spin configuration.
F::NTuple{HFTS, Matrix{T}}
: Fock matrix(s) for one spin configuration.
Eo::NTuple{HFTS, Vector{T}}
: Energies of canonical orbitals.
occu::NTuple{HFTS, NTuple{BN, Int}}
: Occupations of canonical orbitals.
temp::NTuple{HFTS, [HFtempVars](@ref){T, HFT}}
: the intermediate values stored during the Hartree–Fock interactions.
isConverged::Union{Bool, Missing}
: Whether the SCF iteration is converged in the end. When convergence detection is off (see SCFconfig), it is set to missing
.
basis::GTBasis{T, D, BN}
: The basis set used for the Hartree–Fock approximation.
Quiqbox.HFtempVars
— TypeHFtempVars{T, HFT} <: HartreeFockintermediateData{T}
The container to store the intermediate values (only of the one spin configuration) for each iteration during the Hartree–Fock SCF procedure.
≡≡≡ Field(s) ≡≡≡
N::Int
: The number of electrons with the one spin function.
Cs::Vector{Matrix{T}}
: Orbital coefficient matrices.
Ds::Vector{Matrix{T}}
: Density matrices corresponding to only one spin configuration.
Fs::Vector{Matrix{T}}
: Fock matrices.
Es::Vector{T}
: Part of the Hartree–Fock energy corresponding to one spin configuration.
shared.Dtots::Vector{Matrix{T}}
: The total density matrices.
shared.Etots::Vector{T}
: The total Hartree–Fock energy.
NOTE: For unrestricted Hartree–Fock, there are 2 HFtempVars
being updated during the iterations, and changing the field shared.Dtots
or shared.Etots
of one HFtempVars
will affect the other one's.
Quiqbox.SCFconfig
— TypeSCFconfig{T<:Real, L, MS<:NTuple{L, Val}} <: ImmutableParameter{T, SCFconfig}
The struct
for self-consistent field (SCF) iteration configurations.
≡≡≡ Field(s) ≡≡≡
method::MS
: The applied convergence methods. They can be specified as the elements inside an NTuple{L, Symbol}
, which is then input to the constructor of SCFconfig
as the positional argument methods
. The available configuration(s) corresponding to each method in terms of keyword arguments are:
Convergence Method(s) | Configuration(s) | Keyword(s) | Range(s)/Option(s) | Default(s) |
---|---|---|---|---|
:DD | damping strength | dampStrength | [0 , 1 ] | 0.5 |
:DIIS , :EDIIS , :ADIIS | subspace size; DIIS-Method solver; reset threshold¹ | DIISsize ; solver ; resetThreshold | 1 ,2 ...; :LBFGS ...; 1e-14 ... | 10 ; :LBFGS ; N/A |
¹ The reset threshold (resetThreshold::Real
) determines when to clear the memory of the DIIS-based method's subspace and reset the second-to-latest residual vector as the first reference. The reset is executed when the latest computed energy increases an amount above the threshold compared to the second-to-latest computed energy. In default, the threshold is always slightly larger than the machine epsilon of the numerical data type for the SCF computation.
Convergence Methods
:DD
: Direct diagonalization of the Fock matrix.:DIIS
: Direct inversion in the iterative subspace.:EDIIS
: Energy-DIIS.:ADIIS
: DIIS based on the augmented Roothaan–Hall (ARH) energy function.
DIIS-Method Solvers
:LBFGS
: Limited-memory BFGS with box constraints.:LCM
: Lagrange multiplier solver.:SPGB
: Spectral Projected Gradient Method with box constraints.
interval::NTuple{L, T}
: The stopping (or skipping) thresholds for required methods. The last threshold will be the convergence threshold for the SCF procedure. When the last threshold is set to NaN
, there will be no convergence detection.
methodConfig::NTuple{L, Vector{<:Pair}}
: The additional keywords arguments for each method stored as Tuple
s of Pair
s.
secondaryConvRatio::NTuple{2, T}
: The ratios of the secondary convergence criteria to the primary convergence indicator, i.e., the change of the energy (ΔE):
Order | Symbols | Meaning | Default |
---|---|---|---|
1 | RMS(FDS-SDF) | root mean square of the error matrix defined in DIIS | 50 |
2 | RMS(ΔD) | root mean square of the change of the density matrix | 50 |
oscillateThreshold::T
: The threshold for oscillatory convergence.
≡≡≡ Initialization Method(s) ≡≡≡
SCFconfig(methods::NTuple{L, Symbol}, intervals::NTuple{L, T},
config::Dict{Int, <:AbstractVector{<:Pair}}=Dict(1=>Pair[]);
secondaryConvRatio::Union{Real, NTuple{2, Real}}=(50, 50),
oscillateThreshold::Real=1.0e-6) where {L, T} ->
SCFconfig{T, L}
methods
and intervals
are the convergence methods to be applied and their stopping (or skipping) thresholds respectively. config
specifies additional keyword argument(s) for each methods by a Pair
of which the key i::Int
is for i
th method and the pointed AbstractVector{<:Pair}
is the pairs of keyword arguments and their values respectively. If secondaryConvRatio
is Real
, it will be assigned as the value for all the secondary convergence ratios.
SCFconfig(;threshold::AbstractFloat=(0.001, 1.0e-10),
secondaryConvRatio::Union{Real, NTuple{2, Real}}=(50, 50),
oscillateThreshold::Real=defaultOscThreshold) ->
SCFconfig{Float64, 2}
threshold
will update the stopping threshold of the default SCF configuration used in HFconfig() with a new value. In other words, it updates the stopping threshold of :1.0e-10
.
≡≡≡ Example(s) ≡≡≡
julia> SCFconfig((:DD, :ADIIS, :DIIS), (1e-4, 1e-12, 1e-13), Dict(2=>[:solver=>:LCM]));
julia> SCFconfig(threshold=1e-8, oscillateThreshold=1e-5)
SCFconfig{Float64, 2, Tuple{Val{:ADIIS}, Val{:DIIS}}}(method, interval=(0.001, 1.0e-8), methodConfig, secondaryConvRatio, oscillateThreshold)
Quiqbox.GridBox
— TypeGridBox(nGridPerEdge::Int, spacing, center=ntuple(_->eltype(spacing[begin])(0), 3);
canDiff::Bool=true, index::Int=0) ->
GridBox{T, D}
The method of generating a cubic GridBox
. Aside from the common arguments, nGridPerEdge
specifies the number of grids for every dimension. The dimension of the grid box is determined by the dimension of center
.
Quiqbox.GridBox
— TypeGridBox{T, D, NP, GPT<:SpatialPoint{T, D}} <: SpatialStructure{T, D}
A container of multiple D
-dimensional grid points.
≡≡≡ Field(s) ≡≡≡
spacing::NTuple{D, T}
: The distance between adjacent grid points along each dimension.
nPoint::Int
: Total number of the grid points.
point::NTuple{NP, GPT}
: The grid points represented by SpatialPoint
.
param::Tuple{Vararg{ParamBox{T}}}
: All the parameters in the GridBox
.
≡≡≡ Initialization Method(s) ≡≡≡
GridBox(nGrids::NTuple{D, Int}, spacing::NTuple{D, Union{Array{T, 0}, T}},
center::Union{AbstractVector{T}, NTuple{D, T}}=ntuple(_->T(0),Val(D));
canDiff::NTuple{D, Bool}=ntuple(_->true, Val(D)),
index::NTuple{D, Int}=ntuple(_->0, Val(D))) where {T<:AbstractFloat, D} ->
GridBox{T, D}
GridBox(nGrids::NTuple{D, Int}, spacingForAllDim::Union{T, Array{T, 0}},
center::Union{AbstractVector{T}, NTuple{D, T}}=ntuple(_->T(0), Val(D));
canDiff::Bool=true, index::Int=0) where {T<:AbstractFloat, D} ->
GridBox{T, D}
Construct a general D
-dimensional GridBox
.
=== Positional argument(s) ===
nGrids::NTuple{D, Int}
: The numbers of grids along each dimension.
spacing::NTuple{D, Union{Array{T, 0}, T}}
: The spacing between grid points along each dimension.
spacingForAllDim::NTuple{D, Union{Array{T, 0}, T}}
: A single spacing applied for all dimensions.
center::Union{AbstractVector{T}, NTuple{D, T}}
: The coordinate of the geometric center of the grid box.
=== Keyword argument(s) ===
canDiff
::NTuple{D, Bool}: Whether the ParamBox
es of each dimension stored in the constructed GridBox
will be marked as differentiable.
index
::NTuple{D, Int}: The Index(s) that will be assigned to the shared input variable(s) L
of the stored ParamBox
es.
Quiqbox.GDconfig
— TypeGDconfig{T, M, ST<:Union{Array{T, 0}, T}} <: ConfigBox{T, GDconfig, M}
The mutable container of configurations for the default gradient descent method used to optimize parameters.
≡≡≡ Field(s) ≡≡≡
lineSearchMethod::M
: The selected line search method to optimize the step size for each gradient descent iteration. It can be any algorithm defined in the Julia package LineSearches.jl, or itself
so that the step size will be fixed to a constant value during all the iterations.
initialStep::ST
: The value of the initial step size in each iteration. If it's an Array{T, 0}
, then the final step size (optimized by lineSearchMethod
) in current interaction will be set as the initial step size for the next iteration.
stepBound::NTuple{2, T}
: The lower bound and upper bound of the step size to enforce stepBound[begin] ≤ η ≤ stepBound[end]
(where η
is the step size in each iteration). Whenever the size size is out of the boundary, it will be clipped to a bound value.
scaleStepBound::Bool
: Determine whether stepBound
will be scaled so that the norm of the gradient descent in each iteration is constrained within the initially configured stepBound
, i.e., stepBound[begin] ≤ norm(η*∇f) ≤ stepBound[end]
(where ∇f
is the gradient of some function f
in each iteration).
≡≡≡ Initialization Method(s) ≡≡≡
GDconfig(lineSearchMethod::M=BackTracking(),
initialStep::Union{Array{T, 0}, T}=ifelse(M==typeof(itself), 0.1, 1.0);
stepBound::NTuple{2, T}=convert.(eltype(initialStep), (0, Inf)),
scaleStepBound::Bool=ifelse(M==typeof(itself), true, false)) where
{M, T<:AbstractFloat} ->
GDconfig{T, M, typeof(initialStep)}
Quiqbox.POconfig
— TypePOconfig{T, M, CBT<:ConfigBox, TH<:Union{T, NTuple{2, T}},
OM} <: ConfigBox{T, POconfig, M}
The mutable container of parameter optimization configurations.
≡≡≡ Field(s) ≡≡≡
method::Val{M}
: The method that defines the objective function (e.g., HF energy) for the optimization. Available values of M
from Quiqbox are :HFenergy, :DirectRHFenergy.
config::CBT
: The configuration for the selected method
. E.g., for :HFenergy
it's HFconfig
.
target::T
: The target value of the objective function. The difference between the last-step value and the target value will be used for convergence detection. If it's set to NaN
, the gradient of the latest step and the difference between the function values of latest two steps are used instead.
threshold::TH
: The error threshold/thresholds for the function value difference and the gradient both/respectively to determine whether the optimization iteration has converged. When it's (or either of them) set to NaN
, there will be no corresponding convergence detection, and when target
is not NaN
, the threshold for the gradient won't be used because the gradient won't be part of the convergence criteria.
maxStep::Int
: Maximum iteration steps allowed regardless if the iteration converges.
optimizer::F
: Applied parameter optimizer. The default setting is GDconfig
()
. To use a function implemented by the user as the optimizer, it should have the following function signature:
optimizer(f::Function, gf::Function, x0::Vector{T}) where {T} -> Function
where f
is the objective function to be minimized, gf
is a function that returns both the gradient and the returned value of f
given the input value as a vector. x0
is the initial input value. The output of optimizer
, if we name it optimize!
, should have the corresponding function signature:
optimize!(x::Vector{T}, gx::Vector{T}, fx::T) where {T}
where x
, gx
, fx
are the input value, the gradient, and the returned value of f
respectively at one step. In other words, (gx, fx) == (gx, f(x)) == gf(x)
. After accepting those arguments, optimizer
should update (i.e. mutate the elements of) x
so that f(x) will have lower returned value.
saveTrace::NTuple{4, Bool}
: Determine whether saving (by pushing) the intermediate information from all the iterations steps to the output of optimizeParams!
. The types of relevant information are:
Sequence | Corresponding Information |
---|---|
1 | function value of the applied method |
2 | parameter value(s) |
3 | function gradient with respect to the parameter(s) |
4 | complete output of the applied method |
≡≡≡ Initialization Method(s) ≡≡≡
POconfig(;method::Union{Val{M}, Symbol}=HFenergy,
config::ConfigBox=Quiqbox.defaultOFconfigs[method],
target::T=NaN,
threshold::Union{T, NTuple{2, T}}=(5.0e-8, 5.0e-5),
maxStep::Int=200,
optimizer::Function=GDconfig()) where
{T, M} ->
POconfig{T, M}
≡≡≡ Example(s) ≡≡≡
julia> POconfig();
julia> POconfig(maxStep=100);
Quiqbox.CanOrbital
— TypeCanOrbital{T, D, NN} <: AbstractSpinOrbital{T, D}
The spatial part (orbital) of a canonical spin-orbital (the set of which diagonalizes the Fock matrix of a Hartree–Fock state) with its occupation information. This means the maximal occupation number for the mode corresponding to the orbital (namely a canonical orbital) equals 2. Please refer to genCanOrbitals
for the construction of a CanOrbital
.
≡≡≡ Field(s) ≡≡≡
energy::T
: The eigen energy corresponding to the orbital.
index::Int
: The index of the orbital within the same spin configuration.
nuc::NTuple{NN, String}
: The nuclei in the studied system.
nucCoords::NTuple{NN, NTuple{D, T}}
: The coordinates of corresponding nuclei.
occu::NTuple{2, Array{Bool, 0}}
: The occupations of two spin configurations.
orbital::GTBasisFuncs{T, D, 1}
: The spatial orbital part.
Quiqbox.MatterByHF
— TypeMatterByHF{T, D, NN, BN, HFTS} <:MatterData{T, D}
Container of the electronic structure information of a quantum system.
≡≡≡ Field(s) ≡≡≡
Ehf::T
: Hartree–Fock energy of the electronic Hamiltonian.
nuc::NTuple{NN, String}
: The nuclei in the studied system.
nucCoord::NTuple{NN, NTuple{D, T}}
: The coordinates of corresponding nuclei.
Enn::T
: The nuclear repulsion energy.
Ns::NTuple{HFTS, Int}
: The number(s) of electrons with same spin configurations(s). For restricted closed-shell Hartree–Fock (RHF), the single element in .Ns
represents both spin-up electrons and spin-down electrons.
occu::NTuple{HFTS, NTuple{BN, Int}}
: Occupations of canonical orbitals.
occuOrbital::NTuple{HFTS, Tuple{Vararg{CanOrbital{T, D, NN}}}}
: The occupied canonical orbitals.
unocOrbital::NTuple{HFTS, Tuple{Vararg{CanOrbital{T, D, NN}}}}
The unoccupied canonical orbitals.
occuC::NTuple{HFTS, Matrix{T}}
: Coefficient matrix(s) of occupied canonical orbitals.
unocC::NTuple{HFTS, Matrix{T}}
: Coefficient matrix(s) of unoccupied canonical orbitals.
coreHsameSpin::NTuple{HFTS, Matrix{T}}
: Core Hamiltonian(s) (one-body integrals) of the canonical orbitals with same spin configuration(s).
eeIsameSpin::NTuple{HFTS, Array{T, 4}}
: electron-electron interactions (two-body integrals) of the canonical orbitals with same spin configuration(s).
eeIdiffSpin::Matrix{T}
: Coulomb interactions between canonical orbitals with different spins.
basis::GTBasis{T, D, BN}
: The basis set used for the Hartree–Fock approximation.
≡≡≡ Initialization Method(s) ≡≡≡
MatterByHF(HFres::HFfinalVars{T, D, <:Any, NN, BN, HFTS};
roundAtol::Real=getAtolVal(T)) where {T, D, NN, BN, HFTS} ->
MatterByHF{T, D, NN, BN, HFTS}
Construct a MatterByHF
from the result of a Hartree–Fock method HFres
. Each parameter stored in the constructed CanOrbital
s in .occuOrbital
and .unocOrbital
will be rounded to the nearest multiple of roundAtol
; when roundAtol
is set to NaN
, no rounding will be performed.