update grid2d to match grid3d
This commit is contained in:
@ -597,12 +597,12 @@ void MSM::allocate()
|
|||||||
// commgrid using all processors for finest grid level
|
// commgrid using all processors for finest grid level
|
||||||
|
|
||||||
gcall = new Grid3d(lmp,world,1,nx_msm[0],ny_msm[0],nz_msm[0],
|
gcall = new Grid3d(lmp,world,1,nx_msm[0],ny_msm[0],nz_msm[0],
|
||||||
nxlo_in[0],nxhi_in[0],nylo_in[0],
|
nxlo_in[0],nxhi_in[0],nylo_in[0],
|
||||||
nyhi_in[0],nzlo_in[0],nzhi_in[0],
|
nyhi_in[0],nzlo_in[0],nzhi_in[0],
|
||||||
nxlo_out_all,nxhi_out_all,nylo_out_all,
|
nxlo_out_all,nxhi_out_all,nylo_out_all,
|
||||||
nyhi_out_all,nzlo_out_all,nzhi_out_all,
|
nyhi_out_all,nzlo_out_all,nzhi_out_all,
|
||||||
nxlo_out[0],nxhi_out[0],nylo_out[0],
|
nxlo_out[0],nxhi_out[0],nylo_out[0],
|
||||||
nyhi_out[0],nzlo_out[0],nzhi_out[0]);
|
nyhi_out[0],nzlo_out[0],nzhi_out[0]);
|
||||||
|
|
||||||
gcall->setup_comm(ngcall_buf1,ngcall_buf2);
|
gcall->setup_comm(ngcall_buf1,ngcall_buf2);
|
||||||
npergrid = 1;
|
npergrid = 1;
|
||||||
|
|||||||
@ -242,8 +242,9 @@ void *ComputePropertyGrid::get_griddata_by_index(int index)
|
|||||||
void ComputePropertyGrid::allocate_grid()
|
void ComputePropertyGrid::allocate_grid()
|
||||||
{
|
{
|
||||||
if (dimension == 2) {
|
if (dimension == 2) {
|
||||||
grid2d = new Grid2d(lmp, world, nxgrid, nygrid, 0.0, 0, 0.0, nxlo_in, nxhi_in, nylo_in, nyhi_in,
|
grid2d = new Grid2d(lmp, world, nxgrid, nygrid);
|
||||||
nxlo_out, nxhi_out, nylo_out, nyhi_out);
|
grid2d->setup_grid(nxlo_in, nxhi_in, nylo_in, nyhi_in, nxlo_out, nxhi_out, nylo_out, nyhi_out);
|
||||||
|
|
||||||
if (nvalues == 1)
|
if (nvalues == 1)
|
||||||
memory->create2d_offset(vec2d, nylo_out, nyhi_out, nxlo_out, nxhi_out, "property/grid:vec2d");
|
memory->create2d_offset(vec2d, nylo_out, nyhi_out, nxlo_out, nxhi_out, "property/grid:vec2d");
|
||||||
else
|
else
|
||||||
|
|||||||
@ -503,7 +503,7 @@ int DumpGrid::count()
|
|||||||
grid2d = (Grid2d *) compute[field2index[0]]->get_grid_by_index(field2grid[0]);
|
grid2d = (Grid2d *) compute[field2index[0]]->get_grid_by_index(field2grid[0]);
|
||||||
else if (field2source[0] == FIX)
|
else if (field2source[0] == FIX)
|
||||||
grid2d = (Grid2d *) fix[field2index[0]]->get_grid_by_index(field2grid[0]);
|
grid2d = (Grid2d *) fix[field2index[0]]->get_grid_by_index(field2grid[0]);
|
||||||
grid2d->get_bounds(nxlo_in,nxhi_in,nylo_in,nyhi_in);
|
grid2d->get_bounds_owned(nxlo_in,nxhi_in,nylo_in,nyhi_in);
|
||||||
} else {
|
} else {
|
||||||
if (field2source[0] == COMPUTE)
|
if (field2source[0] == COMPUTE)
|
||||||
grid3d = (Grid3d *) compute[field2index[0]]->get_grid_by_index(field2grid[0]);
|
grid3d = (Grid3d *) compute[field2index[0]]->get_grid_by_index(field2grid[0]);
|
||||||
|
|||||||
@ -395,9 +395,10 @@ FixAveGrid::FixAveGrid(LAMMPS *lmp, int narg, char **arg) :
|
|||||||
else if (modegrid) maxdist = 0.0;
|
else if (modegrid) maxdist = 0.0;
|
||||||
|
|
||||||
if (dimension == 2) {
|
if (dimension == 2) {
|
||||||
grid2d = new Grid2d(lmp, world, nxgrid, nygrid, maxdist, 0, 0.5,
|
grid2d = new Grid2d(lmp, world, nxgrid, nygrid);
|
||||||
nxlo_in, nxhi_in, nylo_in, nyhi_in,
|
grid2d->set_distance(maxdist);
|
||||||
nxlo_out, nxhi_out, nylo_out, nyhi_out);
|
grid2d->setup_grid(nxlo_in, nxhi_in, nylo_in, nyhi_in,
|
||||||
|
nxlo_out, nxhi_out, nylo_out, nyhi_out);
|
||||||
|
|
||||||
// ngrid_buf12 converted to nvalues + count
|
// ngrid_buf12 converted to nvalues + count
|
||||||
|
|
||||||
|
|||||||
327
src/grid2d.cpp
327
src/grid2d.cpp
@ -40,118 +40,26 @@ static constexpr int OFFSET = 16384;
|
|||||||
------------------------------------------------------------------------- */
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
/* ----------------------------------------------------------------------
|
/* ----------------------------------------------------------------------
|
||||||
constructor called by all classes except PPPM and MSM
|
constructor to create a 3d distributed grid
|
||||||
comm_caller = caller's communicator
|
gcomm = caller's communicator
|
||||||
nx,ny caller = size of global grid
|
gnx,gny = global grid size
|
||||||
maxdist = max distance outside of proc domain a particle will be
|
|
||||||
extra = additional ghost grid pts needed in each dim, e.g. for stencil
|
|
||||||
shift_caller = 0.0 for grid pt in lower-left corner of grid cell,
|
|
||||||
0.5 for center
|
|
||||||
return:
|
|
||||||
i xy lohi = portion of global grid this proc owns, 0 <= index < N
|
|
||||||
o xy lohi = owned + ghost grid cells needed in all directions
|
|
||||||
for periodic dims, o indices can be < 0 or >= N
|
|
||||||
for non-periodic dims, o indices will be >= 0 and < N
|
|
||||||
since no grid comm is done across non-periodic boundaries
|
|
||||||
NOTE: allow zfactor to be a calling arg for PPPM ?
|
|
||||||
------------------------------------------------------------------------- */
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
Grid2d::Grid2d(LAMMPS *lmp, MPI_Comm comm_caller,
|
Grid2d::Grid2d(LAMMPS *lmp, MPI_Comm gcomm, int gnx, int gny) : Pointers(lmp)
|
||||||
int nx_caller, int ny_caller,
|
|
||||||
double maxdist, int extra, double shift_caller,
|
|
||||||
int &ixlo, int &ixhi, int &iylo, int &iyhi,
|
|
||||||
int &oxlo, int &oxhi, int &oylo, int &oyhi)
|
|
||||||
: Pointers(lmp)
|
|
||||||
{
|
{
|
||||||
// store commnicator and global grid size
|
gridcomm = gcomm;
|
||||||
// set layout mode
|
|
||||||
|
|
||||||
gridcomm = comm_caller;
|
|
||||||
MPI_Comm_rank(gridcomm,&me);
|
MPI_Comm_rank(gridcomm,&me);
|
||||||
MPI_Comm_size(gridcomm,&nprocs);
|
MPI_Comm_size(gridcomm,&nprocs);
|
||||||
|
|
||||||
nx = nx_caller;
|
nx = gnx;
|
||||||
ny = ny_caller;
|
ny = gny;
|
||||||
|
|
||||||
// NOTE: hardwire shift = 0.5 ?
|
// default settings, can be overridden by set() methods
|
||||||
|
|
||||||
shift = shift_caller;
|
|
||||||
|
|
||||||
// define owned grid cells for each proc
|
maxdist = 0.0;
|
||||||
// extend bounds with ghost grid cells in each direction
|
stencil_grid_lo = stencil_grid_hi = 0;
|
||||||
|
stencil_atom_lo = stencil_atom_hi = 0;
|
||||||
double fraclo,frachi;
|
shift_grid = shift_atom = 0.0;
|
||||||
|
|
||||||
if (comm->layout != Comm::LAYOUT_TILED) {
|
|
||||||
fraclo = comm->xsplit[comm->myloc[0]];
|
|
||||||
frachi = comm->xsplit[comm->myloc[0]+1];
|
|
||||||
partition_grid(nx,fraclo,frachi,shift,inxlo,inxhi);
|
|
||||||
fraclo = comm->ysplit[comm->myloc[1]];
|
|
||||||
frachi = comm->ysplit[comm->myloc[1]+1];
|
|
||||||
partition_grid(ny,fraclo,frachi,shift,inylo,inyhi);
|
|
||||||
} else {
|
|
||||||
fraclo = comm->mysplit[0][0];
|
|
||||||
frachi = comm->mysplit[0][1];
|
|
||||||
partition_grid(nx,fraclo,frachi,shift,inxlo,inxhi);
|
|
||||||
fraclo = comm->mysplit[1][0];
|
|
||||||
frachi = comm->mysplit[1][1];
|
|
||||||
partition_grid(ny,fraclo,frachi,shift,inylo,inyhi);
|
|
||||||
}
|
|
||||||
|
|
||||||
ghost_grid(maxdist,extra);
|
|
||||||
|
|
||||||
// error check on size of grid stored by this proc
|
|
||||||
|
|
||||||
bigint total = (bigint) (outxhi - outxlo + 1) * (outyhi - outylo + 1);
|
|
||||||
if (total > MAXSMALLINT)
|
|
||||||
error->one(FLERR, "Too many owned+ghost grid3d points");
|
|
||||||
|
|
||||||
// default = caller grid is allocated to ghost grid
|
|
||||||
// used when computing pack/unpack lists in indices()
|
|
||||||
// NOTE: allow caller to override this
|
|
||||||
|
|
||||||
fullxlo = outxlo;
|
|
||||||
fullxhi = outxhi;
|
|
||||||
fullylo = outylo;
|
|
||||||
fullyhi = outyhi;
|
|
||||||
|
|
||||||
// initialize data structs
|
|
||||||
|
|
||||||
nswap = maxswap = 0;
|
|
||||||
swap = nullptr;
|
|
||||||
|
|
||||||
nsend = nrecv = ncopy = 0;
|
|
||||||
send = nullptr;
|
|
||||||
recv = nullptr;
|
|
||||||
copy = nullptr;
|
|
||||||
requests = nullptr;
|
|
||||||
requests_remap = nullptr;
|
|
||||||
|
|
||||||
xsplit = ysplit = zsplit = nullptr;
|
|
||||||
grid2proc = nullptr;
|
|
||||||
rcbinfo = nullptr;
|
|
||||||
|
|
||||||
nsend_remap = nrecv_remap = self_remap = 0;
|
|
||||||
send_remap = nullptr;
|
|
||||||
recv_remap = nullptr;
|
|
||||||
|
|
||||||
// store info about Comm decomposition needed for remap operation
|
|
||||||
// two Grid instances will exist for duration of remap
|
|
||||||
// each must know Comm decomp at time Grid instance was created
|
|
||||||
|
|
||||||
extract_comm_info();
|
|
||||||
|
|
||||||
// return values
|
|
||||||
|
|
||||||
ixlo = inxlo;
|
|
||||||
ixhi = inxhi;
|
|
||||||
iylo = inylo;
|
|
||||||
iyhi = inyhi;
|
|
||||||
|
|
||||||
oxlo = outxlo;
|
|
||||||
oxhi = outxhi;
|
|
||||||
oylo = outylo;
|
|
||||||
oyhi = outyhi;
|
|
||||||
}
|
}
|
||||||
|
|
||||||
/* ----------------------------------------------------------------------
|
/* ----------------------------------------------------------------------
|
||||||
@ -172,7 +80,6 @@ Grid2d::Grid2d(LAMMPS *lmp, MPI_Comm gcomm,
|
|||||||
: Pointers(lmp)
|
: Pointers(lmp)
|
||||||
{
|
{
|
||||||
// store commnicator and global grid size
|
// store commnicator and global grid size
|
||||||
// set layout mode
|
|
||||||
|
|
||||||
gridcomm = gcomm;
|
gridcomm = gcomm;
|
||||||
MPI_Comm_rank(gridcomm,&me);
|
MPI_Comm_rank(gridcomm,&me);
|
||||||
@ -222,7 +129,6 @@ Grid2d::Grid2d(LAMMPS *lmp, MPI_Comm gcomm, int flag,
|
|||||||
: Pointers(lmp)
|
: Pointers(lmp)
|
||||||
{
|
{
|
||||||
// store commnicator and global grid size
|
// store commnicator and global grid size
|
||||||
// set layout mode
|
|
||||||
|
|
||||||
gridcomm = gcomm;
|
gridcomm = gcomm;
|
||||||
MPI_Comm_rank(gridcomm,&me);
|
MPI_Comm_rank(gridcomm,&me);
|
||||||
@ -239,8 +145,6 @@ Grid2d::Grid2d(LAMMPS *lmp, MPI_Comm gcomm, int flag,
|
|||||||
|
|
||||||
// store grid bounds and proc neighs
|
// store grid bounds and proc neighs
|
||||||
|
|
||||||
// store grid bounds and proc neighs
|
|
||||||
|
|
||||||
if (flag == 1) {
|
if (flag == 1) {
|
||||||
if (layout != Comm::LAYOUT_TILED) {
|
if (layout != Comm::LAYOUT_TILED) {
|
||||||
// this assumes gcomm = world
|
// this assumes gcomm = world
|
||||||
@ -340,7 +244,67 @@ void Grid2d::store(int ixlo, int ixhi, int iylo, int iyhi,
|
|||||||
}
|
}
|
||||||
|
|
||||||
// ----------------------------------------------------------------------
|
// ----------------------------------------------------------------------
|
||||||
// access Grid parameters
|
// set Grid parameters
|
||||||
|
// ----------------------------------------------------------------------
|
||||||
|
|
||||||
|
/* ----------------------------------------------------------------------
|
||||||
|
maxdist = max distance outside proc subdomain a particle can be
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
void Grid2d::set_distance(double distance)
|
||||||
|
{
|
||||||
|
maxdist = distance;
|
||||||
|
}
|
||||||
|
|
||||||
|
/* ----------------------------------------------------------------------
|
||||||
|
# of grid cells beyond an owned grid cell that caller accesses
|
||||||
|
e.g. for a finite different stencial
|
||||||
|
can be different in lo vs hi direction
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
void Grid2d::set_stencil_grid(int lo, int hi)
|
||||||
|
{
|
||||||
|
stencil_grid_lo = lo;
|
||||||
|
stencil_grid_hi = hi;
|
||||||
|
}
|
||||||
|
|
||||||
|
/* ----------------------------------------------------------------------
|
||||||
|
# of grid cells beyond a particle's grid cell that caller updates
|
||||||
|
e.g. for smearing a point charge to the grid
|
||||||
|
can be different in lo vs hi direction
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
void Grid2d::set_stencil_atom(int lo, int hi)
|
||||||
|
{
|
||||||
|
stencil_atom_lo = lo;
|
||||||
|
stencil_atom_hi = hi;
|
||||||
|
}
|
||||||
|
|
||||||
|
/* ----------------------------------------------------------------------
|
||||||
|
shift_grid = offset of position of grid point within grid cell
|
||||||
|
0.5 = cell center, 0.0 = lower-left corner of cell
|
||||||
|
used to determine which proc owns the grid cell within its subdomain
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
void Grid2d::set_shift_grid(double shift)
|
||||||
|
{
|
||||||
|
shift_grid = shift;
|
||||||
|
}
|
||||||
|
|
||||||
|
/* ----------------------------------------------------------------------
|
||||||
|
shift_atom = offset of atoms when caller maps them to grid cells
|
||||||
|
0.5 = half a grid cell, 0.0 = no offset
|
||||||
|
used when computing possible ghost extent
|
||||||
|
PPPM uses 0.5 when order is odd, 0.0 when order is even
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
void Grid2d::set_shift_atom(double shift)
|
||||||
|
{
|
||||||
|
shift_atom = shift;
|
||||||
|
}
|
||||||
|
|
||||||
|
// ----------------------------------------------------------------------
|
||||||
|
// retrieve Grid parameters
|
||||||
// ----------------------------------------------------------------------
|
// ----------------------------------------------------------------------
|
||||||
|
|
||||||
int Grid2d::identical(Grid2d *grid2)
|
int Grid2d::identical(Grid2d *grid2)
|
||||||
@ -348,7 +312,7 @@ int Grid2d::identical(Grid2d *grid2)
|
|||||||
int inxlo2,inxhi2,inylo2,inyhi2;
|
int inxlo2,inxhi2,inylo2,inyhi2;
|
||||||
int outxlo2,outxhi2,outylo2,outyhi2;
|
int outxlo2,outxhi2,outylo2,outyhi2;
|
||||||
|
|
||||||
grid2->get_bounds(inxlo2,inxhi2,inylo2,inyhi2);
|
grid2->get_bounds_owned(inxlo2,inxhi2,inylo2,inyhi2);
|
||||||
grid2->get_bounds_ghost(outxlo2,outxhi2,outylo2,outyhi2);
|
grid2->get_bounds_ghost(outxlo2,outxhi2,outylo2,outyhi2);
|
||||||
|
|
||||||
int flag = 0;
|
int flag = 0;
|
||||||
@ -374,7 +338,7 @@ void Grid2d::get_size(int &nxgrid, int &nygrid)
|
|||||||
|
|
||||||
/* ---------------------------------------------------------------------- */
|
/* ---------------------------------------------------------------------- */
|
||||||
|
|
||||||
void Grid2d::get_bounds(int &xlo, int &xhi, int &ylo, int &yhi)
|
void Grid2d::get_bounds_owned(int &xlo, int &xhi, int &ylo, int &yhi)
|
||||||
{
|
{
|
||||||
xlo = inxlo;
|
xlo = inxlo;
|
||||||
xhi = inxhi;
|
xhi = inxhi;
|
||||||
@ -396,7 +360,99 @@ void Grid2d::get_bounds_ghost(int &xlo, int &xhi, int &ylo, int &yhi)
|
|||||||
// define owned and ghost grid cells
|
// define owned and ghost grid cells
|
||||||
// also store comm and grid partitioning info
|
// also store comm and grid partitioning info
|
||||||
// ----------------------------------------------------------------------
|
// ----------------------------------------------------------------------
|
||||||
|
|
||||||
|
/* ----------------------------------------------------------------------
|
||||||
|
setup grid partition for each proc = owned + ghost cells
|
||||||
|
return:
|
||||||
|
i xy lohi = portion of global grid this proc owns, 0 <= index < N
|
||||||
|
o xy lohi = owned + ghost grid cells in all directions
|
||||||
|
for periodic dims, o indices can be < 0 or >= N
|
||||||
|
for non-periodic dims, o indices will be >= 0 and < N
|
||||||
|
since no grid comm is done across non-periodic boundaries
|
||||||
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
|
void Grid2d::setup_grid(int &ixlo, int &ixhi, int &iylo, int &iyhi,
|
||||||
|
int &oxlo, int &oxhi, int &oylo, int &oyhi)
|
||||||
|
{
|
||||||
|
// owned grid cells = those whose grid point is within proc subdomain
|
||||||
|
// shift_grid = 0.5 for grid point at cell center, 0.0 for lower-left corner
|
||||||
|
|
||||||
|
double fraclo,frachi;
|
||||||
|
|
||||||
|
if (comm->layout != Comm::LAYOUT_TILED) {
|
||||||
|
fraclo = comm->xsplit[comm->myloc[0]];
|
||||||
|
frachi = comm->xsplit[comm->myloc[0]+1];
|
||||||
|
partition_grid(nx,fraclo,frachi,shift_grid,inxlo,inxhi);
|
||||||
|
fraclo = comm->ysplit[comm->myloc[1]];
|
||||||
|
frachi = comm->ysplit[comm->myloc[1]+1];
|
||||||
|
partition_grid(ny,fraclo,frachi,shift_grid,inylo,inyhi);
|
||||||
|
} else {
|
||||||
|
fraclo = comm->mysplit[0][0];
|
||||||
|
frachi = comm->mysplit[0][1];
|
||||||
|
partition_grid(nx,fraclo,frachi,shift_grid,inxlo,inxhi);
|
||||||
|
fraclo = comm->mysplit[1][0];
|
||||||
|
frachi = comm->mysplit[1][1];
|
||||||
|
partition_grid(ny,fraclo,frachi,shift_grid,inylo,inyhi);
|
||||||
|
}
|
||||||
|
|
||||||
|
// extend owned grid bounds with ghost grid cells in each direction
|
||||||
|
|
||||||
|
ghost_grid();
|
||||||
|
|
||||||
|
// error check on size of grid stored by this proc
|
||||||
|
|
||||||
|
bigint total = (bigint) (outxhi - outxlo + 1) * (outyhi - outylo + 1);
|
||||||
|
if (total > MAXSMALLINT)
|
||||||
|
error->one(FLERR, "Too many owned+ghost grid3d points");
|
||||||
|
|
||||||
|
// default = caller grid is allocated to ghost grid
|
||||||
|
// used when computing pack/unpack lists in indices()
|
||||||
|
// NOTE: allow caller to override this
|
||||||
|
|
||||||
|
fullxlo = outxlo;
|
||||||
|
fullxhi = outxhi;
|
||||||
|
fullylo = outylo;
|
||||||
|
fullyhi = outyhi;
|
||||||
|
|
||||||
|
// initialize data structs
|
||||||
|
|
||||||
|
nswap = maxswap = 0;
|
||||||
|
swap = nullptr;
|
||||||
|
|
||||||
|
nsend = nrecv = ncopy = 0;
|
||||||
|
send = nullptr;
|
||||||
|
recv = nullptr;
|
||||||
|
copy = nullptr;
|
||||||
|
requests = nullptr;
|
||||||
|
requests_remap = nullptr;
|
||||||
|
|
||||||
|
xsplit = ysplit = zsplit = nullptr;
|
||||||
|
grid2proc = nullptr;
|
||||||
|
rcbinfo = nullptr;
|
||||||
|
|
||||||
|
nsend_remap = nrecv_remap = self_remap = 0;
|
||||||
|
send_remap = nullptr;
|
||||||
|
recv_remap = nullptr;
|
||||||
|
|
||||||
|
// store info about Comm decomposition needed for remap operation
|
||||||
|
// two Grid instances will exist for duration of remap
|
||||||
|
// each must know Comm decomp at time Grid instance was created
|
||||||
|
|
||||||
|
extract_comm_info();
|
||||||
|
|
||||||
|
// return values
|
||||||
|
|
||||||
|
ixlo = inxlo;
|
||||||
|
ixhi = inxhi;
|
||||||
|
iylo = inylo;
|
||||||
|
iyhi = inyhi;
|
||||||
|
|
||||||
|
oxlo = outxlo;
|
||||||
|
oxhi = outxhi;
|
||||||
|
oylo = outylo;
|
||||||
|
oyhi = outyhi;
|
||||||
|
}
|
||||||
|
|
||||||
/* ----------------------------------------------------------------------
|
/* ----------------------------------------------------------------------
|
||||||
partition a global regular grid into one brick-shaped sub-grid per proc
|
partition a global regular grid into one brick-shaped sub-grid per proc
|
||||||
if grid point is inside my sub-domain I own it,
|
if grid point is inside my sub-domain I own it,
|
||||||
@ -432,7 +488,7 @@ void Grid2d::partition_grid(int ngrid, double fraclo, double frachi,
|
|||||||
ghost cell indices for periodic systems can be < 0 or >= N
|
ghost cell indices for periodic systems can be < 0 or >= N
|
||||||
------------------------------------------------------------------------- */
|
------------------------------------------------------------------------- */
|
||||||
|
|
||||||
void Grid2d::ghost_grid(double maxdist, int extra)
|
void Grid2d::ghost_grid()
|
||||||
{
|
{
|
||||||
double *prd,*boxlo,*sublo,*subhi;
|
double *prd,*boxlo,*sublo,*subhi;
|
||||||
int triclinic = domain->triclinic;
|
int triclinic = domain->triclinic;
|
||||||
@ -455,10 +511,11 @@ void Grid2d::ghost_grid(double maxdist, int extra)
|
|||||||
if (triclinic == 0) dist[0] = dist[1] = dist[2] = maxdist;
|
if (triclinic == 0) dist[0] = dist[1] = dist[2] = maxdist;
|
||||||
else MathExtra::tribbox(domain->h,maxdist,&dist[0]);
|
else MathExtra::tribbox(domain->h,maxdist,&dist[0]);
|
||||||
|
|
||||||
// nlo,nhi = min/max index of global grid cell my owned atoms can be mapped to
|
// lo/hi = min/max index of global grid cells my owned atoms can be mapped to
|
||||||
// including up to maxdist displacement outside my subdomain
|
// includes effects of maxdist and shift_atom settings
|
||||||
// extra = additional ghost layers required by called (e.g. finite diff stencil)
|
// lo/hi can be further extended by stencil_atom and stencil_grid settings
|
||||||
// max of the two quantities = ghost cell layers needed in each dim/dir
|
// all those settings are set by caller
|
||||||
|
// ghost cell layers needed in each dim/dir = max of two extension effects
|
||||||
// OFFSET allows generation of negative indices with static_cast
|
// OFFSET allows generation of negative indices with static_cast
|
||||||
// out xyz lo/hi = index range of owned + ghost cells
|
// out xyz lo/hi = index range of owned + ghost cells
|
||||||
|
|
||||||
@ -466,15 +523,15 @@ void Grid2d::ghost_grid(double maxdist, int extra)
|
|||||||
double dyinv = ny / prd[1];
|
double dyinv = ny / prd[1];
|
||||||
int lo, hi;
|
int lo, hi;
|
||||||
|
|
||||||
lo = static_cast<int>((sublo[0]-dist[0]-boxlo[0]) * dxinv + OFFSET) - OFFSET;
|
lo = static_cast<int>((sublo[0]-dist[0]-boxlo[0]) * dxinv + shift_atom + OFFSET) - OFFSET;
|
||||||
hi = static_cast<int>((subhi[0]+dist[0]-boxlo[0]) * dxinv + OFFSET) - OFFSET;
|
hi = static_cast<int>((subhi[0]+dist[0]-boxlo[0]) * dxinv + shift_atom + OFFSET) - OFFSET;
|
||||||
outxlo = MIN(lo, inxlo - extra);
|
outxlo = MIN(lo - stencil_atom_lo, inxlo - stencil_grid_lo);
|
||||||
outxhi = MAX(hi, inxhi + extra);
|
outxhi = MAX(hi + stencil_atom_hi, inxhi + stencil_grid_hi);
|
||||||
|
|
||||||
lo = static_cast<int>((sublo[1]-dist[1]-boxlo[1]) * dyinv + OFFSET) - OFFSET;
|
lo = static_cast<int>((sublo[1]-dist[1]-boxlo[1]) * dyinv + shift_atom + OFFSET) - OFFSET;
|
||||||
hi = static_cast<int>((subhi[1]+dist[1]-boxlo[1]) * dyinv + OFFSET) - OFFSET;
|
hi = static_cast<int>((subhi[1]+dist[1]-boxlo[1]) * dyinv + shift_atom + OFFSET) - OFFSET;
|
||||||
outylo = MIN(lo, inylo - extra);
|
outxlo = MIN(lo - stencil_atom_lo, inylo - stencil_grid_lo);
|
||||||
outyhi = MAX(hi, inyhi + extra);
|
outxhi = MAX(hi + stencil_atom_hi, inyhi + stencil_grid_hi);
|
||||||
|
|
||||||
// limit out xyz lo/hi indices to global grid for non-periodic dims
|
// limit out xyz lo/hi indices to global grid for non-periodic dims
|
||||||
|
|
||||||
@ -547,7 +604,7 @@ void Grid2d::extract_comm_info()
|
|||||||
// ----------------------------------------------------------------------
|
// ----------------------------------------------------------------------
|
||||||
|
|
||||||
/* ----------------------------------------------------------------------
|
/* ----------------------------------------------------------------------
|
||||||
setup owned/ghost commmunication
|
setup communication of owned/ghost grid cells
|
||||||
either for brick decomp or tiled decomp
|
either for brick decomp or tiled decomp
|
||||||
return sizes of two buffers needed for communication
|
return sizes of two buffers needed for communication
|
||||||
nbuf1 = largest pack or unpack in any Send or Recv or Copy
|
nbuf1 = largest pack or unpack in any Send or Recv or Copy
|
||||||
@ -1268,7 +1325,7 @@ void Grid2d::setup_remap(Grid2d *old, int &nremap_buf1, int &nremap_buf2)
|
|||||||
// overlap_old = vector of overlap info in Overlap data struct
|
// overlap_old = vector of overlap info in Overlap data struct
|
||||||
|
|
||||||
int oldbox[4];
|
int oldbox[4];
|
||||||
old->get_bounds(oldbox[0],oldbox[1],oldbox[2],oldbox[3]);
|
old->get_bounds_owned(oldbox[0],oldbox[1],oldbox[2],oldbox[3]);
|
||||||
pbc[0] = pbc[1] = 0;
|
pbc[0] = pbc[1] = 0;
|
||||||
|
|
||||||
Overlap *overlap_old;
|
Overlap *overlap_old;
|
||||||
@ -1306,7 +1363,7 @@ void Grid2d::setup_remap(Grid2d *old, int &nremap_buf1, int &nremap_buf2)
|
|||||||
// overlap_new = vector of overlap info in Overlap data struct
|
// overlap_new = vector of overlap info in Overlap data struct
|
||||||
|
|
||||||
int newbox[4];
|
int newbox[4];
|
||||||
get_bounds(newbox[0],newbox[1],newbox[2],newbox[3]);
|
get_bounds_owned(newbox[0],newbox[1],newbox[2],newbox[3]);
|
||||||
pbc[0] = pbc[1] = 0;
|
pbc[0] = pbc[1] = 0;
|
||||||
|
|
||||||
Overlap *overlap_new;
|
Overlap *overlap_new;
|
||||||
@ -1585,17 +1642,17 @@ int Grid2d::compute_overlap(int ghostflag, int *box, int *pbc, Overlap *&overlap
|
|||||||
|
|
||||||
// find comm->procgrid indices in each dim for box bounds
|
// find comm->procgrid indices in each dim for box bounds
|
||||||
|
|
||||||
int iproclo = proc_index_uniform(box[0],nx,0,comm->xsplit);
|
int iproclo = proc_index_uniform(box[0],nx,0,xsplit);
|
||||||
int iprochi = proc_index_uniform(box[1],nx,0,comm->xsplit);
|
int iprochi = proc_index_uniform(box[1],nx,0,xsplit);
|
||||||
int jproclo = proc_index_uniform(box[2],ny,1,comm->ysplit);
|
int jproclo = proc_index_uniform(box[2],ny,1,ysplit);
|
||||||
int jprochi = proc_index_uniform(box[3],ny,1,comm->ysplit);
|
int jprochi = proc_index_uniform(box[3],ny,1,ysplit);
|
||||||
|
|
||||||
// compute extent of overlap of box with with each proc's obox
|
// compute extent of overlap of box with with each proc's obox
|
||||||
|
|
||||||
for (int j = jproclo; j <= jprochi; j++)
|
for (int j = jproclo; j <= jprochi; j++)
|
||||||
for (int i = iproclo; i <= iprochi; i++) {
|
for (int i = iproclo; i <= iprochi; i++) {
|
||||||
partition_grid(nx,xsplit[i],xsplit[i+1],shift,obox[0],obox[1]);
|
partition_grid(nx,xsplit[i],xsplit[i+1],shift_grid,obox[0],obox[1]);
|
||||||
partition_grid(ny,ysplit[j],ysplit[j+1],shift,obox[2],obox[3]);
|
partition_grid(ny,ysplit[j],ysplit[j+1],shift_grid,obox[2],obox[3]);
|
||||||
|
|
||||||
if (noverlap_list == maxoverlap_list) grow_overlap();
|
if (noverlap_list == maxoverlap_list) grow_overlap();
|
||||||
overlap_list[noverlap_list].proc = grid2proc[i][j][0];
|
overlap_list[noverlap_list].proc = grid2proc[i][j][0];
|
||||||
@ -1851,7 +1908,7 @@ int Grid2d::proc_index_uniform(int igrid, int n, int dim, double *split)
|
|||||||
}
|
}
|
||||||
|
|
||||||
/* ----------------------------------------------------------------------
|
/* ----------------------------------------------------------------------
|
||||||
compute the grid box for proc within tiled decomposition
|
compute the grid box owned by proc within tiled decomposition
|
||||||
performed recursively until proclower = procupper = proc
|
performed recursively until proclower = procupper = proc
|
||||||
return box = lo/hi bounds of proc's box in 2 dims
|
return box = lo/hi bounds of proc's box in 2 dims
|
||||||
------------------------------------------------------------------------- */
|
------------------------------------------------------------------------- */
|
||||||
|
|||||||
39
src/grid2d.h
39
src/grid2d.h
@ -22,19 +22,27 @@ class Grid2d : protected Pointers {
|
|||||||
public:
|
public:
|
||||||
enum { KSPACE = 0, PAIR = 1, FIX = 2 }; // calling classes
|
enum { KSPACE = 0, PAIR = 1, FIX = 2 }; // calling classes
|
||||||
|
|
||||||
Grid2d(class LAMMPS *, MPI_Comm, int, int, double, int, double,
|
Grid2d(class LAMMPS *, MPI_Comm, int, int);
|
||||||
int &, int &, int &, int &,
|
|
||||||
int &, int &, int &, int &);
|
|
||||||
Grid2d(class LAMMPS *, MPI_Comm, int, int, int, int, int, int, int, int, int, int);
|
Grid2d(class LAMMPS *, MPI_Comm, int, int, int, int, int, int, int, int, int, int);
|
||||||
Grid2d(class LAMMPS *, MPI_Comm, int, int, int, int, int, int, int, int, int, int,
|
Grid2d(class LAMMPS *, MPI_Comm, int, int, int, int, int, int, int, int, int, int,
|
||||||
int, int, int, int, int);
|
int, int, int, int, int);
|
||||||
|
|
||||||
~Grid2d() override;
|
~Grid2d() override;
|
||||||
|
|
||||||
|
void set_distance(double);
|
||||||
|
void set_stencil_grid(int, int);
|
||||||
|
void set_stencil_atom(int, int);
|
||||||
|
void set_shift_grid(double);
|
||||||
|
void set_shift_atom(double);
|
||||||
|
|
||||||
int identical(Grid2d *);
|
int identical(Grid2d *);
|
||||||
void get_size(int &, int &);
|
void get_size(int &, int &);
|
||||||
void get_bounds(int &, int &, int &, int &);
|
void get_bounds_owned(int &, int &, int &, int &);
|
||||||
void get_bounds_ghost(int &, int &, int &, int &);
|
void get_bounds_ghost(int &, int &, int &, int &);
|
||||||
|
|
||||||
|
void setup_grid(int &, int &, int &, int &, int &, int &, int &, int &);
|
||||||
|
|
||||||
void setup_comm(int &, int &);
|
void setup_comm(int &, int &);
|
||||||
int ghost_adjacent();
|
int ghost_adjacent();
|
||||||
void forward_comm(int, void *, int, int, int, void *, void *, MPI_Datatype);
|
void forward_comm(int, void *, int, int, int, void *, void *, MPI_Datatype);
|
||||||
@ -55,6 +63,16 @@ protected:
|
|||||||
// inputs from caller via constructor
|
// inputs from caller via constructor
|
||||||
|
|
||||||
int nx, ny; // size of global grid in both dims
|
int nx, ny; // size of global grid in both dims
|
||||||
|
double maxdist; // distance owned atoms can move outside subdomain
|
||||||
|
int stencil_atom_lo,stencil_atom_hi; // grid cells accessed beyond atom's cell
|
||||||
|
int stencil_grid_lo,stencil_grid_hi; // grid cells accessed beyond owned cell
|
||||||
|
double shift_grid; // location of grid point within grid cell
|
||||||
|
// only affects which proc owns grid cell
|
||||||
|
double shift_atom; // shift applied to atomd when mapped to grid cell by caller
|
||||||
|
// only affects extent of ghost cells
|
||||||
|
|
||||||
|
// extent of my owned and ghost cells
|
||||||
|
|
||||||
int inxlo, inxhi; // inclusive extent of my grid chunk, 0 <= in <= N-1
|
int inxlo, inxhi; // inclusive extent of my grid chunk, 0 <= in <= N-1
|
||||||
int inylo, inyhi;
|
int inylo, inyhi;
|
||||||
int outxlo, outxhi; // inclusive extent of my grid chunk plus
|
int outxlo, outxhi; // inclusive extent of my grid chunk plus
|
||||||
@ -63,9 +81,6 @@ protected:
|
|||||||
int fullxlo, fullxhi; // extent of grid chunk that caller stores
|
int fullxlo, fullxhi; // extent of grid chunk that caller stores
|
||||||
int fullylo, fullyhi; // can be same as out indices or larger
|
int fullylo, fullyhi; // can be same as out indices or larger
|
||||||
|
|
||||||
double shift; // location of grid point within grid cell
|
|
||||||
// only affects which proc owns grid cell
|
|
||||||
|
|
||||||
// -------------------------------------------
|
// -------------------------------------------
|
||||||
// internal variables for BRICK layout
|
// internal variables for BRICK layout
|
||||||
// -------------------------------------------
|
// -------------------------------------------
|
||||||
@ -210,14 +225,14 @@ protected:
|
|||||||
// -------------------------------------------
|
// -------------------------------------------
|
||||||
|
|
||||||
void partition_grid(int, double, double, double, int &, int &);
|
void partition_grid(int, double, double, double, int &, int &);
|
||||||
void ghost_grid(double, int);
|
void ghost_grid();
|
||||||
void extract_comm_info();
|
void extract_comm_info();
|
||||||
|
|
||||||
void store(int, int, int, int, int, int, int, int,
|
void store(int, int, int, int, int, int, int, int,
|
||||||
int, int, int, int, int, int, int, int);
|
int, int, int, int, int, int, int, int);
|
||||||
|
|
||||||
virtual void setup_comm_brick(int &, int &);
|
void setup_comm_brick(int &, int &);
|
||||||
virtual void setup_comm_tiled(int &, int &);
|
void setup_comm_tiled(int &, int &);
|
||||||
int ghost_adjacent_brick();
|
int ghost_adjacent_brick();
|
||||||
int ghost_adjacent_tiled();
|
int ghost_adjacent_tiled();
|
||||||
|
|
||||||
@ -236,7 +251,7 @@ protected:
|
|||||||
void box_drop(int *, int *);
|
void box_drop(int *, int *);
|
||||||
void box_drop_grid(int *, int, int, int &, int *);
|
void box_drop_grid(int *, int, int, int &, int *);
|
||||||
|
|
||||||
virtual void grow_swap();
|
void grow_swap();
|
||||||
void grow_overlap();
|
void grow_overlap();
|
||||||
void deallocate_remap();
|
void deallocate_remap();
|
||||||
|
|
||||||
|
|||||||
@ -83,7 +83,6 @@ Grid3d::Grid3d(LAMMPS *lmp, MPI_Comm gcomm,
|
|||||||
: Pointers(lmp)
|
: Pointers(lmp)
|
||||||
{
|
{
|
||||||
// store commnicator and global grid size
|
// store commnicator and global grid size
|
||||||
// set layout mode
|
|
||||||
|
|
||||||
gridcomm = gcomm;
|
gridcomm = gcomm;
|
||||||
MPI_Comm_rank(gridcomm,&me);
|
MPI_Comm_rank(gridcomm,&me);
|
||||||
@ -134,7 +133,6 @@ Grid3d::Grid3d(LAMMPS *lmp, MPI_Comm gcomm, int flag,
|
|||||||
: Pointers(lmp)
|
: Pointers(lmp)
|
||||||
{
|
{
|
||||||
// store commnicator and global grid size
|
// store commnicator and global grid size
|
||||||
// set layout mode
|
|
||||||
|
|
||||||
gridcomm = gcomm;
|
gridcomm = gcomm;
|
||||||
MPI_Comm_rank(gridcomm,&me);
|
MPI_Comm_rank(gridcomm,&me);
|
||||||
|
|||||||
@ -260,7 +260,7 @@ class Grid3d : protected Pointers {
|
|||||||
void box_drop(int *, int *);
|
void box_drop(int *, int *);
|
||||||
void box_drop_grid(int *, int, int, int &, int *);
|
void box_drop_grid(int *, int, int, int &, int *);
|
||||||
|
|
||||||
virtual void grow_swap();
|
void grow_swap();
|
||||||
void grow_overlap();
|
void grow_overlap();
|
||||||
void deallocate_remap();
|
void deallocate_remap();
|
||||||
|
|
||||||
|
|||||||
Reference in New Issue
Block a user