Simplify PME step workload management
[alexxy/gromacs.git] / src / gromacs / ewald / pme_only.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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, 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 /* IMPORTANT FOR DEVELOPERS:
39  *
40  * Triclinic pme stuff isn't entirely trivial, and we've experienced
41  * some bugs during development (many of them due to me). To avoid
42  * this in the future, please check the following things if you make
43  * changes in this file:
44  *
45  * 1. You should obtain identical (at least to the PME precision)
46  *    energies, forces, and virial for
47  *    a rectangular box and a triclinic one where the z (or y) axis is
48  *    tilted a whole box side. For instance you could use these boxes:
49  *
50  *    rectangular       triclinic
51  *     2  0  0           2  0  0
52  *     0  2  0           0  2  0
53  *     0  0  6           2  2  6
54  *
55  * 2. You should check the energy conservation in a triclinic box.
56  *
57  * It might seem an overkill, but better safe than sorry.
58  * /Erik 001109
59  */
60
61 #include "gmxpre.h"
62
63 #include "pme_only.h"
64
65 #include "config.h"
66
67 #include <cassert>
68 #include <cmath>
69 #include <cstdio>
70 #include <cstdlib>
71 #include <cstring>
72
73 #include <memory>
74 #include <numeric>
75 #include <vector>
76
77 #include "gromacs/domdec/domdec.h"
78 #include "gromacs/ewald/pme.h"
79 #include "gromacs/ewald/pme_coordinate_receiver_gpu.h"
80 #include "gromacs/ewald/pme_force_sender_gpu.h"
81 #include "gromacs/fft/parallel_3dfft.h"
82 #include "gromacs/fileio/pdbio.h"
83 #include "gromacs/gmxlib/network.h"
84 #include "gromacs/gmxlib/nrnb.h"
85 #include "gromacs/gpu_utils/hostallocator.h"
86 #include "gromacs/math/gmxcomplex.h"
87 #include "gromacs/math/units.h"
88 #include "gromacs/math/vec.h"
89 #include "gromacs/mdtypes/commrec.h"
90 #include "gromacs/mdtypes/forceoutput.h"
91 #include "gromacs/mdtypes/inputrec.h"
92 #include "gromacs/mdtypes/simulation_workload.h"
93 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
94 #include "gromacs/timing/cyclecounter.h"
95 #include "gromacs/timing/wallcycle.h"
96 #include "gromacs/utility/fatalerror.h"
97 #include "gromacs/utility/futil.h"
98 #include "gromacs/utility/gmxmpi.h"
99 #include "gromacs/utility/gmxomp.h"
100 #include "gromacs/utility/smalloc.h"
101
102 #include "pme_gpu_internal.h"
103 #include "pme_output.h"
104 #include "pme_pp_communication.h"
105
106 /*! \brief environment variable to enable GPU P2P communication */
107 static const bool c_enableGpuPmePpComms =
108         (getenv("GMX_GPU_PME_PP_COMMS") != nullptr) && GMX_THREAD_MPI && (GMX_GPU == GMX_GPU_CUDA);
109
110 /*! \brief Master PP-PME communication data structure */
111 struct gmx_pme_pp
112 {
113     MPI_Comm             mpi_comm_mysim; /**< MPI communicator for this simulation */
114     std::vector<PpRanks> ppRanks;        /**< The PP partner ranks                 */
115     int                  peerRankId;     /**< The peer PP rank id                  */
116     //@{
117     /**< Vectors of A- and B-state parameters used to transfer vectors to PME ranks  */
118     gmx::PaddedHostVector<real> chargeA;
119     std::vector<real>           chargeB;
120     std::vector<real>           sqrt_c6A;
121     std::vector<real>           sqrt_c6B;
122     std::vector<real>           sigmaA;
123     std::vector<real>           sigmaB;
124     //@}
125     gmx::HostVector<gmx::RVec> x; /**< Vector of atom coordinates to transfer to PME ranks */
126     std::vector<gmx::RVec>     f; /**< Vector of atom forces received from PME ranks */
127     //@{
128     /**< Vectors of MPI objects used in non-blocking communication between multiple PP ranks per PME rank */
129     std::vector<MPI_Request> req;
130     std::vector<MPI_Status>  stat;
131     //@}
132
133     /*! \brief object for receiving coordinates using communications operating on GPU memory space */
134     std::unique_ptr<gmx::PmeCoordinateReceiverGpu> pmeCoordinateReceiverGpu;
135     /*! \brief object for sending PME force using communications operating on GPU memory space */
136     std::unique_ptr<gmx::PmeForceSenderGpu> pmeForceSenderGpu;
137
138     /*! \brief whether GPU direct communications are active for PME-PP transfers */
139     bool useGpuDirectComm = false;
140 };
141
142 /*! \brief Initialize the PME-only side of the PME <-> PP communication */
143 static std::unique_ptr<gmx_pme_pp> gmx_pme_pp_init(const t_commrec* cr)
144 {
145     auto pme_pp = std::make_unique<gmx_pme_pp>();
146
147 #if GMX_MPI
148     int rank;
149
150     pme_pp->mpi_comm_mysim = cr->mpi_comm_mysim;
151     MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
152     auto ppRanks = get_pme_ddranks(cr, rank);
153     pme_pp->ppRanks.reserve(ppRanks.size());
154     for (const auto& ppRankId : ppRanks)
155     {
156         pme_pp->ppRanks.push_back({ ppRankId, 0 });
157     }
158     // The peer PP rank is the last one.
159     pme_pp->peerRankId = pme_pp->ppRanks.back().rankId;
160     pme_pp->req.resize(eCommType_NR * pme_pp->ppRanks.size());
161     pme_pp->stat.resize(eCommType_NR * pme_pp->ppRanks.size());
162 #else
163     GMX_UNUSED_VALUE(cr);
164 #endif
165
166     return pme_pp;
167 }
168
169 static void reset_pmeonly_counters(gmx_wallcycle_t           wcycle,
170                                    gmx_walltime_accounting_t walltime_accounting,
171                                    t_nrnb*                   nrnb,
172                                    int64_t                   step,
173                                    bool                      useGpuForPme)
174 {
175     /* Reset all the counters related to performance over the run */
176     wallcycle_stop(wcycle, ewcRUN);
177     wallcycle_reset_all(wcycle);
178     *nrnb = { 0 };
179     wallcycle_start(wcycle, ewcRUN);
180     walltime_accounting_reset_time(walltime_accounting, step);
181
182     if (useGpuForPme)
183     {
184         resetGpuProfiler();
185     }
186 }
187
188 static gmx_pme_t* gmx_pmeonly_switch(std::vector<gmx_pme_t*>* pmedata,
189                                      const ivec               grid_size,
190                                      real                     ewaldcoeff_q,
191                                      real                     ewaldcoeff_lj,
192                                      const t_commrec*         cr,
193                                      const t_inputrec*        ir)
194 {
195     GMX_ASSERT(pmedata, "Bad PME tuning list pointer");
196     for (auto& pme : *pmedata)
197     {
198         GMX_ASSERT(pme, "Bad PME tuning list element pointer");
199         if (gmx_pme_grid_matches(*pme, grid_size))
200         {
201             /* Here we have found an existing PME data structure that suits us.
202              * However, in the GPU case, we have to reinitialize it - there's only one GPU structure.
203              * This should not cause actual GPU reallocations, at least (the allocated buffers are never shrunk).
204              * So, just some grid size updates in the GPU kernel parameters.
205              * TODO: this should be something like gmx_pme_update_split_params()
206              */
207             gmx_pme_reinit(&pme, cr, pme, ir, grid_size, ewaldcoeff_q, ewaldcoeff_lj);
208             return pme;
209         }
210     }
211
212     const auto& pme          = pmedata->back();
213     gmx_pme_t*  newStructure = nullptr;
214     // Copy last structure with new grid params
215     gmx_pme_reinit(&newStructure, cr, pme, ir, grid_size, ewaldcoeff_q, ewaldcoeff_lj);
216     pmedata->push_back(newStructure);
217     return newStructure;
218 }
219
220 /*! \brief Called by PME-only ranks to receive coefficients and coordinates
221  *
222  * \param[in] pme           PME data structure.
223  * \param[in,out] pme_pp    PME-PP communication structure.
224  * \param[out] natoms       Number of received atoms.
225  * \param[out] box        System box, if received.
226  * \param[out] maxshift_x        Maximum shift in X direction, if received.
227  * \param[out] maxshift_y        Maximum shift in Y direction, if received.
228  * \param[out] lambda_q         Free-energy lambda for electrostatics, if received.
229  * \param[out] lambda_lj         Free-energy lambda for Lennard-Jones, if received.
230  * \param[out] bEnerVir          Set to true if this is an energy/virial calculation step, otherwise
231  * set to false. \param[out] step              MD integration step number. \param[out] grid_size PME
232  * grid size, if received. \param[out] ewaldcoeff_q         Ewald cut-off parameter for
233  * electrostatics, if received. \param[out] ewaldcoeff_lj         Ewald cut-off parameter for
234  * Lennard-Jones, if received. \param[in] useGpuForPme      flag on whether PME is on GPU \param[in]
235  * stateGpu          GPU state propagator object \param[in] runMode           PME run mode
236  *
237  * \retval pmerecvqxX             All parameters were set, chargeA and chargeB can be NULL.
238  * \retval pmerecvqxFINISH        No parameters were set.
239  * \retval pmerecvqxSWITCHGRID    Only grid_size and *ewaldcoeff were set.
240  * \retval pmerecvqxRESETCOUNTERS *step was set.
241  */
242 static int gmx_pme_recv_coeffs_coords(struct gmx_pme_t*            pme,
243                                       gmx_pme_pp*                  pme_pp,
244                                       int*                         natoms,
245                                       matrix                       box,
246                                       int*                         maxshift_x,
247                                       int*                         maxshift_y,
248                                       real*                        lambda_q,
249                                       real*                        lambda_lj,
250                                       gmx_bool*                    bEnerVir,
251                                       int64_t*                     step,
252                                       ivec*                        grid_size,
253                                       real*                        ewaldcoeff_q,
254                                       real*                        ewaldcoeff_lj,
255                                       bool                         useGpuForPme,
256                                       gmx::StatePropagatorDataGpu* stateGpu,
257                                       PmeRunMode gmx_unused runMode)
258 {
259     int status = -1;
260     int nat    = 0;
261
262 #if GMX_MPI
263     unsigned int flags          = 0;
264     int          messages       = 0;
265     bool         atomSetChanged = false;
266
267     do
268     {
269         gmx_pme_comm_n_box_t cnb;
270         cnb.flags = 0;
271
272         /* Receive the send count, box and time step from the peer PP node */
273         MPI_Recv(&cnb, sizeof(cnb), MPI_BYTE, pme_pp->peerRankId, eCommType_CNB,
274                  pme_pp->mpi_comm_mysim, MPI_STATUS_IGNORE);
275
276         /* We accumulate all received flags */
277         flags |= cnb.flags;
278
279         *step = cnb.step;
280
281         if (debug)
282         {
283             fprintf(debug, "PME only rank receiving:%s%s%s%s%s\n",
284                     (cnb.flags & PP_PME_CHARGE) ? " charges" : "",
285                     (cnb.flags & PP_PME_COORD) ? " coordinates" : "",
286                     (cnb.flags & PP_PME_FINISH) ? " finish" : "",
287                     (cnb.flags & PP_PME_SWITCHGRID) ? " switch grid" : "",
288                     (cnb.flags & PP_PME_RESETCOUNTERS) ? " reset counters" : "");
289         }
290
291         pme_pp->useGpuDirectComm = ((cnb.flags & PP_PME_GPUCOMMS) != 0);
292         GMX_ASSERT(!pme_pp->useGpuDirectComm || (pme_pp->pmeForceSenderGpu != nullptr),
293                    "The use of GPU direct communication for PME-PP is enabled, "
294                    "but the PME GPU force reciever object does not exist");
295
296         if (cnb.flags & PP_PME_FINISH)
297         {
298             status = pmerecvqxFINISH;
299         }
300
301         if (cnb.flags & PP_PME_SWITCHGRID)
302         {
303             /* Special case, receive the new parameters and return */
304             copy_ivec(cnb.grid_size, *grid_size);
305             *ewaldcoeff_q  = cnb.ewaldcoeff_q;
306             *ewaldcoeff_lj = cnb.ewaldcoeff_lj;
307
308             status = pmerecvqxSWITCHGRID;
309         }
310
311         if (cnb.flags & PP_PME_RESETCOUNTERS)
312         {
313             /* Special case, receive the step (set above) and return */
314             status = pmerecvqxRESETCOUNTERS;
315         }
316
317         if (cnb.flags & (PP_PME_CHARGE | PP_PME_SQRTC6 | PP_PME_SIGMA))
318         {
319             atomSetChanged = true;
320
321             /* Receive the send counts from the other PP nodes */
322             for (auto& sender : pme_pp->ppRanks)
323             {
324                 if (sender.rankId == pme_pp->peerRankId)
325                 {
326                     sender.numAtoms = cnb.natoms;
327                 }
328                 else
329                 {
330                     MPI_Irecv(&sender.numAtoms, sizeof(sender.numAtoms), MPI_BYTE, sender.rankId,
331                               eCommType_CNB, pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
332                 }
333             }
334             MPI_Waitall(messages, pme_pp->req.data(), pme_pp->stat.data());
335             messages = 0;
336
337             nat = 0;
338             for (const auto& sender : pme_pp->ppRanks)
339             {
340                 nat += sender.numAtoms;
341             }
342
343             if (cnb.flags & PP_PME_CHARGE)
344             {
345                 pme_pp->chargeA.resizeWithPadding(nat);
346             }
347             if (cnb.flags & PP_PME_CHARGEB)
348             {
349                 pme_pp->chargeB.resize(nat);
350             }
351             if (cnb.flags & PP_PME_SQRTC6)
352             {
353                 pme_pp->sqrt_c6A.resize(nat);
354             }
355             if (cnb.flags & PP_PME_SQRTC6B)
356             {
357                 pme_pp->sqrt_c6B.resize(nat);
358             }
359             if (cnb.flags & PP_PME_SIGMA)
360             {
361                 pme_pp->sigmaA.resize(nat);
362             }
363             if (cnb.flags & PP_PME_SIGMAB)
364             {
365                 pme_pp->sigmaB.resize(nat);
366             }
367             pme_pp->x.resize(nat);
368             pme_pp->f.resize(nat);
369
370             /* maxshift is sent when the charges are sent */
371             *maxshift_x = cnb.maxshift_x;
372             *maxshift_y = cnb.maxshift_y;
373
374             /* Receive the charges in place */
375             for (int q = 0; q < eCommType_NR; q++)
376             {
377                 real* bufferPtr;
378
379                 if (!(cnb.flags & (PP_PME_CHARGE << q)))
380                 {
381                     continue;
382                 }
383                 switch (q)
384                 {
385                     case eCommType_ChargeA: bufferPtr = pme_pp->chargeA.data(); break;
386                     case eCommType_ChargeB: bufferPtr = pme_pp->chargeB.data(); break;
387                     case eCommType_SQRTC6A: bufferPtr = pme_pp->sqrt_c6A.data(); break;
388                     case eCommType_SQRTC6B: bufferPtr = pme_pp->sqrt_c6B.data(); break;
389                     case eCommType_SigmaA: bufferPtr = pme_pp->sigmaA.data(); break;
390                     case eCommType_SigmaB: bufferPtr = pme_pp->sigmaB.data(); break;
391                     default: gmx_incons("Wrong eCommType");
392                 }
393                 nat = 0;
394                 for (const auto& sender : pme_pp->ppRanks)
395                 {
396                     if (sender.numAtoms > 0)
397                     {
398                         MPI_Irecv(bufferPtr + nat, sender.numAtoms * sizeof(real), MPI_BYTE,
399                                   sender.rankId, q, pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
400                         nat += sender.numAtoms;
401                         if (debug)
402                         {
403                             fprintf(debug, "Received from PP rank %d: %d %s\n", sender.rankId,
404                                     sender.numAtoms,
405                                     (q == eCommType_ChargeA || q == eCommType_ChargeB) ? "charges"
406                                                                                        : "params");
407                         }
408                     }
409                 }
410             }
411         }
412
413         if (cnb.flags & PP_PME_COORD)
414         {
415             if (atomSetChanged)
416             {
417                 gmx_pme_reinit_atoms(pme, nat, pme_pp->chargeA.data());
418                 if (useGpuForPme)
419                 {
420                     stateGpu->reinit(nat, nat);
421                     pme_gpu_set_device_x(pme, stateGpu->getCoordinates());
422                 }
423                 if (pme_pp->useGpuDirectComm)
424                 {
425                     GMX_ASSERT(runMode == PmeRunMode::GPU,
426                                "GPU Direct PME-PP communication has been enabled, "
427                                "but PME run mode is not PmeRunMode::GPU\n");
428
429                     // This rank will have its data accessed directly by PP rank, so needs to send the remote addresses.
430                     pme_pp->pmeCoordinateReceiverGpu->sendCoordinateBufferAddressToPpRanks(
431                             stateGpu->getCoordinates());
432                     pme_pp->pmeForceSenderGpu->sendForceBufferAddressToPpRanks(
433                             reinterpret_cast<rvec*>(pme_gpu_get_device_f(pme)));
434                 }
435             }
436
437
438             /* The box, FE flag and lambda are sent along with the coordinates
439              *  */
440             copy_mat(cnb.box, box);
441             *lambda_q  = cnb.lambda_q;
442             *lambda_lj = cnb.lambda_lj;
443             *bEnerVir  = ((cnb.flags & PP_PME_ENER_VIR) != 0U);
444             *step      = cnb.step;
445
446             /* Receive the coordinates in place */
447             nat = 0;
448             for (const auto& sender : pme_pp->ppRanks)
449             {
450                 if (sender.numAtoms > 0)
451                 {
452                     if (pme_pp->useGpuDirectComm)
453                     {
454                         pme_pp->pmeCoordinateReceiverGpu->launchReceiveCoordinatesFromPpCudaDirect(
455                                 sender.rankId);
456                     }
457                     else
458                     {
459                         MPI_Irecv(pme_pp->x[nat], sender.numAtoms * sizeof(rvec), MPI_BYTE, sender.rankId,
460                                   eCommType_COORD, pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
461                     }
462                     nat += sender.numAtoms;
463                     if (debug)
464                     {
465                         fprintf(debug,
466                                 "Received from PP rank %d: %d "
467                                 "coordinates\n",
468                                 sender.rankId, sender.numAtoms);
469                     }
470                 }
471             }
472
473             if (pme_pp->useGpuDirectComm)
474             {
475                 pme_pp->pmeCoordinateReceiverGpu->enqueueWaitReceiveCoordinatesFromPpCudaDirect();
476             }
477
478             status = pmerecvqxX;
479         }
480
481         /* Wait for the coordinates and/or charges to arrive */
482         MPI_Waitall(messages, pme_pp->req.data(), pme_pp->stat.data());
483         messages = 0;
484     } while (status == -1);
485 #else
486     GMX_UNUSED_VALUE(pme);
487     GMX_UNUSED_VALUE(pme_pp);
488     GMX_UNUSED_VALUE(box);
489     GMX_UNUSED_VALUE(maxshift_x);
490     GMX_UNUSED_VALUE(maxshift_y);
491     GMX_UNUSED_VALUE(lambda_q);
492     GMX_UNUSED_VALUE(lambda_lj);
493     GMX_UNUSED_VALUE(bEnerVir);
494     GMX_UNUSED_VALUE(step);
495     GMX_UNUSED_VALUE(grid_size);
496     GMX_UNUSED_VALUE(ewaldcoeff_q);
497     GMX_UNUSED_VALUE(ewaldcoeff_lj);
498     GMX_UNUSED_VALUE(useGpuForPme);
499     GMX_UNUSED_VALUE(stateGpu);
500
501     status = pmerecvqxX;
502 #endif
503
504     if (status == pmerecvqxX)
505     {
506         *natoms = nat;
507     }
508
509     return status;
510 }
511
512 #if GMX_MPI
513 /*! \brief Send force data to PP ranks */
514 static void sendFToPP(void* sendbuf, PpRanks receiver, gmx_pme_pp* pme_pp, int* messages)
515 {
516
517     if (pme_pp->useGpuDirectComm)
518     {
519         GMX_ASSERT((pme_pp->pmeForceSenderGpu != nullptr),
520                    "The use of GPU direct communication for PME-PP is enabled, "
521                    "but the PME GPU force reciever object does not exist");
522
523         pme_pp->pmeForceSenderGpu->sendFToPpCudaDirect(receiver.rankId);
524     }
525     else
526     {
527         // Send using MPI
528         MPI_Isend(sendbuf, receiver.numAtoms * sizeof(rvec), MPI_BYTE, receiver.rankId, 0,
529                   pme_pp->mpi_comm_mysim, &pme_pp->req[*messages]);
530         *messages = *messages + 1;
531     }
532 }
533 #endif
534
535 /*! \brief Send the PME mesh force, virial and energy to the PP-only ranks. */
536 static void gmx_pme_send_force_vir_ener(const gmx_pme_t& pme,
537                                         gmx_pme_pp*      pme_pp,
538                                         const PmeOutput& output,
539                                         real             dvdlambda_q,
540                                         real             dvdlambda_lj,
541                                         float            cycles)
542 {
543 #if GMX_MPI
544     gmx_pme_comm_vir_ene_t cve;
545     int                    messages, ind_start, ind_end;
546     cve.cycles = cycles;
547
548     /* Now the evaluated forces have to be transferred to the PP nodes */
549     messages = 0;
550     ind_end  = 0;
551     for (const auto& receiver : pme_pp->ppRanks)
552     {
553         ind_start     = ind_end;
554         ind_end       = ind_start + receiver.numAtoms;
555         void* sendbuf = const_cast<void*>(static_cast<const void*>(output.forces_[ind_start]));
556         if (pme_pp->useGpuDirectComm)
557         {
558             // Data will be transferred directly from GPU.
559             rvec* d_f = reinterpret_cast<rvec*>(pme_gpu_get_device_f(&pme));
560             sendbuf   = reinterpret_cast<void*>(&d_f[ind_start]);
561         }
562         sendFToPP(sendbuf, receiver, pme_pp, &messages);
563     }
564
565     /* send virial and energy to our last PP node */
566     copy_mat(output.coulombVirial_, cve.vir_q);
567     copy_mat(output.lennardJonesVirial_, cve.vir_lj);
568     cve.energy_q     = output.coulombEnergy_;
569     cve.energy_lj    = output.lennardJonesEnergy_;
570     cve.dvdlambda_q  = dvdlambda_q;
571     cve.dvdlambda_lj = dvdlambda_lj;
572     /* check for the signals to send back to a PP node */
573     cve.stop_cond = gmx_get_stop_condition();
574
575     cve.cycles = cycles;
576
577     if (debug)
578     {
579         fprintf(debug, "PME rank sending to PP rank %d: virial and energy\n", pme_pp->peerRankId);
580     }
581     MPI_Isend(&cve, sizeof(cve), MPI_BYTE, pme_pp->peerRankId, 1, pme_pp->mpi_comm_mysim,
582               &pme_pp->req[messages++]);
583
584     /* Wait for the forces to arrive */
585     MPI_Waitall(messages, pme_pp->req.data(), pme_pp->stat.data());
586 #else
587     gmx_call("MPI not enabled");
588     GMX_UNUSED_VALUE(pme);
589     GMX_UNUSED_VALUE(pme_pp);
590     GMX_UNUSED_VALUE(output);
591     GMX_UNUSED_VALUE(dvdlambda_q);
592     GMX_UNUSED_VALUE(dvdlambda_lj);
593     GMX_UNUSED_VALUE(cycles);
594 #endif
595 }
596
597 int gmx_pmeonly(struct gmx_pme_t*         pme,
598                 const t_commrec*          cr,
599                 t_nrnb*                   mynrnb,
600                 gmx_wallcycle*            wcycle,
601                 gmx_walltime_accounting_t walltime_accounting,
602                 t_inputrec*               ir,
603                 PmeRunMode                runMode)
604 {
605     int      ret;
606     int      natoms = 0;
607     matrix   box;
608     real     lambda_q   = 0;
609     real     lambda_lj  = 0;
610     int      maxshift_x = 0, maxshift_y = 0;
611     real     dvdlambda_q, dvdlambda_lj;
612     float    cycles;
613     int      count;
614     gmx_bool bEnerVir = FALSE;
615     int64_t  step;
616
617     /* This data will only use with PME tuning, i.e. switching PME grids */
618     std::vector<gmx_pme_t*> pmedata;
619     pmedata.push_back(pme);
620
621     auto pme_pp = gmx_pme_pp_init(cr);
622
623     std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
624     // TODO the variable below should be queried from the task assignment info
625     const bool useGpuForPme = (runMode == PmeRunMode::GPU) || (runMode == PmeRunMode::Mixed);
626     if (useGpuForPme)
627     {
628         const void* commandStream = pme_gpu_get_device_stream(pme);
629         const void* deviceContext = pme_gpu_get_device_context(pme);
630
631         changePinningPolicy(&pme_pp->chargeA, pme_get_pinning_policy());
632         changePinningPolicy(&pme_pp->x, pme_get_pinning_policy());
633         if (c_enableGpuPmePpComms)
634         {
635             pme_pp->pmeCoordinateReceiverGpu = std::make_unique<gmx::PmeCoordinateReceiverGpu>(
636                     commandStream, pme_pp->mpi_comm_mysim, pme_pp->ppRanks);
637             pme_pp->pmeForceSenderGpu = std::make_unique<gmx::PmeForceSenderGpu>(
638                     commandStream, pme_pp->mpi_comm_mysim, pme_pp->ppRanks);
639         }
640         // TODO: Special PME-only constructor is used here. There is no mechanism to prevent from using the other constructor here.
641         //       This should be made safer.
642         stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
643                 commandStream, deviceContext, GpuApiCallBehavior::Async,
644                 pme_gpu_get_padding_size(pme), wcycle);
645     }
646
647     clear_nrnb(mynrnb);
648
649     count = 0;
650     do /****** this is a quasi-loop over time steps! */
651     {
652         /* The reason for having a loop here is PME grid tuning/switching */
653         do
654         {
655             /* Domain decomposition */
656             ivec newGridSize;
657             real ewaldcoeff_q = 0, ewaldcoeff_lj = 0;
658             ret = gmx_pme_recv_coeffs_coords(pme, pme_pp.get(), &natoms, box, &maxshift_x,
659                                              &maxshift_y, &lambda_q, &lambda_lj, &bEnerVir, &step,
660                                              &newGridSize, &ewaldcoeff_q, &ewaldcoeff_lj,
661                                              useGpuForPme, stateGpu.get(), runMode);
662
663             if (ret == pmerecvqxSWITCHGRID)
664             {
665                 /* Switch the PME grid to newGridSize */
666                 pme = gmx_pmeonly_switch(&pmedata, newGridSize, ewaldcoeff_q, ewaldcoeff_lj, cr, ir);
667             }
668
669             if (ret == pmerecvqxRESETCOUNTERS)
670             {
671                 /* Reset the cycle and flop counters */
672                 reset_pmeonly_counters(wcycle, walltime_accounting, mynrnb, step, useGpuForPme);
673             }
674         } while (ret == pmerecvqxSWITCHGRID || ret == pmerecvqxRESETCOUNTERS);
675
676         if (ret == pmerecvqxFINISH)
677         {
678             /* We should stop: break out of the loop */
679             break;
680         }
681
682         if (count == 0)
683         {
684             wallcycle_start(wcycle, ewcRUN);
685             walltime_accounting_start_time(walltime_accounting);
686         }
687
688         wallcycle_start(wcycle, ewcPMEMESH);
689
690         dvdlambda_q  = 0;
691         dvdlambda_lj = 0;
692
693         // TODO Make a struct of array refs onto these per-atom fields
694         // of pme_pp (maybe box, energy and virial, too; and likewise
695         // from mdatoms for the other call to gmx_pme_do), so we have
696         // fewer lines of code and less parameter passing.
697         gmx::StepWorkload stepWork;
698         stepWork.computeVirial = bEnerVir;
699         stepWork.computeEnergy = bEnerVir;
700         stepWork.computeForces = true;
701         PmeOutput output;
702         if (useGpuForPme)
703         {
704             stepWork.haveDynamicBox      = false;
705             stepWork.useGpuPmeFReduction = pme_pp->useGpuDirectComm;
706             // TODO this should be set properly by gmx_pme_recv_coeffs_coords,
707             // or maybe use inputrecDynamicBox(ir), at the very least - change this when this codepath is tested!
708             pme_gpu_prepare_computation(pme, box, wcycle, stepWork);
709             if (!pme_pp->useGpuDirectComm)
710             {
711                 stateGpu->copyCoordinatesToGpu(gmx::ArrayRef<gmx::RVec>(pme_pp->x), gmx::AtomLocality::All);
712             }
713             // On the separate PME rank we do not need a synchronizer as we schedule everything in a single stream
714             // TODO: with pme on GPU the receive should make a list of synchronizers and pass it here #3157
715             auto xReadyOnDevice = nullptr;
716
717             pme_gpu_launch_spread(pme, xReadyOnDevice, wcycle);
718             pme_gpu_launch_complex_transforms(pme, wcycle, stepWork);
719             pme_gpu_launch_gather(pme, wcycle);
720             output = pme_gpu_wait_finish_task(pme, bEnerVir, wcycle);
721             pme_gpu_reinit_computation(pme, wcycle);
722         }
723         else
724         {
725             GMX_ASSERT(pme_pp->x.size() == static_cast<size_t>(natoms),
726                        "The coordinate buffer should have size natoms");
727
728             gmx_pme_do(pme, pme_pp->x, pme_pp->f, pme_pp->chargeA.data(), pme_pp->chargeB.data(),
729                        pme_pp->sqrt_c6A.data(), pme_pp->sqrt_c6B.data(), pme_pp->sigmaA.data(),
730                        pme_pp->sigmaB.data(), box, cr, maxshift_x, maxshift_y, mynrnb, wcycle,
731                        output.coulombVirial_, output.lennardJonesVirial_, &output.coulombEnergy_,
732                        &output.lennardJonesEnergy_, lambda_q, lambda_lj, &dvdlambda_q,
733                        &dvdlambda_lj, stepWork);
734             output.forces_ = pme_pp->f;
735         }
736
737         cycles = wallcycle_stop(wcycle, ewcPMEMESH);
738         gmx_pme_send_force_vir_ener(*pme, pme_pp.get(), output, dvdlambda_q, dvdlambda_lj, cycles);
739
740         count++;
741     } /***** end of quasi-loop, we stop with the break above */
742     while (TRUE);
743
744     walltime_accounting_end_time(walltime_accounting);
745
746     return 0;
747 }