0e67991867eb750d8bb499a27b24887d84236b33
[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,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.
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/device_stream_manager.h"
86 #include "gromacs/gpu_utils/hostallocator.h"
87 #include "gromacs/math/gmxcomplex.h"
88 #include "gromacs/math/units.h"
89 #include "gromacs/math/vec.h"
90 #include "gromacs/mdtypes/commrec.h"
91 #include "gromacs/mdtypes/forceoutput.h"
92 #include "gromacs/mdtypes/inputrec.h"
93 #include "gromacs/mdtypes/simulation_workload.h"
94 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
95 #include "gromacs/timing/cyclecounter.h"
96 #include "gromacs/timing/wallcycle.h"
97 #include "gromacs/utility/fatalerror.h"
98 #include "gromacs/utility/futil.h"
99 #include "gromacs/utility/gmxmpi.h"
100 #include "gromacs/utility/gmxomp.h"
101 #include "gromacs/utility/smalloc.h"
102
103 #include "pme_gpu_internal.h"
104 #include "pme_output.h"
105 #include "pme_pp_communication.h"
106
107 /*! \brief Master PP-PME communication data structure */
108 struct gmx_pme_pp
109 {
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                  */
113     //@{
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;
121     //@}
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 */
124     //@{
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;
128     //@}
129
130     /*! \brief object for receiving coordinates using communications operating on GPU memory space */
131     std::unique_ptr<gmx::PmeCoordinateReceiverGpu> pmeCoordinateReceiverGpu;
132     /*! \brief object for sending PME force using communications operating on GPU memory space */
133     std::unique_ptr<gmx::PmeForceSenderGpu> pmeForceSenderGpu;
134
135     /*! \brief whether GPU direct communications are active for PME-PP transfers */
136     bool useGpuDirectComm = false;
137     /*! \brief whether GPU direct communications should send forces directly to remote GPU memory */
138     bool sendForcesDirectToPpGpu = 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*            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, WallCycleCounter::Run);
176     wallcycle_reset_all(wcycle);
177     *nrnb = { 0 };
178     wallcycle_start(wcycle, WallCycleCounter::Run);
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] computeEnergyAndVirial Set to true if this is an energy/virial calculation
230  *                                    step, otherwise set to false.
231  * \param[out] step                   MD integration step number.
232  * \param[out] grid_size              PME grid size, if received.
233  * \param[out] ewaldcoeff_q           Ewald cut-off parameter for electrostatics, if received.
234  * \param[out] ewaldcoeff_lj          Ewald cut-off parameter for Lennard-Jones, if received.
235  * \param[in]  useGpuForPme           Flag on whether PME is on GPU.
236  * \param[in]  stateGpu               GPU state propagator object.
237  * \param[in]  runMode                PME run mode.
238  *
239  * \retval pmerecvqxX                 All parameters were set, chargeA and chargeB can be NULL.
240  * \retval pmerecvqxFINISH            No parameters were set.
241  * \retval pmerecvqxSWITCHGRID        Only grid_size and *ewaldcoeff were set.
242  * \retval pmerecvqxRESETCOUNTERS     *step was set.
243  */
244 static int gmx_pme_recv_coeffs_coords(struct gmx_pme_t*            pme,
245                                       gmx_pme_pp*                  pme_pp,
246                                       int*                         natoms,
247                                       matrix                       box,
248                                       int*                         maxshift_x,
249                                       int*                         maxshift_y,
250                                       real*                        lambda_q,
251                                       real*                        lambda_lj,
252                                       gmx_bool*                    computeEnergyAndVirial,
253                                       int64_t*                     step,
254                                       ivec*                        grid_size,
255                                       real*                        ewaldcoeff_q,
256                                       real*                        ewaldcoeff_lj,
257                                       bool                         useGpuForPme,
258                                       gmx::StatePropagatorDataGpu* stateGpu,
259                                       PmeRunMode gmx_unused        runMode)
260 {
261     int status = -1;
262     int nat    = 0;
263
264 #if GMX_MPI
265     int  messages       = 0;
266     bool atomSetChanged = false;
267
268     do
269     {
270         gmx_pme_comm_n_box_t cnb;
271         cnb.flags = 0;
272
273         /* Receive the send count, box and time step from the peer PP node */
274         MPI_Recv(&cnb, sizeof(cnb), MPI_BYTE, pme_pp->peerRankId, eCommType_CNB, pme_pp->mpi_comm_mysim, MPI_STATUS_IGNORE);
275
276         *step = cnb.step;
277
278         if (debug)
279         {
280             fprintf(debug,
281                     "PME only rank receiving:%s%s%s%s%s\n",
282                     (cnb.flags & PP_PME_CHARGE) ? " charges" : "",
283                     (cnb.flags & PP_PME_COORD) ? " coordinates" : "",
284                     (cnb.flags & PP_PME_FINISH) ? " finish" : "",
285                     (cnb.flags & PP_PME_SWITCHGRID) ? " switch grid" : "",
286                     (cnb.flags & PP_PME_RESETCOUNTERS) ? " reset counters" : "");
287         }
288
289         pme_pp->useGpuDirectComm = ((cnb.flags & PP_PME_GPUCOMMS) != 0);
290         GMX_ASSERT(!pme_pp->useGpuDirectComm || (pme_pp->pmeForceSenderGpu != nullptr),
291                    "The use of GPU direct communication for PME-PP is enabled, "
292                    "but the PME GPU force reciever object does not exist");
293         pme_pp->sendForcesDirectToPpGpu = ((cnb.flags & PP_PME_RECVFTOGPU) != 0);
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,
330                               sizeof(sender.numAtoms),
331                               MPI_BYTE,
332                               sender.rankId,
333                               eCommType_CNB,
334                               pme_pp->mpi_comm_mysim,
335                               &pme_pp->req[messages++]);
336                 }
337             }
338             MPI_Waitall(messages, pme_pp->req.data(), pme_pp->stat.data());
339             messages = 0;
340
341             nat = 0;
342             for (const auto& sender : pme_pp->ppRanks)
343             {
344                 nat += sender.numAtoms;
345             }
346
347             if (cnb.flags & PP_PME_CHARGE)
348             {
349                 pme_pp->chargeA.resizeWithPadding(nat);
350             }
351             if (cnb.flags & PP_PME_CHARGEB)
352             {
353                 pme_pp->chargeB.resize(nat);
354             }
355             if (cnb.flags & PP_PME_SQRTC6)
356             {
357                 pme_pp->sqrt_c6A.resize(nat);
358             }
359             if (cnb.flags & PP_PME_SQRTC6B)
360             {
361                 pme_pp->sqrt_c6B.resize(nat);
362             }
363             if (cnb.flags & PP_PME_SIGMA)
364             {
365                 pme_pp->sigmaA.resize(nat);
366             }
367             if (cnb.flags & PP_PME_SIGMAB)
368             {
369                 pme_pp->sigmaB.resize(nat);
370             }
371             pme_pp->x.resize(nat);
372             pme_pp->f.resize(nat);
373
374             /* maxshift is sent when the charges are sent */
375             *maxshift_x = cnb.maxshift_x;
376             *maxshift_y = cnb.maxshift_y;
377
378             /* Receive the charges in place */
379             for (int q = 0; q < eCommType_NR; q++)
380             {
381                 real* bufferPtr;
382
383                 if (!(cnb.flags & (PP_PME_CHARGE << q)))
384                 {
385                     continue;
386                 }
387                 switch (q)
388                 {
389                     case eCommType_ChargeA: bufferPtr = pme_pp->chargeA.data(); break;
390                     case eCommType_ChargeB: bufferPtr = pme_pp->chargeB.data(); break;
391                     case eCommType_SQRTC6A: bufferPtr = pme_pp->sqrt_c6A.data(); break;
392                     case eCommType_SQRTC6B: bufferPtr = pme_pp->sqrt_c6B.data(); break;
393                     case eCommType_SigmaA: bufferPtr = pme_pp->sigmaA.data(); break;
394                     case eCommType_SigmaB: bufferPtr = pme_pp->sigmaB.data(); break;
395                     default: gmx_incons("Wrong eCommType");
396                 }
397                 nat = 0;
398                 for (const auto& sender : pme_pp->ppRanks)
399                 {
400                     if (sender.numAtoms > 0)
401                     {
402                         MPI_Irecv(bufferPtr + nat,
403                                   sender.numAtoms * sizeof(real),
404                                   MPI_BYTE,
405                                   sender.rankId,
406                                   q,
407                                   pme_pp->mpi_comm_mysim,
408                                   &pme_pp->req[messages++]);
409                         nat += sender.numAtoms;
410                         if (debug)
411                         {
412                             fprintf(debug,
413                                     "Received from PP rank %d: %d %s\n",
414                                     sender.rankId,
415                                     sender.numAtoms,
416                                     (q == eCommType_ChargeA || q == eCommType_ChargeB) ? "charges"
417                                                                                        : "params");
418                         }
419                     }
420                 }
421             }
422         }
423
424         if (cnb.flags & PP_PME_COORD)
425         {
426             if (atomSetChanged)
427             {
428                 gmx_pme_reinit_atoms(pme, nat, pme_pp->chargeA, pme_pp->chargeB);
429                 if (useGpuForPme)
430                 {
431                     stateGpu->reinit(nat, nat);
432                     pme_gpu_set_device_x(pme, stateGpu->getCoordinates());
433                 }
434                 if (pme_pp->useGpuDirectComm)
435                 {
436                     GMX_ASSERT(runMode == PmeRunMode::GPU,
437                                "GPU Direct PME-PP communication has been enabled, "
438                                "but PME run mode is not PmeRunMode::GPU\n");
439
440                     // This rank will have its data accessed directly by PP rank, so needs to send the remote addresses.
441                     pme_pp->pmeCoordinateReceiverGpu->sendCoordinateBufferAddressToPpRanks(
442                             stateGpu->getCoordinates());
443                     pme_pp->pmeForceSenderGpu->setForceSendBuffer(pme_gpu_get_device_f(pme));
444                 }
445             }
446
447
448             /* The box, FE flag and lambda are sent along with the coordinates
449              *  */
450             copy_mat(cnb.box, box);
451             *lambda_q               = cnb.lambda_q;
452             *lambda_lj              = cnb.lambda_lj;
453             *computeEnergyAndVirial = ((cnb.flags & PP_PME_ENER_VIR) != 0U);
454             *step                   = cnb.step;
455
456             /* Receive the coordinates in place */
457             nat = 0;
458             for (const auto& sender : pme_pp->ppRanks)
459             {
460                 if (sender.numAtoms > 0)
461                 {
462                     if (pme_pp->useGpuDirectComm)
463                     {
464                         if (GMX_THREAD_MPI)
465                         {
466                             pme_pp->pmeCoordinateReceiverGpu->receiveCoordinatesSynchronizerFromPpCudaDirect(
467                                     sender.rankId);
468                         }
469                         else
470                         {
471                             pme_pp->pmeCoordinateReceiverGpu->launchReceiveCoordinatesFromPpCudaMpi(
472                                     stateGpu->getCoordinates(), nat, sender.numAtoms * sizeof(rvec), sender.rankId);
473                         }
474                     }
475                     else
476                     {
477                         MPI_Irecv(pme_pp->x[nat],
478                                   sender.numAtoms * sizeof(rvec),
479                                   MPI_BYTE,
480                                   sender.rankId,
481                                   eCommType_COORD,
482                                   pme_pp->mpi_comm_mysim,
483                                   &pme_pp->req[messages++]);
484                     }
485                     nat += sender.numAtoms;
486                     if (debug)
487                     {
488                         fprintf(debug,
489                                 "Received from PP rank %d: %d "
490                                 "coordinates\n",
491                                 sender.rankId,
492                                 sender.numAtoms);
493                     }
494                 }
495             }
496
497             if (pme_pp->useGpuDirectComm)
498             {
499                 pme_pp->pmeCoordinateReceiverGpu->synchronizeOnCoordinatesFromPpRanks();
500             }
501
502             status = pmerecvqxX;
503         }
504
505         /* Wait for the coordinates and/or charges to arrive */
506         MPI_Waitall(messages, pme_pp->req.data(), pme_pp->stat.data());
507         messages = 0;
508     } while (status == -1);
509 #else
510     GMX_UNUSED_VALUE(pme);
511     GMX_UNUSED_VALUE(pme_pp);
512     GMX_UNUSED_VALUE(box);
513     GMX_UNUSED_VALUE(maxshift_x);
514     GMX_UNUSED_VALUE(maxshift_y);
515     GMX_UNUSED_VALUE(lambda_q);
516     GMX_UNUSED_VALUE(lambda_lj);
517     GMX_UNUSED_VALUE(computeEnergyAndVirial);
518     GMX_UNUSED_VALUE(step);
519     GMX_UNUSED_VALUE(grid_size);
520     GMX_UNUSED_VALUE(ewaldcoeff_q);
521     GMX_UNUSED_VALUE(ewaldcoeff_lj);
522     GMX_UNUSED_VALUE(useGpuForPme);
523     GMX_UNUSED_VALUE(stateGpu);
524
525     status = pmerecvqxX;
526 #endif
527
528     if (status == pmerecvqxX)
529     {
530         *natoms = nat;
531     }
532
533     return status;
534 }
535
536 /*! \brief Send the PME mesh force, virial and energy to the PP-only ranks. */
537 static void gmx_pme_send_force_vir_ener(const gmx_pme_t& pme,
538                                         gmx_pme_pp*      pme_pp,
539                                         const PmeOutput& output,
540                                         real             dvdlambda_q,
541                                         real             dvdlambda_lj,
542                                         float            cycles)
543 {
544 #if GMX_MPI
545     gmx_pme_comm_vir_ene_t cve;
546     int                    messages, ind_start, ind_end;
547     cve.cycles = cycles;
548
549     /* Now the evaluated forces have to be transferred to the PP nodes */
550     messages = 0;
551     ind_end  = 0;
552     for (const auto& receiver : pme_pp->ppRanks)
553     {
554         ind_start = ind_end;
555         ind_end   = ind_start + receiver.numAtoms;
556         if (pme_pp->useGpuDirectComm)
557         {
558             GMX_ASSERT((pme_pp->pmeForceSenderGpu != nullptr),
559                        "The use of GPU direct communication for PME-PP is enabled, "
560                        "but the PME GPU force reciever object does not exist");
561
562             if (GMX_THREAD_MPI)
563             {
564                 pme_pp->pmeForceSenderGpu->sendFToPpCudaDirect(
565                         receiver.rankId, receiver.numAtoms, pme_pp->sendForcesDirectToPpGpu);
566             }
567             else
568             {
569                 pme_pp->pmeForceSenderGpu->sendFToPpCudaMpi(pme_gpu_get_device_f(&pme),
570                                                             ind_start,
571                                                             receiver.numAtoms * sizeof(rvec),
572                                                             receiver.rankId,
573                                                             &pme_pp->req[messages]);
574
575                 messages++;
576             }
577         }
578         else
579         {
580             void* sendbuf = const_cast<void*>(static_cast<const void*>(output.forces_[ind_start]));
581             // Send using MPI
582             MPI_Isend(sendbuf,
583                       receiver.numAtoms * sizeof(rvec),
584                       MPI_BYTE,
585                       receiver.rankId,
586                       0,
587                       pme_pp->mpi_comm_mysim,
588                       &pme_pp->req[messages]);
589             messages++;
590         }
591     }
592
593     /* send virial and energy to our last PP node */
594     copy_mat(output.coulombVirial_, cve.vir_q);
595     copy_mat(output.lennardJonesVirial_, cve.vir_lj);
596     cve.energy_q     = output.coulombEnergy_;
597     cve.energy_lj    = output.lennardJonesEnergy_;
598     cve.dvdlambda_q  = dvdlambda_q;
599     cve.dvdlambda_lj = dvdlambda_lj;
600     /* check for the signals to send back to a PP node */
601     cve.stop_cond = gmx_get_stop_condition();
602
603     cve.cycles = cycles;
604
605     if (debug)
606     {
607         fprintf(debug, "PME rank sending to PP rank %d: virial and energy\n", pme_pp->peerRankId);
608     }
609     MPI_Isend(&cve, sizeof(cve), MPI_BYTE, pme_pp->peerRankId, 1, pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
610
611     /* Wait for the forces to arrive */
612     MPI_Waitall(messages, pme_pp->req.data(), pme_pp->stat.data());
613 #else
614     GMX_RELEASE_ASSERT(false, "Invalid call to gmx_pme_send_force_vir_ener");
615     GMX_UNUSED_VALUE(pme);
616     GMX_UNUSED_VALUE(pme_pp);
617     GMX_UNUSED_VALUE(output);
618     GMX_UNUSED_VALUE(dvdlambda_q);
619     GMX_UNUSED_VALUE(dvdlambda_lj);
620     GMX_UNUSED_VALUE(cycles);
621 #endif
622 }
623
624 int gmx_pmeonly(struct gmx_pme_t*               pme,
625                 const t_commrec*                cr,
626                 t_nrnb*                         mynrnb,
627                 gmx_wallcycle*                  wcycle,
628                 gmx_walltime_accounting_t       walltime_accounting,
629                 t_inputrec*                     ir,
630                 PmeRunMode                      runMode,
631                 bool                            useGpuPmePpCommunication,
632                 const gmx::DeviceStreamManager* deviceStreamManager)
633 {
634     int     ret;
635     int     natoms = 0;
636     matrix  box;
637     real    lambda_q   = 0;
638     real    lambda_lj  = 0;
639     int     maxshift_x = 0, maxshift_y = 0;
640     real    dvdlambda_q, dvdlambda_lj;
641     float   cycles;
642     int     count;
643     bool    computeEnergyAndVirial = false;
644     int64_t step;
645
646     /* This data will only use with PME tuning, i.e. switching PME grids */
647     std::vector<gmx_pme_t*> pmedata;
648     pmedata.push_back(pme);
649
650     auto pme_pp = gmx_pme_pp_init(cr);
651
652     std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
653     // TODO the variable below should be queried from the task assignment info
654     const bool useGpuForPme = (runMode == PmeRunMode::GPU) || (runMode == PmeRunMode::Mixed);
655     if (useGpuForPme)
656     {
657         GMX_RELEASE_ASSERT(
658                 deviceStreamManager != nullptr,
659                 "Device stream manager can not be nullptr when using GPU in PME-only rank.");
660         GMX_RELEASE_ASSERT(deviceStreamManager->streamIsValid(gmx::DeviceStreamType::Pme),
661                            "Device stream can not be nullptr when using GPU in PME-only rank");
662         changePinningPolicy(&pme_pp->chargeA, pme_get_pinning_policy());
663         changePinningPolicy(&pme_pp->x, pme_get_pinning_policy());
664         if (useGpuPmePpCommunication)
665         {
666             pme_pp->pmeCoordinateReceiverGpu = std::make_unique<gmx::PmeCoordinateReceiverGpu>(
667                     deviceStreamManager->stream(gmx::DeviceStreamType::Pme),
668                     pme_pp->mpi_comm_mysim,
669                     pme_pp->ppRanks);
670             pme_pp->pmeForceSenderGpu =
671                     std::make_unique<gmx::PmeForceSenderGpu>(pme_gpu_get_f_ready_synchronizer(pme),
672                                                              pme_pp->mpi_comm_mysim,
673                                                              deviceStreamManager->context(),
674                                                              pme_pp->ppRanks);
675         }
676         // TODO: Special PME-only constructor is used here. There is no mechanism to prevent from using the other constructor here.
677         //       This should be made safer.
678         stateGpu = std::make_unique<gmx::StatePropagatorDataGpu>(
679                 &deviceStreamManager->stream(gmx::DeviceStreamType::Pme),
680                 deviceStreamManager->context(),
681                 GpuApiCallBehavior::Async,
682                 pme_gpu_get_block_size(pme),
683                 wcycle);
684     }
685
686     clear_nrnb(mynrnb);
687
688     count = 0;
689     do /****** this is a quasi-loop over time steps! */
690     {
691         /* The reason for having a loop here is PME grid tuning/switching */
692         do
693         {
694             /* Domain decomposition */
695             ivec newGridSize;
696             real ewaldcoeff_q = 0, ewaldcoeff_lj = 0;
697             ret = gmx_pme_recv_coeffs_coords(pme,
698                                              pme_pp.get(),
699                                              &natoms,
700                                              box,
701                                              &maxshift_x,
702                                              &maxshift_y,
703                                              &lambda_q,
704                                              &lambda_lj,
705                                              &computeEnergyAndVirial,
706                                              &step,
707                                              &newGridSize,
708                                              &ewaldcoeff_q,
709                                              &ewaldcoeff_lj,
710                                              useGpuForPme,
711                                              stateGpu.get(),
712                                              runMode);
713
714             if (ret == pmerecvqxSWITCHGRID)
715             {
716                 /* Switch the PME grid to newGridSize */
717                 pme = gmx_pmeonly_switch(&pmedata, newGridSize, ewaldcoeff_q, ewaldcoeff_lj, cr, ir);
718             }
719
720             if (ret == pmerecvqxRESETCOUNTERS)
721             {
722                 /* Reset the cycle and flop counters */
723                 reset_pmeonly_counters(wcycle, walltime_accounting, mynrnb, step, useGpuForPme);
724             }
725         } while (ret == pmerecvqxSWITCHGRID || ret == pmerecvqxRESETCOUNTERS);
726
727         if (ret == pmerecvqxFINISH)
728         {
729             /* We should stop: break out of the loop */
730             break;
731         }
732
733         if (count == 0)
734         {
735             wallcycle_start(wcycle, WallCycleCounter::Run);
736             walltime_accounting_start_time(walltime_accounting);
737         }
738
739         wallcycle_start(wcycle, WallCycleCounter::PmeMesh);
740
741         dvdlambda_q  = 0;
742         dvdlambda_lj = 0;
743
744         // TODO Make a struct of array refs onto these per-atom fields
745         // of pme_pp (maybe box, energy and virial, too; and likewise
746         // from mdatoms for the other call to gmx_pme_do), so we have
747         // fewer lines of code and less parameter passing.
748         gmx::StepWorkload stepWork;
749         stepWork.computeVirial = computeEnergyAndVirial;
750         stepWork.computeEnergy = computeEnergyAndVirial;
751         stepWork.computeForces = true;
752         PmeOutput output       = { {}, false, 0, { { 0 } }, 0, 0, { { 0 } }, 0 };
753         if (useGpuForPme)
754         {
755             stepWork.haveDynamicBox      = false;
756             stepWork.useGpuPmeFReduction = pme_pp->useGpuDirectComm;
757             // TODO this should be set properly by gmx_pme_recv_coeffs_coords,
758             // or maybe use inputrecDynamicBox(ir), at the very least - change this when this codepath is tested!
759             pme_gpu_prepare_computation(pme, box, wcycle, stepWork);
760             if (!pme_pp->useGpuDirectComm)
761             {
762                 stateGpu->copyCoordinatesToGpu(gmx::ArrayRef<gmx::RVec>(pme_pp->x), gmx::AtomLocality::All);
763             }
764             // On the separate PME rank we do not need a synchronizer as we schedule everything in a single stream
765             // TODO: with pme on GPU the receive should make a list of synchronizers and pass it here #3157
766             auto xReadyOnDevice = nullptr;
767
768             pme_gpu_launch_spread(pme, xReadyOnDevice, wcycle, lambda_q);
769             pme_gpu_launch_complex_transforms(pme, wcycle, stepWork);
770             pme_gpu_launch_gather(pme, wcycle, lambda_q);
771             output = pme_gpu_wait_finish_task(pme, computeEnergyAndVirial, lambda_q, wcycle);
772             pme_gpu_reinit_computation(pme, wcycle);
773         }
774         else
775         {
776             GMX_ASSERT(pme_pp->x.size() == static_cast<size_t>(natoms),
777                        "The coordinate buffer should have size natoms");
778
779             gmx_pme_do(pme,
780                        pme_pp->x,
781                        pme_pp->f,
782                        pme_pp->chargeA,
783                        pme_pp->chargeB,
784                        pme_pp->sqrt_c6A,
785                        pme_pp->sqrt_c6B,
786                        pme_pp->sigmaA,
787                        pme_pp->sigmaB,
788                        box,
789                        cr,
790                        maxshift_x,
791                        maxshift_y,
792                        mynrnb,
793                        wcycle,
794                        output.coulombVirial_,
795                        output.lennardJonesVirial_,
796                        &output.coulombEnergy_,
797                        &output.lennardJonesEnergy_,
798                        lambda_q,
799                        lambda_lj,
800                        &dvdlambda_q,
801                        &dvdlambda_lj,
802                        stepWork);
803             output.forces_ = pme_pp->f;
804         }
805
806         cycles = wallcycle_stop(wcycle, WallCycleCounter::PmeMesh);
807         gmx_pme_send_force_vir_ener(*pme, pme_pp.get(), output, dvdlambda_q, dvdlambda_lj, cycles);
808
809         count++;
810     } /***** end of quasi-loop, we stop with the break above */
811     while (TRUE);
812
813     walltime_accounting_end_time(walltime_accounting);
814
815     return 0;
816 }