Takes as argument an array of group index for each atom, and create a list of atms in each group using the sparse matrix format (saved as Yale format). This is used by cell lists, polarization group etc.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=ip), | intent(in) | :: | a2g(:) |
Index of polarization group for each MM atom |
||
type(yale_sparse), | intent(out) | :: | g2a |
Indices of atoms included in each polarization group; Atom indeces for the n-th group are found at pg2mm%ci(pg2mm%ri(n):pg2mm%ri(n+1)-1) |
||
integer(kind=ip), | intent(in), | optional | :: | ng_in |
Number of groups if it is not provided in input it is assumed that the number of group equals the largest group index, that is no empty groups are present after the one with the largest index. |
subroutine reverse_grp_tab(a2g, g2a, ng_in)
use mod_memory, only: mallocate, mfree
!! Takes as argument an array of group index for each
!! atom, and create a list of atms in each group using the
!! sparse matrix format (saved as Yale format).
!! This is used by cell lists, polarization group etc.
implicit none
integer(ip), intent(in) :: a2g(:)
!! Index of polarization group for each MM atom
type(yale_sparse), intent(out) :: g2a
!! Indices of atoms included in each polarization group;
!! Atom indeces for the n-th group are found at
!! pg2mm%ci(pg2mm%ri(n):pg2mm%ri(n+1)-1)
integer(ip), intent(in), optional :: ng_in
!! Number of groups if it is not provided in input it is
!! assumed that the number of group equals the largest group
!! index, that is no empty groups are present after the one
!! with the largest index.
integer(ip) :: i, j, na, ng, ig
integer(ip), allocatable :: uc_data(:, :), g_dim(:)
na = size(a2g)
if(present(ng_in)) then
ng = ng_in
else
ng = maxval(a2g)
end if
! Find largest group
call mallocate('reverse_grp_tab [g_dim]', ng, g_dim)
g_dim = 0
do i=1, na
g_dim(a2g(i)) = g_dim(a2g(i)) + 1
end do
! Struct for uncompressed data
call mallocate('reverse_grp_tab [uc_data]', maxval(g_dim), ng, uc_data)
! First invert in an uncompressed structure
uc_data = 0
g_dim = 0
do i=1, na
ig = a2g(i)
g_dim(ig) = g_dim(ig) + 1
uc_data(g_dim(ig),ig) = i
end do
! Compress the list
!g2a%n = ng
!if(.not. allocated(g2a%ri)) &
! call mallocate('reverse_grp_tab [ri]', ng+1, g2a%ri)
!if(.not. allocated(g2a%ci)) &
! call mallocate('reverse_grp_tab [ci]', na, g2a%ci)
!g2a%ri(1) = 1
!do i=1, ng
! g2a%ri(i+1) = g2a%ri(i) + g_dim(i) - 1
! g2a%ci(g2a%ri(i):g2a%ri(i+1)-1) = uc_data(1:g_dim(i)-1,i)
!end do
call compress_list(ng, uc_data, g_dim, g2a)
! Free temporary mem
call mfree('reverse_grp_tab [uc_data]', uc_data)
call mfree('reverse_grp_tab [g_dim]', g_dim)
end subroutine reverse_grp_tab