3a929daf1f5437e337fa8f0122ed1990316c6923
[alexxy/gromacs.git] / src / gromacs / ewald / pme_pp.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 /*! \internal \file
39  *
40  * \brief This file contains function definitions necessary for
41  * managing the offload of long-ranged PME work to separate MPI rank,
42  * for computing energies and forces (Coulomb and LJ).
43  *
44  * \author Berk Hess <hess@kth.se>
45  * \ingroup module_ewald
46  */
47
48 #include "gmxpre.h"
49
50 #include "pme_pp.h"
51
52 #include "config.h"
53
54 #include <cstdio>
55 #include <cstring>
56
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/ewald/pme.h"
60 #include "gromacs/ewald/pme_pp_comm_gpu.h"
61 #include "gromacs/gmxlib/network.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/gmx_omp_nthreads.h"
64 #include "gromacs/mdtypes/commrec.h"
65 #include "gromacs/mdtypes/forceoutput.h"
66 #include "gromacs/mdtypes/forcerec.h"
67 #include "gromacs/mdtypes/interaction_const.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
70 #include "gromacs/nbnxm/nbnxm.h"
71 #include "gromacs/timing/wallcycle.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/gmxmpi.h"
74 #include "gromacs/utility/smalloc.h"
75
76 #include "pme_pp_communication.h"
77
78 /*! \brief Block to wait for communication to PME ranks to complete
79  *
80  * This should be faster with a real non-blocking MPI implementation
81  */
82 static constexpr bool c_useDelayedWait = false;
83
84 /*! \brief Wait for the pending data send requests to PME ranks to complete */
85 static void gmx_pme_send_coeffs_coords_wait(gmx_domdec_t* dd)
86 {
87     if (dd->nreq_pme)
88     {
89 #if GMX_MPI
90         MPI_Waitall(dd->nreq_pme, dd->req_pme, MPI_STATUSES_IGNORE);
91 #endif
92         dd->nreq_pme = 0;
93     }
94 }
95
96 /*! \brief Send data to PME ranks */
97 static void gmx_pme_send_coeffs_coords(t_forcerec*         fr,
98                                        const t_commrec*    cr,
99                                        unsigned int        flags,
100                                        gmx::ArrayRef<real> chargeA,
101                                        gmx::ArrayRef<real> chargeB,
102                                        gmx::ArrayRef<real> c6A,
103                                        gmx::ArrayRef<real> c6B,
104                                        gmx::ArrayRef<real> sigmaA,
105                                        gmx::ArrayRef<real> sigmaB,
106                                        const matrix        box,
107                                        const rvec gmx_unused* x,
108                                        real                   lambda_q,
109                                        real                   lambda_lj,
110                                        int                    maxshift_x,
111                                        int                    maxshift_y,
112                                        int64_t                step,
113                                        bool                   useGpuPmePpComms,
114                                        bool                   reinitGpuPmePpComms,
115                                        bool                   sendCoordinatesFromGpu,
116                                        GpuEventSynchronizer*  coordinatesReadyOnDeviceEvent)
117 {
118     gmx_domdec_t*         dd;
119     gmx_pme_comm_n_box_t* cnb;
120     int                   n;
121
122     dd = cr->dd;
123     n  = dd_numHomeAtoms(*dd);
124
125     if (debug)
126     {
127         fprintf(debug,
128                 "PP rank %d sending to PME rank %d: %d%s%s%s%s\n",
129                 cr->sim_nodeid,
130                 dd->pme_nodeid,
131                 n,
132                 (flags & PP_PME_CHARGE) ? " charges" : "",
133                 (flags & PP_PME_SQRTC6) ? " sqrtC6" : "",
134                 (flags & PP_PME_SIGMA) ? " sigma" : "",
135                 (flags & PP_PME_COORD) ? " coordinates" : "");
136     }
137
138     if (useGpuPmePpComms)
139     {
140         flags |= PP_PME_GPUCOMMS;
141     }
142
143     if (c_useDelayedWait)
144     {
145         /* We can not use cnb until pending communication has finished */
146         gmx_pme_send_coeffs_coords_wait(dd);
147     }
148
149     if (dd->pme_receive_vir_ener)
150     {
151         /* Peer PP node: communicate all data */
152         if (dd->cnb == nullptr)
153         {
154             snew(dd->cnb, 1);
155         }
156         cnb = dd->cnb;
157
158         cnb->flags      = flags;
159         cnb->natoms     = n;
160         cnb->maxshift_x = maxshift_x;
161         cnb->maxshift_y = maxshift_y;
162         cnb->lambda_q   = lambda_q;
163         cnb->lambda_lj  = lambda_lj;
164         cnb->step       = step;
165         if (flags & PP_PME_COORD)
166         {
167             copy_mat(box, cnb->box);
168         }
169 #if GMX_MPI
170         MPI_Isend(cnb,
171                   sizeof(*cnb),
172                   MPI_BYTE,
173                   dd->pme_nodeid,
174                   eCommType_CNB,
175                   cr->mpi_comm_mysim,
176                   &dd->req_pme[dd->nreq_pme++]);
177 #endif
178     }
179     else if (flags & (PP_PME_CHARGE | PP_PME_SQRTC6 | PP_PME_SIGMA))
180     {
181 #if GMX_MPI
182         /* Communicate only the number of atoms */
183         MPI_Isend(&n,
184                   sizeof(n),
185                   MPI_BYTE,
186                   dd->pme_nodeid,
187                   eCommType_CNB,
188                   cr->mpi_comm_mysim,
189                   &dd->req_pme[dd->nreq_pme++]);
190 #endif
191     }
192
193 #if GMX_MPI
194     if (n > 0)
195     {
196         if (flags & PP_PME_CHARGE)
197         {
198             MPI_Isend(chargeA.data(),
199                       n * sizeof(real),
200                       MPI_BYTE,
201                       dd->pme_nodeid,
202                       eCommType_ChargeA,
203                       cr->mpi_comm_mysim,
204                       &dd->req_pme[dd->nreq_pme++]);
205         }
206         if (flags & PP_PME_CHARGEB)
207         {
208             MPI_Isend(chargeB.data(),
209                       n * sizeof(real),
210                       MPI_BYTE,
211                       dd->pme_nodeid,
212                       eCommType_ChargeB,
213                       cr->mpi_comm_mysim,
214                       &dd->req_pme[dd->nreq_pme++]);
215         }
216         if (flags & PP_PME_SQRTC6)
217         {
218             MPI_Isend(c6A.data(),
219                       n * sizeof(real),
220                       MPI_BYTE,
221                       dd->pme_nodeid,
222                       eCommType_SQRTC6A,
223                       cr->mpi_comm_mysim,
224                       &dd->req_pme[dd->nreq_pme++]);
225         }
226         if (flags & PP_PME_SQRTC6B)
227         {
228             MPI_Isend(c6B.data(),
229                       n * sizeof(real),
230                       MPI_BYTE,
231                       dd->pme_nodeid,
232                       eCommType_SQRTC6B,
233                       cr->mpi_comm_mysim,
234                       &dd->req_pme[dd->nreq_pme++]);
235         }
236         if (flags & PP_PME_SIGMA)
237         {
238             MPI_Isend(sigmaA.data(),
239                       n * sizeof(real),
240                       MPI_BYTE,
241                       dd->pme_nodeid,
242                       eCommType_SigmaA,
243                       cr->mpi_comm_mysim,
244                       &dd->req_pme[dd->nreq_pme++]);
245         }
246         if (flags & PP_PME_SIGMAB)
247         {
248             MPI_Isend(sigmaB.data(),
249                       n * sizeof(real),
250                       MPI_BYTE,
251                       dd->pme_nodeid,
252                       eCommType_SigmaB,
253                       cr->mpi_comm_mysim,
254                       &dd->req_pme[dd->nreq_pme++]);
255         }
256         if (flags & PP_PME_COORD)
257         {
258             if (reinitGpuPmePpComms)
259             {
260                 fr->pmePpCommGpu->reinit(n);
261             }
262
263
264             /* MPI_Isend does not accept a const buffer pointer */
265             real* xRealPtr = const_cast<real*>(x[0]);
266             if (useGpuPmePpComms && (fr != nullptr))
267             {
268                 void* sendPtr = sendCoordinatesFromGpu
269                                         ? static_cast<void*>(fr->stateGpu->getCoordinates())
270                                         : static_cast<void*>(xRealPtr);
271                 fr->pmePpCommGpu->sendCoordinatesToPmeCudaDirect(
272                         sendPtr, n, sendCoordinatesFromGpu, coordinatesReadyOnDeviceEvent);
273             }
274             else
275             {
276                 MPI_Isend(xRealPtr,
277                           n * sizeof(rvec),
278                           MPI_BYTE,
279                           dd->pme_nodeid,
280                           eCommType_COORD,
281                           cr->mpi_comm_mysim,
282                           &dd->req_pme[dd->nreq_pme++]);
283             }
284         }
285     }
286 #else
287     GMX_UNUSED_VALUE(fr);
288     GMX_UNUSED_VALUE(chargeA);
289     GMX_UNUSED_VALUE(chargeB);
290     GMX_UNUSED_VALUE(c6A);
291     GMX_UNUSED_VALUE(c6B);
292     GMX_UNUSED_VALUE(sigmaA);
293     GMX_UNUSED_VALUE(sigmaB);
294     GMX_UNUSED_VALUE(reinitGpuPmePpComms);
295     GMX_UNUSED_VALUE(sendCoordinatesFromGpu);
296     GMX_UNUSED_VALUE(coordinatesReadyOnDeviceEvent);
297 #endif
298     if (!c_useDelayedWait)
299     {
300         /* Wait for the data to arrive */
301         /* We can skip this wait as we are sure x and q will not be modified
302          * before the next call to gmx_pme_send_x_q or gmx_pme_receive_f.
303          */
304         gmx_pme_send_coeffs_coords_wait(dd);
305     }
306 }
307
308 void gmx_pme_send_parameters(const t_commrec*           cr,
309                              const interaction_const_t& interactionConst,
310                              bool                       bFreeEnergy_q,
311                              bool                       bFreeEnergy_lj,
312                              gmx::ArrayRef<real>        chargeA,
313                              gmx::ArrayRef<real>        chargeB,
314                              gmx::ArrayRef<real>        sqrt_c6A,
315                              gmx::ArrayRef<real>        sqrt_c6B,
316                              gmx::ArrayRef<real>        sigmaA,
317                              gmx::ArrayRef<real>        sigmaB,
318                              int                        maxshift_x,
319                              int                        maxshift_y)
320 {
321     unsigned int flags = 0;
322
323     if (EEL_PME(interactionConst.eeltype))
324     {
325         flags |= PP_PME_CHARGE;
326     }
327     if (EVDW_PME(interactionConst.vdwtype))
328     {
329         flags |= (PP_PME_SQRTC6 | PP_PME_SIGMA);
330     }
331     if (bFreeEnergy_q || bFreeEnergy_lj)
332     {
333         /* Assumes that the B state flags are in the bits just above
334          * the ones for the A state. */
335         flags |= (flags << 1);
336     }
337
338     gmx_pme_send_coeffs_coords(nullptr,
339                                cr,
340                                flags,
341                                chargeA,
342                                chargeB,
343                                sqrt_c6A,
344                                sqrt_c6B,
345                                sigmaA,
346                                sigmaB,
347                                nullptr,
348                                nullptr,
349                                0,
350                                0,
351                                maxshift_x,
352                                maxshift_y,
353                                -1,
354                                false,
355                                false,
356                                false,
357                                nullptr);
358 }
359
360 void gmx_pme_send_coordinates(t_forcerec*           fr,
361                               const t_commrec*      cr,
362                               const matrix          box,
363                               const rvec*           x,
364                               real                  lambda_q,
365                               real                  lambda_lj,
366                               bool                  computeEnergyAndVirial,
367                               int64_t               step,
368                               bool                  useGpuPmePpComms,
369                               bool                  receiveCoordinateAddressFromPme,
370                               bool                  sendCoordinatesFromGpu,
371                               GpuEventSynchronizer* coordinatesReadyOnDeviceEvent,
372                               gmx_wallcycle*        wcycle)
373 {
374     wallcycle_start(wcycle, ewcPP_PMESENDX);
375
376     unsigned int flags = PP_PME_COORD;
377     if (computeEnergyAndVirial)
378     {
379         flags |= PP_PME_ENER_VIR;
380     }
381     gmx_pme_send_coeffs_coords(fr,
382                                cr,
383                                flags,
384                                {},
385                                {},
386                                {},
387                                {},
388                                {},
389                                {},
390                                box,
391                                x,
392                                lambda_q,
393                                lambda_lj,
394                                0,
395                                0,
396                                step,
397                                useGpuPmePpComms,
398                                receiveCoordinateAddressFromPme,
399                                sendCoordinatesFromGpu,
400                                coordinatesReadyOnDeviceEvent);
401
402     wallcycle_stop(wcycle, ewcPP_PMESENDX);
403 }
404
405 void gmx_pme_send_finish(const t_commrec* cr)
406 {
407     unsigned int flags = PP_PME_FINISH;
408
409     gmx_pme_send_coeffs_coords(
410             nullptr, cr, flags, {}, {}, {}, {}, {}, {}, nullptr, nullptr, 0, 0, 0, 0, -1, false, false, false, nullptr);
411 }
412
413 void gmx_pme_send_switchgrid(const t_commrec* cr, ivec grid_size, real ewaldcoeff_q, real ewaldcoeff_lj)
414 {
415 #if GMX_MPI
416     gmx_pme_comm_n_box_t cnb;
417
418     /* Only let one PP node signal each PME node */
419     if (cr->dd->pme_receive_vir_ener)
420     {
421         cnb.flags = PP_PME_SWITCHGRID;
422         copy_ivec(grid_size, cnb.grid_size);
423         cnb.ewaldcoeff_q  = ewaldcoeff_q;
424         cnb.ewaldcoeff_lj = ewaldcoeff_lj;
425
426         /* We send this, uncommon, message blocking to simplify the code */
427         MPI_Send(&cnb, sizeof(cnb), MPI_BYTE, cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
428     }
429 #else
430     GMX_UNUSED_VALUE(cr);
431     GMX_UNUSED_VALUE(grid_size);
432     GMX_UNUSED_VALUE(ewaldcoeff_q);
433     GMX_UNUSED_VALUE(ewaldcoeff_lj);
434 #endif
435 }
436
437 void gmx_pme_send_resetcounters(const t_commrec gmx_unused* cr, int64_t gmx_unused step)
438 {
439 #if GMX_MPI
440     gmx_pme_comm_n_box_t cnb;
441
442     /* Only let one PP node signal each PME node */
443     if (cr->dd->pme_receive_vir_ener)
444     {
445         cnb.flags = PP_PME_RESETCOUNTERS;
446         cnb.step  = step;
447
448         /* We send this, uncommon, message blocking to simplify the code */
449         MPI_Send(&cnb, sizeof(cnb), MPI_BYTE, cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
450     }
451 #endif
452 }
453
454 /*! \brief Receive virial and energy from PME rank */
455 static void receive_virial_energy(const t_commrec*      cr,
456                                   gmx::ForceWithVirial* forceWithVirial,
457                                   real*                 energy_q,
458                                   real*                 energy_lj,
459                                   real*                 dvdlambda_q,
460                                   real*                 dvdlambda_lj,
461                                   float*                pme_cycles)
462 {
463     gmx_pme_comm_vir_ene_t cve;
464
465     if (cr->dd->pme_receive_vir_ener)
466     {
467         if (debug)
468         {
469             fprintf(debug,
470                     "PP rank %d receiving from PME rank %d: virial and energy\n",
471                     cr->sim_nodeid,
472                     cr->dd->pme_nodeid);
473         }
474 #if GMX_MPI
475         MPI_Recv(&cve, sizeof(cve), MPI_BYTE, cr->dd->pme_nodeid, 1, cr->mpi_comm_mysim, MPI_STATUS_IGNORE);
476 #else
477         memset(&cve, 0, sizeof(cve));
478 #endif
479
480         forceWithVirial->addVirialContribution(cve.vir_q);
481         forceWithVirial->addVirialContribution(cve.vir_lj);
482         *energy_q  = cve.energy_q;
483         *energy_lj = cve.energy_lj;
484         *dvdlambda_q += cve.dvdlambda_q;
485         *dvdlambda_lj += cve.dvdlambda_lj;
486         *pme_cycles = cve.cycles;
487
488         if (cve.stop_cond != gmx_stop_cond_none)
489         {
490             gmx_set_stop_condition(cve.stop_cond);
491         }
492     }
493     else
494     {
495         *energy_q   = 0;
496         *energy_lj  = 0;
497         *pme_cycles = 0;
498     }
499 }
500
501 /*! \brief Recieve force data from PME ranks */
502 static void recvFFromPme(gmx::PmePpCommGpu* pmePpCommGpu,
503                          void*              recvptr,
504                          int                n,
505                          const t_commrec*   cr,
506                          bool               useGpuPmePpComms,
507                          bool               receivePmeForceToGpu)
508 {
509     if (useGpuPmePpComms)
510     {
511         GMX_ASSERT(pmePpCommGpu != nullptr, "Need valid pmePpCommGpu");
512         // Receive directly using CUDA memory copy
513         pmePpCommGpu->receiveForceFromPmeCudaDirect(recvptr, n, receivePmeForceToGpu);
514     }
515     else
516     {
517         // Receive data using MPI
518 #if GMX_MPI
519         MPI_Recv(recvptr, n * sizeof(rvec), MPI_BYTE, cr->dd->pme_nodeid, 0, cr->mpi_comm_mysim, MPI_STATUS_IGNORE);
520 #else
521         GMX_UNUSED_VALUE(cr);
522 #endif
523     }
524 }
525
526
527 void gmx_pme_receive_f(gmx::PmePpCommGpu*    pmePpCommGpu,
528                        const t_commrec*      cr,
529                        gmx::ForceWithVirial* forceWithVirial,
530                        real*                 energy_q,
531                        real*                 energy_lj,
532                        real*                 dvdlambda_q,
533                        real*                 dvdlambda_lj,
534                        bool                  useGpuPmePpComms,
535                        bool                  receivePmeForceToGpu,
536                        float*                pme_cycles)
537 {
538     if (c_useDelayedWait)
539     {
540         /* Wait for the x request to finish */
541         gmx_pme_send_coeffs_coords_wait(cr->dd);
542     }
543
544     const int               natoms = dd_numHomeAtoms(*cr->dd);
545     std::vector<gmx::RVec>& buffer = cr->dd->pmeForceReceiveBuffer;
546     buffer.resize(natoms);
547
548     void* recvptr = reinterpret_cast<void*>(buffer.data());
549     recvFFromPme(pmePpCommGpu, recvptr, natoms, cr, useGpuPmePpComms, receivePmeForceToGpu);
550
551     int nt = gmx_omp_nthreads_get_simple_rvec_task(emntDefault, natoms);
552
553     gmx::ArrayRef<gmx::RVec> f = forceWithVirial->force_;
554
555     if (!receivePmeForceToGpu)
556     {
557         /* Note that we would like to avoid this conditional by putting it
558          * into the omp pragma instead, but then we still take the full
559          * omp parallel for overhead (at least with gcc5).
560          */
561         if (nt == 1)
562         {
563             for (int i = 0; i < natoms; i++)
564             {
565                 f[i] += buffer[i];
566             }
567         }
568         else
569         {
570 #pragma omp parallel for num_threads(nt) schedule(static)
571             for (int i = 0; i < natoms; i++)
572             {
573                 f[i] += buffer[i];
574             }
575         }
576     }
577
578     receive_virial_energy(cr, forceWithVirial, energy_q, energy_lj, dvdlambda_q, dvdlambda_lj, pme_cycles);
579 }