Core Functions
Quiqbox.getCharge
β MethodgetCharge(nucs::Union{AbstractVector{String}, Tuple{Vararg{String}}}) -> Int
Return the total electric charge (in π) of the input nuclei.
Quiqbox.orbitalLin
β FunctionorbitalLin(subshell::String, D::Int=3) -> Tuple{Vararg{NTuple{3, Int}}}
Return all the possible angular momentum configuration(s) within the input subshell
of D
dimension.
Quiqbox.changeMapping
β MethodchangeMapping(pb::ParamBox{T, V}, mapFunction::Function=itself, outSym::Symbol=V;
canDiff::Union{Bool, Array{Bool, 0}}=isDiffParam(pb)) where {T, V} ->
ParamBox{T, outSym}
Change the mapping function of pb
. The symbol of the output variable of the returned ParamBox
can be specified by outSym
, and its differentiability is determined by canDiff
.
Quiqbox.dataOf
β MethoddataOf(pb::ParamBox{T}) where {T} -> Pair{Array{T, 0}, Symbol}
Return the Pair
of the input variable and its symbol.
Quiqbox.disableDiff!
β MethoddisableDiff!(pb::ParamBox) -> Bool
Mark pb
as "non-differentiable" and then return false
.
Quiqbox.enableDiff!
β MethodenableDiff!(pb::ParamBox) -> Bool
Mark pb
as "differentiable" and then return true
.
Quiqbox.fullVarCopy
β MethodfullVarCopy(pb::T) where {T<:ParamBox} -> T
A shallow copy of the input ParamBox
.
Quiqbox.inSymOf
β MethodinSymOf(pb::ParamBox) -> Symbol
Return the symbol of the input variable of pb
.
Quiqbox.inValOf
β MethodinValOf(pb::ParamBox{T}) where {T} -> T
Return the value of the input variable of pb
. Equivalent to pb[]
.
Quiqbox.indVarOf
β MethodindVarOf(pb::ParamBox{T}) -> Pair{}
Return (the name and the value of) the independent variable tied to pb
. Specifically, return the input variable stored in pb
when pb
is marked as differentiable; return the output variable of pb
when pb
is marked as non-differentiable. Thus, it is the variable pb
represents to differentiate any (differentiable) function of ParamBox
es.
Quiqbox.isDiffParam
β MethodisDiffParam(pb::ParamBox) -> Bool
Return the Boolean value of whether pb
is differentiable with respect to its input variable.
Quiqbox.isInSymEqual
β MethodisInSymEqual(pb::ParamBox, sym::Symbol) -> Bool
Return the Boolean value of whether the symbol of pb
's input variable equals sym
.
Quiqbox.isOutSymEqual
β MethodisOutSymEqual(::ParamBox, sym::Symbol) -> Bool
Return the Boolean value of whether the symbol of pb
's output variable equals sym
.
Quiqbox.mapOf
β MethodmapOf(pb::ParamBox) -> Function
Return the mapping function of pb
.
Quiqbox.outSymOf
β MethodoutSymOf(pb::ParamBox) -> Symbol
Return the symbol of the output variable of pb
.
Quiqbox.outValCopy
β MethodoutValCopy(pb::ParamBox{T, V}) where {T} -> ParamBox{T, V, typeof(itself)}
Return a new ParamBox
of which the input variable's value is equal to the output variable's value of pb
.
Quiqbox.outValOf
β MethodoutValOf(pb::ParamBox) -> Number
Return the value of the output variable of pb
. Equivalent to pb()
.
Quiqbox.toggleDiff!
β MethodtoggleDiff!(pb::ParamBox) -> Bool
Toggle the differentiability of the input pb
and then return the updated boolean value.
Quiqbox.absorbNormFactor
β MethodabsorbNormFactor(b::BasisFunc{T, 3, π, GN, PT}) where {T, π, GN, PT} ->
FloatingGTBasisFuncs{T, 3, π, GN, PT}
absorbNormFactor(b::BasisFuncs{T, 3, π, GN, PT}) where {T, π, GN, PT} ->
Vector{<:FloatingGTBasisFuncs{T, 3, π, GN, PT}}
If hasNormFactor(
b) == true
, absorb the normalization factor of each Gaussian-type orbital inside b
into the orbital's corresponding contraction coefficient and then set .normalizeGTO
of b
to false
. Otherwise, directly return b
when it's a BasisFunc
, or [b]
when it's a BasisFuncs
.
Quiqbox.absorbNormFactor
β MethodabsorbNormFactor(bfm::BasisFuncMix{T, 3}) where {T} -> GTBasisFuncs{T, 3}
Apply absorbNormFactor
to every one of the BasisFunc
inside bfm
and them sum them over.
Quiqbox.absorbNormFactor
β MethodabsorbNormFactor(bs::AbstractVector{<:GTBasisFuncs{T, 3}}) where {T} ->
AbstractVector{<:GTBasisFuncs{T, 3}}
absorbNormFactor(bs::Tuple{Vararg{GTBasisFuncs{T, 3}}}) where {T} ->
Tuple{Vararg{GTBasisFuncs{T, 3}}}
Apply absorbNormFactor
to every element inside bs
.
Quiqbox.add
β Methodadd(b1::CompositeGTBasisFuncs{T, D, <:Any, 1},
b2::CompositeGTBasisFuncs{T, D, <:Any, 1};
roundAtol::Real=getAtolVal(T)) where {T, D} ->
CompositeGTBasisFuncs{T, D, <:Any, 1}
Addition between two CompositeGTBasisFuncs{T, D, <:Any, 1}
such as BasisFunc
and BasisFuncMix
. roundAtol
specifies the absolute approximation tolerance of comparing parameters stored in each CompositeGTBasisFuncs
to determine whether they are treated as "equal"; each parameter in the returned CompositeGTBasisFuncs
is set to the nearest exact multiple of 0.5atol
. When roundAtol
is set to NaN
, there will be no approximation nor rounding. This function can be called using +
syntax with the keyword argument set to it default value.
NOTE: For the ParamBox
(stored in the input arguments) that are marked as non-differentiable, they will be fused together if possible to generate new ParamBox
(s) no longer linked to the input variable stored in them.
β‘β‘β‘ Example(s) β‘β‘β‘
julia> bf1 = genBasisFunc([1.,1.,1.], (2.,1.));
julia> bf2 = genBasisFunc([1.,1.,1.], (2.,2.));
julia> bf3 = bf1 + bf2;
julia> bf3.gauss[1].con() == bf1.gauss[1].con() + bf2.gauss[1].con()
true
Quiqbox.assignCenInVal!
β MethodassignCenInVal!(b::FloatingGTBasisFuncs{T, D}, center::AbstractVector{<:Real}) ->
SpatialPoint{T, D}
Change the input value of the ParamBox
stored in b.center
(meaning the output value will also change according to the mapping function). Then, return the altered center.
Quiqbox.centerCoordOf
β MethodcenterCoordOf(bf::FloatingGTBasisFuncs{T}) where {T} -> Vector{T}
Return the center coordinate of the input FloatingGTBasisFuncs
.
Quiqbox.centerNumOf
β MethodcenterNumOf(bf::CompositeGTBasisFuncs) -> Int
Return the number of center (coordinates) associated with the input CompositeGTBasisFuncs
.
Quiqbox.centerOf
β MethodcenterOf(bf::FloatingGTBasisFuncs{T, D}) where {T, D} -> SpatialPoint{T, D}
Return the center of the input FloatingGTBasisFuncs
.
Quiqbox.coordOf
β MethodcoordOf(sp::SpatialPoint{T}) where {T} -> Vector{T}
Get the coordinate represented by the input SpatialPoint
.
Quiqbox.copyBasis
β FunctioncopyBasis(b::GaussFunc, copyOutVal::Bool=true) -> GaussFunc
copyBasis(b::CompositeGTBasisFuncs, copyOutVal::Bool=true) -> CompositeGTBasisFuncs
Return a copy of the input basis. If copyOutVal
is set to true
, then only the output value(s) of the ParamBox
(s) stored in b
will be copied, i.e., outValCopy
is used, otherwise fullVarCopy
is used.
β‘β‘β‘ Example(s) β‘β‘β‘
julia> f(x)=x^2; e = genExponent(3.0, f)
ParamBox{Float64, :Ξ±, β¦}{1}[π][x_Ξ±]β¦ββ§[9.0]
julia> c = genContraction(2.0)
ParamBox{Float64, :d, β¦}{0}[β][d]β¦=β§[2.0]
julia> gf1 = GaussFunc(e, c);
julia> gf2 = copyBasis(gf1)
GaussFunc{Float64, β¦}{0, 0}[ββ][{9.0, 2.0}]
julia> gf1.xpn() == gf2.xpn()
true
julia> (gf1.xpn[] |> gf1.xpn.map) == gf2.xpn[]
true
Quiqbox.decompose
β Functiondecompose(bf::CompositeGTBasisFuncs{T, D}, splitGaussFunc::Bool=false) ->
Matrix{<:BasisFunc{T, D}}
Decompose a CompositeGTBasisFuncs
into a Matrix
of BasisFunc
s. The sum of each column represents one orbital of the input basis function(s). If splitGaussFunc
is true
, then each column consists of the BasisFunc
s with only 1 GaussFunc
.
Quiqbox.dimOf
β MethoddimOf(::DimensionalParamContainer) -> Int
Return the spatial dimension of the input parameterized container such as AbstractSpatialPoint
and QuiqboxBasis
.
Quiqbox.gaussCoeffOf
β MethodgaussCoeffOf(b::FloatingGTBasisFuncs{T}) -> Matrix{T}
Return the exponent and contraction coefficients of each GaussFunc
(in each row of the returned Matrix
) inside b
.
Quiqbox.gaussCoeffOf
β MethodgaussCoeffOf(gf::GaussFunc{T}) -> Matrix{T}
Return the exponent and contraction coefficients of gf
.
Quiqbox.genBFuncsFromText
β MethodgenBFuncsFromText(content::String;
adjustContent::Bool=false,
adjustFunction::Function=sciNotReplace,
excludeFirstNlines::Int=0, excludeLastNlines::Int=0,
center::Union{AbstractArray{T},
NTuple{D, T},
NTuple{D, ParamBox{T}},
SpatialPoint{T, D},
Missing}=(NaN, NaN, NaN),
unlinkCenter::Bool=false,
normalizeGTO::Union{Bool, Missing}=
ifelse(adjustContent, true, missing)) where
{D, T<:AbstractFloat} ->
Array{<:FloatingGTBasisFuncs, 1}
Generate a basis set from content
which is either a basis set String
in Gaussian format or the output from genBasisFuncText
. For the former, adjustContent
needs to be set to true
. adjustFunction
is only applied when adjustContent=true
, which in default is a function
used to detect and convert the format of the scientific notation in content
.
excludeFirstNlines
and excludeLastNlines
are used to exclude first or last few lines of content
if intended. center
is used to assign a center coordinate for all the basis functions from content
; when it's set to missing
, it will try to read the center information in content
, and leave the center as [NaN, NaN, Nan]
if one cannot be found for each corresponding FloatingGTBasisFuncs
. If unlinkCenter = true
, the center of each FloatingGTBasisFuncs
is a Base.deepcopy
of the input center
. Otherwise, they share the same underlying data so changing the value of one will affect others. If the center coordinate is included in content
, it should be right above the subshell information for the FloatingGTBasisFuncs
. E.g.:
"""
X 1.0 0.0 0.0
S 1 1.0 false
2.0 1.0
"""
Finally, normalizeGTO
specifies the field .normalizeGTO
of the generated FloatingGTBasisFuncs
. If it's set to missing
(in default), the normalization configuration of each FloatingGTBasisFuncs
will depend on content
, so different basis functions may have different normalization configurations. However, when it's set to a Bool
value, .normalizeGTO
of all the generated basis functions will be set to that value.
Quiqbox.genBasisFunc
β MethodgenBasisFunc(center::Union{AbstractVector{T}, Tuple{Vararg{T}}, SpatialPoint, Missing},
args..., kws...) where {T<:Union{AbstractFloat, ParamBox}} ->
Union{FloatingGTBasisFuncs{T}, Vector{<:FloatingGTBasisFuncs{T}}}
The constructor of FloatingGTBasisFuncs
, but it also returns different kinds of collections (Vector
) of them based on the input arguments. The first argument center
specifies the center coordinate of the generated FloatingGTBasisFuncs
, and can be left as missing
for later assignment.
NOTE: Although not marked in the return types, if any allowed input argument can lead to a FloatingGTBasisFuncs
with no GaussFunc
stored inside of it, e.g., genBasisFunc([1.0, 2.0, 3.0], ())
, a Quiqbox.EmptyBasisFunc
will be generated instead.
β‘β‘β‘ Method 1 β‘β‘β‘
genBasisFunc(center, GsOrCoeffs, Ls; normalizeGTO=false) ->
FloatingGTBasisFuncs
=== Positional argument(s) ===
GsOrCoeffs::Union{ AbstractGaussFunc{T1}, AbstractVector{<:AbstractGaussFunc{T1}}, Tuple{Vararg{AbstractGaussFunc{T1}}}, NTuple{2, Union{T1, Array{T1, 0}, ParamBox{T1}}, NTuple{2, AbstractVector{<:Union{T1, Array{T1, 0}, ParamBox{T1}}}} } where {T1<:AbstractFloat}
: A collection of concentric GaussFunc
that will be used to construct the basis function. To simplify the procedure, it can also be in the form of a NTuple{2}
of the exponent coefficient(s)::Union{AbstractFloat, AbstractVector{<:AbstractFloat}}
and contraction coefficients::Union{AbstractFloat, AbstractVector{<:AbstractFloat}}
of the GaussFunc
(s) to be input.
Ls::Union{ T2, AbstractVector{T2}, NTuple{<:Any, T2} } where {T2<:Union{NTuple{D, Int}, LTuple{D}} where {D}}
: A collection of at least one angular momentum configuration within the same subshell, in the Cartesian coordinate representation. E.g., for p shell it can be set to ((1,0,0), (0,1,0))
. This will determine the number of spatial orbitals and their angular momentum respectively to be stored in the output FloatingGTBasisFuncs
.
=== Keyword argument(s) ===
normalizeGTO::Bool
: Determine whether the inside GaussFunc
(s) will be normalized in the during the calculation.
=== Example(s) ===
julia> genBasisFunc([0.,0.,0.], GaussFunc(2., 1.), (0, 1, 0))
BasisFunc{Float64, 3, 1, 1, β¦}{0, 0, 0}[(0.0, 0.0, 0.0)][Xβ°YΒΉZβ°]
β‘β‘β‘ Method 2 β‘β‘β‘
genBasisFunc(center, GsOrCoeffs, subshell="s"; normalizeGTO=false) ->
FloatingGTBasisFuncs
genBasisFunc(center, GsOrCoeffs, subshell, lFilter; normalizeGTO=false) ->
FloatingGTBasisFuncs
=== Positional argument(s) ===
subshell::Union{String, Symbol}
: The third argument of the constructor can also be the name of a subshell, which will make sure the output is a BasisFuncs
that contains the spatial orbitals that fully occupy the subshell.
lFilter::Tuple{Vararg{Bool}}
: When this 4th argument is provided, it can determine the orbital(s) to be included based on the given subshell
. The order of the corresponding orbital angular momentum(s) can be inspected using function orbitalLin
.
=== Keyword argument(s) ===
normalizeGTO::Bool
: Same as the one defined in method 1.
=== Example(s) ===
julia> genBasisFunc([0.,0.,0.], (2., 1.), "p")
BasisFuncs{Float64, 3, 1, 1, β¦}{0, 0, 0}[(0.0, 0.0, 0.0)][3/3]
julia> genBasisFunc([0.,0.,0.], (2., 1.5), "p", (true, false, true))
BasisFuncs{Float64, 3, 1, 1, β¦}{0, 0, 0}[(0.0, 0.0, 0.0)][2/3]
β‘β‘β‘ Method 3 β‘β‘β‘
genBasisFunc(center, BSkey, atm="H"; unlinkCenter=false) ->
Vector{<:FloatingGTBasisFuncs}
=== Positional argument(s) ===
BSkey::String
: The name of an existed atomic basis set. The supported options are in ["STO-2G", "STO-3G", "STO-6G", "3-21G", "6-31G", "cc-pVDZ", "cc-pVTZ", "cc-pVQZ"]
.
atm::Union{String, Symbol}
: The name of the atom corresponding to the chosen basis set. The supported options are in ["H", "He", "Li", "Be", "B", "C", "N", "O", "F", "Ne", "Na", "Mg", "Al", "Si", "P", "S", "Cl", "Ar", "K", "Ca"]
.
=== Keyword argument(s) ===
unlinkCenter::Bool
: Determine whether the centers of constructed FloatingGTBasisFuncs
are linked to each other. If set to true
, the center of each FloatingGTBasisFuncs
is a Base.deepcopy
of each other. Otherwise, they share the same underlying data so changing the value of one will affect others.
=== Example(s) ===
julia> genBasisFunc([0.,0.,0.], "6-31G");
julia> genBasisFunc([0.,0.,0.], "STO-3G", "Li");
β‘β‘β‘ Method 4 β‘β‘β‘
genBasisFunc(b::FloatingGTBasisFuncs{T, D}, newFieldVal) where {T, D} ->
FloatingGTBasisFuncs{T, D}
=== Positional argument(s) ===
newFieldVal::Union{ SpatialPoint{T, D}, Tuple{AbstractGaussFunc{T}, Vararg{AbstractGaussFunc{T}}}, Tuple{LTuple{D, π}, Vararg{LTuple{D, π}}} where π, Bool } where {T<:AbstractFloat, D}
: Any one of the fields inside a FloatingGTBasisFuncs
except param
.
This method outputs a FloatingGTBasisFuncs
that has identical fields as the input one except the field that can be replaced by newFieldVal
(and param
if the replaced field contains ParamBox
).
=== Example(s) ===
julia> bf1 = genBasisFunc([1., 2., 3.], (2.0, 1.0), (0, 0, 0));
julia> bf2 = genBasisFunc([1., 2., 3.], (2.0, fill(1.0)), (0, 0, 1));
julia> bf1 = genBasisFunc(bf1, bf2.l);
julia> hasEqual(bf1, bf2)
true
Quiqbox.genBasisFuncText
β MethodgenBasisFuncText(bs::Union{AbstractVector{<:FloatingGTBasisFuncs{T, D}},
Tuple{Vararg{FloatingGTBasisFuncs{T, D}}}};
norm::Real=1.0, printCenter::Bool=true,
groupCenters::Bool=true, roundDigits::Int=-1) where {T, D} ->
Vector{String}
Generate the text of input basis set (consisting of FloatingGTBasisFuncs
). norm
is the additional normalization factor. groupCenters
determines whether the function will group the basis functions with same center together.
Quiqbox.genBasisFuncText
β MethodgenBasisFuncText(bf::FloatingGTBasisFuncs;
norm::Real=1.0, printCenter::Bool=true, roundDigits::Int=-1) -> String
Generate the text of input FloatingGTBasisFuncs
. norm
is the additional normalization factor. If printCenter
is true
, the center coordinate will be added to the first line of the output String
. roundDigits
specifies the rounding digits for the parameters inside bf
; when set to negative, no rounding will be performed.
Quiqbox.genContraction
β MethodgenContraction(pb::ParamBox{T}) where {T<:AbstractFloat} -> ParamBox{T, :d}
Convert a ParamBox
to an exponent coefficient parameter.
Quiqbox.genContraction
β MethodgenContraction(d::T, mapFunction::Function=itself;
canDiff::Bool=ifelse(FLevel(mapFunction)==FLevel{0}, false, true),
inSym::Symbol=x_d) where {T<:AbstractFloat} ->
ParamBox{T, :d}
genContraction(d::Array{T, 0}, mapFunction::Function=itself;
canDiff::Bool=ifelse(FLevel(mapFunction)==FLevel{0}, false, true),
inSym::Symbol=x_d) where {T<:AbstractFloat} ->
ParamBox{T, :d}
genContraction(dData::Pair{Array{T, 0}, Symbol},
mapFunction::Function=itself;
canDiff::Bool=ifelse(FLevel(mapFunction)==FLevel{0}, false, true)) where
{T<:AbstractFloat} ->
ParamBox{T, :d}
Construct a contraction coefficient given a value or variable (with its symbol). mapFunction
, canDiff
, and inSym
work the same way as in a general constructor of a ParamBox
.
Quiqbox.genExponent
β MethodgenExponent(pb::ParamBox{T}) where {T<:AbstractFloat} -> ParamBox{T, :Ξ±}
Convert a ParamBox
to the container of an exponent coefficient.
Quiqbox.genExponent
β MethodgenExponent(e::T, mapFunction::Function=itself;
canDiff::Bool=ifelse(FLevel(mapFunction)==FLevel{0}, false, true),
inSym::Symbol=x_Ξ±) where {T<:AbstractFloat} ->
ParamBox{T, :Ξ±}
genExponent(e::Array{T, 0}, mapFunction::Function=itself;
canDiff::Bool=ifelse(FLevel(mapFunction)==FLevel{0}, false, true),
inSym::Symbol=x_Ξ±) where {T<:AbstractFloat} ->
ParamBox{T, :Ξ±}
genExponent(eData::Pair{Array{T, 0}, Symbol},
mapFunction::Function=itself;
canDiff::Bool=ifelse(FLevel(mapFunction)==FLevel{0}, false, true)) where
{T<:AbstractFloat} ->
ParamBox{T, :Ξ±}
Construct an exponent coefficient given a value or variable (with its symbol). mapFunction
, canDiff
, and inSym
work the same way as in a general constructor of a ParamBox
.
Quiqbox.genSpatialPoint
β FunctiongenSpatialPoint(point::Union{NTuple{D, Union{T, Array{T, 0}, ParamBox{T}}},
AbstractVector},
marker::Symbol=:point) where {D, T<:AbstractFloat} ->
SpatialPoint{T, D}
Construct a SpatialPoint
from a collection of coordinate components. β‘β‘β‘ Example(s) β‘β‘β‘
julia> v1 = [1.0, 2.0, 3.0];
julia> p1 = genSpatialPoint(v1, :p1)
SpatialPoint{Float64, β¦}{0, 0, 0}[βββ][(1.0, 2.0, 3.0)]
julia> p1.marker
:p1
julia> v2 = [fill(1.0), 2.0, 3.0];
julia> p2 = genSpatialPoint(v2); p2[1]
ParamBox{Float64, :X, β¦}{0}[β][X]β¦=β§[1.0]
julia> v2[1][] = 1.2
1.2
julia> p2[1]
ParamBox{Float64, :X, β¦}{0}[β][X]β¦=β§[1.2]
Quiqbox.genSpatialPoint
β MethodgenSpatialPoint(point::Union{Tuple{ParamBox{T}, Vararg{ParamBox{T}}},
AbstractVector{<:ParamBox{T}}},
marker::Symbol=:point) where {T} ->
SpatialPoint{T}
Convert a collection of ParamBox
s to a SpatialPoint
.
Quiqbox.genSpatialPoint
β MethodgenSpatialPoint(comp::T, compIndex::Int,
mapFunction::Function=itself;
canDiff::Bool=ifelse(FLevel(mapFunction)==IL, false, true),
inSym::Symbol=(:x_X, :x_Y, :x_Z)[compIndex]) where {T<:AbstractFloat} ->
ParamBox{T}
genSpatialPoint(comp::Array{T, 0}, compIndex::Int,
mapFunction::Function=itself;
canDiff::Bool=ifelse(FLevel(mapFunction)==IL, false, true),
inSym::Symbol=(:x_X, :x_Y, :x_Z)[compIndex]) where {T<:AbstractFloat} ->
ParamBox{T}
genSpatialPoint(compData::Pair{Array{T, 0}, Symbol}, compIndex::Int,
mapFunction::Function=itself;
canDiff::Bool=ifelse(FLevel(mapFunction)==IL, false, true)) ->
ParamBox{T}
Construct a ParamBox
as the compIndex
th component of a SpatialPoint
given a value or variable (with its symbol). mapFunction
, canDiff
, and inSym
work the same way as in a general constructor of a ParamBox
.
genSpatialPoint(point::ParamBox{T}, compIndex::Int) where {T<:AbstractFloat} ->
ParamBox{T}
Convert a ParamBox
to the compIndex
th component of a SpatialPoint
.
β‘β‘β‘ Example(s) β‘β‘β‘
julia> genSpatialPoint(1.2, 1)
ParamBox{Float64, :X, β¦}{0}[β][X]β¦=β§[1.2]
julia> pointY1 = fill(2.0);
julia> Y1 = genSpatialPoint(pointY1, 2)
ParamBox{Float64, :Y, β¦}{0}[β][Y]β¦=β§[2.0]
julia> pointY1[] = 1.5;
julia> Y1
ParamBox{Float64, :Y, β¦}{0}[β][Y]β¦=β§[1.5]
Quiqbox.getNormFactor
β MethodgetNormFactor(b::FloatingGTBasisFuncs{T, 3}) where {T} -> Array{T}
Return the normalization factors of the Gaussian-type orbitals (GTO) inside the input b
. Each column corresponds to one orbital.
Quiqbox.getParams
β FunctiongetParams(pbc::ParamBox, symbol::Union{Symbol, Missing}=missing) ->
Union{ParamBox, Nothing}
getParams(pbc::ParameterizedContainer, symbol::Union{Symbol, Missing}=missing) ->
AbstractVector{<:ParamBox}
getParams(pbc::Union{AbstractArray{<:T}, Tuple{T, Vararg{T}}},
symbol::Union{Symbol, Missing}=missing) where {T<:QuiqboxContainer} ->
AbstractVector{<:ParamBox}
Return the parameter(s) stored in the input container. If symbol
is set to missing
, then return all the parameter(s). If it's set to a Symbol
tied to a parameter, for example :Ξ±β
, the function will match any pb::
ParamBox
such that indVarOf
(pb)[begin] == :Ξ±β
. If it's set to a Symbol
without any subscript, for example :Ξ±
, the function will match it with all the pb
s such that string(indVarOf(pb)[begin])
contains 'Ξ±'
. If the first argument is a collection, its entries must be ParamBox
containers.
Quiqbox.hasNormFactor
β MethodhasNormFactor(b::FloatingGTBasisFuncs) -> Bool
Indicate whether b
' is be treated as having additional normalization factor(s) which its Gaussian-type orbital(s) will be multiplied by during any calculation.
Quiqbox.lOf
β MethodlOf(::FloatingGTBasisFuncs) -> Int
Return the total orbital angular momentum quantum number (in Cartesian coordinate representation).
Quiqbox.markParams!
β MethodmarkParams!(b::Union{AbstractVector{T}, T, Tuple{T, Vararg{T}}},
filterMapping::Bool=false) where {T<:ParameterizedContainer} ->
Vector{<:ParamBox}
Mark the parameters (ParamBox
) in b
. The parameters that hold the same independent variable (in other words, considered identical in the differentiation procedure) will be marked with same index. filterMapping
determines whether filtering out (i.e. not return) the extra ParamBox
es that hold the same independent variables but represent different output variables. The independent variable held by a ParamBox
can be inspected using indVarOf
.
Quiqbox.mergeBasisFuncsIn
β MethodmergeBasisFuncsIn(bs::Union{AbstractVector{<:GTBasisFuncs{T, D}},
Tuple{Vararg{GTBasisFuncs{T, D}}}};
roundAtol::Real=NaN) where {T, D} ->
Vector{<:GTBasisFuncs{T, D}}
Try merging multiple FloatingGTBasisFuncs
(if there's any) in bs
into FloatingGTBasisFuncs{T, D, <:Any, <:Any, <:Any, ON}
where ON > 1
if possible and then return the resulted basis collection. If no merging is performed, then the returned collection is same as (but not necessarily identical to) bs
.
Quiqbox.mul
β Methodmul(a1::Real, a2::CompositeGTBasisFuncs{T, D, <:Any, 1};
normalizeGTO::Union{Bool, Missing}=missing,
roundAtol::Real=getAtolVal(T)) where {T, D} ->
CompositeGTBasisFuncs{T, D, <:Any, 1}
mul(a1::CompositeGTBasisFuncs{T, D, <:Any, 1}, a2::Real;
normalizeGTO::Union{Bool, Missing}=missing,
roundAtol::Real=getAtolVal(T)) where {T, D} ->
CompositeGTBasisFuncs{T, D, <:Any, 1}
mul(a1::CompositeGTBasisFuncs{T, D, <:Any, 1},
a2::CompositeGTBasisFuncs{T, D, <:Any, 1};
normalizeGTO::Union{Bool, Missing}=missing,
roundAtol::Real=getAtolVal(T)) where {T, D} ->
CompositeGTBasisFuncs{T, D, <:Any, 1}
Multiplication between two CompositeGTBasisFuncs{T, D, <:Any, 1}
(e.g., BasisFunc
and BasisFuncMix
), or a Real
number and a CompositeGTBasisFuncs{T, D, <:Any, 1}
. If normalizeGTO
is set to missing
(in default), The GaussFunc
inside the output containers will be normalized only if every input FloatingGTBasisFuncs
(or inside the input CompositeGTBasisFuncs
) holds hasNormFactor(ai) == true
. roundAtol
specifies the absolute approximation tolerance of comparing parameters stored in each CompositeGTBasisFuncs
to determine whether they are treated as "equal"; each parameter in the returned CompositeGTBasisFuncs
is set to the nearest exact multiple of 0.5atol
. When roundAtol
is set to NaN
, there will be no approximation nor rounding. This function can be called using *
syntax with the keyword argument set to it default value.
β‘β‘β‘ Example(s) β‘β‘β‘
julia> bf1 = genBasisFunc([1.0, 1.0, 1.0], ([2.0, 1.0], [0.1, 0.2]));
julia> gaussCoeffOf(bf1)
2Γ2 Matrix{Float64}:
2.0 0.1
1.0 0.2
julia> bf2 = bf1 * 2;
julia> gaussCoeffOf(bf2)
2Γ2 Matrix{Float64}:
2.0 0.2
1.0 0.4
julia> bf3 = bf1 * bf2;
julia> gaussCoeffOf(bf3)
3Γ2 Matrix{Float64}:
4.0 0.02
3.0 0.08
2.0 0.08
Quiqbox.mul
β Methodmul(gf::GaussFunc{T}, coeff::Real; roundAtol::Real=getAtolVal(T)) where {T} ->
GaussFunc
mul(coeff::Real, gf::GaussFunc{T}; roundAtol::Real=getAtolVal(T)) where {T} ->
GaussFunc
mul(gf1::GaussFunc{T}, gf2::GaussFunc{T}; roundAtol::Real=getAtolVal(T)) where {T} ->
GaussFunc
Multiplication between a Real
number and a GaussFunc
or two GaussFunc
s. roundAtol
specifies the absolute approximation tolerance of comparing parameters stored in each GaussFunc
to determine whether they are treated as "equal"; each parameter in the returned GaussFunc
is set to the nearest exact multiple of 0.5atol
. When roundAtol
is set to NaN
, there will be no approximation nor rounding. This function can be called using *
syntax with the keyword argument set to it default value.
NOTE: For the ParamBox
(stored in the input arguments) that are marked as non-differentiable, they will be fused together if possible to generate new ParamBox
(s) no longer linked to the data (input variable) stored in them.
β‘β‘β‘ Example(s) β‘β‘β‘
julia> gf1 = GaussFunc(3.0, 1.0)
GaussFunc{Float64, β¦}{0, 0}[ββ][{3.0, 1.0}]
julia> gf1 * 2
GaussFunc{Float64, β¦}{0, 0}[ββ][{3.0, 2.0}]
julia> gf1 * gf1
GaussFunc{Float64, β¦}{0, 0}[ββ][{6.0, 1.0}]
julia> gf1 * 2 * gf1
GaussFunc{Float64, β¦}{0, 0}[ββ][{6.0, 2.0}]
Quiqbox.normalizeBasis
β MethodnormalizeBasis(b::BasisFuncs{T, D}) where {T, D} -> Vector{<:FloatingGTBasisFuncs{T, D}}
Normalize each BasisFunc
inside b
and try to merge them back to one BasisFuncs
. If the all the BasisFunc
(s) can be merged, the returned result will be a 1-element Vector
.
Quiqbox.normalizeBasis
β MethodnormalizeBasis(b::GTBasisFuncs{T, D, 1}) where {T, D} -> GTBasisFuncs{T, D, 1}
Multiply the contraction coefficient(s) inside b
by constant coefficient(s) to normalizeBasis the b
, and then return the normalized basis.
Quiqbox.orbitalNumOf
β FunctionorbitalNumOf(subshell::Union{String, Symbol}, D::Integer=3) -> Int
Return the size (number of orbitals) of each subshell in D
dimensional real space.
Quiqbox.orbitalNumOf
β MethodorbitalNumOf(b::QuiqboxBasis) -> Int
Return the numbers of orbitals of the input basis.
Quiqbox.shift
β Methodshift(bf::FloatingGTBasisFuncs{T, D, π, GN, PT, 1},
dl::Union{Vector{Int}, NTuple{D, Int}, LTuple{D}}, op::Function=+) where
{T, D, π, GN, PT} ->
BasisFunc{T, D, <:Any, GN, PT}
Shift (+
as the default binary operator op
) the angular momentum (in Cartesian representation) of the input FloatingGTBasisFuncs
given dl
that specifies the change of each component.
Quiqbox.sortBasis
β MethodsortBasis(bs::Union{AbstractArray{<:CompositeGTBasisFuncs{T, D}},
Tuple{Vararg{CompositeGTBasisFuncs{T, D}}}};
roundAtol::Real=getAtolVal(T)) where {T, D} ->
Vector{<:CompositeGTBasisFuncs{T, D}}
Sort basis functions. roundAtol
specifies the absolute approximation tolerance of comparing parameters stored in each CompositeGTBasisFuncs
to determine whether they are treated as "equal"; when set to NaN
, no approximation will be made during the comparison.
Quiqbox.sortBasis
β MethodsortBasis(b::GTBasis{T, D}; roundAtol::Real=getAtolVal(T)) where {T, D} ->
GTBasis{T, D}
Reconstruct a GTBasis
by sorting the GTBasisFuncs
stored in the input one.
Quiqbox.sortBasis
β MethodsortBasis(bs::Tuple{Vararg{CompositeGTBasisFuncs{T, D}}};
roundAtol::Real=getAtolVal(T)) where {T, D} ->
Tuple{Vararg{CompositeGTBasisFuncs{T, D}}}
Quiqbox.sortBasisFuncs
β MethodsortBasisFuncs(bs::AbstractArray{<:FloatingGTBasisFuncs{T, D}},
groupCenters::Bool=false; roundAtol::Real=getAtolVal(T)) where {T, D} ->
Vector
Sort FloatingGTBasisFuncs
. If groupCenters = true
, Then the function will return an Vector{<:Vector{<:FloatingGTBasisFuncs}}
in which the elements are grouped basis functions with same center coordinates. roundAtol
specifies the absolute approximation tolerance of comparing the center coordinates to determine whether they are treated as "equal"; when set to NaN
, no approximation will be made during the comparison.
Quiqbox.sortBasisFuncs
β MethodsortBasisFuncs(bs::Tuple{FloatingGTBasisFuncs{T, D},
Vararg{FloatingGTBasisFuncs{T, D}}},
groupCenters::Bool=false;
roundAtol::Real=getAtolVal(T)) where {T, D} ->
Tuple
Quiqbox.sortPermBasis
β MethodsortPermBasis(bs::AbstractArray{<:CompositeGTBasisFuncs{T, D}};
roundAtol::Real=getAtolVal(T)) where {T, D} ->
Vector{Int}
Return a Vector
of indices I
such that bs[I] ==
sortBasis
(bs; roundAtol)[I]
.
Quiqbox.sortPermBasisFuncs
β MethodsortPermBasisFuncs(bs::Union{AbstractArray{<:FloatingGTBasisFuncs{T, D}},
Tuple{FloatingGTBasisFuncs{T, D},
Vararg{FloatingGTBasisFuncs{T, D}}}};
roundAtol::Real=getAtolVal(T)) where {T, D} ->
Vector{Int}
Return a Vector
of indices I
such that bs[I] ==
sortBasisFuncs
(bs; roundAtol)[I]
.
Quiqbox.subshellOf
β MethodsubshellOf(::FloatingGTBasisFuncs) -> String
Return the corresponding subshell of the input FloatingGTBasisFuncs
.
Quiqbox.unpackBasis
β MethodunpackBasis(b::GTBasisFuncs{T, D}) -> Vector{<:BasisFunc{T, D}}
Unpack b
to return all the BasisFunc
inside it.
Quiqbox.runHF
β MethodrunHF(bs, nuc, nucCoords, config=HFconfig(), N=getCharge(nuc);
printInfo=true, infoLevel=2) ->
HFfinalVars
runHF(bs, nuc, nucCoords, N=getCharge(nuc), config=HFconfig();
printInfo=true, infoLevel=2) ->
HFfinalVars
Main function to run a HartreeβFock method in Quiqbox. The returned result and relevant information is stored in a HFfinalVars
.
runHFcore(args...; printInfo=false, infoLevel=2) ->
Tuple{Tuple{Vararg{HFtempVars}}, Bool}
The core function of runHF
that accept the same positional arguments as runHF
, except it returns the data (HFtempVars
) collected during the iteration and the boolean result of whether the SCF procedure is converged.
β‘β‘β‘ Positional argument(s) β‘β‘β‘
bs::Union{ BasisSetData{T, D}, AbstractVector{<:AbstractGTBasisFuncs{T, D}}, Tuple{Vararg{AbstractGTBasisFuncs{T, D}}} } where {T, D}
: The basis set used for the HartreeβFock approximation.
nuc::Union{ Tuple{String, Vararg{String, NNMO}} where NNMO, AbstractVector{String} }
: The nuclei in the studied system.
nucCoords::Union{AbstractArray{Tuple{Vararg{T, D}}, 1}, Tuple{Tuple{Vararg{T, D}}, Vararg{Tuple{Vararg{T, D}}, NNMO}}, AbstractVector{<:AbstractVector{<:T}}, Tuple{TL, Vararg{TL, NNMO}} where TL<:(AbstractVector{<:T})} where {T, D, NNMO}
: The coordinates of corresponding nuclei.
config::HFconfig
: The Configuration of selected HartreeβFock method. For more information please refer to HFconfig
.
N::Union{Int, Tuple{Int}, NTuple{2, Int}}
: Total number of electrons, or the number(s) of electrons with same spin configurations(s). NOTE: N::NTuple{2, Int}
is only supported by unrestricted HartreeβFock (UHF).
β‘β‘β‘ Keyword argument(s) β‘β‘β‘
printInfo::Bool
: Whether print out the information of iteration steps and result.
infoLevel::Int
: Printed info's level of details when printInfo=true
. The higher (the absolute value of) it is, more intermediate steps will be printed. Once infoLevel
achieve 5
, every step will be printed.
Quiqbox.runHFcore
β MethodrunHFcore(HTtype, scfConfig, Ns, Hcore, HeeI, S, X, C0, maxStep, earlyStop, saveTrace,
printInfo=false, infoLevel=2) ->
Tuple{Tuple{Vararg{HFtempVars}}, Bool}
Another method of runHFcore
that has the same return value, but takes more underlying data as arguments.
=== Positional argument(s) ===
HTtype::Val{HFT} where HFT
: HartreeβFock method type. Available values of HFT
are :RHF, :UHF.
scfConfig::SCFconfig
: The SCF iteration configuration.
Ns::NTuple{HFTS, Int} where HFTS
: The numbers of electrons with same spin configurations.
Hcore::AbstractMatrix{T} where T
: The core Hamiltonian of the electronic Hamiltonian.
HeeI::AbstractArray{T, 4} where T
: The electron-electron interaction tensor (in the chemists' notation) which includes both the Coulomb interactions and the Exchange Correlations.
S::AbstractMatrix{T} where T
: The overlap matrix of the used basis set.
X::AbstractMatrix{T} where T
: The transformation matrix of S
.
C0::NTuple{HFTS, AbstractMatrix{T}} where {HFTS, T}
: Initial guess of the coefficient matrix(s) of the canonical spin-orbitals.
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 output HFtempVars
of runHFcore
. Its definition is the same as the field .saveTrace
inside a HFconfig
.
printInfo::Bool
: Whether print out the information of iteration steps and result.
infoLevel::Int
: Printed info's level of details when printInfo=true
. The higher (the absolute value of) it is, more intermediate steps and other information will be printed. Once infoLevel
achieve 5
, every step and all available information will be printed.
Quiqbox.gradOfHFenergy
β MethodgradOfHFenergy(par, basis, C, nuc, nucCoords, N=getCharge(nuc)) -> AbstractVector
gradOfHFenergy(par, bs, C, nuc, nucCoords, N=getCharge(nuc)) -> AbstractVector
Two methods of gradOfHFenergy
.
β‘β‘β‘ Positional argument(s) β‘β‘β‘
par::AbstractVector{<:ParamBox}
: The parameters for differentiation.
basis::
GTBasis
{T, D} where {T, D}
: Basis set information.
C::NTuple{<:Any, AbstractMatrix{T}} where T
: The coefficient matrix(s) of the canonical orbitals with respect to the selected basis set.
nuc::Union{ Tuple{String, Vararg{String, NNMO}} where NNMO, AbstractVector{String} }
: The nuclei in the studied system.
nucCoords::Union{AbstractArray{Tuple{Vararg{T, D}}, 1}, Tuple{Tuple{Vararg{T, D}}, Vararg{Tuple{Vararg{T, D}}, NNMO}}, AbstractVector{<:AbstractVector{<:T}}, Tuple{TL, Vararg{TL, NNMO}} where TL<:(AbstractVector{<:T})} where {T, D, NNMO}
: The coordinates of corresponding nuclei.
N::Union{Int, Tuple{Int}, NTuple{2, Int}}
: Total number of electrons, or the number(s) of electrons with same spin configurations(s).
bs::Union{ NTuple{BN, GTBasisFuncs{T, D, 1}}, AbstractVector{<:GTBasisFuncs{T, D, 1}} } where {T, D}
: A collection of basis functions.
NOTE: If any of these two methods is applied, the user needs to make sure the row orders as well as the colum orders of C
are consistent with the element order of bs
(basis.basis
). ``
Quiqbox.gridCoordOf
β MethodgridCoordOf(gb::GridBox{T}) where {T} -> Tuple{Vararg{Vector{T}}}
Return the coordinates of the grid points stored in gb
.
Quiqbox.optimizeParams!
β MethodoptimizeParams!(pbs, bs, nuc, nucCoords,
config=POconfig(), N=getCharge(nuc);
printInfo=true, infoLevel=1) ->
Vector{Any}
optimizeParams!(pbs, bs, nuc, nucCoords,
N=getCharge(nuc), config=POconfig();
printInfo=true, infoLevel=1) ->
Vector{Any}
The main function to optimize the parameters of a given basis set. It returns a Vector
of the relevant information. The first element is the indicator of whether the optimization is converged if the convergence detection is on (i.e., config.threshold
is not NaN
), or else it's set to missing
. More elements may be pushed to the returned result in order depending on config.method
.
=== Positional argument(s) ===
pbs::AbstractVector{<:ParamBox{T}}
: The parameters to be optimized that are extracted from the basis set. If the parameter is marked as "differentiable", the value of its input variable will be optimized.
bs::AbstractVector{<:GTBasisFuncs{T, D}}
: The basis set to be optimized.
nuc::Union{ Tuple{String, Vararg{String, NNMO}} where NNMO, AbstractVector{String} }
: The nuclei in the studied system.
nucCoords::Union{AbstractArray{Tuple{Vararg{T, D}}, 1}, Tuple{Tuple{Vararg{T, D}}, Vararg{Tuple{Vararg{T, D}}, NNMO}}, AbstractVector{<:AbstractVector{<:T}}, Tuple{TL, Vararg{TL, NNMO}} where TL<:(AbstractVector{<:T})} where {T, D, NNMO}
: The coordinates of corresponding nuclei.
config::POconfig
: The Configuration of selected parameter optimization method. For more information please refer to POconfig
.
N::Union{Int, Tuple{Int}, NTuple{2, Int}}
: Total number of electrons, or the number(s) of electrons with same spin configurations(s).
=== Keyword argument(s) ===
printInfo::Bool
: Whether print out the information of iteration steps.
infoLevel::Int
: Printed info's level of details when printInfo=true
. The higher (the absolute value of) it is, more intermediate steps and other information will be printed. Once infoLevel
achieve 5
, every step and all available information will be printed.
Quiqbox.changeHbasis
β MethodchangeHbasis(HFres::HFfinalVars) -> NTuple{2, Any}
Return the one-body and two-body integrals on the basis of the canonical orbitals using the result of a HartreeβFock method HFres
.
Quiqbox.changeHbasis
β MethodchangeHbasis(b::GTBasis{T, D},
nuc::Tuple{String, Vararg{String, NNMO}},
nucCoords::Tuple{NTuple{D, T}, Vararg{NTuple{D, T}, NNMO}},
C::Union{AbstractMatrix{T}, NTuple{2, AbstractMatrix{T}}}) where
{T, D, NNMO} ->
NTuple{2, Any}
Return the one-body and two-body integrals after a change of basis based on the input C
, given the basis set information b
. The type of each element in the returned Tuple
is consistent with the cases where the first argument of changeHbasis
is an AbstractArray
.
Quiqbox.changeHbasis
β MethodchangeHbasis(twoBodyInt::AbstractArray{T, 4},
C1::AbstractMatrix{T}, C2::AbstractMatrix{T}) where {T} ->
AbstractArray{T, 4}
Change the basis of the input two-body integrals twoBodyInt
based on two orbital coefficient matrices C1
and C2
for different spin configurations (e.g., the unrestricted case). The output is a 3-element Tuple
of which the first 2 elements are the spatial integrals of each spin configurations respectively, while the last element is the Coulomb interactions between orbitals with different spins.
Quiqbox.changeHbasis
β MethodchangeHbasis(DbodyInt::AbstractArray{T, D}, C::AbstractMatrix{T}) where {T} ->
AbstractArray{T, D}
Change the basis of the input one-body / two-body integrals DbodyInt
based on the orbital coefficient matrix C
.
Quiqbox.genCanOrbitals
β MethodgenCanOrbitals(fVars::HFfinalVars{T, D, <:Any, NN};
roundAtol::Real=getAtolVal(T)) where {T, D, NN} ->
NTuple{2, Vector{CanOrbital{T, D, NN}}}
Generate the occupied and unoccupied canonical orbitals from the result of a HartreeβFock approximation fVars
. Each parameter stored in the constructed CanOrbital
will be rounded to the nearest multiple of roundAtol
; when roundAtol
is set to NaN
, no rounding will be performed.
Quiqbox.nnRepulsions
β MethodnnRepulsions(nuc::Union{Tuple{String, Vararg{String, NNMO}}, AbstractVector{String}},
nucCoords::Union{AbstractArray{Tuple{Vararg{T, D}}, 1}, Tuple{Tuple{Vararg{T, D}}, Vararg{Tuple{Vararg{T, D}}, NNMO}}, AbstractVector{<:AbstractVector{<:T}}, Tuple{TL, Vararg{TL, NNMO}} where TL<:(AbstractVector{<:T})}) where {NNMO, D, T} ->
T
Return the nuclear repulsion energy given nuclei nuc
and their coordinates nucCoords
.
Quiqbox.coreH
β MethodcoreH(bs::Union{GTBasis, Tuple{Vararg{AbstractGTBasisFuncs{T, D}}},
AbstractVector{<:AbstractGTBasisFuncs{T, D}}},
nuc::Union{Tuple{String, Vararg{String, NNMO}}, AbstractVector{String}},
nucCoords::Union{AbstractArray{Tuple{Vararg{T, D}}, 1}, Tuple{Tuple{Vararg{T, D}}, Vararg{Tuple{Vararg{T, D}}, NNMO}}, AbstractVector{<:AbstractVector{<:T}}, Tuple{TL, Vararg{TL, NNMO}} where TL<:(AbstractVector{<:T})}) where {T, D, NNMO} ->
Matrix{T}
Return the core Hamiltonian given a basis set and the corresponding nuclei with their coordinates (in atomic units).
Quiqbox.coreHij
β MethodcoreHij(bf1::AbstractGTBasisFuncs{T, D, 1}, bf2::AbstractGTBasisFuncs{T, D, 1},
nuc::Union{Tuple{String, Vararg{String, NNMO}}, AbstractVector{String}},
nucCoords::Union{AbstractArray{Tuple{Vararg{T, D}}, 1}, Tuple{Tuple{Vararg{T, D}}, Vararg{Tuple{Vararg{T, D}}, NNMO}}, AbstractVector{<:AbstractVector{<:T}}, Tuple{TL, Vararg{TL, NNMO}} where TL<:(AbstractVector{<:T})}) where {T, D, NNMO} ->
T
Return a matrix element of the core Hamiltonian given two basis functions.
Quiqbox.eKinetic
β MethodeKinetic(bf1::AbstractGTBasisFuncs{T, D, 1}, bf2::AbstractGTBasisFuncs{T, D, 1}) where
{T, D} ->
T
Return the electron kinetic energy between two basis functions.
Quiqbox.eKinetics
β MethodeKinetics(bs::Union{Tuple{Vararg{AbstractGTBasisFuncs{T, D}}},
AbstractVector{<:AbstractGTBasisFuncs{T, D}}}) where {T, D} ->
Matrix{T}
Return the electron kinetic energy matrix given a basis set.
Quiqbox.neAttraction
β MethodneAttraction(bf1::AbstractGTBasisFuncs{T, D, 1}, bf2::AbstractGTBasisFuncs{T, D, 1},
nuc::Union{Tuple{String, Vararg{String, NNMO}}, AbstractVector{String}},
nucCoords::Union{AbstractArray{Tuple{Vararg{T, D}}, 1}, Tuple{Tuple{Vararg{T, D}}, Vararg{Tuple{Vararg{T, D}}, NNMO}}, AbstractVector{<:AbstractVector{<:T}}, Tuple{TL, Vararg{TL, NNMO}} where TL<:(AbstractVector{<:T})}) where {T, D, NNMO} ->
T
Return the nuclear attraction between two basis functions, provided with the nuclei and their coordinates (in the atomic units).
Quiqbox.neAttractions
β MethodneAttractions(bs::Union{Tuple{Vararg{AbstractGTBasisFuncs{T, D}}},
AbstractVector{<:AbstractGTBasisFuncs{T, D}}},
nuc::Union{Tuple{String, Vararg{String, NNMO}}, AbstractVector{String}},
nucCoords::Union{AbstractArray{Tuple{Vararg{T, D}}, 1}, Tuple{Tuple{Vararg{T, D}}, Vararg{Tuple{Vararg{T, D}}, NNMO}}, AbstractVector{<:AbstractVector{<:T}}, Tuple{TL, Vararg{TL, NNMO}} where TL<:(AbstractVector{<:T})}) where {T, D, NNMO} ->
Matrix{T}
Return the nuclear attraction matrix given a basis set and the corresponding nuclei with their coordinates (in atomic units).
Quiqbox.overlap
β Methodoverlap(bf1::AbstractGTBasisFuncs{T, D, 1}, bf2::AbstractGTBasisFuncs{T, D, 1}) where
{T, D, 1} ->
T
Return the orbital overlap between two basis functions.
Quiqbox.overlaps
β Methodoverlaps(bs::Union{Tuple{Vararg{AbstractGTBasisFuncs{T, D}}},
AbstractVector{<:AbstractGTBasisFuncs{T, D}}}) where {T, D} ->
Matrix{T}
Return the orbital overlap matrix given a basis set.
Quiqbox.eeInteraction
β MethodeeInteraction(bf1::AbstractGTBasisFuncs{T, D, 1},
bf2::AbstractGTBasisFuncs{T, D, 1},
bf3::AbstractGTBasisFuncs{T, D, 1},
bf4::AbstractGTBasisFuncs{T, D, 1}) where {T, D} ->
T
Return an electron-electron interaction tensor element given four basis functions (ordered in the chemists' notation).
Quiqbox.eeInteractions
β MethodeeInteractions(bs::Union{Tuple{Vararg{AbstractGTBasisFuncs{T, D}}},
AbstractVector{<:AbstractGTBasisFuncs{T, D}}}) ->
Array{T, 4}
Return the tensor of electron-electron interactions (in the chemists' notation) given a basis set.
Quiqbox.eeIuniqueIndicesOf
β MethodeeIuniqueIndicesOf(basisSetSize::Int) -> Vector{Vector{Int}}
Return the unique matrix element indices (in the chemists' notation) of electron-electron interactions given the size of a basis set.