cell_list_m Module

Sorts atoms using a cell list.

The algorithm to build the cell list partially follows the techniques in Watanabe et al, 2011, “Efficient Implementations of Molecular Dynamics Simulations for Lennard-Jones Systems,” Prog. Theor. Phys. 126, 203–235.

The pairlist is not explicitly built, rather the cells are directly looped over during force calculation.


Uses

  • module~~cell_list_m~~UsesGraph module~cell_list_m cell_list_m module~constants_m constants_m module~cell_list_m->module~constants_m module~simbox_m simbox_m module~cell_list_m->module~simbox_m module~vector_m vector_m module~cell_list_m->module~vector_m iso_fortran_env iso_fortran_env module~constants_m->iso_fortran_env module~simbox_m->module~constants_m module~random_m random_m module~simbox_m->module~random_m module~vector_m->module~constants_m module~qsort_m qsort_m module~vector_m->module~qsort_m module~qsort_m->module~constants_m module~random_m->module~constants_m mkl_vsl_type mkl_vsl_type module~random_m->mkl_vsl_type mkl_vsl mkl_vsl module~random_m->mkl_vsl

Used by

  • module~~cell_list_m~~UsedByGraph module~cell_list_m cell_list_m module~pairtab_m pairtab_m module~pairtab_m->module~cell_list_m module~interaction_m interaction_m module~interaction_m->module~pairtab_m module~bd_solver_m bd_solver_m module~bd_solver_m->module~interaction_m module~setup_m setup_m module~setup_m->module~interaction_m module~setup_m->module~bd_solver_m program~main main program~main->module~setup_m

Contents


Variables

TypeVisibility AttributesNameInitial
real(kind=rp), private, dimension(3):: cell_size =0.0_rp

Cell size along x, y, & z.

integer, private, dimension(3):: nc_max =0

Maximum number of cells along x, y, & z.

integer, private, dimension(3):: nc =0

Number of cells along x, y, & z.

integer, private :: nct_max =0

Maximum total number of cells.

integer, private :: nct =0

Total of cells.

integer, private, dimension(:), allocatable, target:: cells

(na_max,) array. Listing atoms in each cell.

integer, private, dimension(:), allocatable:: cells_pos

(0:nct_max,) index array. Note: 0-based indexing.

type(ivector_t), private :: cell_nbrs

Lists neighbor cells for each cell.

integer, private, dimension(:), allocatable:: cell_nbrs_pos

(0:nct_max,) index array. Note: 0-based indexing.

integer, private, dimension(:), allocatable:: host_cells

(na_max,) array. host_cells(i) stores the linear index of the cell containing atom i. na_max is the total number of atoms under consideration.

integer, private, dimension(:), allocatable:: cell_pop

(0:nct_max-1,) array storing population of each cell. Note: 0-based indexing.

integer, private, parameter, dimension(3,13):: d =reshape([1, 0, 0, 1, 1, 0, -1, 1, 0, 0, 1, 0, 0, 0, 1, -1, 0, 1, 1, 0, 1, -1, -1, 1, 0, -1, 1, 1, -1, 1, -1, 1, 1, 0, 1, 1, 1, 1, 1], [3, 13])

Functions

public function cl_get_num_cells() result(res)

Returns the total number of cells

Arguments

None

Return Value integer


Subroutines

public subroutine cl_init(na_max, cs_min, simbox)

Initializes a cell list.

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: na_max

Maximum number of atoms to be handled.

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

Minimum size (i.e. length) of a cell.

type(smbx_t), intent(in) :: simbox

Simulation box

public subroutine cl_set_cell_size(cs, simbox)

Sets the cell size. The actual cell size may be slightly larger.

Arguments

Type IntentOptional AttributesName
real(kind=rp), intent(in) :: cs
type(smbx_t), intent(in) :: simbox

public subroutine cl_build_cell_nbrs()

Makes a table of neighboring cells.

Arguments

None

public subroutine cl_delete()

Deallocates memory allocated in cl_init.

Arguments

None

public subroutine cl_build(coords)

Sorts atoms into cells for calculating short-range interations

Arguments

Type IntentOptional AttributesName
real(kind=rp), intent(in), dimension(:,:):: coords

public subroutine cl_get_contents(ic, res)

Returns a pointer to the entries of cell with linear index ic.

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: ic
integer, intent(out), dimension(:), pointer:: res

public subroutine cl_get_nbr_cells(ic, res)

Returns a pointer to the neighbor cells of cell with linear index ic.

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: ic
integer, intent(out), dimension(:), pointer:: res

public subroutine cl_print()

Prints a cell list

Arguments

None