Skip to content

Interface subroutines

Memory layout

CC-aims makes strong use of MPI and Scalapack, so that a specific distribution of the input matrix- and tensor-quantities is expected.
CC-aims uses a quadratic process grid for all its computations. Hence, if a the CC-aims is launched with N MPI tasks, the first maxn2 ≤ N n2 of them will be used for computations. Furthermore, CC-aims distributes memory-intensive arrays as contiguous blocks over the relevant tasks.
The exact parameters of distribution are calculated by the CC4S_init_interface subroutine at runtime.

CC4S_init_interface

Before launching the actual interface using the CC4S_interface_driver subroutine, CC4S_init_interface has to be called. This routine requires the system- and calculation-specific parameters of the previous SCF calculation. CC4S_init_interface stores these quantities and alculates parameters related to the task-distribution of the matrix and tensor input quantities needed for the CC4S_interface_driver subroutine. CC4S_init_interface returns these parameters, which allows the calling ab-initio code to distribute these high-dimensional arrays appropriately for CC-aims before calling CC4S_interface_driver.

subroutine CC4S_init_interface(
    n_states, n_basbas, n_basis, n_atoms, n_species, n_k_points, n_kq_points,
    n_spin, max_n_basis_sp, max_n_basbas_sp, n_cells, isreal, n_low_state,
    basbas_atom, species, sp2n_basis_sp, atom2basbas_off, sp2n_basbas_sp,
    atom2basis_off, k_eigvec_loc, kpq_point_list, lb_atom, ub_atom,
    comm, fermi_energy, lbb_row, ubb_row, lbb_col, ubb_col,
    loc_bb_row, loc_bb_col, myid_col, npcol, myid, member, n_tasks, 
    interface_comm, lb_ap, ub_ap, n_ap, lbasis_col, ubasis_col,
    bbxbb_desc, bbxbasis_desc, cntxt
)

Input

argument type/shape description
n_states INTEGER Number of MO/Bloch states
n_basis INTEGER Number of AO/Bloch basis functions
n_atoms INTEGER Number of atoms in molecule/unit cell
n_species INTEGER Number of elements
n_k_points INTEGER Number of reciprocal grid points (k-points)
n_kq_points INTEGER Number of distinct differences of k-points (k-q-points)
n_spin INTEGER =2 for spin-polarized calculation, else 1
max_n_basis_sp INTEGER Maximum number of AO basis functions associated with an atom
max_n_basbas_sp INTEGER Maximum number of auxiliary basis functions/RI basis functions associated with an atom
n_cells INTEGER Number of unit cells in the supercell, used for Fourier transformation from real to reciprocal space
isreal LOGICAL if the SCF eigenvectors are of real or complex
n_low_state INTEGER Index of (energetically) lowest MO/Bloch state to be taken into consideration in CC-aims. n_low_state can assume integer values between 1 and nstates-1
n_basbas_atom INTEGER, DIMENSION(n_atoms) n_basbas_atom(i) contains the number of auxiliary/RI basis functions associated with the i-th atom
species INTEGER, DIMENSION(n_atoms) n_basbas_atom(i) contains the number of auxiliary/RI basis functions associated with the i-th atom
sp2n_basis_sp INTEGER, DIMENSION(n_species) sp2n_basis_sp(i) contains the number of AO basis functions centered around an atom of the i-th element
atom2basbas_off INTEGER, DIMENSION(n_atoms) atom2basbas_off(i) contains the index of the first auxiliary/RI basis function associated with i-th atom
sp2n_basbas_sp INTEGER, DIMENSION(n_species) sp2n_basbas_sp(i) contains the number auxiliary/RI basis functions associated with any atom of the i-th element
atom2basis_off INTEGER, DIMENSION(n_atoms) atom2basis_off(i) contains the index of the first auxiliary/RI basis function associated with the i-th atom
k_eigvec_loc INTEGER, DIMENSION(n_k_points, 2) k_eigvec_loc(i,1) contains the rank of the MPI-task on which the SCF eigenvector matrix of the i-th k-point is allocated. k_eigvec_loc(i,2) contains the local index of the i-th k-point on said task
kpq_point_list INTEGER, DIMENSION(n_k_points, n_k_points) kpq_point_list(i,j) contains the index of the k-point, which results from adding the i-th and the j-th k-points (and backfold to the Brillouin zone)
lb_atom INTEGER, DIMENSION(n_atoms) lb_atom(i) contains the index of the first basis function associated with the i-th atom
ub_atom INTEGER, DIMENSION(n_atoms) ub_atom(i) contains the index of the last basis function associated with the i-th atom
comm INTEGER MPI communicator for CC-aims to be used. (Usually MPI_COMM_WORLD)
fermi_energy DOUBLE PRECISION Chemical potential of the system

Output

argument type/shape description
lbb_row, ubb_row and lbb_col, ubb_col INTEGER Local variable. For matrix-/ tensor-quantities, which have an auxiliary/RI basis functions dimension, these pairs of integers specify the first (lbb_row or lbb_col) and last (ubb_row or ubb_col) index of that auxiliary/RI basis dimension, which is allocated at any given MPI task.Depending on the memory layout, the auxliary/RI basis dimension of an input array is distributed along the rows or columns of the quadratic process grid as described in section. If the auxiliary/RI basis function dimension is distributed along the process rows the index range of task-specific array block is described by lbb_row and ubb_row. Otherwise, lbb_col and ubb_col are relevant for the block-distribution of the array
loc_bb_row and loc_bb_col INTEGER Local variable. For matrix-/ tensor-quantities, which have a auxiliary/RI basis functions dimension, loc_bb_row (loc_bb_col) is the number of auxiliary/RI basis functions contained in the MPI task-local array block of that input array, if the auxiliary/RI basis dimension is distributed along the process row (column)
myid_col INTEGER Local variable. MPI rank denoting the column of the quadratic process grid (see interface_comm below) the present MPI task is located at
npcol INTEGER Number of columns in the quadratic process grids (= Number of rows in process grid)
myid INTEGER Local variable. MPI rank of task with respect to the quadratic process grid (see interface_comm below)
member LOGICAL Local variable. Is .true. if task is part of the quadratic process grid, .false. otherwise
n_tasks INTEGER Number of tasks in the quadratic process grid
interface_comm INTEGER MPI communicator for the quadratic process grid which will be used by CC-aims
lb_ap INTEGER Local variable. Index of first atom-pair associated with the current MPI task
ub_ap INTEGER Local variable. Index of last atom-pair associated with the current MPI task
n_ap INTEGER Local variable. Number of atom-pairs associated with the current MPI task
bbxbb_desc INTEGER, DIMENSION(9) Scalapack array descriptor used for distribution of matrix of shape n_basbas x n_basbas, e.g the auxiliary Coulomb matrix
bbxbasis_desc INTEGER, DIMENSION(9) Scalapack array descriptor used for distribution of matrix of shape n_basbas x n_basis, e.g the RI-V coefficients
cntxt INTEGER Scalapack context used in the distribution of all quantities in CC-aims

CC4S_init_interface_nopbc

This subroutine serves the same purpose as CC4S_init_interface, but is the convenient alternative for non-periodic calculations. Therefore, no k-point- or unit cell-related arguments need to be passed to CC4S_init_interface_nopbc.

subroutine CC4S_init_interface_nopbc(
    n_states, n_basbas, n_basis, n_atoms, n_species,
    n_spin, max_n_basis_sp, max_n_basbas_sp, isreal, n_low_state,
    basbas_atom, species, sp2n_basis_sp, atom2basbas_off, sp2n_basbas_sp,
    atom2basis_off, lb_atom, ub_atom,
    comm, fermi_energy, lbb_row, ubb_row, lbb_col, ubb_col,
    loc_bb_row, loc_bb_col, myid_col, npcol, myid, member, n_tasks, 
    interface_comm, lb_ap, ub_ap, n_ap, lbasis_col, ubasis_col
)

Input

argument type/shape description
n_states INTEGER Number of MO/Bloch states
n_basis INTEGER Number of AO/Bloch basis functions
n_atoms INTEGER Number of atoms in molecule/unit cell
n_species INTEGER Number of elements
n_spin INTEGER =2 for spin-polarized calculation, else 1
max_n_basis_sp INTEGER Maximum number of AO basis functions associated with an atom
max_n_basbas_sp INTEGER Maximum number of auxiliary basis functions/RI basis functions associated with an atom
isreal LOGICAL if the SCF eigenvectors are of real or complex
n_low_state INTEGER Index of (energetically) lowest MO/Bloch state to be taken into consideration in CC-aims. n_low_state can assume integer values between 1 and nstates-1
n_basbas_atom INTEGER, DIMENSION(n_atoms) n_basbas_atom(i) contains the number of auxiliary/RI basis functions associated with the i-th atom
species INTEGER, DIMENSION(n_atoms) n_basbas_atom(i) contains the number of auxiliary/RI basis functions associated with the i-th atom
sp2n_basis_sp INTEGER, DIMENSION(n_species) sp2n_basis_sp(i) contains the number of AO basis functions centered around an atom of the i-th element
atom2basbas_off INTEGER, DIMENSION(n_atoms) atom2basbas_off(i) contains the index of the first auxiliary/RI basis function associated with i-th atom
sp2n_basbas_sp INTEGER, DIMENSION(n_species) sp2n_basbas_sp(i) contains the number auxiliary/RI basis functions associated with any atom of the i-th element
atom2basis_off INTEGER, DIMENSION(n_atoms) atom2basis_off(i) contains the index of the first auxiliary/RI basis function associated with the i-th atom
lb_atom INTEGER, DIMENSION(n_atoms) lb_atom(i) contains the index of the first basis function associated with the i-th atom
ub_atom INTEGER, DIMENSION(n_atoms) ub_atom(i) contains the index of the last basis function associated with the i-th atom
comm INTEGER MPI communicator for CC-aims to be used. (Usually MPI_COMM_WORLD)
fermi_energy DOUBLE PRECISION Chemical potential of the system

Output

argument type/shape description
lbb_row, ubb_row and lbb_col, ubb_col INTEGER Local variable. For matrix-/ tensor-quantities, which have an auxiliary/RI basis functions dimension, these pairs of integers specify the first (lbb_row or lbb_col) and last (ubb_row or ubb_col) index of that auxiliary/RI basis dimension, which is allocated at any given MPI task.Depending on the memory layout, the auxliary/RI basis dimension of an input array is distributed along the rows or columns of the quadratic process grid as described in section. If the auxiliary/RI basis function dimension is distributed along the process rows the index range of task-specific array block is described by lbb_row and ubb_row. Otherwise, lbb_col and ubb_col are relevant for the block-distribution of the array
loc_bb_row and loc_bb_col INTEGER Local variable. For matrix-/ tensor-quantities, which have a auxiliary/RI basis functions dimension, loc_bb_row (loc_bb_col) is the number of auxiliary/RI basis functions contained in the MPI task-local array block of that input array, if the auxiliary/RI basis dimension is distributed along the process row (column)
myid_col INTEGER Local variable. MPI rank denoting the column of the quadratic process grid (see interface_comm below) the present MPI task is located at
npcol INTEGER Number of columns in the quadratic process grids (= Number of rows in process grid)
myid INTEGER Local variable. MPI rank of task with respect to the quadratic process grid (see interface_comm below)
member LOGICAL Local variable. Is .true. if task is part of the quadratic process grid, .false. otherwise
n_tasks INTEGER Number of tasks in the quadratic process grid
interface_comm INTEGER MPI communicator for the quadratic process grid which will be used by CC-aims
lb_ap INTEGER Local variable. Index of first atom-pair associated with the current MPI task
ub_ap INTEGER Local variable. Index of last atom-pair associated with the current MPI task
n_ap INTEGER Local variable. Number of atom-pairs associated with the current MPI task
bbxbb_desc INTEGER, DIMENSION(9) Scalapack array descriptor used for distribution of matrix of shape n_basbas x n_basbas, e.g the auxiliary Coulomb matrix
bbxbasis_desc INTEGER, DIMENSION(9) Scalapack array descriptor used for distribution of matrix ofshape n_basbas x n_basis, e.g the RI-V coefficients
cntxt INTEGER Scalapack context used in the distribution of all quantities in CC-aims

CC4S_interface_driver

After setting up CC-aims via a call to CC4S_init_interface and (re-)distributing the relevant matrix-/tensor-quantities, CC4S_interface_driver can be called. This subroutine handles the construction of the Coulomb vertex, the output of the remaining input quantities of CC4S and the generation of the yaml files.

subroutine CC4S_interface_driver(
    coulmat, lvl_trico, v_trico, 
    scf_eigvecs_cmplx, scf_eigvecs_real, 
    scf_eigvals, comm, debug, singval_thr, 
    singval_num, scf_fock_real, scf_fock_cmplx
)

Input

argument type/shape description
coulmat DOUBLE COMPLEX, DIMENSION(lbb_row:ubb_row, lbb_col:ubb_col, n_kq_points) Task-local block of auxiliary Coulomb matrix. The global array of the auxiliary Coulomb matrix is of shape (n_basbas, n_basbas, n_kq_points) and is distributed both along the rows and columns of the quadratic process grid used by CC-aims.
lvl_trico DOUBLE COMPLEX, DIMENSION(max_n_basis_sp, max_n_basis_sp, max_n_basbas_sp, n_k_points, lb_ap:ub_ap) Task-local block of tensor of RI expansion coefficients if using RI-LVL scheme. The global tensor of RI-expansion coefficients is of shape (max_n_basis_sp, max_n_basis_sp, max_n_basbas_sp, n_k_points, N_ap), where N_ap is the total number of atom pairs of the system
v_trico DOUBLE COMPLEX, DIMENSION(lbb_row:ubb_row,lbasis_col:ubasis_col, n_basis, n_k_points, n_k_points) Task-local tensor of RI-expansion coefficients if using RI-V scheme. The global tensor of RI-V expansion coefficients is of shape (n_basbas, n_basis, n_basis, n_k_points, n_k_points). Depending on the RI scheme, which should be used, either lvl_trico or v_trico needs to be allocated as described and the other one needs to be allocated to have size 1 (e.g if RI-LVL is requested v_trico needs to be of shape (1,1,1,1,1))
scf_eigvecs_cmplx DOUBLE COMPLEX, DIMENSION(n_basis, n_states, n_spin, n_k_points_task) HF- or KS-eigenvector matrix if isreal == .false. n_k_points_task refers to the number of k-points assigned to the local MPI-task by k_eigvec_loc (see. CC4S_init_interface)
scf_eigvecs_real DOUBLE PRECISION, DIMENSION(n_basis, n_states, n_spin, n_kpoints_task) HF- or KS-eigenvector matrix if isreal == .true. n_k_points_task refers to the number of k-points assigned to the local MPI-task by k_eigvec_loc (see. CC4S_init_interface)
scf_eigvals DOUBLE PRECISION, DIMENSION(n_states, n_spin, n_k_points) HF- or KS-eigenvalues
comm INTEGER Global MPI Communicator, which should be identical to comm specified in the CC4S_init_interface-subroutine (usually MPI_COMM_WORLD)
debug LOGICAL if debug == .true., debug mode is enabled. Additional debug information and files will be generated
singval_thr DOUBLE PRECISION After the generation of the Coulomb vertex, a principal component analysis (PCA) is performed to reduce the size of the Coulomb vertex. Singular values below singval_thr will be ignored during the PCA
singval_num INTEGER Number of biggest singular values to be kept during PCA
scf_fock_real DOUBLE PRECISION, DIMENSION(n_states, n_states, n_spin, n_k_points_local) n_k_points_local being the number of k-points assigned to the current MPI task via k_eigvec_loc. Fock matrix if isreal == .true.. This is an optional input, which is only necessary if a KS calculation was performed
scf_fock_cmplx DOUBLE COMPLEX, DIMENSION(n_states, n_states, n_spin, n_k_points_local) n_k_points_local being the number of k-points assigned to the current MPI task via k_eigvec_loc. Fock matrix if isreal == .false.. This is an optional input, which is only necessary if a KS calculation was performed