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.
37 /* IMPORTANT FOR DEVELOPERS:
39 * Triclinic pme stuff isn't entirely trivial, and we've experienced
40 * some bugs during development (many of them due to me). To avoid
41 * this in the future, please check the following things if you make
42 * changes in this file:
44 * 1. You should obtain identical (at least to the PME precision)
45 * energies, forces, and virial for
46 * a rectangular box and a triclinic one where the z (or y) axis is
47 * tilted a whole box side. For instance you could use these boxes:
49 * rectangular triclinic
54 * 2. You should check the energy conservation in a triclinic box.
56 * It might seem an overkill, but better safe than sorry.
74 #include "gromacs/domdec/domdec.h"
75 #include "gromacs/ewald/pme.h"
76 #include "gromacs/fft/parallel_3dfft.h"
77 #include "gromacs/fileio/pdbio.h"
78 #include "gromacs/gmxlib/network.h"
79 #include "gromacs/gmxlib/nrnb.h"
80 #include "gromacs/gpu_utils/hostallocator.h"
81 #include "gromacs/math/gmxcomplex.h"
82 #include "gromacs/math/units.h"
83 #include "gromacs/math/vec.h"
84 #include "gromacs/mdtypes/commrec.h"
85 #include "gromacs/mdtypes/forceoutput.h"
86 #include "gromacs/mdtypes/inputrec.h"
87 #include "gromacs/timing/cyclecounter.h"
88 #include "gromacs/timing/wallcycle.h"
89 #include "gromacs/utility/fatalerror.h"
90 #include "gromacs/utility/futil.h"
91 #include "gromacs/utility/gmxmpi.h"
92 #include "gromacs/utility/gmxomp.h"
93 #include "gromacs/utility/smalloc.h"
95 #include "pme_gpu_internal.h"
96 #include "pme_internal.h"
97 #include "pme_pp_communication.h"
99 //! Contains information about the PP ranks that partner this PME rank.
102 //! The MPI rank ID of this partner PP rank.
104 //! The number of atoms to communicate with this partner PP rank.
108 /*! \brief Master PP-PME communication data structure */
110 MPI_Comm mpi_comm_mysim; /**< MPI communicator for this simulation */
111 std::vector<PpRanks> ppRanks; /**< The PP partner ranks */
112 int peerRankId; /**< The peer PP rank id */
114 /**< Vectors of A- and B-state parameters used to transfer vectors to PME ranks */
115 gmx::PaddedHostVector<real> chargeA;
116 std::vector<real> chargeB;
117 std::vector<real> sqrt_c6A;
118 std::vector<real> sqrt_c6B;
119 std::vector<real> sigmaA;
120 std::vector<real> sigmaB;
122 gmx::HostVector<gmx::RVec> x; /**< Vector of atom coordinates to transfer to PME ranks */
123 std::vector<gmx::RVec> f; /**< Vector of atom forces received from PME ranks */
125 /**< Vectors of MPI objects used in non-blocking communication between multiple PP ranks per PME rank */
126 std::vector<MPI_Request> req;
127 std::vector<MPI_Status> stat;
131 /*! \brief Initialize the PME-only side of the PME <-> PP communication */
132 static std::unique_ptr<gmx_pme_pp> gmx_pme_pp_init(const t_commrec *cr)
134 auto pme_pp = std::make_unique<gmx_pme_pp>();
139 pme_pp->mpi_comm_mysim = cr->mpi_comm_mysim;
140 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
141 auto ppRanks = get_pme_ddranks(cr, rank);
142 pme_pp->ppRanks.reserve(ppRanks.size());
143 for (const auto &ppRankId : ppRanks)
145 pme_pp->ppRanks.push_back({ppRankId, 0});
147 // The peer PP rank is the last one.
148 pme_pp->peerRankId = pme_pp->ppRanks.back().rankId;
149 pme_pp->req.resize(eCommType_NR*pme_pp->ppRanks.size());
150 pme_pp->stat.resize(eCommType_NR*pme_pp->ppRanks.size());
152 GMX_UNUSED_VALUE(cr);
158 static void reset_pmeonly_counters(gmx_wallcycle_t wcycle,
159 gmx_walltime_accounting_t walltime_accounting,
164 /* Reset all the counters related to performance over the run */
165 wallcycle_stop(wcycle, ewcRUN);
166 wallcycle_reset_all(wcycle);
168 wallcycle_start(wcycle, ewcRUN);
169 walltime_accounting_reset_time(walltime_accounting, step);
177 static gmx_pme_t *gmx_pmeonly_switch(std::vector<gmx_pme_t *> *pmedata,
178 const ivec grid_size,
179 real ewaldcoeff_q, real ewaldcoeff_lj,
180 const t_commrec *cr, const t_inputrec *ir)
182 GMX_ASSERT(pmedata, "Bad PME tuning list pointer");
183 for (auto &pme : *pmedata)
185 GMX_ASSERT(pme, "Bad PME tuning list element pointer");
186 if (pme->nkx == grid_size[XX] &&
187 pme->nky == grid_size[YY] &&
188 pme->nkz == grid_size[ZZ])
190 /* Here we have found an existing PME data structure that suits us.
191 * However, in the GPU case, we have to reinitialize it - there's only one GPU structure.
192 * This should not cause actual GPU reallocations, at least (the allocated buffers are never shrunk).
193 * So, just some grid size updates in the GPU kernel parameters.
194 * TODO: this should be something like gmx_pme_update_split_params()
196 gmx_pme_reinit(&pme, cr, pme, ir, grid_size, ewaldcoeff_q, ewaldcoeff_lj);
201 const auto &pme = pmedata->back();
202 gmx_pme_t *newStructure = nullptr;
203 // Copy last structure with new grid params
204 gmx_pme_reinit(&newStructure, cr, pme, ir, grid_size, ewaldcoeff_q, ewaldcoeff_lj);
205 pmedata->push_back(newStructure);
209 /*! \brief Called by PME-only ranks to receive coefficients and coordinates
211 * \param[in,out] pme_pp PME-PP communication structure.
212 * \param[out] natoms Number of received atoms.
213 * \param[out] box System box, if received.
214 * \param[out] maxshift_x Maximum shift in X direction, if received.
215 * \param[out] maxshift_y Maximum shift in Y direction, if received.
216 * \param[out] lambda_q Free-energy lambda for electrostatics, if received.
217 * \param[out] lambda_lj Free-energy lambda for Lennard-Jones, if received.
218 * \param[out] bEnerVir Set to true if this is an energy/virial calculation step, otherwise set to false.
219 * \param[out] step MD integration step number.
220 * \param[out] grid_size PME grid size, if received.
221 * \param[out] ewaldcoeff_q Ewald cut-off parameter for electrostatics, if received.
222 * \param[out] ewaldcoeff_lj Ewald cut-off parameter for Lennard-Jones, if received.
223 * \param[out] atomSetChanged Set to true only if the local domain atom data (charges/coefficients)
224 * has been received (after DD) and should be reinitialized. Otherwise not changed.
226 * \retval pmerecvqxX All parameters were set, chargeA and chargeB can be NULL.
227 * \retval pmerecvqxFINISH No parameters were set.
228 * \retval pmerecvqxSWITCHGRID Only grid_size and *ewaldcoeff were set.
229 * \retval pmerecvqxRESETCOUNTERS *step was set.
231 static int gmx_pme_recv_coeffs_coords(gmx_pme_pp *pme_pp,
243 bool *atomSetChanged)
249 unsigned int flags = 0;
254 gmx_pme_comm_n_box_t cnb;
257 /* Receive the send count, box and time step from the peer PP node */
258 MPI_Recv(&cnb, sizeof(cnb), MPI_BYTE,
259 pme_pp->peerRankId, eCommType_CNB,
260 pme_pp->mpi_comm_mysim, MPI_STATUS_IGNORE);
262 /* We accumulate all received flags */
269 fprintf(debug, "PME only rank receiving:%s%s%s%s%s\n",
270 (cnb.flags & PP_PME_CHARGE) ? " charges" : "",
271 (cnb.flags & PP_PME_COORD ) ? " coordinates" : "",
272 (cnb.flags & PP_PME_FINISH) ? " finish" : "",
273 (cnb.flags & PP_PME_SWITCHGRID) ? " switch grid" : "",
274 (cnb.flags & PP_PME_RESETCOUNTERS) ? " reset counters" : "");
277 if (cnb.flags & PP_PME_FINISH)
279 status = pmerecvqxFINISH;
282 if (cnb.flags & PP_PME_SWITCHGRID)
284 /* Special case, receive the new parameters and return */
285 copy_ivec(cnb.grid_size, *grid_size);
286 *ewaldcoeff_q = cnb.ewaldcoeff_q;
287 *ewaldcoeff_lj = cnb.ewaldcoeff_lj;
289 status = pmerecvqxSWITCHGRID;
292 if (cnb.flags & PP_PME_RESETCOUNTERS)
294 /* Special case, receive the step (set above) and return */
295 status = pmerecvqxRESETCOUNTERS;
298 if (cnb.flags & (PP_PME_CHARGE | PP_PME_SQRTC6 | PP_PME_SIGMA))
300 *atomSetChanged = true;
302 /* Receive the send counts from the other PP nodes */
303 for (auto &sender : pme_pp->ppRanks)
305 if (sender.rankId == pme_pp->peerRankId)
307 sender.numAtoms = cnb.natoms;
311 MPI_Irecv(&sender.numAtoms, sizeof(sender.numAtoms),
313 sender.rankId, eCommType_CNB,
314 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
317 MPI_Waitall(messages, pme_pp->req.data(), pme_pp->stat.data());
321 for (const auto &sender : pme_pp->ppRanks)
323 nat += sender.numAtoms;
326 if (cnb.flags & PP_PME_CHARGE)
328 pme_pp->chargeA.resizeWithPadding(nat);
330 if (cnb.flags & PP_PME_CHARGEB)
332 pme_pp->chargeB.resize(nat);
334 if (cnb.flags & PP_PME_SQRTC6)
336 pme_pp->sqrt_c6A.resize(nat);
338 if (cnb.flags & PP_PME_SQRTC6B)
340 pme_pp->sqrt_c6B.resize(nat);
342 if (cnb.flags & PP_PME_SIGMA)
344 pme_pp->sigmaA.resize(nat);
346 if (cnb.flags & PP_PME_SIGMAB)
348 pme_pp->sigmaB.resize(nat);
350 pme_pp->x.resize(nat);
351 pme_pp->f.resize(nat);
353 /* maxshift is sent when the charges are sent */
354 *maxshift_x = cnb.maxshift_x;
355 *maxshift_y = cnb.maxshift_y;
357 /* Receive the charges in place */
358 for (int q = 0; q < eCommType_NR; q++)
362 if (!(cnb.flags & (PP_PME_CHARGE<<q)))
368 case eCommType_ChargeA: bufferPtr = pme_pp->chargeA.data(); break;
369 case eCommType_ChargeB: bufferPtr = pme_pp->chargeB.data(); break;
370 case eCommType_SQRTC6A: bufferPtr = pme_pp->sqrt_c6A.data(); break;
371 case eCommType_SQRTC6B: bufferPtr = pme_pp->sqrt_c6B.data(); break;
372 case eCommType_SigmaA: bufferPtr = pme_pp->sigmaA.data(); break;
373 case eCommType_SigmaB: bufferPtr = pme_pp->sigmaB.data(); break;
374 default: gmx_incons("Wrong eCommType");
377 for (const auto &sender : pme_pp->ppRanks)
379 if (sender.numAtoms > 0)
381 MPI_Irecv(bufferPtr+nat,
382 sender.numAtoms*sizeof(real),
385 pme_pp->mpi_comm_mysim,
386 &pme_pp->req[messages++]);
387 nat += sender.numAtoms;
390 fprintf(debug, "Received from PP rank %d: %d %s\n",
391 sender.rankId, sender.numAtoms,
392 (q == eCommType_ChargeA ||
393 q == eCommType_ChargeB) ? "charges" : "params");
400 if (cnb.flags & PP_PME_COORD)
402 /* The box, FE flag and lambda are sent along with the coordinates
404 copy_mat(cnb.box, box);
405 *lambda_q = cnb.lambda_q;
406 *lambda_lj = cnb.lambda_lj;
407 *bEnerVir = ((cnb.flags & PP_PME_ENER_VIR) != 0U);
410 /* Receive the coordinates in place */
412 for (const auto &sender : pme_pp->ppRanks)
414 if (sender.numAtoms > 0)
416 MPI_Irecv(pme_pp->x[nat],
417 sender.numAtoms*sizeof(rvec),
419 sender.rankId, eCommType_COORD,
420 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
421 nat += sender.numAtoms;
424 fprintf(debug, "Received from PP rank %d: %d "
426 sender.rankId, sender.numAtoms);
434 /* Wait for the coordinates and/or charges to arrive */
435 MPI_Waitall(messages, pme_pp->req.data(), pme_pp->stat.data());
438 while (status == -1);
440 GMX_UNUSED_VALUE(pme_pp);
441 GMX_UNUSED_VALUE(box);
442 GMX_UNUSED_VALUE(maxshift_x);
443 GMX_UNUSED_VALUE(maxshift_y);
444 GMX_UNUSED_VALUE(lambda_q);
445 GMX_UNUSED_VALUE(lambda_lj);
446 GMX_UNUSED_VALUE(bEnerVir);
447 GMX_UNUSED_VALUE(step);
448 GMX_UNUSED_VALUE(grid_size);
449 GMX_UNUSED_VALUE(ewaldcoeff_q);
450 GMX_UNUSED_VALUE(ewaldcoeff_lj);
451 GMX_UNUSED_VALUE(atomSetChanged);
456 if (status == pmerecvqxX)
464 /*! \brief Send the PME mesh force, virial and energy to the PP-only ranks. */
465 static void gmx_pme_send_force_vir_ener(gmx_pme_pp *pme_pp,
466 const PmeOutput &output,
467 real dvdlambda_q, real dvdlambda_lj,
471 gmx_pme_comm_vir_ene_t cve;
472 int messages, ind_start, ind_end;
475 /* Now the evaluated forces have to be transferred to the PP nodes */
478 for (const auto &receiver : pme_pp->ppRanks)
481 ind_end = ind_start + receiver.numAtoms;
482 if (MPI_Isend(const_cast<void *>(static_cast<const void *>(output.forces_[ind_start])),
483 (ind_end-ind_start)*sizeof(rvec), MPI_BYTE,
485 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]) != 0)
487 gmx_comm("MPI_Isend failed in do_pmeonly");
491 /* send virial and energy to our last PP node */
492 copy_mat(output.coulombVirial_, cve.vir_q);
493 copy_mat(output.lennardJonesVirial_, cve.vir_lj);
494 cve.energy_q = output.coulombEnergy_;
495 cve.energy_lj = output.lennardJonesEnergy_;
496 cve.dvdlambda_q = dvdlambda_q;
497 cve.dvdlambda_lj = dvdlambda_lj;
498 /* check for the signals to send back to a PP node */
499 cve.stop_cond = gmx_get_stop_condition();
505 fprintf(debug, "PME rank sending to PP rank %d: virial and energy\n",
508 MPI_Isend(&cve, sizeof(cve), MPI_BYTE,
509 pme_pp->peerRankId, 1,
510 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
512 /* Wait for the forces to arrive */
513 MPI_Waitall(messages, pme_pp->req.data(), pme_pp->stat.data());
515 gmx_call("MPI not enabled");
516 GMX_UNUSED_VALUE(pme_pp);
517 GMX_UNUSED_VALUE(output);
518 GMX_UNUSED_VALUE(dvdlambda_q);
519 GMX_UNUSED_VALUE(dvdlambda_lj);
520 GMX_UNUSED_VALUE(cycles);
524 int gmx_pmeonly(struct gmx_pme_t *pme,
525 const t_commrec *cr, t_nrnb *mynrnb,
526 gmx_wallcycle *wcycle,
527 gmx_walltime_accounting_t walltime_accounting,
528 t_inputrec *ir, PmeRunMode runMode)
535 int maxshift_x = 0, maxshift_y = 0;
536 real dvdlambda_q, dvdlambda_lj;
539 gmx_bool bEnerVir = FALSE;
542 /* This data will only use with PME tuning, i.e. switching PME grids */
543 std::vector<gmx_pme_t *> pmedata;
544 pmedata.push_back(pme);
546 auto pme_pp = gmx_pme_pp_init(cr);
547 //TODO the variable below should be queried from the task assignment info
548 const bool useGpuForPme = (runMode == PmeRunMode::GPU) || (runMode == PmeRunMode::Mixed);
551 changePinningPolicy(&pme_pp->chargeA, pme_get_pinning_policy());
552 changePinningPolicy(&pme_pp->x, pme_get_pinning_policy());
558 do /****** this is a quasi-loop over time steps! */
560 /* The reason for having a loop here is PME grid tuning/switching */
563 /* Domain decomposition */
565 bool atomSetChanged = false;
566 real ewaldcoeff_q = 0, ewaldcoeff_lj = 0;
567 ret = gmx_pme_recv_coeffs_coords(pme_pp.get(),
570 &maxshift_x, &maxshift_y,
571 &lambda_q, &lambda_lj,
579 if (ret == pmerecvqxSWITCHGRID)
581 /* Switch the PME grid to newGridSize */
582 pme = gmx_pmeonly_switch(&pmedata, newGridSize, ewaldcoeff_q, ewaldcoeff_lj, cr, ir);
587 gmx_pme_reinit_atoms(pme, natoms, pme_pp->chargeA.data());
590 if (ret == pmerecvqxRESETCOUNTERS)
592 /* Reset the cycle and flop counters */
593 reset_pmeonly_counters(wcycle, walltime_accounting, mynrnb, step, useGpuForPme);
596 while (ret == pmerecvqxSWITCHGRID || ret == pmerecvqxRESETCOUNTERS);
598 if (ret == pmerecvqxFINISH)
600 /* We should stop: break out of the loop */
606 wallcycle_start(wcycle, ewcRUN);
607 walltime_accounting_start_time(walltime_accounting);
610 wallcycle_start(wcycle, ewcPMEMESH);
615 // TODO Make a struct of array refs onto these per-atom fields
616 // of pme_pp (maybe box, energy and virial, too; and likewise
617 // from mdatoms for the other call to gmx_pme_do), so we have
618 // fewer lines of code and less parameter passing.
619 const int pmeFlags = GMX_PME_DO_ALL_F | (bEnerVir ? GMX_PME_CALC_ENER_VIR : 0);
620 PmeOutput output = {{}, false, 0, {{0}}, 0, {{0}}};
623 const bool boxChanged = false;
624 const bool useGpuPmeForceReduction = false;
625 //TODO this should be set properly by gmx_pme_recv_coeffs_coords,
626 // or maybe use inputrecDynamicBox(ir), at the very least - change this when this codepath is tested!
627 pme_gpu_prepare_computation(pme, boxChanged, box, wcycle, pmeFlags, useGpuPmeForceReduction);
628 pme_gpu_copy_coordinates_to_gpu(pme, as_rvec_array(pme_pp->x.data()), wcycle);
629 pme_gpu_launch_spread(pme, wcycle);
630 pme_gpu_launch_complex_transforms(pme, wcycle);
631 pme_gpu_launch_gather(pme, wcycle, PmeForceOutputHandling::Set);
632 output = pme_gpu_wait_finish_task(pme, pmeFlags, wcycle);
633 pme_gpu_reinit_computation(pme, wcycle);
637 GMX_ASSERT(pme_pp->x.size() == static_cast<size_t>(natoms), "The coordinate buffer should have size natoms");
639 gmx_pme_do(pme, pme_pp->x, pme_pp->f,
640 pme_pp->chargeA.data(), pme_pp->chargeB.data(),
641 pme_pp->sqrt_c6A.data(), pme_pp->sqrt_c6B.data(),
642 pme_pp->sigmaA.data(), pme_pp->sigmaB.data(), box,
643 cr, maxshift_x, maxshift_y, mynrnb, wcycle,
644 output.coulombVirial_, output.lennardJonesVirial_,
645 &output.coulombEnergy_, &output.lennardJonesEnergy_,
646 lambda_q, lambda_lj, &dvdlambda_q, &dvdlambda_lj,
648 output.forces_ = pme_pp->f;
651 cycles = wallcycle_stop(wcycle, ewcPMEMESH);
653 gmx_pme_send_force_vir_ener(pme_pp.get(), output,
654 dvdlambda_q, dvdlambda_lj, cycles);
657 } /***** end of quasi-loop, we stop with the break above */
660 walltime_accounting_end_time(walltime_accounting);