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