reverse_grp_tab Subroutine

public subroutine reverse_grp_tab(a2g, g2a, ng_in)

Uses

  • proc~~reverse_grp_tab~~UsesGraph proc~reverse_grp_tab reverse_grp_tab module~mod_memory mod_memory proc~reverse_grp_tab->module~mod_memory module~mod_io mod_io module~mod_memory->module~mod_io module~mod_constants mod_constants module~mod_memory->module~mod_constants iso_c_binding iso_c_binding module~mod_memory->iso_c_binding module~mod_io->module~mod_constants module~mod_constants->iso_c_binding

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.

Arguments

Type IntentOptional 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.


Calls

proc~~reverse_grp_tab~~CallsGraph proc~reverse_grp_tab reverse_grp_tab interface~mallocate mallocate proc~reverse_grp_tab->interface~mallocate proc~compress_list compress_list proc~reverse_grp_tab->proc~compress_list interface~mfree mfree proc~reverse_grp_tab->interface~mfree proc~r_alloc1 r_alloc1 interface~mallocate->proc~r_alloc1 proc~r_alloc3 r_alloc3 interface~mallocate->proc~r_alloc3 proc~i_alloc2 i_alloc2 interface~mallocate->proc~i_alloc2 proc~i_alloc1 i_alloc1 interface~mallocate->proc~i_alloc1 proc~r_alloc2 r_alloc2 interface~mallocate->proc~r_alloc2 proc~i_alloc3 i_alloc3 interface~mallocate->proc~i_alloc3 proc~l_alloc1 l_alloc1 interface~mallocate->proc~l_alloc1 proc~l_alloc2 l_alloc2 interface~mallocate->proc~l_alloc2 proc~compress_list->interface~mallocate proc~compress_list->interface~mfree proc~i_free1 i_free1 interface~mfree->proc~i_free1 proc~l_free1 l_free1 interface~mfree->proc~l_free1 proc~r_free3 r_free3 interface~mfree->proc~r_free3 proc~r_free1 r_free1 interface~mfree->proc~r_free1 proc~i_free2 i_free2 interface~mfree->proc~i_free2 proc~r_free2 r_free2 interface~mfree->proc~r_free2 proc~l_free2 l_free2 interface~mfree->proc~l_free2 proc~i_free3 i_free3 interface~mfree->proc~i_free3 proc~chk_free chk_free proc~i_free1->proc~chk_free proc~l_free1->proc~chk_free proc~r_free3->proc~chk_free proc~chk_alloc chk_alloc proc~r_alloc1->proc~chk_alloc proc~memory_init memory_init proc~r_alloc1->proc~memory_init proc~r_alloc3->proc~chk_alloc proc~r_alloc3->proc~memory_init proc~i_alloc2->proc~chk_alloc proc~i_alloc2->proc~memory_init proc~i_alloc1->proc~chk_alloc proc~i_alloc1->proc~memory_init proc~r_alloc2->proc~chk_alloc proc~r_alloc2->proc~memory_init proc~i_alloc3->proc~chk_alloc proc~i_alloc3->proc~memory_init proc~l_alloc1->proc~chk_alloc proc~l_alloc1->proc~memory_init proc~l_alloc2->proc~chk_alloc proc~l_alloc2->proc~memory_init proc~r_free1->proc~chk_free proc~i_free2->proc~chk_free proc~r_free2->proc~chk_free proc~l_free2->proc~chk_free proc~i_free3->proc~chk_free proc~fatal_error fatal_error proc~chk_free->proc~fatal_error proc~chk_alloc->proc~fatal_error proc~ommp_message ommp_message proc~fatal_error->proc~ommp_message proc~close_output close_output proc~fatal_error->proc~close_output proc~close_output->proc~ommp_message

Called by

proc~~reverse_grp_tab~~CalledByGraph proc~reverse_grp_tab reverse_grp_tab proc~nl_update nl_update proc~nl_update->proc~reverse_grp_tab proc~mmpol_prepare mmpol_prepare proc~mmpol_prepare->proc~reverse_grp_tab proc~nl_init nl_init proc~nl_init->proc~nl_update proc~mmpol_init_from_mmp mmpol_init_from_mmp proc~mmpol_init_from_mmp->proc~mmpol_prepare proc~ommp_system_from_qm_helper ommp_system_from_qm_helper proc~ommp_system_from_qm_helper->proc~mmpol_prepare proc~assign_vdw assign_vdw proc~ommp_system_from_qm_helper->proc~assign_vdw proc~mmpol_init_from_xyz mmpol_init_from_xyz proc~mmpol_init_from_xyz->proc~mmpol_prepare proc~mmpol_init_from_xyz->proc~assign_vdw proc~vdw_init vdw_init proc~vdw_init->proc~nl_init proc~vdw_set_cutoff vdw_set_cutoff proc~vdw_set_cutoff->proc~nl_init proc~c_ommp_system_from_qm_helper C_ommp_system_from_qm_helper proc~c_ommp_system_from_qm_helper->proc~ommp_system_from_qm_helper proc~ommp_init_xyz ommp_init_xyz proc~ommp_init_xyz->proc~mmpol_init_from_xyz proc~ommp_init_mmp ommp_init_mmp proc~ommp_init_mmp->proc~mmpol_init_from_mmp proc~assign_vdw->proc~vdw_init proc~c_ommp_init_xyz C_ommp_init_xyz proc~c_ommp_init_xyz->proc~ommp_init_xyz proc~qm_helper_init_vdw qm_helper_init_vdw proc~qm_helper_init_vdw->proc~vdw_init proc~ommp_set_vdw_cutoff ommp_set_vdw_cutoff proc~ommp_set_vdw_cutoff->proc~vdw_set_cutoff proc~c_ommp_init_mmp C_ommp_init_mmp proc~c_ommp_init_mmp->proc~ommp_init_mmp proc~qm_helper_init_vdw_prm qm_helper_init_vdw_prm proc~qm_helper_init_vdw_prm->proc~assign_vdw proc~c_ommp_set_vdw_cutoff C_ommp_set_vdw_cutoff proc~c_ommp_set_vdw_cutoff->proc~ommp_set_vdw_cutoff proc~c_ommp_qm_helper_init_vdw C_ommp_qm_helper_init_vdw proc~c_ommp_qm_helper_init_vdw->proc~qm_helper_init_vdw proc~ommp_create_link_atom ommp_create_link_atom proc~ommp_create_link_atom->proc~qm_helper_init_vdw_prm proc~c_ommp_qm_helper_init_vdw_prm C_ommp_qm_helper_init_vdw_prm proc~c_ommp_qm_helper_init_vdw_prm->proc~qm_helper_init_vdw_prm proc~c_ommp_create_link_atom C_ommp_create_link_atom proc~c_ommp_create_link_atom->proc~ommp_create_link_atom

Contents

Source Code


Source Code

        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