Expand description
This project contains libcint (C language) FFI bindings, wrapper and build-from-source.
§Quick Links to Important Parts
CInt
basic struct for defining molecule for electron integral.CInt::integrate
: integrate function (partly PySCF’smol.intor
counterpart, column-major which is the same to convention of libcint C library, and the same memory layout to PySCF’s 2/3-center integrals).CInt::integrate_row_major
: integrate function with row-major output (the same shape convention to PySCF, and same memory layout to PySCF 4-center integrals).CInt::integrate_cross
: integrate function with multipleCInt
s, useful for integral with auxiliary basis sets.CInt::integrate_with_args
,CInt::integrate_args_builder
: integrate function with arguments builder (advanced options).- cint_wrapper and cecp_wrapper: supported integrators.
§Introduction
libcint is a C library for GTO (gaussian-type orbital) electronic integral, can be applied in computational chemistry, and has already been applied extensively in PySCF.
Resources | Badges |
---|---|
Crate | |
API Document | |
FFI Binding (libcint) | v6.1.2 |
FFI Binding (qcint) | v6.1.2 |
ECP Source Code (PySCF) | v2.9.0 |
This crate is not official bindgen project, nither libcint, PySCF, nor REST. It is originally intended to be some pioneer work for possible future development of rest_libcint wrapper.
§Minimal Example
The most important function is CInt::integrate
. It is somehow similar
to PySCF’s mol.intor(intor, aosym, shls_slice)
.
use libcint::prelude::*;
// This is for testing and development purpose only.
// For actual usage, you should initialize `CInt` with your own molecule data.
let cint_data: CInt = init_h2o_def2_tzvp();
// int1e_ipkin: the integral of kinetic energy operator with derivative on the first orbital
// [mu, nu, comp], column-major, same data memory layout to PySCF's 2/3-center intor
let (out, shape): (Vec<f64>, Vec<usize>) = cint_data.integrate("int1e_ipkin", None, None).into();
assert_eq!(shape, [43, 43, 3]);
// [comp, mu, nu], row-major, same shape to PySCF intor
let (out, shape): (Vec<f64>, Vec<usize>) = cint_data.integrate_row_major("int1e_ipkin", None, None).into();
assert_eq!(shape, [3, 43, 43]);
§For users from PySCF
This library corresponds to some parts of pyscf.gto
module in PySCF.
Similar parts are:
CInt::integrate
corresponds tomol.intor(intor, aosym, shls_slice)
;- Various get/set/with-clauses methods;
CInt
corresponds to(mol._atm, mol._bas, mol._env)
in PySCF (addingmol._ecpbas
for ECP).
Differences are:
- Recall that PySCF’s
Mole
class handles three parts: Python input from user (mol.atom
,mol.basis
, etc.), Python intermediates (mol._basis
, etc.), data for C-FFI (mol._atm
,mol._bas
,mol._env
). In rust’sCInt
, it only handles the last part. - This library does not handle molecule input (coords, spin, charge, etc.) and basis set parse. Molecule initialization should be performed by user.
§Installation and Cargo Features
§Install with pre-compiled libcint.so
(recommended)
If you have already compiled libcint.so
, then put path of this shared
library in CINT_DIR
or LD_LIBRARY_PATH
(or REST_EXT_DIR
). Then you
just use this library in your Cargo.toml
file by
[dependencies]
libcint = { version = "0.1" }
§Install and also build-from-source
If you have not compiled libcint.so
or libcint.a
, then you are suggested
to use this library by specifying some cargo features:
[dependencies]
libcint = { version = "0.1", features = ["build_from_source", "static"] }
The source code will be automatically downloaded from github, and cargo will handle the building process.
If access to github is not available, you can use environment variable
CINT_SRC
to specify source mirror of
sunqm/libcint
or sunqm/qcint.
§Cargo features
- Default features: None of any listed below (use library provided by system or user, using sunqm/libcint, dynamic linking, without F12 and 4c1e support).
build_from_source
: Trigger of C language library libcint building. This performs by CMake; source code will be automatically downloaded from github (if environment variableCINT_SRC
not specified).static
: Use static library for linking. This will require static linklibcint.a
, and dynamic linklibquadmath.so
.qcint
: Use sunqm/qcint instead of sunqm/libcint. Some integrals will not be available ifqcint
does not supports that. This will also change URL source if cargo featurebuild_from_source
specified.with_f12
: Whether F12 integrals (int2e_stg
,int2e_yp
, etc.) are supported.with_4c1e
: Whether 4c1e integrals (int4c1e
, etc.) are supported.
§Shell environment variables
CINT_DIR
,LD_LIBRARY_PATH
,REST_EXT_DIR
: Your compiled library path oflibcint.so
andlibcint.a
. This crate will try to find if this library is in directory root, ordirectory_root/lib
. May override the library built by cargo featurebuild_from_source
.CINT_SRC
: Source of libcint or qcint (must be a git repository). Only works with cargo featurebuild_from_source
.CINT_VIR
: Version of libcint or qcint (e.g.v6.1.2
, must starts withv
prefix). Only works with cargo featurebuild_from_source
.
§50 lines RHF with Rust
This is an example code to compute the RHF energy of H2O molecule using
libcint
as electronic integral library (corresponding to some parts ofpyscf.gto
);rstsr
as tensor library (corresponding to NumPy and SciPy).
You will see that except for import and preparation, the code with core algorithms is 27 lines (blank and printing lines included). For comparison, the same code in Python is 23 lines.
use libcint::prelude::*;
use rstsr::prelude::*;
pub type Tsr = Tensor<f64, DeviceBLAS>;
pub type TsrView<'a> = TensorView<'a, f64, DeviceBLAS>;
/// Obtain integrals (in row-major, same to PySCF but reverse of libcint)
pub fn intor_row_major(cint_data: &CInt, intor: &str) -> Tsr {
// use up all rayon available threads for tensor operations
let device = DeviceBLAS::default();
// intor, "s1", full_shls_slice
let (out, shape) = cint_data.integrate(intor, None, None).into();
// row-major by transposition of col-major shape
rt::asarray((out, shape.f(), &device)).into_reverse_axes()
}
fn main() {
let device = DeviceBLAS::default();
// Assuming H2O/def2-TZVP data for `CInt` has been prepared
let cint_data = init_h2o_def2_tzvp();
/* #region rhf total energy algorithms */
let atom_coords = {
let coords = cint_data.atom_coords();
let coords = coords.into_iter().flatten().collect::<Vec<f64>>();
rt::asarray((coords, &device)).into_shape((-1, 3))
};
let atom_charges = rt::asarray((cint_data.atom_charges(), &device));
let mut dist = rt::sci::cdist((atom_coords.view(), atom_coords.view()));
dist.diagonal_mut(None).fill(f64::INFINITY);
let eng_nuc = 0.5 * (&atom_charges * atom_charges.i((.., None)) / dist).sum();
println!("Nuclear repulsion energy: {eng_nuc}");
let hcore = intor_row_major(&cint_data, "int1e_kin") + intor_row_major(&cint_data, "int1e_nuc");
let ovlp = intor_row_major(&cint_data, "int1e_ovlp");
let int2e = intor_row_major(&cint_data, "int2e");
let nocc = 5; // hardcoded for H2O, 5 occupied orbitals
let mut dm = ovlp.zeros_like();
for _ in 0..40 {
// hardcoded SCF iterations
let fock = &hcore + ((1.0_f64 * &int2e - 0.5_f64 * int2e.swapaxes(1, 2)) * &dm).sum_axes([-1, -2]);
let (_, c) = rt::linalg::eigh((&fock, &ovlp)).into();
dm = 2.0_f64 * c.i((.., ..nocc)) % c.i((.., ..nocc)).t();
}
let eng_scratch = &hcore + ((0.5_f64 * &int2e - 0.25_f64 * int2e.swapaxes(1, 2)) * &dm).sum_axes([-1, -2]);
let eng_elec = (dm * &eng_scratch).sum();
println!("Total elec energy: {eng_elec}");
println!("Total RHF energy: {}", eng_nuc + eng_elec);
/* #endregion */
}
Modules§
- cint
- Main
CInt
with related struct definition, and impl ofCInt::integrate
. - cint_
change - Getter/setter and with-clause temporary changes to
CInt
instance. - cint_
crafter - Integral crafter for
CInt
instance (the main integral functions). - cint_
fill_ col_ major - Implementation of integral at lower API level (crafting, col-major).
- cint_
fill_ grids - Implementation of grids inplementation at lower API level.
- cint_
fill_ row_ major - Implementation of integral at lower API level (crafting, row-major).
- cint_
initialize - Initializers for the CINT library.
- cint_
prop - Properties of the
CInt
instance (compute and memory cost is not essential in this module). - ffi
- FFI binding and basic wrapper.
- prelude
- Use
libcint::prelude::*
to import this crate. - test_
mol - Data and functions only for doc and unit testing.
- util
- Utility functions and structs (not for users).