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