Julia implementation of various electron-correlation methods (main focus on coupled cluster methods).
The integrals are obtained from a FCIDUMP file or calculated using the GaussianBasis
package.
Canonical | DF | Only closed-shell | |
---|---|---|---|
RHF | ❌ | ✔️ | ❗ |
UHF | ❌ | ✔️ | |
BO-HF | ✔️ | ||
MCSCF | ❌ | 🔧 | |
MP2 | ✔️ | ||
CCSD | ✔️ | ✔️ | |
RCCSD | ✔️ | ||
UCCSD | ✔️ | ||
ΛCCSD | ✔️ | ✔️ | |
ΛUCCSD | ✔️ | ||
CCSD(T) | ✔️ | ✔️ | |
UCCSD(T) | ✔️ | ||
ΛCCSD(T) | ✔️ | ✔️ | |
ΛUCCSD(T) | ✔️ | ||
FR-CCSD | ✔️ | ||
2D-CCSD | ✔️ | ||
DCSD | ✔️ | ✔️ | |
RDCSD | ✔️ | ||
UDCSD | ✔️ | ||
ΛDCSD | ✔️ | ✔️ | |
ΛUDCSD | ✔️ | ||
FR-DCSD | ✔️ | ||
2D-DCSD | ✔️ | ||
SVD-DCSD | ✔️ | ✔️ | ❗ |
SVD-DC-CCSDT | ✔️ | ✔️ | ❗ |
Requirements: julia (>1.8)
Packages: LinearAlgebra, NPZ, Mmap, TensorOperations, Printf, IterativeSolvers, GaussianBasis, DocStringExtensions, MKL(optional)
For a development version of ElemCo.jl
, clone the repository and create an alias to set the project to the ElemCo.jl
directory,
alias jlm='julia --project=<path_to_ElemCo.jl>'
Now the command jlm
can be used to start the calculations,
jlm input.jl
Default scratch directory path on Windows is the first environment variable found in the ordered list TMP
, TEMP
, USERPROFILE
.
On all other operating systems TMPDIR
, TMP
, TEMP
, and TEMPDIR
. If none of these are found, the path /tmp
is used.
Default scratch folder name is elemcojlscr
.
Variable names fcidump
, geometry
and basis
are reserved for the file name of FCIDUMP, geometry specification and basis sets, respectively.
The ground state energy can be calculated using the DCSD method with the following script:
using ElemCo
@print_input
fcidump = "../test/H2O.FCIDUMP"
@cc dcsd
In order to calculate the ground state energy of the water molecule using the DCSD method, the following script can be used:
using ElemCo
@print_input
geometry="bohr
O 0.000000000 0.000000000 -0.130186067
H1 0.000000000 1.489124508 1.033245507
H2 0.000000000 -1.489124508 1.033245507"
basis = Dict("ao"=>"cc-pVDZ",
"jkfit"=>"cc-pvtz-jkfit",
"mp2fit"=>"cc-pvdz-rifit")
@dfhf
@cc dcsd
The @dfhf
macro calculates the density-fitted Hartree-Fock energy and orbitals
and then DCSD calculation is performed using density-fitted integrals.
Further example scripts are provided in the examples
directory.
Documentation is available at https://elem.co.il.
Electron coil
A poem by Bing
In the heart of an atom lies a tiny core
Where protons and neutrons are tightly bound
But around this nucleus, there's so much more
A cloud of electrons that swirls around
They don't orbit in circles like planets do
But jump and spin in quantum states
They can be here and there and everywhere too
And sometimes they even change their mates
This is the electron coil, the turmoil of the shell
The source of light and heat and power
The force that makes the atoms repel or gel
The spark that ignites the cosmic flower