Basis set
The basis set is defined as a dictionary, where the keys are the types of the basis sets, and the values are the basis set:
basis = Dict("ao"=>"cc-pVTZ",
"jkfit"=>"cc-pvtz-jkfit",
"mpfit"=>"cc-pvtz-mpfit")The basis set dictionary contains three keys: ao, jkfit, and mpfit. The ao key contains the basis set for the AO integrals, the jkfit key contains the basis set for the density fitting integrals in the Hartree-Fock calculations, and the mpfit key contains the fitting basis set for the correlated calculations.
Alternatively, you can define the basis set using a string that defines the AO basis. In this case, the jkfit and mpfit basis names will be generated automatically. Here's an example of how you can define the basis set using a string:
basis = "cc-pVDZ"Common acronyms are also supported for the basis set names, e.g., cc-pVDZ can be written as vdz, and aug-cc-pVTZ can be written as avtz.
ElemCo.BasisSets — ModuleBasisSetsModule for working with basis sets.
The basis set is stored in the BasisSet instance. The basis set can be generated using the generate_basis function.
Exported functions and types
ElemCo.BasisSets.AngularShell — TypeAngularShellType for angular shells, i.e, subshells with the same angular momentum. For general contracted basis sets, the angular shell is a collection of all subshells with the same l quantum number. For some other basis sets (e.g., the def2-family), the angular shell can be a single subshell with a specific l quantum number. id is the index of the angular shell in the basis set.
element::String: element symbol (e.g., "H")l::Int64: angular momentumexponents::Vector{Float64}: array of exponentssubshells::Vector{ElemCo.BasisSets.BasisContraction}: array of subshells (contractions)id::Int64: index of the angular shell in the basis set
ElemCo.BasisSets.BasisCenter — TypeBasisCenterA basis center (atom) with basis functions.
name::String: basis center name (e.g., "H1")position::StaticArraysCore.SVector{3, Float64}: atomic position in Bohr (3D vector)atomic_number::Int64: atomic numberbasis::String: name of the basis set (e.g., "cc-pVDZ")shells::Vector{ElemCo.BasisSets.AngularShell}: array of angular shells
ElemCo.BasisSets.BasisContraction — TypeBasisContractionA basis contraction. exprange is the range of primitives (from exponents in the angular shell).
ElemCo.BasisSets.BasisSet — TypeBasisSetA basis set with basis centers (atoms) and basis functions.
centers::Vector{ElemCo.BasisSets.BasisCenter}: array of basis centers (atoms) with basis functions.shell_indices::Vector{CartesianIndex{2}}: indices for angular shells.center_ranges::Vector{UnitRange{Int64}}: center ranges for each basis set in a combined set.shell_ranges::Vector{UnitRange{Int64}}: angular shell ranges for each basis sets in a combined set.cartesian::Bool: cartesian basis setlib::ElemCo.BasisSets.ILibcint5: infos for integral library (at the moment only libcint5 is possible).
ElemCo.BasisSets.ILibcint5 — TypeILibcint5Infos for Libcint5 integral library.
ElemCo.BasisSets.ILibcint5 — MethodILibcint5(atoms::Vector{BasisCenter}, cartesian::Bool)Prepare the infos for Libcint5 integral library.
ElemCo.BasisSets.ao_list — Functionao_list(basis::BasisSet, ibas=1)Return the list of atomic orbitals in the basis set.
For a combined basis set, use ibas to select the basis set.
ElemCo.BasisSets.basis_name — Functionbasis_name(atoms, type="ao")Return the name of the basis set (or unknown if not found). atoms can be a single atom ::Atom or a system ::AbstractSystem.
ElemCo.BasisSets.center_range — Functioncenter_range(bs::BasisSet, i::Int=1)Return the range of centers for the ith basis set.
The range is used to access the centers in the basis set, e.g., bs.centers[i] for i in center_range(bs, 1) gives the centers of the first basis set.
ElemCo.BasisSets.coefficients_1mat — Methodcoefficients_1mat(ashell::AngularShell, cartesian::Bool)Return a single contraction matrix of the coefficients in the angular shell (nprimitives × nsubshells). The contractions are normalized.
The missing coefficients are set to zero in the matrix.
ElemCo.BasisSets.combine — Methodcombine(bs1::BasisSet, bs2::BasisSet)Combine two basis sets.
The centers are concatenated. The center/shell ranges (center_ranges/shell_ranges) corresponding to the centers/shells for each basis set can be used to access the centers/shells in the combined basis set, e.g. bs.centers[i] for i in bs.center_ranges[1] gives the centers of the first basis set in the combined set.
ElemCo.BasisSets.generate_basis — Functiongenerate_basis(EC::ECInfo, type="ao"; basisset::AbstractString="")Generate basis sets for integral calculations.
The basis set is stored in BasisSet object. type can be "ao", "mpfit" or "jkfit". If basisset is provided, it is used as the basis set.
ElemCo.BasisSets.generate_basis — Functiongenerate_basis(ms::AbstractSystem, type="ao"; cartesian=false, basisset::AbstractString="")Generate basis sets for integral calculations.
The basis set is stored in BasisSet object. type can be "ao", "mpfit" or "jkfit". If basisset is provided, it is used as the basis set.
ElemCo.BasisSets.guess_norb — Methodguess_norb(EC::AbstractECInfo)Guess the number of orbitals in the system.
ElemCo.BasisSets.is_cartesian — Methodis_cartesian(bs::BasisSet)Check if the basis set is Cartesian.
ElemCo.BasisSets.max_l — Methodmax_l(basis::BasisSet)Return the maximum angular momentum in the basis set.
ElemCo.BasisSets.n_angularshells — Methodn_angularshells(atom::BasisCenter)Return the number of angular shells in the basis set for atom.
ElemCo.BasisSets.n_angularshells — Methodn_angularshells(atoms::BasisSet)Return the number of angular shells in the basis set.
ElemCo.BasisSets.n_ao — Methodn_ao(ashell::AngularShell, cartesian::Bool)Return the number of atomic orbitals in the angular shell.
ElemCo.BasisSets.n_ao — Methodn_ao(atom::BasisCenter, cartesian::Bool)Return the number of atomic orbitals in the basis set for atom.
ElemCo.BasisSets.n_ao — Methodn_ao(atoms::BasisSet)Return the number of atomic orbitals in the basis set.
ElemCo.BasisSets.n_coefficients — Methodn_coefficients(ashell::AngularShell)Return the total number of coefficients in the angular shell.
ElemCo.BasisSets.n_coefficients — Methodn_coefficients(subshell::BasisContraction)Return the number of coefficients for the subshell.
ElemCo.BasisSets.n_coefficients_1mat — Methodn_coefficients_1mat(ashell::AngularShell)Return the number of coefficients in the angular shell for a single contraction matrix.
The missing coefficients will be set to zero.
ElemCo.BasisSets.n_primitives — Methodn_primitives(ashell::AngularShell)Return the total number of primitives in the angular shell.
ElemCo.BasisSets.n_primitives — Methodn_primitives(atom::BasisCenter)Return the number of primitives in the basis set for atom.
ElemCo.BasisSets.n_primitives — Methodn_primitives(subshell::BasisContraction)Return the number of primitives for the subshell.
ElemCo.BasisSets.n_primitives — Methodn_primitives(atoms::BasisSet)Return the number of primitives in the basis set.
ElemCo.BasisSets.n_subshells — Methodn_subshells(ashell::AngularShell)Return the number of subshells in the angular shell.
ElemCo.BasisSets.n_subshells — Methodn_subshells(atom::BasisCenter)Return the number of subshells in the basis set for atom.
ElemCo.BasisSets.n_subshells — Methodn_subshells(atoms::BasisSet)Return the number of subshells in the basis set.
ElemCo.BasisSets.normalize_contraction — Methodnormalize_contraction(subshell::BasisContraction, ashell::AngularShell, cartesian::Bool)Normalize the contraction coefficients in subshell.
The subshell has to be part of the angular shell ashell. Return the normalized contraction.
ElemCo.BasisSets.print_ao — Methodprint_ao(ao::AbstractAtomicOrbital, basis::BasisSet)
Print the atomic orbital.
ElemCo.BasisSets.shell_range — Functionshell_range(bs::BasisSet, i::Int=1)Return the range of angular shells for the ith basis set.
The range is used to access the angular shells in the basis set, e.g., bs[i] for i in shell_range(bs, 1) gives the angular shells of the first basis set.
ElemCo.BasisSets.subshell_char — Methodsubshell_char(l)Return the character for the subshell with angular momentum l.
Internal functions and types
ElemCo.BasisSets.AbstractILib — TypeAbstractILibAbstract type for infos for integral libraries.
ElemCo.BasisSets.CartesianAtomicOrbital — TypeCartesianAtomicOrbital
Represents an atomic orbital in a cartesian basis set.
icenter::UInt16: index of the center in the basis set objectiangularshell::UInt16: index of the angular shell in the basis set objectisubshell::UInt8: index of the subshell in the angular shelln::UInt8: principal quantum numberl::UInt8: orbital angular momentum quantum numberml::Int8: magnetic "quantum number" (ml = -l:l(l+1)/2)
ElemCo.BasisSets.SphericalAtomicOrbital — TypeSphericalAtomicOrbital
Represents an atomic orbital in a spherical basis set.
icenter::UInt16: index of the center in the basis set objectiangularshell::UInt16: index of the angular shell in the basis set objectisubshell::UInt8: index of the contraction in the angular shelln::UInt8: principal quantum numberl::UInt8: orbital angular momentum quantum numberml::Int8: magnetic quantum number (ml = -l:l)
Base.getindex — MethodBase.getindex(bs::BasisSet, i::CartesianIndex{2})Return the the angular shell.
Base.getindex — MethodBase.getindex(bs::BasisSet, i::Int, j::Int)Return the the angular shell.
Base.getindex — MethodBase.getindex(bs::BasisSet, i::Int)Return the the angular shell.
Base.iterate — FunctionBase.iterate(bs::BasisSet, state=1)Iterate over the basis angular shells in the basis set.
Base.length — MethodBase.length(bs::BasisSet)Return the number of angular shells in the basis set.
ElemCo.BasisSets.add_subshell! — Methodadd_subshell!(ashell::AngularShell, exprange, contraction)Add a subshell to the angular shell.
ElemCo.BasisSets.basis_file — Methodbasis_file(basis_name::AbstractString)Return the full path to the basis set file.
ElemCo.BasisSets.full_basis_name — Methodfull_basis_name(basis_name::AbstractString)Return the full basis name and version number (if given as *.v[0-2], otherwise -1 is returned).
I.e,
[a][wc/c]vXz*->[aug-]cc-p[wc/c]vXz*svp*->def2-svp*[tq]zvp*->def2-[tq]zvp*
Additionally check for version number (e.g., vdz.v2)
ElemCo.BasisSets.generate_angularshell — Methodgenerate_angularshell(elem, l, exponents)Generate an angular shell with angular momentum l and exponents. The contractions have to be added later. Return an angular shell of type AngularShell.
ElemCo.BasisSets.guess_basis_name — Methodguess_basis_name(atom::Atom, type)Guess the name of the basis set. type can be "ao", "mpfit" or "jkfit".
ElemCo.BasisSets.n_ao4subshell — Methodn_ao4subshell(ashell::AngularShell, cartesian::Bool)Return the number of atomic orbitals for the subshell.
ElemCo.BasisSets.normalize_cartesian_contraction — Methodnormalize_cartesian_contraction(contraction, exponents, l)Normalize the subshell.
Return the normalized contraction.
ElemCo.BasisSets.normalize_spherical_contraction — Methodnormalize_spherical_contraction(contraction, exponents, l)Normalize the spherical subshell.
Return the normalized contraction.
ElemCo.BasisSets.parse_basis — Methodparse_basis(basis_name::String, atom::Atom)Search and parse the basis set for a given atom.
Return a list of angular shells AngularShell.
ElemCo.BasisSets.parse_basis_block — Methodparse_basis_block(basisblock::AbstractString, atom::Atom)Parse the basis block for a given atom.
Return a list of angular shells AngularShell. The basis block is in the Molpro format:
!commentss,p,d,f,g,hangular momentumc, <from>.<to>contraction coefficients for primitives
Example cc-pVDZ for H atom:
s, H , 13.0100000, 1.9620000, 0.4446000, 0.1220000
c, 1.4, 0.0196850, 0.1379770, 0.4781480, 0.5012400
c, 4.4, 1.0000000
p, H , 0.7270000
c, 1.1, 1.0000000For generally-contracted basis sets (like the one above), one angular shell is created for each angular momentum type s,p,d,f,g,h with the corresponding exponents and contraction coefficients. For other basis sets, like the def2-SVP, each contraction is a separate angular shell:
! hydrogen (4s,1p) -> [2s,1p]
s, H , 13.0107010, 1.9622572, 0.44453796, 0.12194962
c, 1.3, 0.19682158E-01, 0.13796524, 0.47831935
c, 4.4, 1.0000000
p, H , 0.8000000
c, 1.1, 1.0000000ElemCo.BasisSets.parse_contraction — Methodparse_contraction(conline::AbstractString)Parse contraction coefficients from a line in the basis block.
Return the range of exponents and the contraction coefficients as a tuple. The line is in the Molpro format: c, 1.4, 0.0196850, 0.1379770, 0.4781480, 0.5012400 where c is the contraction, 1.4 is the exponent range, and the rest are the coefficients.
ElemCo.BasisSets.parse_exponents — Methodparse_exponents(expline::AbstractString)Parse exponents from a line in the basis block.
Return the angular momentum and exponents as a tuple. The line is in the Molpro format: s, H , 13.0100000, 1.9620000, 0.4446000, 0.1220000 where s is the angular momentum, H is the element symbol, and the rest are the exponents.
ElemCo.BasisSets.read_basis_block — Methodread_basis_block(basisfile::AbstractString, atom::Atom)Read the basis block for a given atom.
The basis library is in the Molpro format:
!comments- basis block starts with
! <elementname> .... - basis block ends with
!or} - basis block contains:
s,p,d,f,g,hangular momentumc, <from>.<to>contraction coefficients for primitives
Example cc-pVDZ for H atom:
!
! hydrogen (4s,1p) -> [2s,1p]
s, H , 13.0100000, 1.9620000, 0.4446000, 0.1220000
c, 1.4, 0.0196850, 0.1379770, 0.4781480, 0.5012400
c, 4.4, 1.0000000
p, H , 0.7270000
c, 1.1, 1.0000000
!ElemCo.BasisSets.set_id! — Methodset_id!(ashells::AbstractArray{AngularShell}, start_id)Set the id for each angular shell in the array. Return the next id.
ElemCo.BasisSets.set_id! — Methodset_id!(centers::AbstractArray{BasisCenter}, start_id)Set the id for each angular shell in the array of centers. Return the next id.
ElemCo.BasisSets.split_angular_shell — Methodsplit_angular_shell(ashell::AngularShell)If the ranges of exponents do not overlap, split the angular shell into separate angular shells for each subshell. The shells are kept together only if one is a subset of the other.