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.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
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]) |
Returns the total number of cells
Initializes a cell list.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
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 |
Sets the cell size. The actual cell size may be slightly larger.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=rp), | intent(in) | :: | cs | |||
type(smbx_t), | intent(in) | :: | simbox |
Makes a table of neighboring cells.
Deallocates memory allocated in cl_init
.
Sorts atoms into cells for calculating short-range interations
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=rp), | intent(in), | dimension(:,:) | :: | coords |
Returns a pointer to the entries of cell with linear index ic.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | ic | |||
integer, | intent(out), | dimension(:), pointer | :: | res |
Returns a pointer to the neighbor cells of cell with linear index ic.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | ic | |||
integer, | intent(out), | dimension(:), pointer | :: | res |
Prints a cell list