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