Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fix: Update function calls in pyabacus to align with new function signature in hpsi_func #5176

Merged
merged 3 commits into from
Sep 26, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
13 changes: 8 additions & 5 deletions python/pyabacus/src/py_diago_dav_subspace.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -110,11 +110,14 @@ class PyDiagoDavSubspace
bool scf_type,
hsolver::diag_comm_info comm_info
) {
auto hpsi_func = [mm_op] (std::complex<double> *hpsi_out,
std::complex<double> *psi_in, const int nband_in,
const int nbasis_in, const int band_index1,
const int band_index2)
{
auto hpsi_func = [mm_op] (
std::complex<double> *psi_in,
std::complex<double> *hpsi_out,
const int nband_in,
const int nbasis_in,
const int band_index1,
const int band_index2
) {
// Note: numpy's py::array_t is row-major, but
// our raw pointer-array is column-major
py::array_t<std::complex<double>, py::array::f_style> psi({nbasis_in, band_index2 - band_index1 + 1});
Expand Down
4 changes: 2 additions & 2 deletions python/pyabacus/src/py_diago_david.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -109,8 +109,8 @@ class PyDiagoDavid
hsolver::diag_comm_info comm_info
) {
auto hpsi_func = [mm_op] (
std::complex<double> *hpsi_out,
std::complex<double> *psi_in,
std::complex<double> *psi_in,
std::complex<double> *hpsi_out,
const int nband_in,
const int nbasis_in,
const int band_index1,
Expand Down
30 changes: 16 additions & 14 deletions python/pyabacus/src/pyabacus/hsolver/_hsolver.py
Original file line number Diff line number Diff line change
Expand Up @@ -16,7 +16,7 @@ def rank(self) -> int: ...
def nproc(self) -> int: ...

def dav_subspace(
mm_op: Callable[[NDArray[np.complex128]], NDArray[np.complex128]],
mvv_op: Callable[[NDArray[np.complex128]], NDArray[np.complex128]],
init_v: NDArray[np.complex128],
dim: int,
num_eigs: int,
Expand All @@ -32,9 +32,10 @@ def dav_subspace(

Parameters
----------
mm_op : Callable[[NDArray[np.complex128]], NDArray[np.complex128]],
The operator to be diagonalized, which is a function that takes a matrix as input
and returns a matrix mv_op(X) = H * X as output.
mvv_op : Callable[[NDArray[np.complex128]], NDArray[np.complex128]],
The operator to be diagonalized, which is a function that takes a set of
vectors X = [x1, ..., xN] as input and returns a matrix(vector block)
mvv_op(X) = H * X ([Hx1, ..., HxN]) as output.
init_v : NDArray[np.complex128]
The initial guess for the eigenvectors.
dim : int
Expand Down Expand Up @@ -68,8 +69,8 @@ def dav_subspace(
v : NDArray[np.complex128]
The eigenvectors corresponding to the eigenvalues.
"""
if not callable(mm_op):
raise TypeError("mm_op must be a callable object.")
if not callable(mvv_op):
raise TypeError("mvv_op must be a callable object.")

if is_occupied is None:
is_occupied = [True] * num_eigs
Expand All @@ -86,7 +87,7 @@ def dav_subspace(
assert dav_ndim * num_eigs < dim * comm_info.nproc, "dav_ndim * num_eigs must be less than dim * comm_info.nproc."

_ = _diago_obj_dav_subspace.diag(
mm_op,
mvv_op,
pre_condition,
dav_ndim,
tol,
Expand All @@ -103,7 +104,7 @@ def dav_subspace(
return e, v

def davidson(
mm_op: Callable[[NDArray[np.complex128]], NDArray[np.complex128]],
mvv_op: Callable[[NDArray[np.complex128]], NDArray[np.complex128]],
init_v: NDArray[np.complex128],
dim: int,
num_eigs: int,
Expand All @@ -119,9 +120,10 @@ def davidson(

Parameters
----------
mm_op : Callable[[NDArray[np.complex128]], NDArray[np.complex128]],
The operator to be diagonalized, which is a function that takes a matrix as input
and returns a matrix mv_op(X) = H * X as output.
mvv_op : Callable[[NDArray[np.complex128]], NDArray[np.complex128]],
The operator to be diagonalized, which is a function that takes a set of
vectors X = [x1, ..., xN] as input and returns a matrix(vector block)
mvv_op(X) = H * X ([Hx1, ..., HxN]) as output.
init_v : NDArray[np.complex128]
The initial guess for the eigenvectors.
dim : int
Expand All @@ -146,8 +148,8 @@ def davidson(
v : NDArray[np.complex128]
The eigenvectors corresponding to the eigenvalues.
"""
if not callable(mm_op):
raise TypeError("mm_op must be a callable object.")
if not callable(mvv_op):
raise TypeError("mvv_op must be a callable object.")

if init_v.ndim != 1 or init_v.dtype != np.complex128:
init_v = init_v.flatten().astype(np.complex128, order='C')
Expand All @@ -159,7 +161,7 @@ def davidson(
comm_info = hsolver.diag_comm_info(0, 1)

_ = _diago_obj_dav_subspace.diag(
mm_op,
mvv_op,
pre_condition,
dav_ndim,
tol,
Expand Down
Loading