mod_neighbor_list Module


Uses

  • module~~mod_neighbor_list~~UsesGraph module~mod_neighbor_list mod_neighbor_list module~mod_memory mod_memory module~mod_neighbor_list->module~mod_memory module~mod_adjacency_mat mod_adjacency_mat module~mod_neighbor_list->module~mod_adjacency_mat module~mod_io mod_io module~mod_neighbor_list->module~mod_io module~mod_memory->module~mod_io iso_c_binding iso_c_binding module~mod_memory->iso_c_binding module~mod_constants mod_constants module~mod_memory->module~mod_constants module~mod_adjacency_mat->module~mod_memory module~mod_io->module~mod_constants module~mod_constants->iso_c_binding

Used by


Contents


Derived Types

type, public ::  ommp_neigh_list

Components

Type Visibility Attributes Name Initial
integer(kind=ip), public :: n

Number of particles in the system

real(kind=rp), public :: cutoff

Cutoff distance for the neighbor list

integer(kind=ip), public :: cellf = 1

Number of subdivision of 'standard' cell. If you create smaller cells (cutoff/cellf) you can sensibly increase the performanc of the algorithm. See 10.1016/S0010-4655(98)00203-3

real(kind=rp), public :: celld

Dimension of each cell

integer(kind=ip), public :: ncell(3)

Number of cell along each dimension

integer(kind=ip), public :: ncells

Total number of cells

real(kind=rp), public :: offset(3)

Coordinates offset along each dimension (that is the minimum value of each coordinate.

integer(kind=ip), public :: nneigh

Number of neighbor cells

integer(kind=ip), public, allocatable :: neigh_offset(:)

Offsets of neighbor cells

integer(kind=ip), public, allocatable :: p2c(:)

Cell of each particle

type(yale_sparse), public :: c2p

Particles contained in each cell, sparse matrix format


Subroutines

public subroutine nl_init(nl, c, cutoff, f)

Arguments

Type IntentOptional Attributes Name
type(ommp_neigh_list), intent(inout) :: nl

Neigh list object to initialize

real(kind=rp), intent(in) :: c(:,:)

Coordinates in input

real(kind=rp), intent(in) :: cutoff

Cut off distance

integer(kind=ip), intent(in) :: f

Subdivision required for each cell

public subroutine nl_terminate(nl)

Arguments

Type IntentOptional Attributes Name
type(ommp_neigh_list), intent(inout) :: nl

public subroutine nl_update(nl, c)

TODO this should be improved

Arguments

Type IntentOptional Attributes Name
type(ommp_neigh_list), intent(inout) :: nl

Neigh list object to initialize

real(kind=rp), intent(in) :: c(3,nl%n)

Coordinates in input

public subroutine get_ith_nl(nl, i, c, neigh, dist, nn)

Once that the neighbor list have been initialized and updated, this function provide a logical array for atom i with all interactions that should be computed and corresponding distances.

Arguments

Type IntentOptional Attributes Name
type(ommp_neigh_list), intent(in) :: nl

Neigh list object

integer(kind=ip), intent(in) :: i

Index of atom for which the neigbor list is required

real(kind=rp), intent(in) :: c(3,nl%n)

Coordinates in input

integer(kind=ip), intent(out) :: neigh(nl%n)

Integer array with neighbors' indexes. Only the first nn elements are valid

real(kind=rp), intent(out) :: dist(nl%n)

Array for returning distances. Only the first nn elements are valid

integer(kind=ip), intent(out) :: nn

Number of neighbors