Octopus
isdf_serial_oct_m Module Reference

Serial prototype for benchmarking and validating ISDF implementation. More...

Detailed Description

Serial prototype for benchmarking and validating ISDF implementation.

Functions/Subroutines

subroutine, public isdf_serial_interpolation_vectors (isdf, namespace, mesh, st, indices, phi_mu, P_r_mu, isdf_vectors)
 Construct interpolative separable density fitting (ISDF) vectors and other intermediate quantities relevant for density-fitting adaptively- compressed exchange. More...
 
subroutine collate_batches_get_state (mesh, st, max_state, psi)
 Loop over states per block, which makes applying the maximum state limit much simpler Use this to compare to the output of collate_batches Not distributed in states or domain. More...
 
subroutine sample_phi_at_centroids (phi_r, indices, phi_mu)
 Sample KS states at centroid points. More...
 
subroutine construct_zct (zct)
 Construct the product of Z and C matrices from the element-wise product of the quasi-density matrix. More...
 
subroutine construct_cct (indices, zct, cct)
 Construct the product \(CC^T\) from \(ZC^T\) by masking the first dimension of \(ZC^T\). More...
 
subroutine construct_density_matrix_packed (phi, phi_mu, P_r_mu)
 @ brief Construct the density matrix with shape (np, n_int). Denoted packed, because it expects phi in state-major format. More...
 
subroutine construct_density_matrix_with_occ_packed (st, phi, phi_mu, P_r_mu)
 
subroutine, public quantify_error_and_visualise (isdf, namespace, st, space, mesh, ions, indices, isdf_vectors, output_cubes)
 Wrapper for quantifying the error in the expansion of the product basis. More...
 
subroutine approximate_pair_products (psi_mu, zeta, product_basis)
 Construct a set of approximate pair products using the ISDF interpolation vectors. More...
 
subroutine error_in_product_basis (mesh, product_basis, approx_product_basis, error, mean_error)
 Quantify the error in the product basis expansion. More...
 
subroutine generate_product_state_cubes (namespace, space, mesh, ions, file_prefix, data, limits)
 Helper function to output a set of pair product states. More...
 
subroutine, public isdf_serial_ace_compute_potentials (exxop, namespace, space, mesh, st, Vx_on_st, kpoints)
 ISDF wrapper computing interpolation points and vectors, which are used to build the potential \(W_j=\left(V_x\left[\left\{\varphi_i\right\}\right] \psi_j\right)(\mathbf{r})\) used in adaptively-compressed exchange. More...
 
subroutine isdf_serial_ace_w_unpacked (namespace, P_r_mu, isdf_vectors, psi_mu, poisson_solver, cam, st, W_ace)
 Compute the action of the exchange potential on KS states for adaptively-compressed exchange. More...
 
subroutine isdf_serial_ace_batch_w (isdf, st, W_ace, Vx_on_st)
 Put the bare array representation of W into a batch. More...
 

Variables

integer, parameter ik = 1
 Hard-coded for Gamma-point, spin-unpolarised calculations. More...
 

Function/Subroutine Documentation

◆ isdf_serial_interpolation_vectors()

subroutine, public isdf_serial_oct_m::isdf_serial_interpolation_vectors ( type(isdf_options_t), intent(in)  isdf,
type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(in)  mesh,
type(states_elec_t), intent(in)  st,
integer(int64), dimension(:), intent(in), contiguous  indices,
real(real64), dimension(:, :), intent(out), allocatable  phi_mu,
real(real64), dimension(:, :), intent(out), allocatable  P_r_mu,
real(real64), dimension(:, :), intent(out), allocatable  isdf_vectors 
)

Construct interpolative separable density fitting (ISDF) vectors and other intermediate quantities relevant for density-fitting adaptively- compressed exchange.

Parameters
[in]stKS states
[in]indicesInterpolation point indices
[out]phi_muBatched array of wave functions
[out]p_r_muQuasi-density matrix,
[out]isdf_vectorsInterpolation vectors

Definition at line 169 of file isdf_serial.F90.

◆ collate_batches_get_state()

subroutine isdf_serial_oct_m::collate_batches_get_state ( class(mesh_t), intent(in)  mesh,
type(states_elec_t), intent(in)  st,
integer, intent(in)  max_state,
real(real64), dimension(:, :), intent(out), allocatable  psi 
)
private

Loop over states per block, which makes applying the maximum state limit much simpler Use this to compare to the output of collate_batches Not distributed in states or domain.

Parameters
[in]stStates
[in]max_statemax state in ISDF expansion
[out]psiPacked states, \(\varphi_i(\mathbf{r})\)

Definition at line 290 of file isdf_serial.F90.

◆ sample_phi_at_centroids()

subroutine isdf_serial_oct_m::sample_phi_at_centroids ( real(real64), dimension(:, :), intent(in), contiguous  phi_r,
integer(int64), dimension(:), intent(in), contiguous  indices,
real(real64), dimension(:, :), intent(out), contiguous  phi_mu 
)
private

Sample KS states at centroid points.

Parameters
[in]phi_rPhi at all grid points
[in]indicescentroid to grid map
[out]phi_muPhi at centroid points

Definition at line 323 of file isdf_serial.F90.

◆ construct_zct()

subroutine isdf_serial_oct_m::construct_zct ( real(real64), dimension(:, :), intent(inout), contiguous  zct)
private

Construct the product of Z and C matrices from the element-wise product of the quasi-density matrix.

In this instance, the routine only accepts one set of Kohn-sham states:

\[ (ZC^T)_{\mathbf{r}, \mathbf{r}_{\mu}} = P^{\varphi}(\mathbf{r}, \mathbf{r}_{\mu}) \circ P^{\varphi}(\mathbf{r}, \mathbf{r}_{\mu}) \quad \text{for all elements} \; \mathbf{r} \; \text{and} \; \mathbf{r}_{\mu} \]

Parameters
[in,out]zctIn: Quasi-density matrix

Definition at line 364 of file isdf_serial.F90.

◆ construct_cct()

subroutine isdf_serial_oct_m::construct_cct ( integer(int64), dimension(:), intent(in), contiguous  indices,
real(real64), dimension(:, :), intent(in), contiguous  zct,
real(real64), dimension(:, :), intent(out), contiguous  cct 
)
private

Construct the product \(CC^T\) from \(ZC^T\) by masking the first dimension of \(ZC^T\).

Parameters
[in]indicesInterpolation point indices on global grid
[in]zctProduct ZC^T
[out]cctProduct CC^T

Definition at line 390 of file isdf_serial.F90.

◆ construct_density_matrix_packed()

subroutine isdf_serial_oct_m::construct_density_matrix_packed ( real(real64), dimension(:, :), intent(in), contiguous  phi,
real(real64), dimension(:, :), intent(in), contiguous  phi_mu,
real(real64), dimension(:, :), intent(out), contiguous  P_r_mu 
)
private

@ brief Construct the density matrix with shape (np, n_int). Denoted packed, because it expects phi in state-major format.

Construct the density matrix:

\[ P^{\varphi}\left(\mathbf{r}, \mathbf{r}_{\mu}\right) = \sum_{i=1}^m \varphi_i(\mathbf{r}) \varphi_i\left(\mathbf{r}_{\mu}\right). \]

Parameters
[in]phiA set of states defined on real-space grid
[in]phi_muKS states, only defined at interpolation points
[out]p_r_muDensity matrix of shape (np, n_int)

Definition at line 430 of file isdf_serial.F90.

◆ construct_density_matrix_with_occ_packed()

subroutine isdf_serial_oct_m::construct_density_matrix_with_occ_packed ( type(states_elec_t), intent(in)  st,
real(real64), dimension(:, :), intent(in)  phi,
real(real64), dimension(:, :), intent(in)  phi_mu,
real(real64), dimension(:, :), intent(out)  P_r_mu 
)
private
Parameters
[in]stCurrent Kohn-Sham states
[in]phiA set of states defined on real-space grid
[in]phi_muKS states, only defined at interpolation points
[out]p_r_muDensity matrix of shape (np, n_int)

Definition at line 462 of file isdf_serial.F90.

◆ quantify_error_and_visualise()

subroutine, public isdf_serial_oct_m::quantify_error_and_visualise ( type(isdf_options_t), intent(in)  isdf,
type(namespace_t), intent(in)  namespace,
type(states_elec_t), intent(in)  st,
class(space_t), intent(in)  space,
class(mesh_t), intent(in)  mesh,
class(ions_t), intent(in), pointer  ions,
integer(int64), dimension(:), intent(in), contiguous  indices,
real(real64), dimension(:, :), intent(inout), allocatable  isdf_vectors,
logical, intent(in)  output_cubes 
)

Wrapper for quantifying the error in the expansion of the product basis.

Parameters
[in]indicesInterpolation indices
[in,out]isdf_vectorsISDF vectors

Definition at line 504 of file isdf_serial.F90.

◆ approximate_pair_products()

subroutine isdf_serial_oct_m::approximate_pair_products ( real(real64), dimension(:, :), intent(in), contiguous  psi_mu,
real(real64), dimension(:, :), intent(in), contiguous  zeta,
real(real64), dimension(:, :), intent(out), contiguous  product_basis 
)
private

Construct a set of approximate pair products using the ISDF interpolation vectors.

\[ \varphi_i(\mathbf{r}) \varphi_j(\mathbf{r}) \approx \sum_{\mu =1} \varphi_i(\mathbf{r}_\mu) \varphi_j(\mathbf{r}_\mu) \zeta_\mu(\mathbf{r}) \]

which is formulated in terms of a Khatri-Rao product to form, \(\varphi_i(\mathbf{r}_\mu) \varphi_j(\mathbf{r}_\mu)\), followed by the matrix-matrix product, \( \left[ \varphi_{ij, \mu} \right] \left[\zeta_{\mathbf{r},\mu} \right]^T \) to contract over the interpolation index.

Parameters
[in]psi_muStates sampled at interpolation points
[in]zetaInterpolation vectors
[out]product_basisApprox product basis

Definition at line 597 of file isdf_serial.F90.

◆ error_in_product_basis()

subroutine isdf_serial_oct_m::error_in_product_basis ( class(mesh_t), intent(in)  mesh,
real(real64), dimension(:, :), intent(in), contiguous  product_basis,
real(real64), dimension(:, :), intent(in), contiguous  approx_product_basis,
real(real64), dimension(:), intent(out), contiguous  error,
real(real64), intent(out)  mean_error 
)
private

Quantify the error in the product basis expansion.

\[ \left[ \frac{dV}{N_p}\sum_{\mathbf{r}} \left[ \varphi_i(\mathbf{r})\varphi_j(\mathbf{r}) - \tilde{\varphi}_i(\mathbf{r})\tilde{\varphi}_j(\mathbf{r}) \right]^2 \right]^{1/2} \]

Parameters
[in]meshMesh instance
[in]product_basisProduct basis with shape (mn_states, np)
[in]approx_product_basisApprox product basis with shape (mn_states, np)

Definition at line 634 of file isdf_serial.F90.

◆ generate_product_state_cubes()

subroutine isdf_serial_oct_m::generate_product_state_cubes ( type(namespace_t), intent(in)  namespace,
class(space_t), intent(in)  space,
class(mesh_t), intent(in)  mesh,
class(ions_t), intent(in), pointer  ions,
character(len=*), intent(in)  file_prefix,
real(real64), dimension(:, :), intent(in), contiguous  data,
integer, dimension(2), intent(in), optional  limits 
)
private

Helper function to output a set of pair product states.

Parameters
[in]limitsOrdered inner loop to outer loop

Definition at line 683 of file isdf_serial.F90.

◆ isdf_serial_ace_compute_potentials()

subroutine, public isdf_serial_oct_m::isdf_serial_ace_compute_potentials ( type(exchange_operator_t), intent(in)  exxop,
type(namespace_t), intent(in)  namespace,
class(space_t), intent(in)  space,
class(mesh_t), intent(in)  mesh,
type(states_elec_t), intent(in)  st,
type(states_elec_t), intent(out)  Vx_on_st,
type(kpoints_t), intent(in)  kpoints 
)

ISDF wrapper computing interpolation points and vectors, which are used to build the potential \(W_j=\left(V_x\left[\left\{\varphi_i\right\}\right] \psi_j\right)(\mathbf{r})\) used in adaptively-compressed exchange.

The routine performs the steps:

  1. Compute ISDF vectors, and also return \(\psi_i(\mathbf{r}_\mu)\) and \(P_{\mu, \mathbf{r}\), required for the construction of \(W_j\).
  2. Compute \(W_j\) as a bare array, then put it into a batched container
Note
If one wishes to also compute \(M_{ij}\) in the ISDF basis:
  • P_r_mu can be masked using the interpolation point indices defined in centroids obj
  • A new routine is required to compute \(M_{\mu,\nu}\) in the ISDF basis
  • A new routine is required to rotate M from the ISDF basis to the KS basis.
Parameters
[in]exxopExchange instance, containing Poisson solver settings,
[in]spaceOnly included so API is consistent
[in]meshReal-space grid
[in]stCurrent Kohn-Sham states
[in]kpointsOnly included so API is consistent
[out]vx_on_stV_X operating on all states (W in ACE)

Definition at line 736 of file isdf_serial.F90.

◆ isdf_serial_ace_w_unpacked()

subroutine isdf_serial_oct_m::isdf_serial_ace_w_unpacked ( type(namespace_t), intent(in)  namespace,
real(real64), dimension(:, :), intent(in), contiguous  P_r_mu,
real(real64), dimension(:, :), intent(in), contiguous  isdf_vectors,
real(real64), dimension(:, :), intent(in), contiguous  psi_mu,
type(poisson_t), intent(in)  poisson_solver,
type(xc_cam_t), intent(in)  cam,
type(states_elec_t), intent(in)  st,
real(real64), dimension(:, :), intent(out), contiguous  W_ace 
)
private

Compute the action of the exchange potential on KS states for adaptively-compressed exchange.

\[ W_j(\mathbf{r}) = - \sum_{\mu} f P_{\mathbf{r}, \mu} V_{\mathbf{r}, \mu} \phi_j(\mathbf{r}_\mu). \]

Returns W_ace with shape(np, nst), as this is optimal for putting into a batch.

Parameters
[in]p_r_muQuasi-density matrix, shape(np, n_int)
[in]isdf_vectorsISDF vectors \( \{ \zeta \} \), with shape(np, n_int)
[in]psi_mu\( \phi_j(\mathbf{r}_\mu) \), with shape(nocc, n_int)
[in]poisson_solverPoisson solver settings
[in]camHybrid CAM parameters
[in]stCurrent Kohn-Sham states
[out]w_aceACE W-matrix, expected with shape(np, nocc)

Definition at line 786 of file isdf_serial.F90.

◆ isdf_serial_ace_batch_w()

subroutine isdf_serial_oct_m::isdf_serial_ace_batch_w ( type(isdf_options_t), intent(in)  isdf,
type(states_elec_t), intent(in)  st,
real(real64), dimension(:, :), intent(in), contiguous  W_ace,
type(states_elec_t), intent(out)  Vx_on_st 
)
private

Put the bare array representation of W into a batch.

REQUIRES that W_ace has shape (np, nst).

Parameters
[in]stKS states
[in]w_aceACE W-matrix, returned with shape(np, nst)
[out]vx_on_stKS states having had V_X act on them

Definition at line 860 of file isdf_serial.F90.

Variable Documentation

◆ ik

integer, parameter isdf_serial_oct_m::ik = 1
private

Hard-coded for Gamma-point, spin-unpolarised calculations.

Definition at line 162 of file isdf_serial.F90.