2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2005,2006,2007,2008,2009 by the GROMACS development team.
5 * Copyright (c) 2010,2011,2012,2013,2014 by the GROMACS development team.
6 * Copyright (c) 2015,2016,2017,2018,2019 by the GROMACS development team.
7 * Copyright (c) 2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 * \brief Implements atom redistribution functions.
42 * \author Berk Hess <hess@kth.se>
43 * \ingroup module_domdec
48 #include "redistribute.h"
52 #include "gromacs/domdec/domdec_network.h"
53 #include "gromacs/domdec/ga2la.h"
54 #include "gromacs/gmxlib/nrnb.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/mdlib/gmx_omp_nthreads.h"
57 #include "gromacs/mdtypes/forcerec.h"
58 #include "gromacs/mdtypes/nblist.h"
59 #include "gromacs/mdtypes/state.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/gmxomp.h"
64 #include "domdec_internal.h"
67 /* The size per charge group of the cggl_flag buffer in gmx_domdec_comm_t */
68 static constexpr int DD_CGIBS = 2;
70 /* The flags for the cggl_flag buffer in gmx_domdec_comm_t */
72 /* The lower 16 bits are reserved for the charge group size */
73 static constexpr int DD_FLAG_NRCG = 65535;
75 /* Returns which bit tells whether to move a group forward along dimension d */
76 static inline int DD_FLAG_FW(int d)
78 return 1 << (16 + d * 2);
81 /* Returns which bit tells whether to move a group backward along dimension d */
82 static inline int DD_FLAG_BW(int d)
84 return 1 << (16 + d * 2 + 1);
87 static void copyMovedAtomsToBufferPerAtom(gmx::ArrayRef<const int> move,
91 gmx_domdec_comm_t* comm)
93 int pos_vec[DIM * 2] = { 0 };
95 for (gmx::index i = 0; i < move.ssize(); i++)
97 /* Skip moved atoms */
98 const int m = move[i];
101 /* Copy to the communication buffer */
102 pos_vec[m] += 1 + vec;
103 copy_rvec(src[i], comm->cgcm_state[m][pos_vec[m]++]);
104 pos_vec[m] += nvec - vec - 1;
109 static void copyMovedUpdateGroupCogs(gmx::ArrayRef<const int> move,
111 gmx::ArrayRef<const gmx::RVec> coordinates,
112 gmx_domdec_comm_t* comm)
114 int pos_vec[DIM * 2] = { 0 };
116 for (gmx::index g = 0; g < move.ssize(); g++)
118 /* Skip moved atoms */
119 const int m = move[g];
122 /* Copy to the communication buffer */
123 const gmx::RVec& cog =
124 (comm->systemInfo.useUpdateGroups ? comm->updateGroupsCog->cogForAtom(g)
126 copy_rvec(cog, comm->cgcm_state[m][pos_vec[m]]);
127 pos_vec[m] += 1 + nvec;
132 static void clear_and_mark_ind(gmx::ArrayRef<const int> move,
133 gmx::ArrayRef<const int> globalAtomIndices,
137 for (gmx::index a = 0; a < move.ssize(); a++)
141 /* Clear the global indices */
142 ga2la->erase(globalAtomIndices[a]);
143 /* Signal that this atom has moved using the ns cell index.
144 * Here we set it to -1. fill_grid will change it
145 * from -1 to NSGRID_SIGNAL_MOVED_FAC*grid->ncells.
152 static void print_cg_move(FILE* fplog,
153 const gmx_domdec_t* dd,
158 gmx_bool bHaveCgcmOld,
164 const gmx_domdec_comm_t* comm = dd->comm;
167 fprintf(fplog, "\nStep %" PRId64 ":\n", step);
169 if (comm->systemInfo.useUpdateGroups)
171 mesg += "The update group starting at atom";
177 mesg += gmx::formatString(
178 " %d moved more than the distance allowed by the domain decomposition", ddglatnr(dd, cg));
181 mesg += gmx::formatString(" (%f)", limitd);
183 mesg += gmx::formatString(" in direction %c\n", dim2char(dim));
184 fprintf(fplog, "%s", mesg.c_str());
187 "distance out of cell %f\n",
188 dir == 1 ? pos_d - comm->cell_x1[dim] : pos_d - comm->cell_x0[dim]);
191 fprintf(fplog, "Old coordinates: %8.3f %8.3f %8.3f\n", cm_old[XX], cm_old[YY], cm_old[ZZ]);
193 fprintf(fplog, "New coordinates: %8.3f %8.3f %8.3f\n", cm_new[XX], cm_new[YY], cm_new[ZZ]);
195 "Old cell boundaries in direction %c: %8.3f %8.3f\n",
197 comm->old_cell_x0[dim],
198 comm->old_cell_x1[dim]);
200 "New cell boundaries in direction %c: %8.3f %8.3f\n",
206 [[noreturn]] static void cg_move_error(FILE* fplog,
207 const gmx_domdec_t* dd,
212 gmx_bool bHaveCgcmOld,
220 print_cg_move(fplog, dd, step, cg, dim, dir, bHaveCgcmOld, limitd, cm_old, cm_new, pos_d);
222 print_cg_move(stderr, dd, step, cg, dim, dir, bHaveCgcmOld, limitd, cm_old, cm_new, pos_d);
224 "One or more atoms moved too far between two domain decomposition steps.\n"
225 "This usually means that your system is not well equilibrated");
228 static void rotate_state_atom(t_state* state, int a)
230 if (state->flags & enumValueToBitMask(StateEntry::X))
232 auto x = makeArrayRef(state->x);
233 /* Rotate the complete state; for a rectangular box only */
234 x[a][YY] = state->box[YY][YY] - x[a][YY];
235 x[a][ZZ] = state->box[ZZ][ZZ] - x[a][ZZ];
237 if (state->flags & enumValueToBitMask(StateEntry::V))
239 auto v = makeArrayRef(state->v);
240 v[a][YY] = -v[a][YY];
241 v[a][ZZ] = -v[a][ZZ];
243 if (state->flags & enumValueToBitMask(StateEntry::Cgp))
245 auto cg_p = makeArrayRef(state->cg_p);
246 cg_p[a][YY] = -cg_p[a][YY];
247 cg_p[a][ZZ] = -cg_p[a][ZZ];
251 /* Returns a pointer to a buffer of size \p numAtomsNew with contents 0 from numAtomsOld upwards
253 * Note: numAtomsOld should either be 0 or match the current buffer size.
255 static int* getMovedBuffer(gmx_domdec_comm_t* comm, size_t numAtomsOld, size_t numAtomsNew)
257 std::vector<int>& movedBuffer = comm->movedBuffer;
259 GMX_RELEASE_ASSERT(numAtomsOld == 0 || movedBuffer.size() == numAtomsOld,
260 "numAtomsOld should either be 0 or match the current size");
262 /* Contents up to numAtomsOld should be preserved, clear from numAtomsOld */
263 if (numAtomsOld == 0)
267 movedBuffer.resize(numAtomsNew);
269 return movedBuffer.data();
272 /* Bounds to determine whether an atom group moved to far between DD steps */
275 rvec distance; /* Maximum distance over which a group can move */
276 rvec lower; /* Lower bound for group location */
277 rvec upper; /* Upper bound for group location */
280 /* Compute flag that tells where to move an atom group
282 * The input dev tells whether the group moves fw,remains,bw (-1,0,1) along
284 * The return value is -1 when the group does not move.
285 * The return has move flags set when the group does move and the lower 4 bits
286 * are (mis)used to tell which is the first dimension (bit 1,2,3) the group
287 * needs to be moved along and in which direction (bit 0 not set for fw
288 * and set for bw), however only values in [0,6) are used.
290 static int computeMoveFlag(const gmx_domdec_t& dd, const ivec& dev)
293 int firstMoveDimValue = -1;
294 for (int d = 0; d < dd.ndim; d++)
296 const int dim = dd.dim[d];
299 flag |= DD_FLAG_FW(d);
300 if (firstMoveDimValue == -1)
302 firstMoveDimValue = d * 2;
305 else if (dev[dim] == -1)
307 flag |= DD_FLAG_BW(d);
308 if (firstMoveDimValue == -1)
310 if (dd.numCells[dim] > 2)
312 firstMoveDimValue = d * 2 + 1;
316 firstMoveDimValue = d * 2;
322 return firstMoveDimValue + flag;
325 /* Determine to which atoms in the range \p cg_start, \p cg_end should go.
327 * Returns in the move array where the atoms should go.
329 static void calc_cg_move(FILE* fplog,
337 const MoveLimits& moveLimits,
340 gmx::ArrayRef<int> move)
342 const int npbcdim = dd->unitCellInfo.npbcdim;
343 auto x = makeArrayRef(state->x);
345 for (int a = cg_start; a < cg_end; a++)
347 // TODO: Rename this center of geometry variable to cogNew
349 copy_rvec(x[a], cm_new);
352 /* Do pbc and check DD cell boundary crossings */
353 for (int d = DIM - 1; d >= 0; d--)
355 if (dd->numCells[d] > 1)
357 bool bScrew = (dd->unitCellInfo.haveScrewPBC && d == XX);
358 /* Determine the location of this cg in lattice coordinates */
359 real pos_d = cm_new[d];
362 for (int d2 = d + 1; d2 < DIM; d2++)
364 pos_d += cm_new[d2] * tcm[d2][d];
367 /* Put the charge group in the triclinic unit-cell */
368 if (pos_d >= cell_x1[d])
370 if (pos_d >= moveLimits.upper[d])
373 fplog, dd, step, a, d, 1, false, moveLimits.distance[d], cm_new, cm_new, pos_d);
376 if (dd->ci[d] == dd->numCells[d] - 1)
378 rvec_dec(cm_new, state->box[d]);
381 cm_new[YY] = state->box[YY][YY] - cm_new[YY];
382 cm_new[ZZ] = state->box[ZZ][ZZ] - cm_new[ZZ];
384 rvec_dec(x[a], state->box[d]);
387 rotate_state_atom(state, a);
391 else if (pos_d < cell_x0[d])
393 if (pos_d < moveLimits.lower[d])
396 fplog, dd, step, a, d, -1, false, moveLimits.distance[d], cm_new, cm_new, pos_d);
401 rvec_inc(cm_new, state->box[d]);
404 cm_new[YY] = state->box[YY][YY] - cm_new[YY];
405 cm_new[ZZ] = state->box[ZZ][ZZ] - cm_new[ZZ];
407 rvec_inc(x[a], state->box[d]);
410 rotate_state_atom(state, a);
415 else if (d < npbcdim)
417 /* Put the charge group in the rectangular unit-cell */
418 while (cm_new[d] >= state->box[d][d])
420 rvec_dec(cm_new, state->box[d]);
421 rvec_dec(x[a], state->box[d]);
423 while (cm_new[d] < 0)
425 rvec_inc(cm_new, state->box[d]);
426 rvec_inc(x[a], state->box[d]);
431 /* Temporarily store the flag in move */
432 move[a] = computeMoveFlag(*dd, dev);
438 /* Constructor that purposely does not initialize anything */
445 /* Determine to which domains update groups in the range \p groupBegin, \p groupEnd should go.
447 * Returns in the move array where the groups should go.
448 * Also updates the COGs and coordinates for jumps over periodic boundaries.
450 static void calcGroupMove(FILE* fplog,
452 const gmx_domdec_t* dd,
453 const t_state* state,
458 const MoveLimits& moveLimits,
461 gmx::ArrayRef<PbcAndFlag> pbcAndFlags)
463 GMX_RELEASE_ASSERT(!dd->unitCellInfo.haveScrewPBC, "Screw PBC is not supported here");
465 const int npbcdim = dd->unitCellInfo.npbcdim;
467 gmx::UpdateGroupsCog* updateGroupsCog = dd->comm->updateGroupsCog.get();
469 for (int g = groupBegin; g < groupEnd; g++)
472 gmx::RVec& cog = updateGroupsCog->cog(g);
473 gmx::RVec cogOld = cog;
476 /* Do pbc and check DD cell boundary crossings */
477 for (int d = DIM - 1; d >= 0; d--)
479 if (dd->numCells[d] > 1)
481 /* Determine the location of this COG in lattice coordinates */
485 for (int d2 = d + 1; d2 < DIM; d2++)
487 pos_d += cog[d2] * tcm[d2][d];
490 /* Put the COG in the triclinic unit-cell */
491 if (pos_d >= cell_x1[d])
493 if (pos_d >= moveLimits.upper[d])
496 fplog, dd, step, g, d, 1, true, moveLimits.distance[d], cogOld, cog, pos_d);
499 if (dd->ci[d] == dd->numCells[d] - 1)
501 rvec_dec(cog, state->box[d]);
504 else if (pos_d < cell_x0[d])
506 if (pos_d < moveLimits.lower[d])
509 fplog, dd, step, g, d, -1, true, moveLimits.distance[d], cogOld, cog, pos_d);
514 rvec_inc(cog, state->box[d]);
518 else if (d < npbcdim)
520 /* Put the COG in the rectangular unit-cell */
521 while (cog[d] >= state->box[d][d])
523 rvec_dec(cog, state->box[d]);
527 rvec_inc(cog, state->box[d]);
532 /* Store the PBC and move flag, so we can later apply them to the atoms */
533 PbcAndFlag& pbcAndFlag = pbcAndFlags[g];
535 rvec_sub(cog, cogOld, pbcAndFlag.pbcShift);
536 pbcAndFlag.moveFlag = computeMoveFlag(*dd, dev);
540 static void applyPbcAndSetMoveFlags(const gmx::UpdateGroupsCog& updateGroupsCog,
541 gmx::ArrayRef<const PbcAndFlag> pbcAndFlags,
544 gmx::ArrayRef<gmx::RVec> atomCoords,
545 gmx::ArrayRef<int> move)
547 for (int a = atomBegin; a < atomEnd; a++)
549 const PbcAndFlag& pbcAndFlag = pbcAndFlags[updateGroupsCog.cogIndex(a)];
550 rvec_inc(atomCoords[a], pbcAndFlag.pbcShift);
551 /* Temporarily store the flag in move */
552 move[a] = pbcAndFlag.moveFlag;
556 void dd_redistribute_cg(FILE* fplog,
565 gmx_domdec_comm_t* comm = dd->comm;
567 if (dd->unitCellInfo.haveScrewPBC)
569 check_screw_box(state->box);
572 // Positions are always present, so there's nothing to flag
573 bool bV = (state->flags & enumValueToBitMask(StateEntry::V)) != 0;
574 bool bCGP = (state->flags & enumValueToBitMask(StateEntry::Cgp)) != 0;
576 DDBufferAccess<int> moveBuffer(comm->intBuffer, dd->numHomeAtoms);
577 gmx::ArrayRef<int> move = moveBuffer.buffer;
579 const int npbcdim = dd->unitCellInfo.npbcdim;
581 rvec cell_x0, cell_x1;
582 MoveLimits moveLimits;
583 for (int d = 0; (d < DIM); d++)
585 moveLimits.distance[d] = dd->comm->cellsize_min[d];
586 if (d >= npbcdim && dd->ci[d] == 0)
588 cell_x0[d] = -GMX_FLOAT_MAX;
592 cell_x0[d] = comm->cell_x0[d];
594 if (d >= npbcdim && dd->ci[d] == dd->numCells[d] - 1)
596 cell_x1[d] = GMX_FLOAT_MAX;
600 cell_x1[d] = comm->cell_x1[d];
604 moveLimits.lower[d] = comm->old_cell_x0[d] - moveLimits.distance[d];
605 moveLimits.upper[d] = comm->old_cell_x1[d] + moveLimits.distance[d];
609 /* We check after communication if a charge group moved
610 * more than one cell. Set the pre-comm check limit to float_max.
612 moveLimits.lower[d] = -GMX_FLOAT_MAX;
613 moveLimits.upper[d] = GMX_FLOAT_MAX;
618 make_tric_corr_matrix(npbcdim, state->box, tcm);
620 const int nthread = gmx_omp_nthreads_get(ModuleMultiThread::Domdec);
622 /* Compute the center of geometry for all home charge groups
623 * and put them in the box and determine where they should go.
625 std::vector<PbcAndFlag> pbcAndFlags(
626 comm->systemInfo.useUpdateGroups ? comm->updateGroupsCog->numCogs() : 0);
628 #pragma omp parallel num_threads(nthread)
632 const int thread = gmx_omp_get_thread_num();
634 if (comm->systemInfo.useUpdateGroups)
636 const auto& updateGroupsCog = *comm->updateGroupsCog;
637 const int numGroups = updateGroupsCog.numCogs();
647 (thread * numGroups) / nthread,
648 ((thread + 1) * numGroups) / nthread,
650 /* We need a barrier as atoms below can be in a COG of a different thread */
652 const int numHomeAtoms = comm->atomRanges.numHomeAtoms();
653 applyPbcAndSetMoveFlags(updateGroupsCog,
655 (thread * numHomeAtoms) / nthread,
656 ((thread + 1) * numHomeAtoms) / nthread,
662 /* Here we handle single atoms or charge groups */
672 (thread * dd->numHomeAtoms) / nthread,
673 ((thread + 1) * dd->numHomeAtoms) / nthread,
677 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
680 // The counts of atoms to move, forward or backward, over the
681 // possible DIM dimensions.
682 int nat[DIM * 2] = { 0 };
683 for (int cg = 0; cg < dd->numHomeAtoms; cg++)
687 // The value in move[cg] was computed by computeMoveFlags
688 // and describes how this atom should move between domains.
689 const int flag = move[cg] & ~DD_FLAG_NRCG;
690 // mc contains 4 bits that tell which is the first
691 // dimension (bit 1,2,3) that the group needs to be moved
692 // along, and in which direction (bit 0; not set for fw
693 // and set for bw). However the value is always in
695 const int mc = move[cg] & DD_FLAG_NRCG;
698 std::vector<int>& cggl_flag = comm->cggl_flag[mc];
700 /* TODO: See if we can use push_back instead */
701 if ((nat[mc] + 1) * DD_CGIBS > gmx::index(cggl_flag.size()))
703 cggl_flag.resize((nat[mc] + 1) * DD_CGIBS);
705 cggl_flag[nat[mc] * DD_CGIBS] = dd->globalAtomGroupIndices[cg];
706 cggl_flag[nat[mc] * DD_CGIBS + 1] = flag;
711 inc_nrnb(nrnb, eNR_CGCM, comm->atomRanges.numHomeAtoms());
712 inc_nrnb(nrnb, eNR_RESETX, dd->numHomeAtoms);
715 for (int i = 0; i < dd->ndim * 2; i++)
717 *ncg_moved += nat[i];
730 /* Make sure the communication buffers are large enough */
731 for (int mc = 0; mc < dd->ndim * 2; mc++)
733 size_t nvr = nat[mc] * (1 + nvec);
734 if (nvr > comm->cgcm_state[mc].size())
736 comm->cgcm_state[mc].resize(nvr);
740 /* With update groups we send over their COGs.
741 * Without update groups we send the moved atom coordinates
742 * over twice. This is so the code further down can be used
743 * without many conditionals both with and without update groups.
745 copyMovedUpdateGroupCogs(move, nvec, state->x, comm);
748 copyMovedAtomsToBufferPerAtom(move, nvec, vectorIndex++, state->x.rvec_array(), comm);
751 copyMovedAtomsToBufferPerAtom(move, nvec, vectorIndex++, state->v.rvec_array(), comm);
755 copyMovedAtomsToBufferPerAtom(move, nvec, vectorIndex++, state->cg_p.rvec_array(), comm);
758 int* moved = getMovedBuffer(comm, 0, dd->numHomeAtoms);
760 clear_and_mark_ind(move, dd->globalAtomIndices, dd->ga2la, moved);
762 /* Now we can remove the excess global atom-group indices from the list */
763 dd->globalAtomGroupIndices.resize(dd->numHomeAtoms);
765 /* We reuse the intBuffer without reacquiring since we are in the same scope */
766 DDBufferAccess<int>& flagBuffer = moveBuffer;
768 gmx::ArrayRef<const gmx::AtomInfoWithinMoleculeBlock> atomInfoForEachMoleculeBlock =
769 fr->atomInfoForEachMoleculeBlock;
771 /* Temporarily store atoms passed to our rank at the end of the range */
772 int home_pos_at = dd->numHomeAtoms;
773 for (int d = 0; d < dd->ndim; d++)
775 DDBufferAccess<gmx::RVec> rvecBuffer(comm->rvecBuffer, 0);
777 const int dim = dd->dim[d];
778 int totalAtomsReceived = 0;
780 for (int dir = 0; dir < (dd->numCells[dim] == 2 ? 1 : 2); dir++)
782 const int cdd = d * 2 + dir;
783 /* Communicate the atom counts */
786 fprintf(debug, "Sending ddim %d dir %d: nat %d\n", d, dir, nat[cdd]);
788 int atomCountToReceive;
789 ddSendrecv(dd, d, dir, &nat[cdd], 1, &atomCountToReceive, 1);
791 flagBuffer.resize((totalAtomsReceived + atomCountToReceive) * DD_CGIBS);
793 /* Communicate the charge group indices, sizes and flags */
797 comm->cggl_flag[cdd].data(),
799 flagBuffer.buffer.data() + totalAtomsReceived * DD_CGIBS,
800 atomCountToReceive * DD_CGIBS);
802 const int nvs = nat[cdd] * (1 + nvec);
803 const int i = atomCountToReceive * (1 + nvec);
804 rvecBuffer.resize(nvr + i);
806 /* Communicate cgcm and state */
810 as_rvec_array(comm->cgcm_state[cdd].data()),
812 as_rvec_array(rvecBuffer.buffer.data()) + nvr,
814 totalAtomsReceived += atomCountToReceive;
818 dd_resize_atominfo_and_state(fr, state, home_pos_at + totalAtomsReceived);
820 /* Process the received charge or update groups */
822 for (int cg = 0; cg < totalAtomsReceived; cg++)
824 /* Extract the move flags and COG for the charge or update group */
825 int flag = flagBuffer.buffer[cg * DD_CGIBS + 1];
826 const gmx::RVec& cog = rvecBuffer.buffer[buf_pos];
828 if (dim >= npbcdim && dd->numCells[dim] > 2)
830 /* No pbc in this dim and more than one domain boundary.
831 * We do a separate check if a charge group didn't move too far.
833 if (((flag & DD_FLAG_FW(d)) && cog[dim] > cell_x1[dim])
834 || ((flag & DD_FLAG_BW(d)) && cog[dim] < cell_x0[dim]))
836 rvec pos = { cog[0], cog[1], cog[2] };
838 fplog, dd, step, cg, dim, (flag & DD_FLAG_FW(d)) ? 1 : 0, false, 0, pos, pos, pos[dim]);
843 if (d < dd->ndim - 1)
845 /* Check which direction this cg should go */
846 for (int d2 = d + 1; (d2 < dd->ndim && mc == -1); d2++)
848 if (isDlbOn(dd->comm))
850 /* The cell boundaries for dimension d2 are not equal
851 * for each cell row of the lower dimension(s),
852 * therefore we might need to redetermine where
855 const int dim2 = dd->dim[d2];
856 /* The DD grid is not staggered at the box boundaries,
857 * so we do not need to handle boundary crossings.
858 * This also means we do not have to handle PBC here.
860 if (!((dd->ci[dim2] == dd->numCells[dim2] - 1 && (flag & DD_FLAG_FW(d2)))
861 || (dd->ci[dim2] == 0 && (flag & DD_FLAG_BW(d2)))))
863 /* Clear the two flags for this dimension */
864 flag &= ~(DD_FLAG_FW(d2) | DD_FLAG_BW(d2));
865 /* Determine the location of this cg
866 * in lattice coordinates
868 real pos_d = cog[dim2];
871 for (int d3 = dim2 + 1; d3 < DIM; d3++)
873 pos_d += cog[d3] * tcm[d3][dim2];
877 GMX_ASSERT(dim2 >= 0 && dim2 < DIM, "Keep the static analyzer happy");
879 /* Check if we need to move this group
880 * to an adjacent cell because of the
883 if (pos_d >= cell_x1[dim2] && dd->ci[dim2] != dd->numCells[dim2] - 1)
885 flag |= DD_FLAG_FW(d2);
887 else if (pos_d < cell_x0[dim2] && dd->ci[dim2] != 0)
889 flag |= DD_FLAG_BW(d2);
892 flagBuffer.buffer[cg * DD_CGIBS + 1] = flag;
895 /* Set to which neighboring cell this cg should go */
896 if (flag & DD_FLAG_FW(d2))
900 else if (flag & DD_FLAG_BW(d2))
902 if (dd->numCells[dd->dim[d2]] > 2)
916 /* Set the global charge group index and size */
917 const int globalAtomGroupIndex = flagBuffer.buffer[cg * DD_CGIBS];
918 dd->globalAtomGroupIndices.push_back(globalAtomGroupIndex);
919 /* Skip the COG entry in the buffer */
923 fr->atomInfo[home_pos_at] =
924 ddGetAtomInfo(atomInfoForEachMoleculeBlock, globalAtomGroupIndex);
926 auto x = makeArrayRef(state->x);
927 auto v = makeArrayRef(state->v);
928 auto cg_p = makeArrayRef(state->cg_p);
929 rvec* rvecPtr = as_rvec_array(rvecBuffer.buffer.data());
930 copy_rvec(rvecPtr[buf_pos++], x[home_pos_at]);
933 copy_rvec(rvecPtr[buf_pos++], v[home_pos_at]);
937 copy_rvec(rvecPtr[buf_pos++], cg_p[home_pos_at]);
943 /* Reallocate the buffers if necessary */
944 if ((nat[mc] + 1) * DD_CGIBS > gmx::index(comm->cggl_flag[mc].size()))
946 comm->cggl_flag[mc].resize((nat[mc] + 1) * DD_CGIBS);
948 size_t nvr = nat[mc] * (1 + nvec);
949 if (nvr + 1 + nvec > comm->cgcm_state[mc].size())
951 comm->cgcm_state[mc].resize(nvr + 1 + nvec);
953 /* Copy from the receive to the send buffers */
954 memcpy(comm->cggl_flag[mc].data() + nat[mc] * DD_CGIBS,
955 flagBuffer.buffer.data() + cg * DD_CGIBS,
956 DD_CGIBS * sizeof(int));
957 memcpy(comm->cgcm_state[mc][nvr],
958 rvecBuffer.buffer.data() + buf_pos,
959 (1 + nvec) * sizeof(rvec));
966 /* Note that the indices are now only partially up to date
967 * and numHomeAtoms and nat_home are not the real count, since there are
968 * "holes" in the arrays for the charge groups that moved to neighbors.
971 /* We need to clear the moved flags for the received atoms,
972 * because the moved buffer will be passed to the nbnxm gridding call.
974 moved = getMovedBuffer(comm, dd->numHomeAtoms, home_pos_at);
976 for (int i = dd->numHomeAtoms; i < home_pos_at; i++)
981 dd->numHomeAtoms = home_pos_at;
982 comm->atomRanges.setEnd(DDAtomRanges::Type::Home, home_pos_at);
987 "Finished repartitioning: cgs moved out %d, new home %d\n",
989 dd->numHomeAtoms - *ncg_moved);