Expand description
This project contains libcint (C language) FFI bindings, wrapper and build-from-source.
§Quick Links to Important Parts
CIntbasic struct for defining molecule for electron integral.CInt::integrate: integrate function (partly PySCF’smol.intorcounterpart, 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::eval_gto: evaluate GTO values on grids, useful for DFT computation (PySCF’smol.eval_aoormol.eval_gtocounterpart, with same memory layout to PySCF, but note the returned shape is in column-major).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 multipleCInts, useful for integral with auxiliary basis sets.CInt::integrate_with_args,CInt::integrate_args_builder: integrate function with arguments builder (advanced options).CIntMol(featurebse): molecule builder from TOML/JSON, a weaker version of PySCF’sMoleclass (handling only atom coordinates and basis sets).- 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, neither 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]);§Building Molecule from TOML/JSON (feature bse)
With the bse feature, you can build a CIntMol object from TOML/JSON:
use libcint::prelude::*;
let toml = r#"
atom = "O 0 0 0; H 0 0 0.9572; H 0 0.9266 -0.239987"
basis = "STO-3G"
"#;
let mol = CIntMol::from_toml(toml).unwrap();
let ovlp = mol.cint.integrate("int1e_ovlp", None, None).out.unwrap();See parse::mole::CIntMol and parse::serde_build for details.
Full TOML specification is documented in
libcint/src/parse/from_toml_docs.md.
§For users from PySCF
This library corresponds to some parts of pyscf.gto module in PySCF.
Similar parts are:
CInt::integratecorresponds tomol.intor(intor, aosym, shls_slice);- Various get/set/with-clauses methods;
CIntcorresponds to(mol._atm, mol._bas, mol._env)in PySCF (addingmol._ecpbasfor ECP).CIntMol(featurebse) corresponds to PySCF’sMoleclass for molecule input parsing (atom coordinates, basis sets, ECP).
Differences are:
CIntonly handles C-FFI data (_atm,_bas,_env,_ecpbas), similar to PySCF’s internal arrays.CIntMolprovides molecule input parsing (atom coordinates, basis/ECP from BSE, unit conversion, ghost atoms), comparable to PySCF’sMole.build()flow.
§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:
bse(molecule builder from TOML/JSON with Basis Set Exchange support). bse: EnableCIntMol::from_tomlandCIntMol::from_jsonfor building molecules from text input. Requiresbsecrate for basis set data.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_SRCnot 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 ifqcintdoes not supports that. This will also change URL source if cargo featurebuild_from_sourcespecified.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.soandlibcint.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 withvprefix). 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
libcintas electronic integral library (corresponding to some parts ofpyscf.gto);rstsras 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.
//! H2O RHF calculation example using libcint.
//!
//! This example demonstrates a simple RHF calculation on H2O/def2-TZVP.
use libcint::prelude::*;
use rstsr::prelude::*;
pub type Tsr = Tensor<f64, DeviceCpu>;
pub type TsrView<'a> = TensorView<'a, f64, DeviceCpu>;
/// Obtain integrals (in row-major, same to PySCF but reverse of libcint)
pub fn intor_row_major(cint: &CInt, intor: &str) -> Tsr {
// use up all rayon available threads for tensor operations
let device = DeviceCpu::default();
// intor, "s1", full_shls_slice
let (out, shape) = cint.integrate_row_major(intor, None, None).into();
// row-major by transposition of col-major shape
rt::asarray((out, shape.c(), &device))
}
// This attribute line is for CI to only run this example when "bse" feature is
// enabled. API User is not supposed to add this line.
#[cfg(feature = "bse")]
fn main() {
let device = DeviceCpu::default();
// Assuming H2O/def2-TZVP data for `CInt` has been prepared
let mol_input = r#"
atom = "O; H 1 0.94; H 1 0.94 2 104.5"
basis = "def2-TZVP"
"#;
let cint_mol = CIntMol::from_toml(mol_input);
let cint = cint_mol.cint;
/* #region rhf total energy algorithms */
let atom_coords = {
let coords = cint.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.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, "int1e_kin") + intor_row_major(&cint, "int1e_nuc");
let ovlp = intor_row_major(&cint, "int1e_ovlp");
let int2e = intor_row_major(&cint, "int2e");
let nocc = 5; // hardcoded for H2O, 5 occupied orbitals
let mut dm = ovlp.zeros_like();
for _idx_iter 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 (_mo_energy, mo_coeff) = rt::linalg::eigh((&fock, &ovlp)).into();
dm = 2.0_f64 * mo_coeff.i((.., ..nocc)) % mo_coeff.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 */
assert!((eng_nuc + eng_elec - -76.05945519696209).abs() < 1e-8)
}
// ignore the following lines
#[cfg(not(feature = "bse"))]
fn main() {
println!("This example requires the 'bse' feature to be enabled.");
}Modules§
- cint
- Main
CIntwith related struct definition, and impl ofCInt::integrate. - cint_
change - Getter/setter and with-clause temporary changes to
CIntinstance. - cint_
crafter - Integral crafter for
CIntinstance (the main integral functions). - cint_
fill_ col_ major - Implementation of integral at lower API level (crafting, col-major).
- cint_
fill_ grids - Implementation of grids implementation 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
CIntinstance (compute and memory cost is not essential in this module). - cint_
result - ffi
- FFI binding and basic wrapper.
- gto
- GTO (gaussian-type atomic orbitals) values on coordinate grids.
- parse
- Generation of molecule object from input coordinates and basis sets.
- 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).