2020-08-18 14:27:28 +02:00
|
|
|
/*
|
2022-09-05 10:39:42 +02:00
|
|
|
* Copyright (C) 2022 NHR@FAU, University Erlangen-Nuremberg.
|
|
|
|
* All rights reserved. This file is part of MD-Bench.
|
|
|
|
* Use of this source code is governed by a LGPL-3.0
|
|
|
|
* license that can be found in the LICENSE file.
|
2020-08-18 14:27:28 +02:00
|
|
|
*/
|
|
|
|
#include <stdlib.h>
|
|
|
|
#include <stdio.h>
|
|
|
|
#include <math.h>
|
|
|
|
|
|
|
|
#include <neighbor.h>
|
|
|
|
#include <parameter.h>
|
|
|
|
#include <atom.h>
|
2022-01-25 00:43:10 +01:00
|
|
|
#include <util.h>
|
2020-08-18 14:27:28 +02:00
|
|
|
|
|
|
|
#define SMALL 1.0e-6
|
|
|
|
#define FACTOR 0.999
|
|
|
|
|
2022-01-27 03:07:31 +01:00
|
|
|
static MD_FLOAT xprd, yprd, zprd;
|
2022-01-25 00:43:10 +01:00
|
|
|
static MD_FLOAT bininvx, bininvy;
|
|
|
|
static int mbinxlo, mbinylo;
|
|
|
|
static int nbinx, nbiny;
|
|
|
|
static int mbinx, mbiny; // n bins in x, y
|
2020-08-18 14:27:28 +02:00
|
|
|
static int *bincount;
|
|
|
|
static int *bins;
|
2022-03-09 02:25:39 +01:00
|
|
|
static int *bin_nclusters;
|
|
|
|
static int *bin_clusters;
|
2020-08-18 14:27:28 +02:00
|
|
|
static int mbins; //total number of bins
|
|
|
|
static int atoms_per_bin; // max atoms per bin
|
2022-01-25 00:43:10 +01:00
|
|
|
static int clusters_per_bin; // max clusters per bin
|
2020-11-05 12:41:44 +01:00
|
|
|
static MD_FLOAT cutneigh;
|
|
|
|
static MD_FLOAT cutneighsq; // neighbor cutoff squared
|
2020-08-18 14:27:28 +02:00
|
|
|
static int nmax;
|
2020-11-05 12:41:44 +01:00
|
|
|
static int nstencil; // # of bins in stencil
|
|
|
|
static int* stencil; // stencil list of bin offsets
|
2022-01-25 00:43:10 +01:00
|
|
|
static MD_FLOAT binsizex, binsizey;
|
2020-08-18 14:27:28 +02:00
|
|
|
|
2022-01-25 00:43:10 +01:00
|
|
|
static int coord2bin(MD_FLOAT, MD_FLOAT);
|
|
|
|
static MD_FLOAT bindist(int, int);
|
|
|
|
|
2020-08-18 14:27:28 +02:00
|
|
|
/* exported subroutines */
|
2022-01-27 03:07:31 +01:00
|
|
|
void initNeighbor(Neighbor *neighbor, Parameter *param) {
|
2020-11-05 12:41:44 +01:00
|
|
|
MD_FLOAT neighscale = 5.0 / 6.0;
|
2021-03-30 01:54:56 +02:00
|
|
|
xprd = param->nx * param->lattice;
|
|
|
|
yprd = param->ny * param->lattice;
|
2022-01-27 03:07:31 +01:00
|
|
|
zprd = param->nz * param->lattice;
|
2020-08-18 14:27:28 +02:00
|
|
|
cutneigh = param->cutneigh;
|
|
|
|
nmax = 0;
|
|
|
|
atoms_per_bin = 8;
|
2022-03-10 01:31:50 +01:00
|
|
|
clusters_per_bin = (atoms_per_bin / CLUSTER_M) + 10;
|
2020-08-18 14:27:28 +02:00
|
|
|
stencil = NULL;
|
|
|
|
bins = NULL;
|
|
|
|
bincount = NULL;
|
2022-03-09 02:25:39 +01:00
|
|
|
bin_clusters = NULL;
|
|
|
|
bin_nclusters = NULL;
|
2022-03-22 23:47:05 +01:00
|
|
|
neighbor->half_neigh = param->half_neigh;
|
2020-08-18 14:27:28 +02:00
|
|
|
neighbor->maxneighs = 100;
|
|
|
|
neighbor->numneigh = NULL;
|
|
|
|
neighbor->neighbors = NULL;
|
|
|
|
}
|
|
|
|
|
2022-01-27 03:07:31 +01:00
|
|
|
void setupNeighbor(Parameter *param, Atom *atom) {
|
2020-11-05 12:41:44 +01:00
|
|
|
MD_FLOAT coord;
|
2022-01-25 00:43:10 +01:00
|
|
|
int mbinxhi, mbinyhi;
|
2020-08-18 14:27:28 +02:00
|
|
|
int nextx, nexty, nextz;
|
2021-11-30 01:33:55 +01:00
|
|
|
|
|
|
|
if(param->input_file != NULL) {
|
|
|
|
xprd = param->xprd;
|
|
|
|
yprd = param->yprd;
|
2022-01-27 03:07:31 +01:00
|
|
|
zprd = param->zprd;
|
2021-11-30 01:33:55 +01:00
|
|
|
}
|
|
|
|
|
|
|
|
// TODO: update lo and hi for standard case and use them here instead
|
2020-11-05 12:41:44 +01:00
|
|
|
MD_FLOAT xlo = 0.0; MD_FLOAT xhi = xprd;
|
|
|
|
MD_FLOAT ylo = 0.0; MD_FLOAT yhi = yprd;
|
2022-01-27 03:07:31 +01:00
|
|
|
MD_FLOAT zlo = 0.0; MD_FLOAT zhi = zprd;
|
2020-08-18 14:27:28 +02:00
|
|
|
|
2022-01-27 03:07:31 +01:00
|
|
|
MD_FLOAT atom_density = ((MD_FLOAT)(atom->Nlocal)) / ((xhi - xlo) * (yhi - ylo) * (zhi - zlo));
|
2022-03-09 02:25:39 +01:00
|
|
|
MD_FLOAT atoms_in_cell = MAX(CLUSTER_M, CLUSTER_N);
|
2022-02-24 16:42:58 +01:00
|
|
|
MD_FLOAT targetsizex = cbrt(atoms_in_cell / atom_density);
|
|
|
|
MD_FLOAT targetsizey = cbrt(atoms_in_cell / atom_density);
|
|
|
|
nbinx = MAX(1, (int)ceil((xhi - xlo) / targetsizex));
|
|
|
|
nbiny = MAX(1, (int)ceil((yhi - ylo) / targetsizey));
|
|
|
|
binsizex = (xhi - xlo) / nbinx;
|
|
|
|
binsizey = (yhi - ylo) / nbiny;
|
2022-01-27 03:07:31 +01:00
|
|
|
bininvx = 1.0 / binsizex;
|
|
|
|
bininvy = 1.0 / binsizey;
|
2022-02-24 16:42:58 +01:00
|
|
|
cutneighsq = cutneigh * cutneigh;
|
2022-01-27 03:07:31 +01:00
|
|
|
|
2020-08-18 14:27:28 +02:00
|
|
|
coord = xlo - cutneigh - SMALL * xprd;
|
2022-03-09 02:25:39 +01:00
|
|
|
mbinxlo = (int)(coord * bininvx);
|
|
|
|
if(coord < 0.0) { mbinxlo = mbinxlo - 1; }
|
2020-08-18 14:27:28 +02:00
|
|
|
coord = xhi + cutneigh + SMALL * xprd;
|
2022-03-09 02:25:39 +01:00
|
|
|
mbinxhi = (int)(coord * bininvx);
|
2020-08-18 14:27:28 +02:00
|
|
|
|
|
|
|
coord = ylo - cutneigh - SMALL * yprd;
|
2022-03-09 02:25:39 +01:00
|
|
|
mbinylo = (int)(coord * bininvy);
|
|
|
|
if(coord < 0.0) { mbinylo = mbinylo - 1; }
|
2020-08-18 14:27:28 +02:00
|
|
|
coord = yhi + cutneigh + SMALL * yprd;
|
2022-03-09 02:25:39 +01:00
|
|
|
mbinyhi = (int)(coord * bininvy);
|
2020-08-18 14:27:28 +02:00
|
|
|
|
|
|
|
mbinxlo = mbinxlo - 1;
|
|
|
|
mbinxhi = mbinxhi + 1;
|
|
|
|
mbinx = mbinxhi - mbinxlo + 1;
|
|
|
|
|
|
|
|
mbinylo = mbinylo - 1;
|
|
|
|
mbinyhi = mbinyhi + 1;
|
|
|
|
mbiny = mbinyhi - mbinylo + 1;
|
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
nextx = (int)(cutneigh * bininvx);
|
|
|
|
nexty = (int)(cutneigh * bininvy);
|
2020-08-18 14:27:28 +02:00
|
|
|
if(nextx * binsizex < FACTOR * cutneigh) nextx++;
|
|
|
|
if(nexty * binsizey < FACTOR * cutneigh) nexty++;
|
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
if (stencil) { free(stencil); }
|
2022-01-25 00:43:10 +01:00
|
|
|
stencil = (int *) malloc((2 * nexty + 1) * (2 * nextx + 1) * sizeof(int));
|
2020-08-18 14:27:28 +02:00
|
|
|
nstencil = 0;
|
2022-01-25 00:43:10 +01:00
|
|
|
|
|
|
|
for(int j = -nexty; j <= nexty; j++) {
|
|
|
|
for(int i = -nextx; i <= nextx; i++) {
|
|
|
|
if(bindist(i, j) < cutneighsq) {
|
|
|
|
stencil[nstencil++] = j * mbinx + i;
|
2020-08-18 14:27:28 +02:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
if(bincount) { free(bincount); }
|
|
|
|
if(bins) { free(bins); }
|
|
|
|
if(bin_nclusters) { free(bin_nclusters); }
|
|
|
|
if(bin_clusters) { free(bin_clusters); }
|
2022-01-25 00:43:10 +01:00
|
|
|
mbins = mbinx * mbiny;
|
2020-08-18 14:27:28 +02:00
|
|
|
bincount = (int*) malloc(mbins * sizeof(int));
|
|
|
|
bins = (int*) malloc(mbins * atoms_per_bin * sizeof(int));
|
2022-03-09 02:25:39 +01:00
|
|
|
bin_nclusters = (int*) malloc(mbins * sizeof(int));
|
|
|
|
bin_clusters = (int*) malloc(mbins * clusters_per_bin * sizeof(int));
|
2022-02-24 16:42:58 +01:00
|
|
|
|
|
|
|
/*
|
|
|
|
DEBUG_MESSAGE("lo, hi = (%e, %e, %e), (%e, %e, %e)\n", xlo, ylo, zlo, xhi, yhi, zhi);
|
|
|
|
DEBUG_MESSAGE("binsize = %e, %e\n", binsizex, binsizey);
|
|
|
|
DEBUG_MESSAGE("mbin lo, hi = (%d, %d), (%d, %d)\n", mbinxlo, mbinylo, mbinxhi, mbinyhi);
|
|
|
|
DEBUG_MESSAGE("mbins = %d (%d x %d)\n", mbins, mbinx, mbiny);
|
|
|
|
DEBUG_MESSAGE("nextx = %d, nexty = %d\n", nextx, nexty);
|
|
|
|
*/
|
2020-08-18 14:27:28 +02:00
|
|
|
}
|
|
|
|
|
2022-01-25 00:43:10 +01:00
|
|
|
MD_FLOAT getBoundingBoxDistanceSq(Atom *atom, int ci, int cj) {
|
2022-03-09 02:25:39 +01:00
|
|
|
MD_FLOAT dl = atom->iclusters[ci].bbminx - atom->jclusters[cj].bbmaxx;
|
|
|
|
MD_FLOAT dh = atom->jclusters[cj].bbminx - atom->iclusters[ci].bbmaxx;
|
2022-01-25 00:43:10 +01:00
|
|
|
MD_FLOAT dm = MAX(dl, dh);
|
|
|
|
MD_FLOAT dm0 = MAX(dm, 0.0);
|
|
|
|
MD_FLOAT d2 = dm0 * dm0;
|
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
dl = atom->iclusters[ci].bbminy - atom->jclusters[cj].bbmaxy;
|
|
|
|
dh = atom->jclusters[cj].bbminy - atom->iclusters[ci].bbmaxy;
|
2022-01-25 00:43:10 +01:00
|
|
|
dm = MAX(dl, dh);
|
|
|
|
dm0 = MAX(dm, 0.0);
|
|
|
|
d2 += dm0 * dm0;
|
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
dl = atom->iclusters[ci].bbminz - atom->jclusters[cj].bbmaxz;
|
|
|
|
dh = atom->jclusters[cj].bbminz - atom->iclusters[ci].bbmaxz;
|
2022-01-25 00:43:10 +01:00
|
|
|
dm = MAX(dl, dh);
|
|
|
|
dm0 = MAX(dm, 0.0);
|
|
|
|
d2 += dm0 * dm0;
|
|
|
|
return d2;
|
|
|
|
}
|
|
|
|
|
|
|
|
int atomDistanceInRange(Atom *atom, int ci, int cj, MD_FLOAT rsq) {
|
2022-03-09 02:25:39 +01:00
|
|
|
int ci_vec_base = CI_VECTOR_BASE_INDEX(ci);
|
|
|
|
int cj_vec_base = CJ_VECTOR_BASE_INDEX(cj);
|
|
|
|
MD_FLOAT *ci_x = &atom->cl_x[ci_vec_base];
|
|
|
|
MD_FLOAT *cj_x = &atom->cl_x[cj_vec_base];
|
|
|
|
|
|
|
|
for(int cii = 0; cii < atom->iclusters[ci].natoms; cii++) {
|
|
|
|
for(int cjj = 0; cjj < atom->jclusters[cj].natoms; cjj++) {
|
|
|
|
MD_FLOAT delx = ci_x[CL_X_OFFSET + cii] - cj_x[CL_X_OFFSET + cjj];
|
|
|
|
MD_FLOAT dely = ci_x[CL_Y_OFFSET + cii] - cj_x[CL_Y_OFFSET + cjj];
|
|
|
|
MD_FLOAT delz = ci_x[CL_Z_OFFSET + cii] - cj_x[CL_Z_OFFSET + cjj];
|
2022-01-25 00:43:10 +01:00
|
|
|
if(delx * delx + dely * dely + delz * delz < rsq) {
|
|
|
|
return 1;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
return 0;
|
|
|
|
}
|
|
|
|
|
2023-03-23 00:58:25 +01:00
|
|
|
/* Returns a diagonal or off-diagonal interaction mask for plain C lists */
|
|
|
|
static unsigned int get_imask(int rdiag, int ci, int cj) {
|
|
|
|
return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
|
|
|
|
}
|
|
|
|
|
|
|
|
/* Returns a diagonal or off-diagonal interaction mask for cj-size=2 */
|
|
|
|
static unsigned int get_imask_simd_j2(int rdiag, int ci, int cj) {
|
|
|
|
return (rdiag && ci * 2 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_0
|
|
|
|
: (rdiag && ci * 2 + 1 == cj ? NBNXN_INTERACTION_MASK_DIAG_J2_1
|
|
|
|
: NBNXN_INTERACTION_MASK_ALL));
|
|
|
|
}
|
|
|
|
|
|
|
|
/* Returns a diagonal or off-diagonal interaction mask for cj-size=4 */
|
|
|
|
static unsigned int get_imask_simd_j4(int rdiag, int ci, int cj) {
|
|
|
|
return (rdiag && ci == cj ? NBNXN_INTERACTION_MASK_DIAG : NBNXN_INTERACTION_MASK_ALL);
|
|
|
|
}
|
|
|
|
|
|
|
|
/* Returns a diagonal or off-diagonal interaction mask for cj-size=8 */
|
|
|
|
static unsigned int get_imask_simd_j8(int rdiag, int ci, int cj) {
|
|
|
|
return (rdiag && ci == cj * 2 ? NBNXN_INTERACTION_MASK_DIAG_J8_0
|
|
|
|
: (rdiag && ci == cj * 2 + 1 ? NBNXN_INTERACTION_MASK_DIAG_J8_1
|
|
|
|
: NBNXN_INTERACTION_MASK_ALL));
|
|
|
|
}
|
|
|
|
|
|
|
|
#if VECTOR_WIDTH == 2
|
|
|
|
# define get_imask_simd_4xn get_imask_simd_j2
|
|
|
|
#elif VECTOR_WIDTH== 4
|
|
|
|
# define get_imask_simd_4xn get_imask_simd_j4
|
|
|
|
#elif VECTOR_WIDTH == 8
|
|
|
|
# define get_imask_simd_4xn get_imask_simd_j8
|
|
|
|
# define get_imask_simd_2xnn get_imask_simd_j4
|
|
|
|
#elif VECTOR_WIDTH == 16
|
|
|
|
# define get_imask_simd_2xnn get_imask_simd_j8
|
|
|
|
#else
|
|
|
|
# error "Invalid cluster configuration"
|
|
|
|
#endif
|
|
|
|
|
2022-01-27 03:07:31 +01:00
|
|
|
void buildNeighbor(Atom *atom, Neighbor *neighbor) {
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("buildNeighbor start\n");
|
2020-08-18 14:27:28 +02:00
|
|
|
|
|
|
|
/* extend atom arrays if necessary */
|
2022-03-10 22:33:41 +01:00
|
|
|
if(atom->Nclusters_local > nmax) {
|
|
|
|
nmax = atom->Nclusters_local;
|
2020-08-18 14:27:28 +02:00
|
|
|
if(neighbor->numneigh) free(neighbor->numneigh);
|
|
|
|
if(neighbor->neighbors) free(neighbor->neighbors);
|
|
|
|
neighbor->numneigh = (int*) malloc(nmax * sizeof(int));
|
2023-03-23 00:58:25 +01:00
|
|
|
neighbor->neighbors = (NeighborCluster*) malloc(nmax * neighbor->maxneighs * sizeof(NeighborCluster));
|
2020-08-18 14:27:28 +02:00
|
|
|
}
|
|
|
|
|
2022-01-31 17:49:22 +01:00
|
|
|
MD_FLOAT bbx = 0.5 * (binsizex + binsizex);
|
|
|
|
MD_FLOAT bby = 0.5 * (binsizey + binsizey);
|
|
|
|
MD_FLOAT rbb_sq = MAX(0.0, cutneigh - 0.5 * sqrt(bbx * bbx + bby * bby));
|
|
|
|
rbb_sq = rbb_sq * rbb_sq;
|
2020-08-18 14:27:28 +02:00
|
|
|
int resize = 1;
|
|
|
|
|
|
|
|
/* loop over each atom, storing neighbors */
|
|
|
|
while(resize) {
|
|
|
|
int new_maxneighs = neighbor->maxneighs;
|
|
|
|
resize = 0;
|
|
|
|
|
2022-01-25 00:43:10 +01:00
|
|
|
for(int ci = 0; ci < atom->Nclusters_local; ci++) {
|
2022-03-22 23:47:05 +01:00
|
|
|
int ci_cj1 = CJ1_FROM_CI(ci);
|
2023-03-23 00:58:25 +01:00
|
|
|
NeighborCluster *neighptr = &(neighbor->neighbors[ci * neighbor->maxneighs]);
|
2020-08-18 14:27:28 +02:00
|
|
|
int n = 0;
|
2022-03-10 01:31:50 +01:00
|
|
|
int ibin = atom->icluster_bin[ci];
|
2022-03-09 02:25:39 +01:00
|
|
|
MD_FLOAT ibb_xmin = atom->iclusters[ci].bbminx;
|
|
|
|
MD_FLOAT ibb_xmax = atom->iclusters[ci].bbmaxx;
|
|
|
|
MD_FLOAT ibb_ymin = atom->iclusters[ci].bbminy;
|
|
|
|
MD_FLOAT ibb_ymax = atom->iclusters[ci].bbmaxy;
|
|
|
|
MD_FLOAT ibb_zmin = atom->iclusters[ci].bbminz;
|
|
|
|
MD_FLOAT ibb_zmax = atom->iclusters[ci].bbmaxz;
|
2022-01-25 00:43:10 +01:00
|
|
|
|
2020-08-18 14:27:28 +02:00
|
|
|
for(int k = 0; k < nstencil; k++) {
|
|
|
|
int jbin = ibin + stencil[k];
|
2022-03-09 02:25:39 +01:00
|
|
|
int *loc_bin = &bin_clusters[jbin * clusters_per_bin];
|
2022-01-25 00:43:10 +01:00
|
|
|
int cj, m = -1;
|
2022-02-07 18:00:21 +01:00
|
|
|
MD_FLOAT jbb_xmin, jbb_xmax, jbb_ymin, jbb_ymax, jbb_zmin, jbb_zmax;
|
2022-03-09 02:25:39 +01:00
|
|
|
const int c = bin_nclusters[jbin];
|
2022-01-27 03:07:31 +01:00
|
|
|
|
|
|
|
if(c > 0) {
|
2022-02-07 18:28:53 +01:00
|
|
|
MD_FLOAT dl, dh, dm, dm0, d_bb_sq;
|
|
|
|
|
|
|
|
do {
|
2022-01-27 03:07:31 +01:00
|
|
|
m++;
|
|
|
|
cj = loc_bin[m];
|
2022-03-22 23:47:05 +01:00
|
|
|
if(neighbor->half_neigh && ci_cj1 > cj) {
|
|
|
|
continue;
|
|
|
|
}
|
2022-03-09 02:25:39 +01:00
|
|
|
jbb_zmin = atom->jclusters[cj].bbminz;
|
|
|
|
jbb_zmax = atom->jclusters[cj].bbmaxz;
|
2022-02-07 18:28:53 +01:00
|
|
|
dl = ibb_zmin - jbb_zmax;
|
|
|
|
dh = jbb_zmin - ibb_zmax;
|
|
|
|
dm = MAX(dl, dh);
|
|
|
|
dm0 = MAX(dm, 0.0);
|
|
|
|
d_bb_sq = dm0 * dm0;
|
|
|
|
} while(m + 1 < c && d_bb_sq > cutneighsq);
|
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
jbb_xmin = atom->jclusters[cj].bbminx;
|
|
|
|
jbb_xmax = atom->jclusters[cj].bbmaxx;
|
|
|
|
jbb_ymin = atom->jclusters[cj].bbminy;
|
|
|
|
jbb_ymax = atom->jclusters[cj].bbmaxy;
|
2022-01-27 03:07:31 +01:00
|
|
|
|
2022-02-01 00:46:12 +01:00
|
|
|
while(m < c) {
|
2022-03-22 23:47:05 +01:00
|
|
|
if(!neighbor->half_neigh || ci_cj1 <= cj) {
|
|
|
|
dl = ibb_zmin - jbb_zmax;
|
|
|
|
dh = jbb_zmin - ibb_zmax;
|
|
|
|
dm = MAX(dl, dh);
|
|
|
|
dm0 = MAX(dm, 0.0);
|
|
|
|
d_bb_sq = dm0 * dm0;
|
|
|
|
|
|
|
|
/*if(d_bb_sq > cutneighsq) {
|
|
|
|
break;
|
|
|
|
}*/
|
|
|
|
|
|
|
|
dl = ibb_ymin - jbb_ymax;
|
|
|
|
dh = jbb_ymin - ibb_ymax;
|
|
|
|
dm = MAX(dl, dh);
|
|
|
|
dm0 = MAX(dm, 0.0);
|
|
|
|
d_bb_sq += dm0 * dm0;
|
|
|
|
|
|
|
|
dl = ibb_xmin - jbb_xmax;
|
|
|
|
dh = jbb_xmin - ibb_xmax;
|
|
|
|
dm = MAX(dl, dh);
|
|
|
|
dm0 = MAX(dm, 0.0);
|
|
|
|
d_bb_sq += dm0 * dm0;
|
|
|
|
|
|
|
|
if(d_bb_sq < cutneighsq) {
|
|
|
|
if(d_bb_sq < rbb_sq || atomDistanceInRange(atom, ci, cj, cutneighsq)) {
|
2023-03-23 00:58:25 +01:00
|
|
|
unsigned int imask;
|
|
|
|
#if CLUSTER_N == (VECTOR_WIDTH / 2) // 2xnn
|
|
|
|
imask = get_imask_simd_2xnn(neighbor->half_neigh, ci, cj);
|
|
|
|
#else // 4xn
|
|
|
|
imask = get_imask_simd_4xn(neighbor->half_neigh, ci, cj);
|
|
|
|
#endif
|
|
|
|
|
|
|
|
neighptr[n].cj = cj;
|
|
|
|
neighptr[n].imask = imask;
|
|
|
|
n++;
|
2022-03-22 23:47:05 +01:00
|
|
|
}
|
2022-01-25 00:43:10 +01:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2022-01-27 03:07:31 +01:00
|
|
|
m++;
|
|
|
|
if(m < c) {
|
|
|
|
cj = loc_bin[m];
|
2022-03-09 02:25:39 +01:00
|
|
|
jbb_xmin = atom->jclusters[cj].bbminx;
|
|
|
|
jbb_xmax = atom->jclusters[cj].bbmaxx;
|
|
|
|
jbb_ymin = atom->jclusters[cj].bbminy;
|
|
|
|
jbb_ymax = atom->jclusters[cj].bbmaxy;
|
|
|
|
jbb_zmin = atom->jclusters[cj].bbminz;
|
|
|
|
jbb_zmax = atom->jclusters[cj].bbmaxz;
|
2022-01-27 03:07:31 +01:00
|
|
|
}
|
|
|
|
}
|
2020-08-18 14:27:28 +02:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
// Fill neighbor list with dummy values to fit vector width
|
|
|
|
if(CLUSTER_N < VECTOR_WIDTH) {
|
|
|
|
while(n % (VECTOR_WIDTH / CLUSTER_N)) {
|
2023-03-23 00:58:25 +01:00
|
|
|
neighptr[n].cj = atom->dummy_cj; // Last cluster is always a dummy cluster
|
|
|
|
neighptr[n].imask = 0;
|
|
|
|
n++;
|
2022-02-02 18:00:44 +01:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2022-01-25 00:43:10 +01:00
|
|
|
neighbor->numneigh[ci] = n;
|
2020-08-18 14:27:28 +02:00
|
|
|
if(n >= neighbor->maxneighs) {
|
|
|
|
resize = 1;
|
|
|
|
|
|
|
|
if(n >= new_maxneighs) {
|
|
|
|
new_maxneighs = n;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
if(resize) {
|
2022-02-01 20:16:04 +01:00
|
|
|
fprintf(stdout, "RESIZE %d\n", neighbor->maxneighs);
|
2020-08-18 14:27:28 +02:00
|
|
|
neighbor->maxneighs = new_maxneighs * 1.2;
|
|
|
|
free(neighbor->neighbors);
|
2023-03-23 00:58:25 +01:00
|
|
|
neighbor->neighbors = (NeighborCluster*) malloc(atom->Nmax * neighbor->maxneighs * sizeof(int));
|
2020-08-18 14:27:28 +02:00
|
|
|
}
|
|
|
|
}
|
2022-01-31 17:49:22 +01:00
|
|
|
|
|
|
|
/*
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("\ncutneighsq = %f, rbb_sq = %f\n", cutneighsq, rbb_sq);
|
2022-01-31 23:46:20 +01:00
|
|
|
for(int ci = 0; ci < 6; ci++) {
|
2022-03-10 22:33:41 +01:00
|
|
|
int ci_vec_base = CI_VECTOR_BASE_INDEX(ci);
|
|
|
|
MD_FLOAT *ci_x = &atom->cl_x[ci_vec_base];
|
2022-01-31 17:49:22 +01:00
|
|
|
int* neighptr = &(neighbor->neighbors[ci * neighbor->maxneighs]);
|
2022-02-01 20:16:04 +01:00
|
|
|
|
|
|
|
DEBUG_MESSAGE("Cluster %d, bbx = {%f, %f}, bby = {%f, %f}, bbz = {%f, %f}\n",
|
|
|
|
ci,
|
2022-03-10 22:33:41 +01:00
|
|
|
atom->iclusters[ci].bbminx,
|
|
|
|
atom->iclusters[ci].bbmaxx,
|
|
|
|
atom->iclusters[ci].bbminy,
|
|
|
|
atom->iclusters[ci].bbmaxy,
|
|
|
|
atom->iclusters[ci].bbminz,
|
|
|
|
atom->iclusters[ci].bbmaxz);
|
2022-02-01 20:16:04 +01:00
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
for(int cii = 0; cii < CLUSTER_M; cii++) {
|
2022-03-10 22:33:41 +01:00
|
|
|
DEBUG_MESSAGE("%f, %f, %f\n", ci_x[CL_X_OFFSET + cii], ci_x[CL_Y_OFFSET + cii], ci_x[CL_Z_OFFSET + cii]);
|
2022-01-31 17:49:22 +01:00
|
|
|
}
|
|
|
|
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("Neighbors:\n");
|
2022-01-31 17:49:22 +01:00
|
|
|
for(int k = 0; k < neighbor->numneigh[ci]; k++) {
|
2022-03-10 22:33:41 +01:00
|
|
|
int cj = neighptr[k];
|
|
|
|
int cj_vec_base = CJ_VECTOR_BASE_INDEX(cj);
|
|
|
|
MD_FLOAT *cj_x = &atom->cl_x[cj_vec_base];
|
2022-02-01 20:16:04 +01:00
|
|
|
|
|
|
|
DEBUG_MESSAGE(" Cluster %d, bbx = {%f, %f}, bby = {%f, %f}, bbz = {%f, %f}\n",
|
|
|
|
cj,
|
2022-03-10 22:33:41 +01:00
|
|
|
atom->jclusters[cj].bbminx,
|
|
|
|
atom->jclusters[cj].bbmaxx,
|
|
|
|
atom->jclusters[cj].bbminy,
|
|
|
|
atom->jclusters[cj].bbmaxy,
|
|
|
|
atom->jclusters[cj].bbminz,
|
|
|
|
atom->jclusters[cj].bbmaxz);
|
|
|
|
|
|
|
|
for(int cjj = 0; cjj < CLUSTER_N; cjj++) {
|
|
|
|
DEBUG_MESSAGE(" %f, %f, %f\n", cj_x[CL_X_OFFSET + cjj], cj_x[CL_Y_OFFSET + cjj], cj_x[CL_Z_OFFSET + cjj]);
|
2022-01-31 17:49:22 +01:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
*/
|
|
|
|
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("buildNeighbor end\n");
|
2020-08-18 14:27:28 +02:00
|
|
|
}
|
|
|
|
|
2022-02-28 17:20:39 +01:00
|
|
|
void pruneNeighbor(Parameter *param, Atom *atom, Neighbor *neighbor) {
|
|
|
|
DEBUG_MESSAGE("pruneNeighbor start\n");
|
|
|
|
//MD_FLOAT cutsq = param->cutforce * param->cutforce;
|
|
|
|
MD_FLOAT cutsq = cutneighsq;
|
|
|
|
|
|
|
|
for(int ci = 0; ci < atom->Nclusters_local; ci++) {
|
2023-03-23 00:58:25 +01:00
|
|
|
NeighborCluster *neighs = &neighbor->neighbors[ci * neighbor->maxneighs];
|
2022-02-28 17:20:39 +01:00
|
|
|
int numneighs = neighbor->numneigh[ci];
|
|
|
|
int k = 0;
|
|
|
|
|
|
|
|
// Remove dummy clusters if necessary
|
2022-03-09 02:25:39 +01:00
|
|
|
if(CLUSTER_N < VECTOR_WIDTH) {
|
2023-03-23 00:58:25 +01:00
|
|
|
while(neighs[numneighs - 1].cj == atom->dummy_cj) {
|
2022-02-28 17:20:39 +01:00
|
|
|
numneighs--;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
while(k < numneighs) {
|
2023-03-23 00:58:25 +01:00
|
|
|
int cj = neighs[k].cj;
|
2022-02-28 17:20:39 +01:00
|
|
|
if(atomDistanceInRange(atom, ci, cj, cutsq)) {
|
|
|
|
k++;
|
|
|
|
} else {
|
|
|
|
numneighs--;
|
|
|
|
neighs[k] = neighs[numneighs];
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
// Readd dummy clusters if necessary
|
2022-03-09 02:25:39 +01:00
|
|
|
if(CLUSTER_N < VECTOR_WIDTH) {
|
|
|
|
while(numneighs % (VECTOR_WIDTH / CLUSTER_N)) {
|
2023-03-23 00:58:25 +01:00
|
|
|
neighs[numneighs].cj = atom->dummy_cj; // Last cluster is always a dummy cluster
|
|
|
|
neighs[numneighs].imask = 0;
|
|
|
|
numneighs++;
|
2022-02-28 17:20:39 +01:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
neighbor->numneigh[ci] = numneighs;
|
|
|
|
}
|
|
|
|
|
|
|
|
DEBUG_MESSAGE("pruneNeighbor end\n");
|
|
|
|
}
|
|
|
|
|
2020-08-18 14:27:28 +02:00
|
|
|
/* internal subroutines */
|
2022-01-25 00:43:10 +01:00
|
|
|
MD_FLOAT bindist(int i, int j) {
|
2020-11-05 12:41:44 +01:00
|
|
|
MD_FLOAT delx, dely, delz;
|
2020-08-18 14:27:28 +02:00
|
|
|
|
|
|
|
if(i > 0) {
|
|
|
|
delx = (i - 1) * binsizex;
|
|
|
|
} else if(i == 0) {
|
|
|
|
delx = 0.0;
|
|
|
|
} else {
|
|
|
|
delx = (i + 1) * binsizex;
|
|
|
|
}
|
|
|
|
|
|
|
|
if(j > 0) {
|
|
|
|
dely = (j - 1) * binsizey;
|
|
|
|
} else if(j == 0) {
|
|
|
|
dely = 0.0;
|
|
|
|
} else {
|
|
|
|
dely = (j + 1) * binsizey;
|
|
|
|
}
|
|
|
|
|
2022-01-25 00:43:10 +01:00
|
|
|
return (delx * delx + dely * dely);
|
2020-08-18 14:27:28 +02:00
|
|
|
}
|
|
|
|
|
2022-01-25 00:43:10 +01:00
|
|
|
int coord2bin(MD_FLOAT xin, MD_FLOAT yin) {
|
|
|
|
int ix, iy;
|
2020-08-18 14:27:28 +02:00
|
|
|
|
|
|
|
if(xin >= xprd) {
|
|
|
|
ix = (int)((xin - xprd) * bininvx) + nbinx - mbinxlo;
|
|
|
|
} else if(xin >= 0.0) {
|
|
|
|
ix = (int)(xin * bininvx) - mbinxlo;
|
|
|
|
} else {
|
|
|
|
ix = (int)(xin * bininvx) - mbinxlo - 1;
|
|
|
|
}
|
|
|
|
|
|
|
|
if(yin >= yprd) {
|
|
|
|
iy = (int)((yin - yprd) * bininvy) + nbiny - mbinylo;
|
|
|
|
} else if(yin >= 0.0) {
|
|
|
|
iy = (int)(yin * bininvy) - mbinylo;
|
|
|
|
} else {
|
|
|
|
iy = (int)(yin * bininvy) - mbinylo - 1;
|
|
|
|
}
|
|
|
|
|
2022-01-25 00:43:10 +01:00
|
|
|
return (iy * mbinx + ix + 1);
|
|
|
|
}
|
|
|
|
|
|
|
|
void coord2bin2D(MD_FLOAT xin, MD_FLOAT yin, int *ix, int *iy) {
|
|
|
|
if(xin >= xprd) {
|
|
|
|
*ix = (int)((xin - xprd) * bininvx) + nbinx - mbinxlo;
|
|
|
|
} else if(xin >= 0.0) {
|
|
|
|
*ix = (int)(xin * bininvx) - mbinxlo;
|
2020-08-18 14:27:28 +02:00
|
|
|
} else {
|
2022-01-25 00:43:10 +01:00
|
|
|
*ix = (int)(xin * bininvx) - mbinxlo - 1;
|
2020-08-18 14:27:28 +02:00
|
|
|
}
|
|
|
|
|
2022-01-25 00:43:10 +01:00
|
|
|
if(yin >= yprd) {
|
|
|
|
*iy = (int)((yin - yprd) * bininvy) + nbiny - mbinylo;
|
|
|
|
} else if(yin >= 0.0) {
|
|
|
|
*iy = (int)(yin * bininvy) - mbinylo;
|
|
|
|
} else {
|
|
|
|
*iy = (int)(yin * bininvy) - mbinylo - 1;
|
|
|
|
}
|
2020-08-18 14:27:28 +02:00
|
|
|
}
|
|
|
|
|
2022-01-27 03:07:31 +01:00
|
|
|
void binAtoms(Atom *atom) {
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("binAtoms start\n");
|
2020-08-18 14:27:28 +02:00
|
|
|
int resize = 1;
|
|
|
|
|
|
|
|
while(resize > 0) {
|
|
|
|
resize = 0;
|
|
|
|
|
|
|
|
for(int i = 0; i < mbins; i++) {
|
|
|
|
bincount[i] = 0;
|
|
|
|
}
|
|
|
|
|
2022-01-31 23:46:20 +01:00
|
|
|
for(int i = 0; i < atom->Nlocal; i++) {
|
2022-01-25 00:43:10 +01:00
|
|
|
int ibin = coord2bin(atom_x(i), atom_y(i));
|
2020-08-18 14:27:28 +02:00
|
|
|
if(bincount[ibin] < atoms_per_bin) {
|
|
|
|
int ac = bincount[ibin]++;
|
|
|
|
bins[ibin * atoms_per_bin + ac] = i;
|
|
|
|
} else {
|
|
|
|
resize = 1;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
if(resize) {
|
|
|
|
free(bins);
|
|
|
|
atoms_per_bin *= 2;
|
|
|
|
bins = (int*) malloc(mbins * atoms_per_bin * sizeof(int));
|
|
|
|
}
|
|
|
|
}
|
2022-02-01 20:16:04 +01:00
|
|
|
|
|
|
|
DEBUG_MESSAGE("binAtoms end\n");
|
2020-08-18 14:27:28 +02:00
|
|
|
}
|
2021-11-03 00:57:24 +01:00
|
|
|
|
2022-01-25 00:43:10 +01:00
|
|
|
// TODO: Use pigeonhole sorting
|
2022-01-27 03:07:31 +01:00
|
|
|
void sortAtomsByZCoord(Atom *atom) {
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("sortAtomsByZCoord start\n");
|
2022-01-27 03:07:31 +01:00
|
|
|
for(int bin = 0; bin < mbins; bin++) {
|
|
|
|
int c = bincount[bin];
|
|
|
|
int *bin_ptr = &bins[bin * atoms_per_bin];
|
|
|
|
|
|
|
|
for(int ac_i = 0; ac_i < c; ac_i++) {
|
|
|
|
int i = bin_ptr[ac_i];
|
|
|
|
int min_ac = ac_i;
|
|
|
|
int min_idx = i;
|
|
|
|
MD_FLOAT min_z = atom_z(i);
|
|
|
|
|
|
|
|
for(int ac_j = ac_i + 1; ac_j < c; ac_j++) {
|
|
|
|
int j = bin_ptr[ac_j];
|
|
|
|
MD_FLOAT zj = atom_z(j);
|
|
|
|
if(zj < min_z) {
|
2022-01-31 17:49:22 +01:00
|
|
|
min_ac = ac_j;
|
2022-01-27 03:07:31 +01:00
|
|
|
min_idx = j;
|
|
|
|
min_z = zj;
|
|
|
|
}
|
2022-01-25 00:43:10 +01:00
|
|
|
}
|
2021-11-03 00:57:24 +01:00
|
|
|
|
2022-01-27 03:07:31 +01:00
|
|
|
bin_ptr[ac_i] = min_idx;
|
|
|
|
bin_ptr[min_ac] = i;
|
|
|
|
}
|
2021-11-03 00:57:24 +01:00
|
|
|
}
|
2022-02-01 20:16:04 +01:00
|
|
|
|
|
|
|
DEBUG_MESSAGE("sortAtomsByZCoord end\n");
|
2022-01-25 00:43:10 +01:00
|
|
|
}
|
|
|
|
|
2022-01-27 03:07:31 +01:00
|
|
|
void buildClusters(Atom *atom) {
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("buildClusters start\n");
|
2022-01-27 03:07:31 +01:00
|
|
|
atom->Nclusters_local = 0;
|
|
|
|
|
2022-01-25 21:00:11 +01:00
|
|
|
/* bin local atoms */
|
2022-01-27 03:07:31 +01:00
|
|
|
binAtoms(atom);
|
|
|
|
sortAtomsByZCoord(atom);
|
2022-01-25 21:00:11 +01:00
|
|
|
|
2022-01-25 00:43:10 +01:00
|
|
|
for(int bin = 0; bin < mbins; bin++) {
|
2022-01-27 03:07:31 +01:00
|
|
|
int c = bincount[bin];
|
|
|
|
int ac = 0;
|
2022-03-09 02:25:39 +01:00
|
|
|
int nclusters = ((c + CLUSTER_M - 1) / CLUSTER_M);
|
|
|
|
if(CLUSTER_N > CLUSTER_M && nclusters % 2) { nclusters++; }
|
2022-01-27 03:07:31 +01:00
|
|
|
for(int cl = 0; cl < nclusters; cl++) {
|
|
|
|
const int ci = atom->Nclusters_local;
|
|
|
|
if(ci >= atom->Nclusters_max) {
|
|
|
|
growClusters(atom);
|
|
|
|
}
|
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
int ci_sca_base = CI_SCALAR_BASE_INDEX(ci);
|
|
|
|
int ci_vec_base = CI_VECTOR_BASE_INDEX(ci);
|
|
|
|
MD_FLOAT *ci_x = &atom->cl_x[ci_vec_base];
|
|
|
|
MD_FLOAT *ci_v = &atom->cl_v[ci_vec_base];
|
|
|
|
int *ci_type = &atom->cl_type[ci_sca_base];
|
2022-01-27 03:07:31 +01:00
|
|
|
MD_FLOAT bbminx = INFINITY, bbmaxx = -INFINITY;
|
|
|
|
MD_FLOAT bbminy = INFINITY, bbmaxy = -INFINITY;
|
|
|
|
MD_FLOAT bbminz = INFINITY, bbmaxz = -INFINITY;
|
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
atom->iclusters[ci].natoms = 0;
|
|
|
|
for(int cii = 0; cii < CLUSTER_M; cii++) {
|
2022-01-27 03:07:31 +01:00
|
|
|
if(ac < c) {
|
|
|
|
int i = bins[bin * atoms_per_bin + ac];
|
|
|
|
MD_FLOAT xtmp = atom_x(i);
|
|
|
|
MD_FLOAT ytmp = atom_y(i);
|
|
|
|
MD_FLOAT ztmp = atom_z(i);
|
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
ci_x[CL_X_OFFSET + cii] = xtmp;
|
|
|
|
ci_x[CL_Y_OFFSET + cii] = ytmp;
|
|
|
|
ci_x[CL_Z_OFFSET + cii] = ztmp;
|
|
|
|
ci_v[CL_X_OFFSET + cii] = atom->vx[i];
|
|
|
|
ci_v[CL_Y_OFFSET + cii] = atom->vy[i];
|
|
|
|
ci_v[CL_Z_OFFSET + cii] = atom->vz[i];
|
2022-01-27 03:07:31 +01:00
|
|
|
|
|
|
|
// TODO: To create the bounding boxes faster, we can use SIMD operations
|
|
|
|
if(bbminx > xtmp) { bbminx = xtmp; }
|
|
|
|
if(bbmaxx < xtmp) { bbmaxx = xtmp; }
|
|
|
|
if(bbminy > ytmp) { bbminy = ytmp; }
|
|
|
|
if(bbmaxy < ytmp) { bbmaxy = ytmp; }
|
2022-01-31 17:49:22 +01:00
|
|
|
if(bbminz > ztmp) { bbminz = ztmp; }
|
|
|
|
if(bbmaxz < ztmp) { bbmaxz = ztmp; }
|
2022-01-27 03:07:31 +01:00
|
|
|
|
2022-03-10 01:31:50 +01:00
|
|
|
ci_type[cii] = atom->type[i];
|
2022-03-09 02:25:39 +01:00
|
|
|
atom->iclusters[ci].natoms++;
|
2022-01-27 03:07:31 +01:00
|
|
|
} else {
|
2022-03-09 02:25:39 +01:00
|
|
|
ci_x[CL_X_OFFSET + cii] = INFINITY;
|
|
|
|
ci_x[CL_Y_OFFSET + cii] = INFINITY;
|
|
|
|
ci_x[CL_Z_OFFSET + cii] = INFINITY;
|
2022-01-27 03:07:31 +01:00
|
|
|
}
|
|
|
|
|
|
|
|
ac++;
|
|
|
|
}
|
|
|
|
|
2022-03-10 01:31:50 +01:00
|
|
|
atom->icluster_bin[ci] = bin;
|
2022-03-09 02:25:39 +01:00
|
|
|
atom->iclusters[ci].bbminx = bbminx;
|
|
|
|
atom->iclusters[ci].bbmaxx = bbmaxx;
|
|
|
|
atom->iclusters[ci].bbminy = bbminy;
|
|
|
|
atom->iclusters[ci].bbmaxy = bbmaxy;
|
|
|
|
atom->iclusters[ci].bbminz = bbminz;
|
|
|
|
atom->iclusters[ci].bbmaxz = bbmaxz;
|
2022-01-27 03:07:31 +01:00
|
|
|
atom->Nclusters_local++;
|
|
|
|
}
|
2022-01-25 00:43:10 +01:00
|
|
|
}
|
2022-01-31 17:49:22 +01:00
|
|
|
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("buildClusters end\n");
|
2022-01-31 23:46:20 +01:00
|
|
|
}
|
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
void defineJClusters(Atom *atom) {
|
|
|
|
DEBUG_MESSAGE("defineJClusters start\n");
|
|
|
|
|
2022-03-10 22:33:41 +01:00
|
|
|
for(int ci = 0; ci < atom->Nclusters_local; ci++) {
|
2022-03-09 02:25:39 +01:00
|
|
|
int cj0 = CJ0_FROM_CI(ci);
|
|
|
|
|
|
|
|
if(CLUSTER_M == CLUSTER_N) {
|
|
|
|
atom->jclusters[cj0].bbminx = atom->iclusters[ci].bbminx;
|
|
|
|
atom->jclusters[cj0].bbmaxx = atom->iclusters[ci].bbmaxx;
|
|
|
|
atom->jclusters[cj0].bbminy = atom->iclusters[ci].bbminy;
|
|
|
|
atom->jclusters[cj0].bbmaxy = atom->iclusters[ci].bbmaxy;
|
|
|
|
atom->jclusters[cj0].bbminz = atom->iclusters[ci].bbminz;
|
|
|
|
atom->jclusters[cj0].bbmaxz = atom->iclusters[ci].bbmaxz;
|
|
|
|
atom->jclusters[cj0].natoms = atom->iclusters[ci].natoms;
|
|
|
|
|
|
|
|
} else if(CLUSTER_M > CLUSTER_N) {
|
|
|
|
int cj1 = CJ1_FROM_CI(ci);
|
|
|
|
int ci_vec_base = CI_VECTOR_BASE_INDEX(ci);
|
|
|
|
MD_FLOAT *ci_x = &atom->cl_x[ci_vec_base];
|
|
|
|
MD_FLOAT bbminx = INFINITY, bbmaxx = -INFINITY;
|
|
|
|
MD_FLOAT bbminy = INFINITY, bbmaxy = -INFINITY;
|
|
|
|
MD_FLOAT bbminz = INFINITY, bbmaxz = -INFINITY;
|
|
|
|
|
|
|
|
for(int cii = 0; cii < MAX(atom->iclusters[ci].natoms, CLUSTER_N); cii++) {
|
|
|
|
MD_FLOAT xtmp = ci_x[CL_X_OFFSET + cii];
|
|
|
|
MD_FLOAT ytmp = ci_x[CL_Y_OFFSET + cii];
|
|
|
|
MD_FLOAT ztmp = ci_x[CL_Z_OFFSET + cii];
|
|
|
|
|
|
|
|
// TODO: To create the bounding boxes faster, we can use SIMD operations
|
|
|
|
if(bbminx > xtmp) { bbminx = xtmp; }
|
|
|
|
if(bbmaxx < xtmp) { bbmaxx = xtmp; }
|
|
|
|
if(bbminy > ytmp) { bbminy = ytmp; }
|
|
|
|
if(bbmaxy < ytmp) { bbmaxy = ytmp; }
|
|
|
|
if(bbminz > ztmp) { bbminz = ztmp; }
|
|
|
|
if(bbmaxz < ztmp) { bbmaxz = ztmp; }
|
|
|
|
}
|
|
|
|
|
|
|
|
atom->jclusters[cj0].bbminx = bbminx;
|
|
|
|
atom->jclusters[cj0].bbmaxx = bbmaxx;
|
|
|
|
atom->jclusters[cj0].bbminy = bbminy;
|
|
|
|
atom->jclusters[cj0].bbmaxy = bbmaxy;
|
|
|
|
atom->jclusters[cj0].bbminz = bbminz;
|
|
|
|
atom->jclusters[cj0].bbmaxz = bbmaxz;
|
|
|
|
atom->jclusters[cj0].natoms = MAX(atom->iclusters[ci].natoms, CLUSTER_N);
|
|
|
|
|
|
|
|
bbminx = INFINITY, bbmaxx = -INFINITY;
|
|
|
|
bbminy = INFINITY, bbmaxy = -INFINITY;
|
|
|
|
bbminz = INFINITY, bbmaxz = -INFINITY;
|
|
|
|
|
|
|
|
for(int cii = CLUSTER_N; cii < atom->iclusters[ci].natoms; cii++) {
|
|
|
|
MD_FLOAT xtmp = ci_x[CL_X_OFFSET + cii];
|
|
|
|
MD_FLOAT ytmp = ci_x[CL_Y_OFFSET + cii];
|
|
|
|
MD_FLOAT ztmp = ci_x[CL_Z_OFFSET + cii];
|
|
|
|
|
|
|
|
// TODO: To create the bounding boxes faster, we can use SIMD operations
|
|
|
|
if(bbminx > xtmp) { bbminx = xtmp; }
|
|
|
|
if(bbmaxx < xtmp) { bbmaxx = xtmp; }
|
|
|
|
if(bbminy > ytmp) { bbminy = ytmp; }
|
|
|
|
if(bbmaxy < ytmp) { bbmaxy = ytmp; }
|
|
|
|
if(bbminz > ztmp) { bbminz = ztmp; }
|
|
|
|
if(bbmaxz < ztmp) { bbmaxz = ztmp; }
|
|
|
|
}
|
|
|
|
|
|
|
|
atom->jclusters[cj1].bbminx = bbminx;
|
|
|
|
atom->jclusters[cj1].bbmaxx = bbmaxx;
|
|
|
|
atom->jclusters[cj1].bbminy = bbminy;
|
|
|
|
atom->jclusters[cj1].bbmaxy = bbmaxy;
|
|
|
|
atom->jclusters[cj1].bbminz = bbminz;
|
|
|
|
atom->jclusters[cj1].bbmaxz = bbmaxz;
|
|
|
|
atom->jclusters[cj1].natoms = MIN(0, atom->iclusters[ci].natoms - CLUSTER_N);
|
|
|
|
|
|
|
|
} else {
|
|
|
|
if(ci % 2 == 0) {
|
2022-03-10 22:33:41 +01:00
|
|
|
const int ci1 = ci + 1;
|
|
|
|
atom->jclusters[cj0].bbminx = MIN(atom->iclusters[ci].bbminx, atom->iclusters[ci1].bbminx);
|
|
|
|
atom->jclusters[cj0].bbmaxx = MAX(atom->iclusters[ci].bbmaxx, atom->iclusters[ci1].bbmaxx);
|
|
|
|
atom->jclusters[cj0].bbminy = MIN(atom->iclusters[ci].bbminy, atom->iclusters[ci1].bbminy);
|
|
|
|
atom->jclusters[cj0].bbmaxy = MAX(atom->iclusters[ci].bbmaxy, atom->iclusters[ci1].bbmaxy);
|
|
|
|
atom->jclusters[cj0].bbminz = MIN(atom->iclusters[ci].bbminz, atom->iclusters[ci1].bbminz);
|
|
|
|
atom->jclusters[cj0].bbmaxz = MAX(atom->iclusters[ci].bbmaxz, atom->iclusters[ci1].bbmaxz);
|
|
|
|
atom->jclusters[cj0].natoms = atom->iclusters[ci].natoms + atom->iclusters[ci1].natoms;
|
2022-03-09 02:25:39 +01:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
DEBUG_MESSAGE("defineJClusters end\n");
|
|
|
|
}
|
|
|
|
|
2022-01-31 23:46:20 +01:00
|
|
|
void binClusters(Atom *atom) {
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("binClusters start\n");
|
2022-01-31 23:46:20 +01:00
|
|
|
|
2022-01-31 17:49:22 +01:00
|
|
|
/*
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("Nghost = %d\n", atom->Nclusters_ghost);
|
2022-01-31 23:46:20 +01:00
|
|
|
for(int ci = atom->Nclusters_local; ci < atom->Nclusters_local + 4; ci++) {
|
2022-01-31 17:49:22 +01:00
|
|
|
MD_FLOAT *cptr = cluster_pos_ptr(ci);
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("Cluster %d:\n", ci);
|
|
|
|
DEBUG_MESSAGE("bin=%d, Natoms=%d, bbox={%f,%f},{%f,%f},{%f,%f}\n",
|
2022-03-10 01:31:50 +01:00
|
|
|
atom->icluster_bin[ci],
|
2022-02-01 20:16:04 +01:00
|
|
|
atom->clusters[ci].natoms,
|
|
|
|
atom->clusters[ci].bbminx,
|
|
|
|
atom->clusters[ci].bbmaxx,
|
|
|
|
atom->clusters[ci].bbminy,
|
|
|
|
atom->clusters[ci].bbmaxy,
|
|
|
|
atom->clusters[ci].bbminz,
|
|
|
|
atom->clusters[ci].bbmaxz);
|
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
for(int cii = 0; cii < CLUSTER_M; cii++) {
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("%f, %f, %f\n", cluster_x(cptr, cii), cluster_y(cptr, cii), cluster_z(cptr, cii));
|
2022-01-31 17:49:22 +01:00
|
|
|
}
|
|
|
|
}
|
|
|
|
*/
|
|
|
|
|
2022-01-27 03:07:31 +01:00
|
|
|
const int nlocal = atom->Nclusters_local;
|
2022-03-10 22:33:41 +01:00
|
|
|
const int jfac = MAX(1, CLUSTER_N / CLUSTER_M);
|
|
|
|
const int ncj = atom->Nclusters_local / jfac;
|
2022-01-27 03:07:31 +01:00
|
|
|
|
2022-03-09 02:25:39 +01:00
|
|
|
int resize = 1;
|
2022-01-25 00:43:10 +01:00
|
|
|
while(resize > 0) {
|
|
|
|
resize = 0;
|
2022-01-27 03:07:31 +01:00
|
|
|
|
2022-01-25 00:43:10 +01:00
|
|
|
for(int bin = 0; bin < mbins; bin++) {
|
2022-03-09 02:25:39 +01:00
|
|
|
bin_nclusters[bin] = 0;
|
2022-01-27 03:07:31 +01:00
|
|
|
}
|
2022-01-25 00:43:10 +01:00
|
|
|
|
2022-01-27 03:07:31 +01:00
|
|
|
for(int ci = 0; ci < nlocal && !resize; ci++) {
|
2022-03-09 02:25:39 +01:00
|
|
|
// Assure we add this j-cluster only once in the bin
|
|
|
|
if(!(CLUSTER_M < CLUSTER_N && ci % 2)) {
|
2022-03-10 01:31:50 +01:00
|
|
|
int bin = atom->icluster_bin[ci];
|
2022-03-09 02:25:39 +01:00
|
|
|
int c = bin_nclusters[bin];
|
|
|
|
if(c + 1 < clusters_per_bin) {
|
|
|
|
bin_clusters[bin * clusters_per_bin + c] = CJ0_FROM_CI(ci);
|
|
|
|
bin_nclusters[bin]++;
|
|
|
|
|
|
|
|
if(CLUSTER_M > CLUSTER_N) {
|
|
|
|
int cj1 = CJ1_FROM_CI(ci);
|
2022-03-09 17:23:49 +01:00
|
|
|
if(atom->jclusters[cj1].natoms > 0) {
|
2022-03-09 02:25:39 +01:00
|
|
|
bin_clusters[bin * clusters_per_bin + c + 1] = cj1;
|
|
|
|
bin_nclusters[bin]++;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
} else {
|
|
|
|
resize = 1;
|
|
|
|
}
|
2022-01-27 03:07:31 +01:00
|
|
|
}
|
|
|
|
}
|
2022-01-25 00:43:10 +01:00
|
|
|
|
2022-01-31 17:49:22 +01:00
|
|
|
for(int cg = 0; cg < atom->Nclusters_ghost && !resize; cg++) {
|
2022-03-10 22:33:41 +01:00
|
|
|
const int cj = ncj + cg;
|
|
|
|
int ix = -1, iy = -1;
|
|
|
|
MD_FLOAT xtmp, ytmp;
|
|
|
|
|
|
|
|
if(atom->jclusters[cj].natoms > 0) {
|
|
|
|
int cj_vec_base = CJ_VECTOR_BASE_INDEX(cj);
|
|
|
|
MD_FLOAT *cj_x = &atom->cl_x[cj_vec_base];
|
|
|
|
MD_FLOAT cj_minz = atom->jclusters[cj].bbminz;
|
|
|
|
|
|
|
|
xtmp = cj_x[CL_X_OFFSET + 0];
|
|
|
|
ytmp = cj_x[CL_Y_OFFSET + 0];
|
|
|
|
coord2bin2D(xtmp, ytmp, &ix, &iy);
|
|
|
|
ix = MAX(MIN(ix, mbinx - 1), 0);
|
|
|
|
iy = MAX(MIN(iy, mbiny - 1), 0);
|
|
|
|
for(int cjj = 1; cjj < atom->jclusters[cj].natoms; cjj++) {
|
|
|
|
int nix, niy;
|
|
|
|
xtmp = cj_x[CL_X_OFFSET + cjj];
|
|
|
|
ytmp = cj_x[CL_Y_OFFSET + cjj];
|
|
|
|
coord2bin2D(xtmp, ytmp, &nix, &niy);
|
|
|
|
nix = MAX(MIN(nix, mbinx - 1), 0);
|
|
|
|
niy = MAX(MIN(niy, mbiny - 1), 0);
|
|
|
|
|
|
|
|
// Always put the cluster on the bin of its innermost atom so
|
|
|
|
// the cluster should be closer to local clusters
|
|
|
|
if(atom->PBCx[cg] > 0 && ix > nix) { ix = nix; }
|
|
|
|
if(atom->PBCx[cg] < 0 && ix < nix) { ix = nix; }
|
|
|
|
if(atom->PBCy[cg] > 0 && iy > niy) { iy = niy; }
|
|
|
|
if(atom->PBCy[cg] < 0 && iy < niy) { iy = niy; }
|
|
|
|
}
|
2022-02-07 18:00:21 +01:00
|
|
|
|
2022-03-10 22:33:41 +01:00
|
|
|
int bin = iy * mbinx + ix + 1;
|
|
|
|
int c = bin_nclusters[bin];
|
|
|
|
if(c < clusters_per_bin) {
|
|
|
|
// Insert the current ghost cluster in the bin keeping clusters
|
|
|
|
// sorted by z coordinate
|
|
|
|
int inserted = 0;
|
|
|
|
for(int i = 0; i < c; i++) {
|
|
|
|
int last_cl = bin_clusters[bin * clusters_per_bin + i];
|
|
|
|
if(atom->jclusters[last_cl].bbminz > cj_minz) {
|
|
|
|
bin_clusters[bin * clusters_per_bin + i] = cj;
|
|
|
|
|
|
|
|
for(int j = i + 1; j <= c; j++) {
|
|
|
|
int tmp = bin_clusters[bin * clusters_per_bin + j];
|
|
|
|
bin_clusters[bin * clusters_per_bin + j] = last_cl;
|
|
|
|
last_cl = tmp;
|
2022-03-09 02:25:39 +01:00
|
|
|
}
|
2022-02-07 18:00:21 +01:00
|
|
|
|
2022-03-10 22:33:41 +01:00
|
|
|
inserted = 1;
|
|
|
|
break;
|
2022-03-09 02:25:39 +01:00
|
|
|
}
|
|
|
|
}
|
2022-03-10 22:33:41 +01:00
|
|
|
|
|
|
|
if(!inserted) {
|
|
|
|
bin_clusters[bin * clusters_per_bin + c] = cj;
|
|
|
|
}
|
|
|
|
|
|
|
|
bin_nclusters[bin]++;
|
|
|
|
} else {
|
|
|
|
resize = 1;
|
2022-03-09 02:25:39 +01:00
|
|
|
}
|
2022-01-25 00:43:10 +01:00
|
|
|
}
|
|
|
|
}
|
2021-11-03 00:57:24 +01:00
|
|
|
|
2022-01-25 00:43:10 +01:00
|
|
|
if(resize) {
|
2022-03-09 02:25:39 +01:00
|
|
|
free(bin_clusters);
|
2022-01-25 00:43:10 +01:00
|
|
|
clusters_per_bin *= 2;
|
2022-03-09 02:25:39 +01:00
|
|
|
bin_clusters = (int*) malloc(mbins * clusters_per_bin * sizeof(int));
|
2021-11-03 00:57:24 +01:00
|
|
|
}
|
|
|
|
}
|
2022-01-31 17:49:22 +01:00
|
|
|
|
|
|
|
/*
|
2022-03-09 02:25:39 +01:00
|
|
|
DEBUG_MESSAGE("bin_nclusters\n");
|
|
|
|
for(int i = 0; i < mbins; i++) { DEBUG_MESSAGE("%d, ", bin_nclusters[i]); }
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("\n");
|
2022-01-31 17:49:22 +01:00
|
|
|
*/
|
|
|
|
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("binClusters stop\n");
|
2022-01-25 00:43:10 +01:00
|
|
|
}
|
|
|
|
|
2022-01-27 03:07:31 +01:00
|
|
|
void updateSingleAtoms(Atom *atom) {
|
2022-02-01 20:16:04 +01:00
|
|
|
DEBUG_MESSAGE("updateSingleAtoms start\n");
|
2022-01-27 03:07:31 +01:00
|
|
|
int Natom = 0;
|
|
|
|
|
|
|
|
for(int ci = 0; ci < atom->Nclusters_local; ci++) {
|
2022-03-09 02:25:39 +01:00
|
|
|
int ci_vec_base = CI_VECTOR_BASE_INDEX(ci);
|
|
|
|
MD_FLOAT *ci_x = &atom->cl_x[ci_vec_base];
|
|
|
|
MD_FLOAT *ci_v = &atom->cl_v[ci_vec_base];
|
2022-01-27 03:07:31 +01:00
|
|
|
|
2022-03-09 17:23:49 +01:00
|
|
|
for(int cii = 0; cii < atom->iclusters[ci].natoms; cii++) {
|
2022-03-09 02:25:39 +01:00
|
|
|
atom_x(Natom) = ci_x[CL_X_OFFSET + cii];
|
|
|
|
atom_y(Natom) = ci_x[CL_Y_OFFSET + cii];
|
|
|
|
atom_z(Natom) = ci_x[CL_Z_OFFSET + cii];
|
|
|
|
atom->vx[Natom] = ci_v[CL_X_OFFSET + cii];
|
|
|
|
atom->vy[Natom] = ci_v[CL_Y_OFFSET + cii];
|
|
|
|
atom->vz[Natom] = ci_v[CL_Z_OFFSET + cii];
|
2022-01-27 03:07:31 +01:00
|
|
|
Natom++;
|
2022-01-25 00:43:10 +01:00
|
|
|
}
|
2022-01-27 03:07:31 +01:00
|
|
|
}
|
2021-11-03 00:57:24 +01:00
|
|
|
|
2022-01-27 03:07:31 +01:00
|
|
|
if(Natom != atom->Nlocal) {
|
|
|
|
fprintf(stderr, "updateSingleAtoms(): Number of atoms changed!\n");
|
2022-01-25 00:43:10 +01:00
|
|
|
}
|
2022-02-01 20:16:04 +01:00
|
|
|
|
|
|
|
DEBUG_MESSAGE("updateSingleAtoms stop\n");
|
2021-11-03 00:57:24 +01:00
|
|
|
}
|