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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,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 This file contains function definitions necessary for
41 * managing the offload of long-ranged PME work to separate MPI rank,
42 * for computing energies and forces (Coulomb and LJ).
44 * \author Berk Hess <hess@kth.se>
45 * \ingroup module_ewald
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/ewald/pme.h"
60 #include "gromacs/ewald/pme_pp_comm_gpu.h"
61 #include "gromacs/gmxlib/network.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/gmx_omp_nthreads.h"
64 #include "gromacs/mdtypes/commrec.h"
65 #include "gromacs/mdtypes/forceoutput.h"
66 #include "gromacs/mdtypes/forcerec.h"
67 #include "gromacs/mdtypes/interaction_const.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
70 #include "gromacs/nbnxm/nbnxm.h"
71 #include "gromacs/timing/wallcycle.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/gmxmpi.h"
74 #include "gromacs/utility/smalloc.h"
76 #include "pme_pp_communication.h"
78 /*! \brief Block to wait for communication to PME ranks to complete
80 * This should be faster with a real non-blocking MPI implementation
82 static constexpr bool c_useDelayedWait = false;
84 /*! \brief Wait for the pending data send requests to PME ranks to complete */
85 static void gmx_pme_send_coeffs_coords_wait(gmx_domdec_t* dd)
90 MPI_Waitall(dd->nreq_pme, dd->req_pme, MPI_STATUSES_IGNORE);
96 /*! \brief Send data to PME ranks */
97 static void gmx_pme_send_coeffs_coords(t_forcerec* fr,
100 gmx::ArrayRef<const real> chargeA,
101 gmx::ArrayRef<const real> chargeB,
102 gmx::ArrayRef<const real> c6A,
103 gmx::ArrayRef<const real> c6B,
104 gmx::ArrayRef<const real> sigmaA,
105 gmx::ArrayRef<const real> sigmaB,
107 gmx::ArrayRef<const gmx::RVec> x,
113 bool useGpuPmePpComms,
114 bool reinitGpuPmePpComms,
115 bool sendCoordinatesFromGpu,
116 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
119 gmx_pme_comm_n_box_t* cnb;
123 n = dd_numHomeAtoms(*dd);
128 "PP rank %d sending to PME rank %d: %d%s%s%s%s\n",
132 (flags & PP_PME_CHARGE) ? " charges" : "",
133 (flags & PP_PME_SQRTC6) ? " sqrtC6" : "",
134 (flags & PP_PME_SIGMA) ? " sigma" : "",
135 (flags & PP_PME_COORD) ? " coordinates" : "");
138 if (useGpuPmePpComms)
140 flags |= PP_PME_GPUCOMMS;
143 if (c_useDelayedWait)
145 /* We can not use cnb until pending communication has finished */
146 gmx_pme_send_coeffs_coords_wait(dd);
149 if (dd->pme_receive_vir_ener)
151 /* Peer PP node: communicate all data */
152 if (dd->cnb == nullptr)
160 cnb->maxshift_x = maxshift_x;
161 cnb->maxshift_y = maxshift_y;
162 cnb->lambda_q = lambda_q;
163 cnb->lambda_lj = lambda_lj;
165 if (flags & PP_PME_COORD)
167 copy_mat(box, cnb->box);
176 &dd->req_pme[dd->nreq_pme++]);
179 else if (flags & (PP_PME_CHARGE | PP_PME_SQRTC6 | PP_PME_SIGMA))
182 /* Communicate only the number of atoms */
189 &dd->req_pme[dd->nreq_pme++]);
196 if (flags & PP_PME_CHARGE)
198 MPI_Isend(chargeA.data(),
204 &dd->req_pme[dd->nreq_pme++]);
206 if (flags & PP_PME_CHARGEB)
208 MPI_Isend(chargeB.data(),
214 &dd->req_pme[dd->nreq_pme++]);
216 if (flags & PP_PME_SQRTC6)
218 MPI_Isend(c6A.data(),
224 &dd->req_pme[dd->nreq_pme++]);
226 if (flags & PP_PME_SQRTC6B)
228 MPI_Isend(c6B.data(),
234 &dd->req_pme[dd->nreq_pme++]);
236 if (flags & PP_PME_SIGMA)
238 MPI_Isend(sigmaA.data(),
244 &dd->req_pme[dd->nreq_pme++]);
246 if (flags & PP_PME_SIGMAB)
248 MPI_Isend(sigmaB.data(),
254 &dd->req_pme[dd->nreq_pme++]);
256 if (flags & PP_PME_COORD)
258 if (reinitGpuPmePpComms)
260 fr->pmePpCommGpu->reinit(n);
263 if (useGpuPmePpComms && (fr != nullptr))
265 if (sendCoordinatesFromGpu)
267 fr->pmePpCommGpu->sendCoordinatesToPmeFromGpu(
268 fr->stateGpu->getCoordinates(), n, coordinatesReadyOnDeviceEvent);
272 fr->pmePpCommGpu->sendCoordinatesToPmeFromCpu(
273 const_cast<gmx::RVec*>(x.data()), n, coordinatesReadyOnDeviceEvent);
284 &dd->req_pme[dd->nreq_pme++]);
289 GMX_UNUSED_VALUE(fr);
290 GMX_UNUSED_VALUE(chargeA);
291 GMX_UNUSED_VALUE(chargeB);
292 GMX_UNUSED_VALUE(c6A);
293 GMX_UNUSED_VALUE(c6B);
294 GMX_UNUSED_VALUE(sigmaA);
295 GMX_UNUSED_VALUE(sigmaB);
297 GMX_UNUSED_VALUE(reinitGpuPmePpComms);
298 GMX_UNUSED_VALUE(sendCoordinatesFromGpu);
299 GMX_UNUSED_VALUE(coordinatesReadyOnDeviceEvent);
301 if (!c_useDelayedWait)
303 /* Wait for the data to arrive */
304 /* We can skip this wait as we are sure x and q will not be modified
305 * before the next call to gmx_pme_send_x_q or gmx_pme_receive_f.
307 gmx_pme_send_coeffs_coords_wait(dd);
311 void gmx_pme_send_parameters(const t_commrec* cr,
312 const interaction_const_t& interactionConst,
315 gmx::ArrayRef<const real> chargeA,
316 gmx::ArrayRef<const real> chargeB,
317 gmx::ArrayRef<const real> sqrt_c6A,
318 gmx::ArrayRef<const real> sqrt_c6B,
319 gmx::ArrayRef<const real> sigmaA,
320 gmx::ArrayRef<const real> sigmaB,
324 unsigned int flags = 0;
326 if (EEL_PME(interactionConst.eeltype))
328 flags |= PP_PME_CHARGE;
330 if (EVDW_PME(interactionConst.vdwtype))
332 flags |= (PP_PME_SQRTC6 | PP_PME_SIGMA);
334 if (bFreeEnergy_q || bFreeEnergy_lj)
336 /* Assumes that the B state flags are in the bits just above
337 * the ones for the A state. */
338 flags |= (flags << 1);
341 gmx_pme_send_coeffs_coords(nullptr,
351 gmx::ArrayRef<gmx::RVec>(),
363 void gmx_pme_send_coordinates(t_forcerec* fr,
366 gmx::ArrayRef<const gmx::RVec> x,
369 bool computeEnergyAndVirial,
371 bool useGpuPmePpComms,
372 bool receiveCoordinateAddressFromPme,
373 bool sendCoordinatesFromGpu,
374 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent,
375 gmx_wallcycle* wcycle)
377 wallcycle_start(wcycle, WallCycleCounter::PpPmeSendX);
379 unsigned int flags = PP_PME_COORD;
380 if (computeEnergyAndVirial)
382 flags |= PP_PME_ENER_VIR;
384 gmx_pme_send_coeffs_coords(fr,
401 receiveCoordinateAddressFromPme,
402 sendCoordinatesFromGpu,
403 coordinatesReadyOnDeviceEvent);
405 wallcycle_stop(wcycle, WallCycleCounter::PpPmeSendX);
408 void gmx_pme_send_finish(const t_commrec* cr)
410 unsigned int flags = PP_PME_FINISH;
412 gmx_pme_send_coeffs_coords(
413 nullptr, cr, flags, {}, {}, {}, {}, {}, {}, nullptr, gmx::ArrayRef<gmx::RVec>(), 0, 0, 0, 0, -1, false, false, false, nullptr);
416 void gmx_pme_send_switchgrid(const t_commrec* cr, ivec grid_size, real ewaldcoeff_q, real ewaldcoeff_lj)
419 gmx_pme_comm_n_box_t cnb;
421 /* Only let one PP node signal each PME node */
422 if (cr->dd->pme_receive_vir_ener)
424 cnb.flags = PP_PME_SWITCHGRID;
425 copy_ivec(grid_size, cnb.grid_size);
426 cnb.ewaldcoeff_q = ewaldcoeff_q;
427 cnb.ewaldcoeff_lj = ewaldcoeff_lj;
429 /* We send this, uncommon, message blocking to simplify the code */
430 MPI_Send(&cnb, sizeof(cnb), MPI_BYTE, cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
433 GMX_UNUSED_VALUE(cr);
434 GMX_UNUSED_VALUE(grid_size);
435 GMX_UNUSED_VALUE(ewaldcoeff_q);
436 GMX_UNUSED_VALUE(ewaldcoeff_lj);
440 void gmx_pme_send_resetcounters(const t_commrec gmx_unused* cr, int64_t gmx_unused step)
443 gmx_pme_comm_n_box_t cnb;
445 /* Only let one PP node signal each PME node */
446 if (cr->dd->pme_receive_vir_ener)
448 cnb.flags = PP_PME_RESETCOUNTERS;
451 /* We send this, uncommon, message blocking to simplify the code */
452 MPI_Send(&cnb, sizeof(cnb), MPI_BYTE, cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
457 /*! \brief Receive virial and energy from PME rank */
458 static void receive_virial_energy(const t_commrec* cr,
459 gmx::ForceWithVirial* forceWithVirial,
466 gmx_pme_comm_vir_ene_t cve;
468 if (cr->dd->pme_receive_vir_ener)
473 "PP rank %d receiving from PME rank %d: virial and energy\n",
478 MPI_Recv(&cve, sizeof(cve), MPI_BYTE, cr->dd->pme_nodeid, 1, cr->mpi_comm_mysim, MPI_STATUS_IGNORE);
480 memset(&cve, 0, sizeof(cve));
483 forceWithVirial->addVirialContribution(cve.vir_q);
484 forceWithVirial->addVirialContribution(cve.vir_lj);
485 *energy_q = cve.energy_q;
486 *energy_lj = cve.energy_lj;
487 *dvdlambda_q += cve.dvdlambda_q;
488 *dvdlambda_lj += cve.dvdlambda_lj;
489 *pme_cycles = cve.cycles;
491 if (cve.stop_cond != StopCondition::None)
493 gmx_set_stop_condition(cve.stop_cond);
504 /*! \brief Recieve force data from PME ranks */
505 static void recvFFromPme(gmx::PmePpCommGpu* pmePpCommGpu,
509 bool useGpuPmePpComms,
510 bool receivePmeForceToGpu)
512 if (useGpuPmePpComms)
514 GMX_ASSERT(pmePpCommGpu != nullptr, "Need valid pmePpCommGpu");
515 // Receive forces from PME rank
516 pmePpCommGpu->receiveForceFromPme(static_cast<gmx::RVec*>(recvptr), n, receivePmeForceToGpu);
520 // Receive data using MPI
522 MPI_Recv(recvptr, n * sizeof(rvec), MPI_BYTE, cr->dd->pme_nodeid, 0, cr->mpi_comm_mysim, MPI_STATUS_IGNORE);
524 GMX_UNUSED_VALUE(cr);
531 void gmx_pme_receive_f(gmx::PmePpCommGpu* pmePpCommGpu,
533 gmx::ForceWithVirial* forceWithVirial,
538 bool useGpuPmePpComms,
539 bool receivePmeForceToGpu,
542 if (c_useDelayedWait)
544 /* Wait for the x request to finish */
545 gmx_pme_send_coeffs_coords_wait(cr->dd);
548 const int natoms = dd_numHomeAtoms(*cr->dd);
549 std::vector<gmx::RVec>& buffer = cr->dd->pmeForceReceiveBuffer;
550 buffer.resize(natoms);
552 void* recvptr = reinterpret_cast<void*>(buffer.data());
553 recvFFromPme(pmePpCommGpu, recvptr, natoms, cr, useGpuPmePpComms, receivePmeForceToGpu);
555 int nt = gmx_omp_nthreads_get_simple_rvec_task(ModuleMultiThread::Default, natoms);
557 gmx::ArrayRef<gmx::RVec> f = forceWithVirial->force_;
559 if (!receivePmeForceToGpu)
561 /* Note that we would like to avoid this conditional by putting it
562 * into the omp pragma instead, but then we still take the full
563 * omp parallel for overhead (at least with gcc5).
567 for (int i = 0; i < natoms; i++)
574 #pragma omp parallel for num_threads(nt) schedule(static)
575 for (int i = 0; i < natoms; i++)
582 receive_virial_energy(cr, forceWithVirial, energy_q, energy_lj, dvdlambda_q, dvdlambda_lj, pme_cycles);