Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / domdec / collect.cpp
blob00fd425223a0feabac62d442b4bb6d2f0065185b
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 /* \internal \file
37 * \brief Implements functions to collect state data to the master rank.
39 * \author Berk Hess <hess@kth.se>
40 * \ingroup module_domdec
43 #include "gmxpre.h"
45 #include "collect.h"
47 #include "config.h"
49 #include "gromacs/domdec/domdec_network.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/mdtypes/state.h"
52 #include "gromacs/utility/fatalerror.h"
54 #include "atomdistribution.h"
55 #include "distribute.h"
56 #include "domdec_internal.h"
58 static void dd_collect_cg(gmx_domdec_t *dd,
59 const t_state *state_local)
61 if (state_local->ddp_count == dd->comm->master_cg_ddp_count)
63 /* The master has the correct distribution */
64 return;
67 gmx::ArrayRef<const int> atomGroups;
68 int nat_home = 0;
70 if (state_local->ddp_count == dd->ddp_count)
72 /* The local state and DD are in sync, use the DD indices */
73 atomGroups = gmx::constArrayRefFromArray(dd->globalAtomGroupIndices.data(), dd->ncg_home);
74 nat_home = dd->comm->atomRanges.numHomeAtoms();
76 else if (state_local->ddp_count_cg_gl == state_local->ddp_count)
78 /* The DD is out of sync with the local state, but we have stored
79 * the cg indices with the local state, so we can use those.
81 const t_block &cgs_gl = dd->comm->cgs_gl;
83 atomGroups = state_local->cg_gl;
84 nat_home = 0;
85 for (const int &i : atomGroups)
87 nat_home += cgs_gl.index[i + 1] - cgs_gl.index[i];
90 else
92 gmx_incons("Attempted to collect a vector for a state for which the charge group distribution is unknown");
95 AtomDistribution *ma = dd->ma.get();
97 /* Collect the charge group and atom counts on the master */
98 int localBuffer[2] = { static_cast<int>(atomGroups.size()), nat_home };
99 dd_gather(dd, 2*sizeof(int), localBuffer,
100 DDMASTER(dd) ? ma->intBuffer.data() : nullptr);
102 if (DDMASTER(dd))
104 int groupOffset = 0;
105 for (int rank = 0; rank < dd->nnodes; rank++)
107 auto &domainGroups = ma->domainGroups[rank];
108 int numGroups = ma->intBuffer[2*rank];
110 domainGroups.atomGroups = gmx::constArrayRefFromArray(ma->atomGroups.data() + groupOffset, numGroups);
112 domainGroups.numAtoms = ma->intBuffer[2*rank + 1];
114 groupOffset += numGroups;
117 if (debug)
119 fprintf(debug, "Initial charge group distribution: ");
120 for (int rank = 0; rank < dd->nnodes; rank++)
122 fprintf(debug, " %td", ma->domainGroups[rank].atomGroups.ssize());
124 fprintf(debug, "\n");
127 /* Make byte counts and indices */
128 int offset = 0;
129 for (int rank = 0; rank < dd->nnodes; rank++)
131 int numGroups = ma->domainGroups[rank].atomGroups.size();
132 ma->intBuffer[rank] = numGroups*sizeof(int);
133 ma->intBuffer[dd->nnodes + rank] = offset*sizeof(int);
134 offset += numGroups;
138 /* Collect the charge group indices on the master */
139 dd_gatherv(dd,
140 atomGroups.size()*sizeof(int), atomGroups.data(),
141 DDMASTER(dd) ? ma->intBuffer.data() : nullptr,
142 DDMASTER(dd) ? ma->intBuffer.data() + dd->nnodes : nullptr,
143 DDMASTER(dd) ? ma->atomGroups.data() : nullptr);
145 dd->comm->master_cg_ddp_count = state_local->ddp_count;
148 static void dd_collect_vec_sendrecv(gmx_domdec_t *dd,
149 gmx::ArrayRef<const gmx::RVec> lv,
150 gmx::ArrayRef<gmx::RVec> v)
152 if (!DDMASTER(dd))
154 #if GMX_MPI
155 const int numHomeAtoms = dd->comm->atomRanges.numHomeAtoms();
156 MPI_Send(const_cast<void *>(static_cast<const void *>(lv.data())), numHomeAtoms*sizeof(rvec), MPI_BYTE,
157 dd->masterrank, dd->rank, dd->mpi_comm_all);
158 #endif
160 else
162 AtomDistribution &ma = *dd->ma;
164 /* Copy the master coordinates to the global array */
165 const t_block &cgs_gl = dd->comm->cgs_gl;
167 int rank = dd->masterrank;
168 int localAtom = 0;
169 for (const int &i : ma.domainGroups[rank].atomGroups)
171 for (int globalAtom = cgs_gl.index[i]; globalAtom < cgs_gl.index[i + 1]; globalAtom++)
173 copy_rvec(lv[localAtom++], v[globalAtom]);
177 for (int rank = 0; rank < dd->nnodes; rank++)
179 if (rank != dd->rank)
181 const auto &domainGroups = ma.domainGroups[rank];
183 GMX_RELEASE_ASSERT(v.data() != ma.rvecBuffer.data(), "We need different communication and return buffers");
185 /* When we send/recv instead of scatter/gather, we might need
186 * to increase the communication buffer size here.
188 if (static_cast<size_t>(domainGroups.numAtoms) > ma.rvecBuffer.size())
190 ma.rvecBuffer.resize(domainGroups.numAtoms);
193 #if GMX_MPI
194 MPI_Recv(ma.rvecBuffer.data(), domainGroups.numAtoms*sizeof(rvec), MPI_BYTE, rank,
195 rank, dd->mpi_comm_all, MPI_STATUS_IGNORE);
196 #endif
197 int localAtom = 0;
198 for (const int &cg : domainGroups.atomGroups)
200 for (int globalAtom = cgs_gl.index[cg]; globalAtom < cgs_gl.index[cg + 1]; globalAtom++)
202 copy_rvec(ma.rvecBuffer[localAtom++], v[globalAtom]);
210 static void dd_collect_vec_gatherv(gmx_domdec_t *dd,
211 gmx::ArrayRef<const gmx::RVec> lv,
212 gmx::ArrayRef<gmx::RVec> v)
214 int *recvCounts = nullptr;
215 int *displacements = nullptr;
217 if (DDMASTER(dd))
219 get_commbuffer_counts(dd->ma.get(), &recvCounts, &displacements);
222 const int numHomeAtoms = dd->comm->atomRanges.numHomeAtoms();
223 dd_gatherv(dd, numHomeAtoms*sizeof(rvec), lv.data(), recvCounts, displacements,
224 DDMASTER(dd) ? dd->ma->rvecBuffer.data() : nullptr);
226 if (DDMASTER(dd))
228 const AtomDistribution &ma = *dd->ma;
230 const t_block &cgs_gl = dd->comm->cgs_gl;
232 int bufferAtom = 0;
233 for (int rank = 0; rank < dd->nnodes; rank++)
235 const auto &domainGroups = ma.domainGroups[rank];
236 for (const int &cg : domainGroups.atomGroups)
238 for (int globalAtom = cgs_gl.index[cg]; globalAtom < cgs_gl.index[cg + 1]; globalAtom++)
240 copy_rvec(ma.rvecBuffer[bufferAtom++], v[globalAtom]);
247 void dd_collect_vec(gmx_domdec_t *dd,
248 const t_state *state_local,
249 gmx::ArrayRef<const gmx::RVec> lv,
250 gmx::ArrayRef<gmx::RVec> v)
252 dd_collect_cg(dd, state_local);
254 if (dd->nnodes <= c_maxNumRanksUseSendRecvForScatterAndGather)
256 dd_collect_vec_sendrecv(dd, lv, v);
258 else
260 dd_collect_vec_gatherv(dd, lv, v);
265 void dd_collect_state(gmx_domdec_t *dd,
266 const t_state *state_local, t_state *state)
268 int nh = state_local->nhchainlength;
270 if (DDMASTER(dd))
272 GMX_RELEASE_ASSERT(state->nhchainlength == nh, "The global and local Nose-Hoover chain lengths should match");
274 for (int i = 0; i < efptNR; i++)
276 state->lambda[i] = state_local->lambda[i];
278 state->fep_state = state_local->fep_state;
279 state->veta = state_local->veta;
280 state->vol0 = state_local->vol0;
281 copy_mat(state_local->box, state->box);
282 copy_mat(state_local->boxv, state->boxv);
283 copy_mat(state_local->svir_prev, state->svir_prev);
284 copy_mat(state_local->fvir_prev, state->fvir_prev);
285 copy_mat(state_local->pres_prev, state->pres_prev);
287 for (int i = 0; i < state_local->ngtc; i++)
289 for (int j = 0; j < nh; j++)
291 state->nosehoover_xi[i*nh+j] = state_local->nosehoover_xi[i*nh+j];
292 state->nosehoover_vxi[i*nh+j] = state_local->nosehoover_vxi[i*nh+j];
294 state->therm_integral[i] = state_local->therm_integral[i];
296 for (int i = 0; i < state_local->nnhpres; i++)
298 for (int j = 0; j < nh; j++)
300 state->nhpres_xi[i*nh+j] = state_local->nhpres_xi[i*nh+j];
301 state->nhpres_vxi[i*nh+j] = state_local->nhpres_vxi[i*nh+j];
304 state->baros_integral = state_local->baros_integral;
305 state->pull_com_prev_step = state_local->pull_com_prev_step;
307 if (state_local->flags & (1 << estX))
309 auto globalXRef = state ? state->x : gmx::ArrayRef<gmx::RVec>();
310 dd_collect_vec(dd, state_local, state_local->x, globalXRef);
312 if (state_local->flags & (1 << estV))
314 auto globalVRef = state ? state->v : gmx::ArrayRef<gmx::RVec>();
315 dd_collect_vec(dd, state_local, state_local->v, globalVRef);
317 if (state_local->flags & (1 << estCGP))
319 auto globalCgpRef = state ? state->cg_p : gmx::ArrayRef<gmx::RVec>();
320 dd_collect_vec(dd, state_local, state_local->cg_p, globalCgpRef);