Increase const correctness in gmx_pme_send_parameters
[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<const real>      chargeA,
101                                        gmx::ArrayRef<const real>      chargeB,
102                                        gmx::ArrayRef<const real>      c6A,
103                                        gmx::ArrayRef<const real>      c6B,
104                                        gmx::ArrayRef<const real>      sigmaA,
105                                        gmx::ArrayRef<const real>      sigmaB,
106                                        const matrix                   box,
107                                        gmx::ArrayRef<const gmx::RVec> 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             if (useGpuPmePpComms && (fr != nullptr))
264             {
265                 if (sendCoordinatesFromGpu)
266                 {
267                     fr->pmePpCommGpu->sendCoordinatesToPmeFromGpu(
268                             fr->stateGpu->getCoordinates(), n, coordinatesReadyOnDeviceEvent);
269                 }
270                 else
271                 {
272                     fr->pmePpCommGpu->sendCoordinatesToPmeFromCpu(
273                             const_cast<gmx::RVec*>(x.data()), n, coordinatesReadyOnDeviceEvent);
274                 }
275             }
276             else
277             {
278                 MPI_Isend(x.data(),
279                           n * sizeof(rvec),
280                           MPI_BYTE,
281                           dd->pme_nodeid,
282                           eCommType_COORD,
283                           cr->mpi_comm_mysim,
284                           &dd->req_pme[dd->nreq_pme++]);
285             }
286         }
287     }
288 #else
289     GMX_UNUSED_VALUE(fr);
290     GMX_UNUSED_VALUE(chargeA);
291     GMX_UNUSED_VALUE(chargeB);
292     GMX_UNUSED_VALUE(c6A);
293     GMX_UNUSED_VALUE(c6B);
294     GMX_UNUSED_VALUE(sigmaA);
295     GMX_UNUSED_VALUE(sigmaB);
296     GMX_UNUSED_VALUE(reinitGpuPmePpComms);
297     GMX_UNUSED_VALUE(sendCoordinatesFromGpu);
298     GMX_UNUSED_VALUE(coordinatesReadyOnDeviceEvent);
299 #endif
300     if (!c_useDelayedWait)
301     {
302         /* Wait for the data to arrive */
303         /* We can skip this wait as we are sure x and q will not be modified
304          * before the next call to gmx_pme_send_x_q or gmx_pme_receive_f.
305          */
306         gmx_pme_send_coeffs_coords_wait(dd);
307     }
308 }
309
310 void gmx_pme_send_parameters(const t_commrec*           cr,
311                              const interaction_const_t& interactionConst,
312                              bool                       bFreeEnergy_q,
313                              bool                       bFreeEnergy_lj,
314                              gmx::ArrayRef<const real>  chargeA,
315                              gmx::ArrayRef<const real>  chargeB,
316                              gmx::ArrayRef<const real>  sqrt_c6A,
317                              gmx::ArrayRef<const real>  sqrt_c6B,
318                              gmx::ArrayRef<const real>  sigmaA,
319                              gmx::ArrayRef<const real>  sigmaB,
320                              int                        maxshift_x,
321                              int                        maxshift_y)
322 {
323     unsigned int flags = 0;
324
325     if (EEL_PME(interactionConst.eeltype))
326     {
327         flags |= PP_PME_CHARGE;
328     }
329     if (EVDW_PME(interactionConst.vdwtype))
330     {
331         flags |= (PP_PME_SQRTC6 | PP_PME_SIGMA);
332     }
333     if (bFreeEnergy_q || bFreeEnergy_lj)
334     {
335         /* Assumes that the B state flags are in the bits just above
336          * the ones for the A state. */
337         flags |= (flags << 1);
338     }
339
340     gmx_pme_send_coeffs_coords(nullptr,
341                                cr,
342                                flags,
343                                chargeA,
344                                chargeB,
345                                sqrt_c6A,
346                                sqrt_c6B,
347                                sigmaA,
348                                sigmaB,
349                                nullptr,
350                                gmx::ArrayRef<gmx::RVec>(),
351                                0,
352                                0,
353                                maxshift_x,
354                                maxshift_y,
355                                -1,
356                                false,
357                                false,
358                                false,
359                                nullptr);
360 }
361
362 void gmx_pme_send_coordinates(t_forcerec*                    fr,
363                               const t_commrec*               cr,
364                               const matrix                   box,
365                               gmx::ArrayRef<const gmx::RVec> x,
366                               real                           lambda_q,
367                               real                           lambda_lj,
368                               bool                           computeEnergyAndVirial,
369                               int64_t                        step,
370                               bool                           useGpuPmePpComms,
371                               bool                           receiveCoordinateAddressFromPme,
372                               bool                           sendCoordinatesFromGpu,
373                               GpuEventSynchronizer*          coordinatesReadyOnDeviceEvent,
374                               gmx_wallcycle*                 wcycle)
375 {
376     wallcycle_start(wcycle, WallCycleCounter::PpPmeSendX);
377
378     unsigned int flags = PP_PME_COORD;
379     if (computeEnergyAndVirial)
380     {
381         flags |= PP_PME_ENER_VIR;
382     }
383     gmx_pme_send_coeffs_coords(fr,
384                                cr,
385                                flags,
386                                {},
387                                {},
388                                {},
389                                {},
390                                {},
391                                {},
392                                box,
393                                x,
394                                lambda_q,
395                                lambda_lj,
396                                0,
397                                0,
398                                step,
399                                useGpuPmePpComms,
400                                receiveCoordinateAddressFromPme,
401                                sendCoordinatesFromGpu,
402                                coordinatesReadyOnDeviceEvent);
403
404     wallcycle_stop(wcycle, WallCycleCounter::PpPmeSendX);
405 }
406
407 void gmx_pme_send_finish(const t_commrec* cr)
408 {
409     unsigned int flags = PP_PME_FINISH;
410
411     gmx_pme_send_coeffs_coords(
412             nullptr, cr, flags, {}, {}, {}, {}, {}, {}, nullptr, gmx::ArrayRef<gmx::RVec>(), 0, 0, 0, 0, -1, false, false, false, nullptr);
413 }
414
415 void gmx_pme_send_switchgrid(const t_commrec* cr, ivec grid_size, real ewaldcoeff_q, real ewaldcoeff_lj)
416 {
417 #if GMX_MPI
418     gmx_pme_comm_n_box_t cnb;
419
420     /* Only let one PP node signal each PME node */
421     if (cr->dd->pme_receive_vir_ener)
422     {
423         cnb.flags = PP_PME_SWITCHGRID;
424         copy_ivec(grid_size, cnb.grid_size);
425         cnb.ewaldcoeff_q  = ewaldcoeff_q;
426         cnb.ewaldcoeff_lj = ewaldcoeff_lj;
427
428         /* We send this, uncommon, message blocking to simplify the code */
429         MPI_Send(&cnb, sizeof(cnb), MPI_BYTE, cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
430     }
431 #else
432     GMX_UNUSED_VALUE(cr);
433     GMX_UNUSED_VALUE(grid_size);
434     GMX_UNUSED_VALUE(ewaldcoeff_q);
435     GMX_UNUSED_VALUE(ewaldcoeff_lj);
436 #endif
437 }
438
439 void gmx_pme_send_resetcounters(const t_commrec gmx_unused* cr, int64_t gmx_unused step)
440 {
441 #if GMX_MPI
442     gmx_pme_comm_n_box_t cnb;
443
444     /* Only let one PP node signal each PME node */
445     if (cr->dd->pme_receive_vir_ener)
446     {
447         cnb.flags = PP_PME_RESETCOUNTERS;
448         cnb.step  = step;
449
450         /* We send this, uncommon, message blocking to simplify the code */
451         MPI_Send(&cnb, sizeof(cnb), MPI_BYTE, cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
452     }
453 #endif
454 }
455
456 /*! \brief Receive virial and energy from PME rank */
457 static void receive_virial_energy(const t_commrec*      cr,
458                                   gmx::ForceWithVirial* forceWithVirial,
459                                   real*                 energy_q,
460                                   real*                 energy_lj,
461                                   real*                 dvdlambda_q,
462                                   real*                 dvdlambda_lj,
463                                   float*                pme_cycles)
464 {
465     gmx_pme_comm_vir_ene_t cve;
466
467     if (cr->dd->pme_receive_vir_ener)
468     {
469         if (debug)
470         {
471             fprintf(debug,
472                     "PP rank %d receiving from PME rank %d: virial and energy\n",
473                     cr->sim_nodeid,
474                     cr->dd->pme_nodeid);
475         }
476 #if GMX_MPI
477         MPI_Recv(&cve, sizeof(cve), MPI_BYTE, cr->dd->pme_nodeid, 1, cr->mpi_comm_mysim, MPI_STATUS_IGNORE);
478 #else
479         memset(&cve, 0, sizeof(cve));
480 #endif
481
482         forceWithVirial->addVirialContribution(cve.vir_q);
483         forceWithVirial->addVirialContribution(cve.vir_lj);
484         *energy_q  = cve.energy_q;
485         *energy_lj = cve.energy_lj;
486         *dvdlambda_q += cve.dvdlambda_q;
487         *dvdlambda_lj += cve.dvdlambda_lj;
488         *pme_cycles = cve.cycles;
489
490         if (cve.stop_cond != StopCondition::None)
491         {
492             gmx_set_stop_condition(cve.stop_cond);
493         }
494     }
495     else
496     {
497         *energy_q   = 0;
498         *energy_lj  = 0;
499         *pme_cycles = 0;
500     }
501 }
502
503 /*! \brief Recieve force data from PME ranks */
504 static void recvFFromPme(gmx::PmePpCommGpu* pmePpCommGpu,
505                          void*              recvptr,
506                          int                n,
507                          const t_commrec*   cr,
508                          bool               useGpuPmePpComms,
509                          bool               receivePmeForceToGpu)
510 {
511     if (useGpuPmePpComms)
512     {
513         GMX_ASSERT(pmePpCommGpu != nullptr, "Need valid pmePpCommGpu");
514         // Receive forces from PME rank
515         pmePpCommGpu->receiveForceFromPme(static_cast<gmx::RVec*>(recvptr), n, receivePmeForceToGpu);
516     }
517     else
518     {
519         // Receive data using MPI
520 #if GMX_MPI
521         MPI_Recv(recvptr, n * sizeof(rvec), MPI_BYTE, cr->dd->pme_nodeid, 0, cr->mpi_comm_mysim, MPI_STATUS_IGNORE);
522 #else
523         GMX_UNUSED_VALUE(cr);
524 #endif
525     }
526 }
527
528
529 void gmx_pme_receive_f(gmx::PmePpCommGpu*    pmePpCommGpu,
530                        const t_commrec*      cr,
531                        gmx::ForceWithVirial* forceWithVirial,
532                        real*                 energy_q,
533                        real*                 energy_lj,
534                        real*                 dvdlambda_q,
535                        real*                 dvdlambda_lj,
536                        bool                  useGpuPmePpComms,
537                        bool                  receivePmeForceToGpu,
538                        float*                pme_cycles)
539 {
540     if (c_useDelayedWait)
541     {
542         /* Wait for the x request to finish */
543         gmx_pme_send_coeffs_coords_wait(cr->dd);
544     }
545
546     const int               natoms = dd_numHomeAtoms(*cr->dd);
547     std::vector<gmx::RVec>& buffer = cr->dd->pmeForceReceiveBuffer;
548     buffer.resize(natoms);
549
550     void* recvptr = reinterpret_cast<void*>(buffer.data());
551     recvFFromPme(pmePpCommGpu, recvptr, natoms, cr, useGpuPmePpComms, receivePmeForceToGpu);
552
553     int nt = gmx_omp_nthreads_get_simple_rvec_task(ModuleMultiThread::Default, natoms);
554
555     gmx::ArrayRef<gmx::RVec> f = forceWithVirial->force_;
556
557     if (!receivePmeForceToGpu)
558     {
559         /* Note that we would like to avoid this conditional by putting it
560          * into the omp pragma instead, but then we still take the full
561          * omp parallel for overhead (at least with gcc5).
562          */
563         if (nt == 1)
564         {
565             for (int i = 0; i < natoms; i++)
566             {
567                 f[i] += buffer[i];
568             }
569         }
570         else
571         {
572 #pragma omp parallel for num_threads(nt) schedule(static)
573             for (int i = 0; i < natoms; i++)
574             {
575                 f[i] += buffer[i];
576             }
577         }
578     }
579
580     receive_virial_energy(cr, forceWithVirial, energy_q, energy_lj, dvdlambda_q, dvdlambda_lj, pme_cycles);
581 }