2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief This file contains function definitions necessary for
40 * managing the offload of long-ranged PME work to separate MPI rank,
41 * for computing energies and forces (Coulomb and LJ).
43 * \author Berk Hess <hess@kth.se>
44 * \ingroup module_ewald
54 #include "gromacs/domdec/domdec.h"
55 #include "gromacs/domdec/domdec_struct.h"
56 #include "gromacs/ewald/pme.h"
57 #include "gromacs/ewald/pme_pp_comm_gpu.h"
58 #include "gromacs/gmxlib/network.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/mdlib/gmx_omp_nthreads.h"
61 #include "gromacs/mdtypes/commrec.h"
62 #include "gromacs/mdtypes/forceoutput.h"
63 #include "gromacs/mdtypes/forcerec.h"
64 #include "gromacs/mdtypes/interaction_const.h"
65 #include "gromacs/mdtypes/md_enums.h"
66 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
67 #include "gromacs/nbnxm/nbnxm.h"
68 #include "gromacs/timing/wallcycle.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/gmxmpi.h"
71 #include "gromacs/utility/smalloc.h"
73 #include "pme_internal.h"
74 #include "pme_pp_communication.h"
76 /*! \brief Block to wait for communication to PME ranks to complete
78 * This should be faster with a real non-blocking MPI implementation
80 static constexpr bool c_useDelayedWait = false;
82 /*! \brief Wait for the pending data send requests to PME ranks to complete */
83 static void gmx_pme_send_coeffs_coords_wait(gmx_domdec_t* dd)
88 MPI_Waitall(dd->nreq_pme, dd->req_pme, MPI_STATUSES_IGNORE);
94 /*! \brief Send data to PME ranks */
95 static void gmx_pme_send_coeffs_coords(t_forcerec* fr,
98 real gmx_unused* chargeA,
99 real gmx_unused* chargeB,
100 real gmx_unused* c6A,
101 real gmx_unused* c6B,
102 real gmx_unused* sigmaA,
103 real gmx_unused* sigmaB,
105 const rvec gmx_unused* x,
111 bool useGpuPmePpComms,
112 bool reinitGpuPmePpComms,
113 bool sendCoordinatesFromGpu,
114 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
117 gmx_pme_comm_n_box_t* cnb;
121 n = dd_numHomeAtoms(*dd);
125 fprintf(debug, "PP rank %d sending to PME rank %d: %d%s%s%s%s\n", cr->sim_nodeid, dd->pme_nodeid,
126 n, (flags & PP_PME_CHARGE) ? " charges" : "", (flags & PP_PME_SQRTC6) ? " sqrtC6" : "",
127 (flags & PP_PME_SIGMA) ? " sigma" : "", (flags & PP_PME_COORD) ? " coordinates" : "");
130 if (useGpuPmePpComms)
132 flags |= PP_PME_GPUCOMMS;
135 if (c_useDelayedWait)
137 /* We can not use cnb until pending communication has finished */
138 gmx_pme_send_coeffs_coords_wait(dd);
141 if (dd->pme_receive_vir_ener)
143 /* Peer PP node: communicate all data */
144 if (dd->cnb == nullptr)
152 cnb->maxshift_x = maxshift_x;
153 cnb->maxshift_y = maxshift_y;
154 cnb->lambda_q = lambda_q;
155 cnb->lambda_lj = lambda_lj;
157 if (flags & PP_PME_COORD)
159 copy_mat(box, cnb->box);
162 MPI_Isend(cnb, sizeof(*cnb), MPI_BYTE, dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim,
163 &dd->req_pme[dd->nreq_pme++]);
166 else if (flags & (PP_PME_CHARGE | PP_PME_SQRTC6 | PP_PME_SIGMA))
169 /* Communicate only the number of atoms */
170 MPI_Isend(&n, sizeof(n), MPI_BYTE, dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim,
171 &dd->req_pme[dd->nreq_pme++]);
178 if (flags & PP_PME_CHARGE)
180 MPI_Isend(chargeA, n * sizeof(real), MPI_BYTE, dd->pme_nodeid, eCommType_ChargeA,
181 cr->mpi_comm_mysim, &dd->req_pme[dd->nreq_pme++]);
183 if (flags & PP_PME_CHARGEB)
185 MPI_Isend(chargeB, n * sizeof(real), MPI_BYTE, dd->pme_nodeid, eCommType_ChargeB,
186 cr->mpi_comm_mysim, &dd->req_pme[dd->nreq_pme++]);
188 if (flags & PP_PME_SQRTC6)
190 MPI_Isend(c6A, n * sizeof(real), MPI_BYTE, dd->pme_nodeid, eCommType_SQRTC6A,
191 cr->mpi_comm_mysim, &dd->req_pme[dd->nreq_pme++]);
193 if (flags & PP_PME_SQRTC6B)
195 MPI_Isend(c6B, n * sizeof(real), MPI_BYTE, dd->pme_nodeid, eCommType_SQRTC6B,
196 cr->mpi_comm_mysim, &dd->req_pme[dd->nreq_pme++]);
198 if (flags & PP_PME_SIGMA)
200 MPI_Isend(sigmaA, n * sizeof(real), MPI_BYTE, dd->pme_nodeid, eCommType_SigmaA,
201 cr->mpi_comm_mysim, &dd->req_pme[dd->nreq_pme++]);
203 if (flags & PP_PME_SIGMAB)
205 MPI_Isend(sigmaB, n * sizeof(real), MPI_BYTE, dd->pme_nodeid, eCommType_SigmaB,
206 cr->mpi_comm_mysim, &dd->req_pme[dd->nreq_pme++]);
208 if (flags & PP_PME_COORD)
210 if (reinitGpuPmePpComms)
212 fr->pmePpCommGpu->reinit(n);
216 /* MPI_Isend does not accept a const buffer pointer */
217 real* xRealPtr = const_cast<real*>(x[0]);
218 if (useGpuPmePpComms && (fr != nullptr))
220 void* sendPtr = sendCoordinatesFromGpu
221 ? static_cast<void*>(fr->stateGpu->getCoordinates())
222 : static_cast<void*>(xRealPtr);
223 fr->pmePpCommGpu->sendCoordinatesToPmeCudaDirect(sendPtr, n, sendCoordinatesFromGpu,
224 coordinatesReadyOnDeviceEvent);
228 MPI_Isend(xRealPtr, n * sizeof(rvec), MPI_BYTE, dd->pme_nodeid, eCommType_COORD,
229 cr->mpi_comm_mysim, &dd->req_pme[dd->nreq_pme++]);
234 GMX_UNUSED_VALUE(fr);
235 GMX_UNUSED_VALUE(reinitGpuPmePpComms);
236 GMX_UNUSED_VALUE(sendCoordinatesFromGpu);
237 GMX_UNUSED_VALUE(coordinatesReadyOnDeviceEvent);
239 if (!c_useDelayedWait)
241 /* Wait for the data to arrive */
242 /* We can skip this wait as we are sure x and q will not be modified
243 * before the next call to gmx_pme_send_x_q or gmx_pme_receive_f.
245 gmx_pme_send_coeffs_coords_wait(dd);
249 void gmx_pme_send_parameters(const t_commrec* cr,
250 const interaction_const_t* ic,
251 gmx_bool bFreeEnergy_q,
252 gmx_bool bFreeEnergy_lj,
262 unsigned int flags = 0;
264 if (EEL_PME(ic->eeltype))
266 flags |= PP_PME_CHARGE;
268 if (EVDW_PME(ic->vdwtype))
270 flags |= (PP_PME_SQRTC6 | PP_PME_SIGMA);
272 if (bFreeEnergy_q || bFreeEnergy_lj)
274 /* Assumes that the B state flags are in the bits just above
275 * the ones for the A state. */
276 flags |= (flags << 1);
279 gmx_pme_send_coeffs_coords(nullptr, cr, flags, chargeA, chargeB, sqrt_c6A, sqrt_c6B, sigmaA,
280 sigmaB, nullptr, nullptr, 0, 0, maxshift_x, maxshift_y, -1, false,
281 false, false, nullptr);
284 void gmx_pme_send_coordinates(t_forcerec* fr,
292 bool useGpuPmePpComms,
293 bool receiveCoordinateAddressFromPme,
294 bool sendCoordinatesFromGpu,
295 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent,
296 gmx_wallcycle* wcycle)
298 wallcycle_start(wcycle, ewcPP_PMESENDX);
300 unsigned int flags = PP_PME_COORD;
303 flags |= PP_PME_ENER_VIR;
305 gmx_pme_send_coeffs_coords(fr, cr, flags, nullptr, nullptr, nullptr, nullptr, nullptr, nullptr,
306 box, x, lambda_q, lambda_lj, 0, 0, step, useGpuPmePpComms,
307 receiveCoordinateAddressFromPme, sendCoordinatesFromGpu,
308 coordinatesReadyOnDeviceEvent);
310 wallcycle_stop(wcycle, ewcPP_PMESENDX);
313 void gmx_pme_send_finish(const t_commrec* cr)
315 unsigned int flags = PP_PME_FINISH;
317 gmx_pme_send_coeffs_coords(nullptr, cr, flags, nullptr, nullptr, nullptr, nullptr, nullptr, nullptr,
318 nullptr, nullptr, 0, 0, 0, 0, -1, false, false, false, nullptr);
321 void gmx_pme_send_switchgrid(const t_commrec* cr, ivec grid_size, real ewaldcoeff_q, real ewaldcoeff_lj)
324 gmx_pme_comm_n_box_t cnb;
326 /* Only let one PP node signal each PME node */
327 if (cr->dd->pme_receive_vir_ener)
329 cnb.flags = PP_PME_SWITCHGRID;
330 copy_ivec(grid_size, cnb.grid_size);
331 cnb.ewaldcoeff_q = ewaldcoeff_q;
332 cnb.ewaldcoeff_lj = ewaldcoeff_lj;
334 /* We send this, uncommon, message blocking to simplify the code */
335 MPI_Send(&cnb, sizeof(cnb), MPI_BYTE, cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
338 GMX_UNUSED_VALUE(cr);
339 GMX_UNUSED_VALUE(grid_size);
340 GMX_UNUSED_VALUE(ewaldcoeff_q);
341 GMX_UNUSED_VALUE(ewaldcoeff_lj);
345 void gmx_pme_send_resetcounters(const t_commrec gmx_unused* cr, int64_t gmx_unused step)
348 gmx_pme_comm_n_box_t cnb;
350 /* Only let one PP node signal each PME node */
351 if (cr->dd->pme_receive_vir_ener)
353 cnb.flags = PP_PME_RESETCOUNTERS;
356 /* We send this, uncommon, message blocking to simplify the code */
357 MPI_Send(&cnb, sizeof(cnb), MPI_BYTE, cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
362 /*! \brief Receive virial and energy from PME rank */
363 static void receive_virial_energy(const t_commrec* cr,
364 gmx::ForceWithVirial* forceWithVirial,
371 gmx_pme_comm_vir_ene_t cve;
373 if (cr->dd->pme_receive_vir_ener)
377 fprintf(debug, "PP rank %d receiving from PME rank %d: virial and energy\n",
378 cr->sim_nodeid, cr->dd->pme_nodeid);
381 MPI_Recv(&cve, sizeof(cve), MPI_BYTE, cr->dd->pme_nodeid, 1, cr->mpi_comm_mysim, MPI_STATUS_IGNORE);
383 memset(&cve, 0, sizeof(cve));
386 forceWithVirial->addVirialContribution(cve.vir_q);
387 forceWithVirial->addVirialContribution(cve.vir_lj);
388 *energy_q = cve.energy_q;
389 *energy_lj = cve.energy_lj;
390 *dvdlambda_q += cve.dvdlambda_q;
391 *dvdlambda_lj += cve.dvdlambda_lj;
392 *pme_cycles = cve.cycles;
394 if (cve.stop_cond != gmx_stop_cond_none)
396 gmx_set_stop_condition(cve.stop_cond);
407 /*! \brief Recieve force data from PME ranks */
408 static void recvFFromPme(gmx::PmePpCommGpu* pmePpCommGpu,
412 bool useGpuPmePpComms,
413 bool receivePmeForceToGpu)
415 if (useGpuPmePpComms)
417 GMX_ASSERT(pmePpCommGpu != nullptr, "Need valid pmePpCommGpu");
418 // Receive directly using CUDA memory copy
419 pmePpCommGpu->receiveForceFromPmeCudaDirect(recvptr, n, receivePmeForceToGpu);
423 // Receive data using MPI
425 MPI_Recv(recvptr, n * sizeof(rvec), MPI_BYTE, cr->dd->pme_nodeid, 0, cr->mpi_comm_mysim,
428 GMX_UNUSED_VALUE(cr);
434 void gmx_pme_receive_f(gmx::PmePpCommGpu* pmePpCommGpu,
436 gmx::ForceWithVirial* forceWithVirial,
441 bool useGpuPmePpComms,
442 bool receivePmeForceToGpu,
445 if (c_useDelayedWait)
447 /* Wait for the x request to finish */
448 gmx_pme_send_coeffs_coords_wait(cr->dd);
451 const int natoms = dd_numHomeAtoms(*cr->dd);
452 std::vector<gmx::RVec>& buffer = cr->dd->pmeForceReceiveBuffer;
453 buffer.resize(natoms);
455 void* recvptr = reinterpret_cast<void*>(buffer.data());
456 recvFFromPme(pmePpCommGpu, recvptr, natoms, cr, useGpuPmePpComms, receivePmeForceToGpu);
458 int nt = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, natoms);
460 gmx::ArrayRef<gmx::RVec> f = forceWithVirial->force_;
462 if (!receivePmeForceToGpu)
464 /* Note that we would like to avoid this conditional by putting it
465 * into the omp pragma instead, but then we still take the full
466 * omp parallel for overhead (at least with gcc5).
470 for (int i = 0; i < natoms; i++)
477 #pragma omp parallel for num_threads(nt) schedule(static)
478 for (int i = 0; i < natoms; i++)
485 receive_virial_energy(cr, forceWithVirial, energy_q, energy_lj, dvdlambda_q, dvdlambda_lj, pme_cycles);