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.

Data Types

type  isdf_options_t
 ISDF input options. More...
 

Functions/Subroutines

subroutine isdf_options_init (this, namespace, st)
 Initialise isdf_inp_options_t instance. More...
 
subroutine, public serial_interpolative_separable_density_fitting_vectors (namespace, mesh, st, int_indices, phi, isdf_vectors)
 Compute a set of ISDF interpolation vectors in serial, for code validation. More...
 
subroutine isdf_construct_interpolation_vectors_packed (namespace, mesh, phi, indices, isdf_vectors)
 Compute a set of ISDF interpolation vectors, where intermediate quantities such as phi are constructed such that the first index runs over states, and the second runs over grid points. 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 collate_batches (np, st, max_state, psi)
 Put batches into a single 2D array. 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_all_centroids_packed (phi_mu, P_mu_nu)
 @ brief Construct the density matrix with shape (n_int, n_int). Denoted packed, because it expects phi in state-major format. More...
 
subroutine, public quantify_error_and_visualise (namespace, st, space, mesh, ions, phi, 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 output_matrix (namespace, fname, matrix)
 Helper routine to output a 2D matrix. 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...
 
integer function highest_occupied_index (st, ik_index)
 Return the index of highest occupied Kohn-Sham state for k-point ik. More...
 

Variables

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

Function/Subroutine Documentation

◆ isdf_options_init()

subroutine isdf_serial_oct_m::isdf_options_init ( class(isdf_options_t), intent(out)  this,
type(namespace_t), intent(in)  namespace,
type(states_elec_t), intent(in)  st 
)
private

Initialise isdf_inp_options_t instance.

Definition at line 167 of file isdf_serial.F90.

◆ serial_interpolative_separable_density_fitting_vectors()

subroutine, public isdf_serial_oct_m::serial_interpolative_separable_density_fitting_vectors ( type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(in)  mesh,
type(states_elec_t), intent(in)  st,
integer(int64), dimension(:), intent(in)  int_indices,
real(real64), dimension(:, :), intent(out), allocatable  phi,
real(real64), dimension(:, :), intent(out), allocatable  isdf_vectors 
)

Compute a set of ISDF interpolation vectors in serial, for code validation.

  • No domain or state parallelisation
    Parameters
    [in]stKS states
    [in]int_indicesGrid indices of interpolation points
    [out]phiExpanded batch array, \(\varphi_i(\mathbf{r}_\mu)\),
    [out]isdf_vectorsInterpolation vectors

Definition at line 198 of file isdf_serial.F90.

◆ isdf_construct_interpolation_vectors_packed()

subroutine isdf_serial_oct_m::isdf_construct_interpolation_vectors_packed ( type(namespace_t), intent(in)  namespace,
class(mesh_t), intent(in)  mesh,
real(real64), dimension(:, :), intent(in)  phi,
integer(int64), dimension(:), intent(in)  indices,
real(real64), dimension(:, :), intent(out), allocatable  isdf_vectors 
)
private

Compute a set of ISDF interpolation vectors, where intermediate quantities such as phi are constructed such that the first index runs over states, and the second runs over grid points.

Parameters
[in]phiExpanded batches
[in]indicesInterpolation indices, with a max element value, np_global
[out]isdf_vectorsInterpolation vectors

Definition at line 269 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 364 of file isdf_serial.F90.

◆ collate_batches()

subroutine isdf_serial_oct_m::collate_batches ( integer, intent(in)  np,
type(states_elec_t), intent(in)  st,
integer, intent(in)  max_state,
real(real64), dimension(:, :), intent(out), allocatable  psi 
)
private

Put batches into a single 2D array.

Parameters
[in]npNumber of grid points in domain
[in]stStates
[in]max_stateMax number of states to include
[out]psiPacked states, \(\varphi_i(\mathbf{r})\)

Definition at line 402 of file isdf_serial.F90.

◆ sample_phi_at_centroids()

subroutine isdf_serial_oct_m::sample_phi_at_centroids ( real(real64), dimension(:, :), intent(in)  phi_r,
integer(int64), dimension(:), intent(in)  indices,
real(real64), dimension(:, :), intent(out)  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 452 of file isdf_serial.F90.

◆ construct_zct()

subroutine isdf_serial_oct_m::construct_zct ( real(real64), dimension(:, :), intent(inout)  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 493 of file isdf_serial.F90.

◆ construct_cct()

subroutine isdf_serial_oct_m::construct_cct ( integer(int64), dimension(:), intent(in)  indices,
real(real64), dimension(:, :), intent(in)  zct,
real(real64), dimension(:, :), intent(out)  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 "@CMAKE_C_COMPILER@"^T

Definition at line 519 of file isdf_serial.F90.

◆ construct_density_matrix_packed()

subroutine isdf_serial_oct_m::construct_density_matrix_packed ( real(real64), dimension(:, :), intent(in)  phi,
real(real64), dimension(:, :), intent(in)  phi_mu,
real(real64), dimension(:, :), intent(out)  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 559 of file isdf_serial.F90.

◆ construct_density_matrix_all_centroids_packed()

subroutine isdf_serial_oct_m::construct_density_matrix_all_centroids_packed ( real(real64), dimension(:, :), intent(in)  phi_mu,
real(real64), dimension(:, :), intent(out)  P_mu_nu 
)
private

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

Construct the density matrix:

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

Parameters
[in]phi_muKS states, only defined at interpolation points
[out]p_mu_nuDensity matrix of shape (n_int, n_int)

Definition at line 599 of file isdf_serial.F90.

◆ quantify_error_and_visualise()

subroutine, public isdf_serial_oct_m::quantify_error_and_visualise ( 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,
real(real64), dimension(:, :), intent(inout), allocatable  phi,
integer(int64), dimension(:), intent(in)  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,out]phiKS states
[in]indicesInterpolation indices
[in,out]isdf_vectorsISDF vectors

Definition at line 624 of file isdf_serial.F90.

◆ approximate_pair_products()

subroutine isdf_serial_oct_m::approximate_pair_products ( real(real64), dimension(:, :), intent(in)  psi_mu,
real(real64), dimension(:, :), intent(in)  zeta,
real(real64), dimension(:, :), intent(out)  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 725 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)  product_basis,
real(real64), dimension(:, :), intent(in)  approx_product_basis,
real(real64), dimension(:), intent(out)  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 762 of file isdf_serial.F90.

◆ output_matrix()

subroutine isdf_serial_oct_m::output_matrix ( type(namespace_t), intent(in)  namespace,
character(len=*), intent(in)  fname,
real(real64), dimension(:, :), intent(in)  matrix 
)
private

Helper routine to output a 2D matrix.

Definition at line 811 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)  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 841 of file isdf_serial.F90.

◆ highest_occupied_index()

integer function isdf_serial_oct_m::highest_occupied_index ( type(states_elec_t), intent(in)  st,
integer, intent(in)  ik_index 
)
private

Return the index of highest occupied Kohn-Sham state for k-point ik.

This is only valid for smear == SMEAR_SEMICONDUCTOR. This routines also assumes a contiguous distribution of occupied states. It will not, for example return the correct index for the input: ``` Occupations 0 | 1 %

Parameters
[in]ik_indexk-point index
Returns
Index of highest occupied KS state

Definition at line 890 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 150 of file isdf_serial.F90.