Refactor md_enums
[alexxy/gromacs.git] / src / gromacs / ewald / pme.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 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  * computing energies and forces for the PME long-ranged part (Coulomb
42  * and LJ).
43  *
44  * \author Erik Lindahl <erik@kth.se>
45  * \author Berk Hess <hess@kth.se>
46  * \ingroup module_ewald
47  */
48 /* IMPORTANT FOR DEVELOPERS:
49  *
50  * Triclinic pme stuff isn't entirely trivial, and we've experienced
51  * some bugs during development (many of them due to me). To avoid
52  * this in the future, please check the following things if you make
53  * changes in this file:
54  *
55  * 1. You should obtain identical (at least to the PME precision)
56  *    energies, forces, and virial for
57  *    a rectangular box and a triclinic one where the z (or y) axis is
58  *    tilted a whole box side. For instance you could use these boxes:
59  *
60  *    rectangular       triclinic
61  *     2  0  0           2  0  0
62  *     0  2  0           0  2  0
63  *     0  0  6           2  2  6
64  *
65  * 2. You should check the energy conservation in a triclinic box.
66  *
67  * It might seem an overkill, but better safe than sorry.
68  * /Erik 001109
69  */
70
71 #include "gmxpre.h"
72
73 #include "pme.h"
74
75 #include "config.h"
76
77 #include <cassert>
78 #include <cmath>
79 #include <cstdio>
80 #include <cstdlib>
81 #include <cstring>
82
83 #include <algorithm>
84 #include <list>
85
86 #include "gromacs/domdec/domdec.h"
87 #include "gromacs/ewald/ewald_utils.h"
88 #include "gromacs/fft/parallel_3dfft.h"
89 #include "gromacs/fileio/pdbio.h"
90 #include "gromacs/gmxlib/network.h"
91 #include "gromacs/gmxlib/nrnb.h"
92 #include "gromacs/hardware/hw_info.h"
93 #include "gromacs/math/gmxcomplex.h"
94 #include "gromacs/math/invertmatrix.h"
95 #include "gromacs/math/units.h"
96 #include "gromacs/math/vec.h"
97 #include "gromacs/math/vectypes.h"
98 #include "gromacs/mdtypes/commrec.h"
99 #include "gromacs/mdtypes/forcerec.h"
100 #include "gromacs/mdtypes/inputrec.h"
101 #include "gromacs/mdtypes/md_enums.h"
102 #include "gromacs/mdtypes/simulation_workload.h"
103 #include "gromacs/pbcutil/pbc.h"
104 #include "gromacs/timing/cyclecounter.h"
105 #include "gromacs/timing/wallcycle.h"
106 #include "gromacs/timing/walltime_accounting.h"
107 #include "gromacs/topology/topology.h"
108 #include "gromacs/utility/basedefinitions.h"
109 #include "gromacs/utility/cstringutil.h"
110 #include "gromacs/utility/exceptions.h"
111 #include "gromacs/utility/fatalerror.h"
112 #include "gromacs/utility/gmxmpi.h"
113 #include "gromacs/utility/gmxomp.h"
114 #include "gromacs/utility/logger.h"
115 #include "gromacs/utility/real.h"
116 #include "gromacs/utility/smalloc.h"
117 #include "gromacs/utility/stringutil.h"
118 #include "gromacs/utility/unique_cptr.h"
119
120 #include "calculate_spline_moduli.h"
121 #include "pme_gather.h"
122 #include "pme_gpu_internal.h"
123 #include "pme_grid.h"
124 #include "pme_internal.h"
125 #include "pme_redistribute.h"
126 #include "pme_solve.h"
127 #include "pme_spline_work.h"
128 #include "pme_spread.h"
129
130 /*! \brief Help build a descriptive message in \c error if there are
131  * \c errorReasons why PME on GPU is not supported.
132  *
133  * \returns Whether the lack of errorReasons indicate there is support. */
134 static bool addMessageIfNotSupported(const std::list<std::string>& errorReasons, std::string* error)
135 {
136     bool isSupported = errorReasons.empty();
137     if (!isSupported && error)
138     {
139         std::string regressionTestMarker = "PME GPU does not support";
140         // this prefix is tested for in the regression tests script gmxtest.pl
141         *error = regressionTestMarker;
142         if (errorReasons.size() == 1)
143         {
144             *error += " " + errorReasons.back();
145         }
146         else
147         {
148             *error += ": " + gmx::joinStrings(errorReasons, "; ");
149         }
150         *error += ".";
151     }
152     return isSupported;
153 }
154
155 bool pme_gpu_supports_build(std::string* error)
156 {
157     std::list<std::string> errorReasons;
158     if (GMX_DOUBLE)
159     {
160         errorReasons.emplace_back("a double-precision build");
161     }
162     if (!GMX_GPU)
163     {
164         errorReasons.emplace_back("a non-GPU build");
165     }
166     if (GMX_GPU_SYCL)
167     {
168         errorReasons.emplace_back("SYCL build"); // SYCL-TODO
169     }
170     return addMessageIfNotSupported(errorReasons, error);
171 }
172
173 bool pme_gpu_supports_hardware(const gmx_hw_info_t gmx_unused& hwinfo, std::string* error)
174 {
175     std::list<std::string> errorReasons;
176
177     if (GMX_GPU_OPENCL)
178     {
179 #ifdef __APPLE__
180         errorReasons.emplace_back("Apple OS X operating system");
181 #endif
182     }
183     return addMessageIfNotSupported(errorReasons, error);
184 }
185
186 bool pme_gpu_supports_input(const t_inputrec& ir, std::string* error)
187 {
188     std::list<std::string> errorReasons;
189     if (!EEL_PME(ir.coulombtype))
190     {
191         errorReasons.emplace_back("systems that do not use PME for electrostatics");
192     }
193     if (ir.pme_order != 4)
194     {
195         errorReasons.emplace_back("interpolation orders other than 4");
196     }
197     if (EVDW_PME(ir.vdwtype))
198     {
199         errorReasons.emplace_back("Lennard-Jones PME");
200     }
201     if (!EI_DYNAMICS(ir.eI))
202     {
203         errorReasons.emplace_back(
204                 "Cannot compute PME interactions on a GPU, because PME GPU requires a dynamical "
205                 "integrator (md, sd, etc).");
206     }
207     return addMessageIfNotSupported(errorReasons, error);
208 }
209
210 /*! \brief \libinternal
211  * Finds out if PME with given inputs is possible to run on GPU.
212  * This function is an internal final check, validating the whole PME structure on creation,
213  * but it still duplicates the preliminary checks from the above (externally exposed) pme_gpu_supports_input() - just in case.
214  *
215  * \param[in]  pme          The PME structure.
216  * \param[out] error        The error message if the input is not supported on GPU.
217  * \returns                 True if this PME input is possible to run on GPU, false otherwise.
218  */
219 static bool pme_gpu_check_restrictions(const gmx_pme_t* pme, std::string* error)
220 {
221     std::list<std::string> errorReasons;
222     if (pme->nnodes != 1)
223     {
224         errorReasons.emplace_back("PME decomposition");
225     }
226     if (pme->pme_order != 4)
227     {
228         errorReasons.emplace_back("interpolation orders other than 4");
229     }
230     if (pme->doLJ)
231     {
232         errorReasons.emplace_back("Lennard-Jones PME");
233     }
234     if (GMX_DOUBLE)
235     {
236         errorReasons.emplace_back("double precision");
237     }
238     if (!GMX_GPU)
239     {
240         errorReasons.emplace_back("non-GPU build of GROMACS");
241     }
242     if (GMX_GPU_SYCL)
243     {
244         errorReasons.emplace_back("SYCL build of GROMACS"); // SYCL-TODO
245     }
246     return addMessageIfNotSupported(errorReasons, error);
247 }
248
249 PmeRunMode pme_run_mode(const gmx_pme_t* pme)
250 {
251     GMX_ASSERT(pme != nullptr, "Expecting valid PME data pointer");
252     return pme->runMode;
253 }
254
255 gmx::PinningPolicy pme_get_pinning_policy()
256 {
257     return gmx::PinningPolicy::PinnedIfSupported;
258 }
259
260 /*! \brief Number of bytes in a cache line.
261  *
262  * Must also be a multiple of the SIMD and SIMD4 register size, to
263  * preserve alignment.
264  */
265 const int gmxCacheLineSize = 64;
266
267 //! Set up coordinate communication
268 static void setup_coordinate_communication(PmeAtomComm* atc)
269 {
270     int nslab, n, i;
271     int fw, bw;
272
273     nslab = atc->nslab;
274
275     n = 0;
276     for (i = 1; i <= nslab / 2; i++)
277     {
278         fw = (atc->nodeid + i) % nslab;
279         bw = (atc->nodeid - i + nslab) % nslab;
280         if (n < nslab - 1)
281         {
282             atc->slabCommSetup[n].node_dest = fw;
283             atc->slabCommSetup[n].node_src  = bw;
284             n++;
285         }
286         if (n < nslab - 1)
287         {
288             atc->slabCommSetup[n].node_dest = bw;
289             atc->slabCommSetup[n].node_src  = fw;
290             n++;
291         }
292     }
293 }
294
295 /*! \brief Round \p n up to the next multiple of \p f */
296 static int mult_up(int n, int f)
297 {
298     return ((n + f - 1) / f) * f;
299 }
300
301 /*! \brief Return estimate of the load imbalance from the PME grid not being a good match for the number of PME ranks */
302 static double estimate_pme_load_imbalance(struct gmx_pme_t* pme)
303 {
304     int    nma, nmi;
305     double n1, n2, n3;
306
307     nma = pme->nnodes_major;
308     nmi = pme->nnodes_minor;
309
310     n1 = mult_up(pme->nkx, nma) * mult_up(pme->nky, nmi) * pme->nkz;
311     n2 = mult_up(pme->nkx, nma) * mult_up(pme->nkz, nmi) * pme->nky;
312     n3 = mult_up(pme->nky, nma) * mult_up(pme->nkz, nmi) * pme->nkx;
313
314     /* pme_solve is roughly double the cost of an fft */
315
316     return (n1 + n2 + 3 * n3) / static_cast<double>(6 * pme->nkx * pme->nky * pme->nkz);
317 }
318
319 #ifndef DOXYGEN
320
321 PmeAtomComm::PmeAtomComm(MPI_Comm   PmeMpiCommunicator,
322                          const int  numThreads,
323                          const int  pmeOrder,
324                          const int  dimIndex,
325                          const bool doSpread) :
326     dimind(dimIndex),
327     bSpread(doSpread),
328     pme_order(pmeOrder),
329     nthread(numThreads),
330     spline(nthread)
331 {
332     if (PmeMpiCommunicator != MPI_COMM_NULL)
333     {
334         mpi_comm = PmeMpiCommunicator;
335 #    if GMX_MPI
336         MPI_Comm_size(mpi_comm, &nslab);
337         MPI_Comm_rank(mpi_comm, &nodeid);
338 #    endif
339     }
340     if (debug)
341     {
342         fprintf(debug, "For PME atom communication in dimind %d: nslab %d rank %d\n", dimind, nslab, nodeid);
343     }
344
345     if (nslab > 1)
346     {
347         slabCommSetup.resize(nslab);
348         setup_coordinate_communication(this);
349
350         count_thread.resize(nthread);
351         for (auto& countThread : count_thread)
352         {
353             countThread.resize(nslab);
354         }
355     }
356
357     if (nthread > 1)
358     {
359         threadMap.resize(nthread);
360
361 #    pragma omp parallel for num_threads(nthread) schedule(static)
362         for (int thread = 0; thread < nthread; thread++)
363         {
364             try
365             {
366                 /* Allocate buffer with padding to avoid cache polution */
367                 threadMap[thread].nBuffer.resize(nthread + 2 * gmxCacheLineSize);
368                 threadMap[thread].n = threadMap[thread].nBuffer.data() + gmxCacheLineSize;
369             }
370             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
371         }
372     }
373 }
374
375 #endif // !DOXYGEN
376
377 /*! \brief Initialize data structure for communication */
378 static void init_overlap_comm(pme_overlap_t* ol, int norder, MPI_Comm comm, int nnodes, int nodeid, int ndata, int commplainsize)
379 {
380     gmx_bool bCont;
381
382     ol->mpi_comm = comm;
383     ol->nnodes   = nnodes;
384     ol->nodeid   = nodeid;
385
386     /* Linear translation of the PME grid won't affect reciprocal space
387      * calculations, so to optimize we only interpolate "upwards",
388      * which also means we only have to consider overlap in one direction.
389      * I.e., particles on this node might also be spread to grid indices
390      * that belong to higher nodes (modulo nnodes)
391      */
392
393     ol->s2g0.resize(ol->nnodes + 1);
394     ol->s2g1.resize(ol->nnodes);
395     if (debug)
396     {
397         fprintf(debug, "PME slab boundaries:");
398     }
399     for (int i = 0; i < nnodes; i++)
400     {
401         /* s2g0 the local interpolation grid start.
402          * s2g1 the local interpolation grid end.
403          * Since in calc_pidx we divide particles, and not grid lines,
404          * spatially uniform along dimension x or y, we need to round
405          * s2g0 down and s2g1 up.
406          */
407         ol->s2g0[i] = (i * ndata + 0) / nnodes;
408         ol->s2g1[i] = ((i + 1) * ndata + nnodes - 1) / nnodes + norder - 1;
409
410         if (debug)
411         {
412             fprintf(debug, "  %3d %3d", ol->s2g0[i], ol->s2g1[i]);
413         }
414     }
415     ol->s2g0[nnodes] = ndata;
416     if (debug)
417     {
418         fprintf(debug, "\n");
419     }
420
421     /* Determine with how many nodes we need to communicate the grid overlap */
422     int testRankCount = 0;
423     do
424     {
425         testRankCount++;
426         bCont = FALSE;
427         for (int i = 0; i < nnodes; i++)
428         {
429             if ((i + testRankCount < nnodes && ol->s2g1[i] > ol->s2g0[i + testRankCount])
430                 || (i + testRankCount >= nnodes && ol->s2g1[i] > ol->s2g0[i + testRankCount - nnodes] + ndata))
431             {
432                 bCont = TRUE;
433             }
434         }
435     } while (bCont && testRankCount < nnodes);
436
437     ol->comm_data.resize(testRankCount - 1);
438     ol->send_size = 0;
439
440     for (size_t b = 0; b < ol->comm_data.size(); b++)
441     {
442         pme_grid_comm_t* pgc = &ol->comm_data[b];
443
444         /* Send */
445         pgc->send_id  = (ol->nodeid + (b + 1)) % ol->nnodes;
446         int fft_start = ol->s2g0[pgc->send_id];
447         int fft_end   = ol->s2g0[pgc->send_id + 1];
448         if (pgc->send_id < nodeid)
449         {
450             fft_start += ndata;
451             fft_end += ndata;
452         }
453         int send_index1  = ol->s2g1[nodeid];
454         send_index1      = std::min(send_index1, fft_end);
455         pgc->send_index0 = fft_start;
456         pgc->send_nindex = std::max(0, send_index1 - pgc->send_index0);
457         ol->send_size += pgc->send_nindex;
458
459         /* We always start receiving to the first index of our slab */
460         pgc->recv_id    = (ol->nodeid - (b + 1) + ol->nnodes) % ol->nnodes;
461         fft_start       = ol->s2g0[ol->nodeid];
462         fft_end         = ol->s2g0[ol->nodeid + 1];
463         int recv_index1 = ol->s2g1[pgc->recv_id];
464         if (pgc->recv_id > nodeid)
465         {
466             recv_index1 -= ndata;
467         }
468         recv_index1      = std::min(recv_index1, fft_end);
469         pgc->recv_index0 = fft_start;
470         pgc->recv_nindex = std::max(0, recv_index1 - pgc->recv_index0);
471     }
472
473 #if GMX_MPI
474     /* Communicate the buffer sizes to receive */
475     MPI_Status stat;
476     for (size_t b = 0; b < ol->comm_data.size(); b++)
477     {
478         MPI_Sendrecv(&ol->send_size,
479                      1,
480                      MPI_INT,
481                      ol->comm_data[b].send_id,
482                      b,
483                      &ol->comm_data[b].recv_size,
484                      1,
485                      MPI_INT,
486                      ol->comm_data[b].recv_id,
487                      b,
488                      ol->mpi_comm,
489                      &stat);
490     }
491 #endif
492
493     /* For non-divisible grid we need pme_order iso pme_order-1 */
494     ol->sendbuf.resize(norder * commplainsize);
495     ol->recvbuf.resize(norder * commplainsize);
496 }
497
498 int minimalPmeGridSize(int pmeOrder)
499 {
500     /* The actual grid size limitations are:
501      *   serial:        >= pme_order
502      *   DD, no OpenMP: >= 2*(pme_order - 1)
503      *   DD, OpenMP:    >= pme_order + 1
504      * But we use the maximum for simplicity since in practice there is not
505      * much performance difference between pme_order and 2*(pme_order -1).
506      */
507     int minimalSize = 2 * (pmeOrder - 1);
508
509     GMX_RELEASE_ASSERT(pmeOrder >= 3, "pmeOrder has to be >= 3");
510     GMX_RELEASE_ASSERT(minimalSize >= pmeOrder + 1, "The grid size should be >= pmeOrder + 1");
511
512     return minimalSize;
513 }
514
515 bool gmx_pme_check_restrictions(int pme_order, int nkx, int nky, int nkz, int numPmeDomainsAlongX, bool useThreads, bool errorsAreFatal)
516 {
517     if (pme_order > PME_ORDER_MAX)
518     {
519         if (!errorsAreFatal)
520         {
521             return false;
522         }
523
524         std::string message = gmx::formatString(
525                 "pme_order (%d) is larger than the maximum allowed value (%d). Modify and "
526                 "recompile the code if you really need such a high order.",
527                 pme_order,
528                 PME_ORDER_MAX);
529         GMX_THROW(gmx::InconsistentInputError(message));
530     }
531
532     const int minGridSize = minimalPmeGridSize(pme_order);
533     if (nkx < minGridSize || nky < minGridSize || nkz < minGridSize)
534     {
535         if (!errorsAreFatal)
536         {
537             return false;
538         }
539         std::string message =
540                 gmx::formatString("The PME grid sizes need to be >= 2*(pme_order-1) (%d)", minGridSize);
541         GMX_THROW(gmx::InconsistentInputError(message));
542     }
543
544     /* Check for a limitation of the (current) sum_fftgrid_dd code.
545      * We only allow multiple communication pulses in dim 1, not in dim 0.
546      */
547     if (useThreads
548         && (nkx < numPmeDomainsAlongX * pme_order && nkx != numPmeDomainsAlongX * (pme_order - 1)))
549     {
550         if (!errorsAreFatal)
551         {
552             return false;
553         }
554         gmx_fatal(FARGS,
555                   "The number of PME grid lines per rank along x is %g. But when using OpenMP "
556                   "threads, the number of grid lines per rank along x should be >= pme_order (%d) "
557                   "or = pmeorder-1. To resolve this issue, use fewer ranks along x (and possibly "
558                   "more along y and/or z) by specifying -dd manually.",
559                   nkx / static_cast<double>(numPmeDomainsAlongX),
560                   pme_order);
561     }
562
563     return true;
564 }
565
566 /*! \brief Round \p enumerator */
567 static int div_round_up(int enumerator, int denominator)
568 {
569     return (enumerator + denominator - 1) / denominator;
570 }
571
572 gmx_pme_t* gmx_pme_init(const t_commrec*     cr,
573                         const NumPmeDomains& numPmeDomains,
574                         const t_inputrec*    ir,
575                         gmx_bool             bFreeEnergy_q,
576                         gmx_bool             bFreeEnergy_lj,
577                         gmx_bool             bReproducible,
578                         real                 ewaldcoeff_q,
579                         real                 ewaldcoeff_lj,
580                         int                  nthread,
581                         PmeRunMode           runMode,
582                         PmeGpu*              pmeGpu,
583                         const DeviceContext* deviceContext,
584                         const DeviceStream*  deviceStream,
585                         const PmeGpuProgram* pmeGpuProgram,
586                         const gmx::MDLogger& mdlog)
587 {
588     int  use_threads, sum_use_threads, i;
589     ivec ndata;
590
591     if (debug)
592     {
593         fprintf(debug, "Creating PME data structures.\n");
594     }
595
596     gmx::unique_cptr<gmx_pme_t, gmx_pme_destroy> pme(new gmx_pme_t());
597
598     pme->sum_qgrid_tmp    = nullptr;
599     pme->sum_qgrid_dd_tmp = nullptr;
600
601     pme->buf_nalloc = 0;
602
603     pme->nnodes  = 1;
604     pme->bPPnode = TRUE;
605
606     pme->nnodes_major = numPmeDomains.x;
607     pme->nnodes_minor = numPmeDomains.y;
608
609     if (numPmeDomains.x * numPmeDomains.y > 1)
610     {
611         pme->mpi_comm = cr->mpi_comm_mygroup;
612
613 #if GMX_MPI
614         MPI_Comm_rank(pme->mpi_comm, &pme->nodeid);
615         MPI_Comm_size(pme->mpi_comm, &pme->nnodes);
616 #endif
617         if (pme->nnodes != numPmeDomains.x * numPmeDomains.y)
618         {
619             gmx_incons("PME rank count mismatch");
620         }
621     }
622     else
623     {
624         pme->mpi_comm = MPI_COMM_NULL;
625     }
626
627     if (pme->nnodes == 1)
628     {
629         pme->mpi_comm_d[0] = MPI_COMM_NULL;
630         pme->mpi_comm_d[1] = MPI_COMM_NULL;
631         pme->ndecompdim    = 0;
632         pme->nodeid_major  = 0;
633         pme->nodeid_minor  = 0;
634     }
635     else
636     {
637         if (numPmeDomains.y == 1)
638         {
639             pme->mpi_comm_d[0] = pme->mpi_comm;
640             pme->mpi_comm_d[1] = MPI_COMM_NULL;
641             pme->ndecompdim    = 1;
642             pme->nodeid_major  = pme->nodeid;
643             pme->nodeid_minor  = 0;
644         }
645         else if (numPmeDomains.x == 1)
646         {
647             pme->mpi_comm_d[0] = MPI_COMM_NULL;
648             pme->mpi_comm_d[1] = pme->mpi_comm;
649             pme->ndecompdim    = 1;
650             pme->nodeid_major  = 0;
651             pme->nodeid_minor  = pme->nodeid;
652         }
653         else
654         {
655             if (pme->nnodes % numPmeDomains.x != 0)
656             {
657                 gmx_incons(
658                         "For 2D PME decomposition, #PME ranks must be divisible by the number of "
659                         "domains along x");
660             }
661             pme->ndecompdim = 2;
662
663 #if GMX_MPI
664             MPI_Comm_split(pme->mpi_comm,
665                            pme->nodeid % numPmeDomains.y,
666                            pme->nodeid,
667                            &pme->mpi_comm_d[0]); /* My communicator along major dimension */
668             MPI_Comm_split(pme->mpi_comm,
669                            pme->nodeid / numPmeDomains.y,
670                            pme->nodeid,
671                            &pme->mpi_comm_d[1]); /* My communicator along minor dimension */
672
673             MPI_Comm_rank(pme->mpi_comm_d[0], &pme->nodeid_major);
674             MPI_Comm_size(pme->mpi_comm_d[0], &pme->nnodes_major);
675             MPI_Comm_rank(pme->mpi_comm_d[1], &pme->nodeid_minor);
676             MPI_Comm_size(pme->mpi_comm_d[1], &pme->nnodes_minor);
677 #endif
678         }
679     }
680     // cr is always initialized if there is a a PP rank, so we can safely assume
681     // that when it is not, like in ewald tests, we not on a PP rank.
682     pme->bPPnode = ((cr != nullptr && cr->duty != 0) && thisRankHasDuty(cr, DUTY_PP));
683
684     pme->nthread = nthread;
685
686     /* Check if any of the PME MPI ranks uses threads */
687     use_threads = (pme->nthread > 1 ? 1 : 0);
688 #if GMX_MPI
689     if (pme->nnodes > 1)
690     {
691         MPI_Allreduce(&use_threads, &sum_use_threads, 1, MPI_INT, MPI_SUM, pme->mpi_comm);
692     }
693     else
694 #endif
695     {
696         sum_use_threads = use_threads;
697     }
698     pme->bUseThreads = (sum_use_threads > 0);
699
700     if (ir->pbcType == PbcType::Screw)
701     {
702         gmx_fatal(FARGS, "pme does not (yet) work with pbc = screw");
703     }
704
705     /* NOTE:
706      * It is likely that the current gmx_pme_do() routine supports calculating
707      * only Coulomb or LJ while gmx_pme_init() configures for both,
708      * but that has never been tested.
709      * It is likely that the current gmx_pme_do() routine supports calculating,
710      * not calculating free-energy for Coulomb and/or LJ while gmx_pme_init()
711      * configures with free-energy, but that has never been tested.
712      */
713     pme->doCoulomb = EEL_PME(ir->coulombtype);
714     pme->doLJ      = EVDW_PME(ir->vdwtype);
715     pme->bFEP_q    = ((ir->efep != FreeEnergyPerturbationType::No) && bFreeEnergy_q);
716     pme->bFEP_lj   = ((ir->efep != FreeEnergyPerturbationType::No) && bFreeEnergy_lj);
717     pme->bFEP      = (pme->bFEP_q || pme->bFEP_lj);
718     pme->nkx       = ir->nkx;
719     pme->nky       = ir->nky;
720     pme->nkz       = ir->nkz;
721     pme->bP3M = (ir->coulombtype == CoulombInteractionType::P3mAD || getenv("GMX_PME_P3M") != nullptr);
722     pme->pme_order     = ir->pme_order;
723     pme->ewaldcoeff_q  = ewaldcoeff_q;
724     pme->ewaldcoeff_lj = ewaldcoeff_lj;
725
726     /* Always constant electrostatics coefficients */
727     pme->epsilon_r = ir->epsilon_r;
728
729     /* Always constant LJ coefficients */
730     pme->ljpme_combination_rule = ir->ljpme_combination_rule;
731
732     // The box requires scaling with nwalls = 2, we store that condition as well
733     // as the scaling factor
734     delete pme->boxScaler;
735     pme->boxScaler = new EwaldBoxZScaler(*ir);
736
737     /* If we violate restrictions, generate a fatal error here */
738     gmx_pme_check_restrictions(
739             pme->pme_order, pme->nkx, pme->nky, pme->nkz, pme->nnodes_major, pme->bUseThreads, true);
740
741     if (pme->nnodes > 1)
742     {
743         double imbal;
744
745 #if GMX_MPI
746         MPI_Type_contiguous(DIM, GMX_MPI_REAL, &(pme->rvec_mpi));
747         MPI_Type_commit(&(pme->rvec_mpi));
748 #endif
749
750         /* Note that the coefficient spreading and force gathering, which usually
751          * takes about the same amount of time as FFT+solve_pme,
752          * is always fully load balanced
753          * (unless the coefficient distribution is inhomogeneous).
754          */
755
756         imbal = estimate_pme_load_imbalance(pme.get());
757         if (imbal >= 1.2 && pme->nodeid_major == 0 && pme->nodeid_minor == 0)
758         {
759             GMX_LOG(mdlog.warning)
760                     .asParagraph()
761                     .appendTextFormatted(
762                             "NOTE: The load imbalance in PME FFT and solve is %d%%.\n"
763                             "      For optimal PME load balancing\n"
764                             "      PME grid_x (%d) and grid_y (%d) should be divisible by "
765                             "#PME_ranks_x "
766                             "(%d)\n"
767                             "      and PME grid_y (%d) and grid_z (%d) should be divisible by "
768                             "#PME_ranks_y "
769                             "(%d)",
770                             gmx::roundToInt((imbal - 1) * 100),
771                             pme->nkx,
772                             pme->nky,
773                             pme->nnodes_major,
774                             pme->nky,
775                             pme->nkz,
776                             pme->nnodes_minor);
777         }
778     }
779
780     /* For non-divisible grid we need pme_order iso pme_order-1 */
781     /* In sum_qgrid_dd x overlap is copied in place: take padding into account.
782      * y is always copied through a buffer: we don't need padding in z,
783      * but we do need the overlap in x because of the communication order.
784      */
785     init_overlap_comm(&pme->overlap[0],
786                       pme->pme_order,
787                       pme->mpi_comm_d[0],
788                       pme->nnodes_major,
789                       pme->nodeid_major,
790                       pme->nkx,
791                       (div_round_up(pme->nky, pme->nnodes_minor) + pme->pme_order)
792                               * (pme->nkz + pme->pme_order - 1));
793
794     /* Along overlap dim 1 we can send in multiple pulses in sum_fftgrid_dd.
795      * We do this with an offset buffer of equal size, so we need to allocate
796      * extra for the offset. That's what the (+1)*pme->nkz is for.
797      */
798     init_overlap_comm(&pme->overlap[1],
799                       pme->pme_order,
800                       pme->mpi_comm_d[1],
801                       pme->nnodes_minor,
802                       pme->nodeid_minor,
803                       pme->nky,
804                       (div_round_up(pme->nkx, pme->nnodes_major) + pme->pme_order + 1) * pme->nkz);
805
806     /* Double-check for a limitation of the (current) sum_fftgrid_dd code.
807      * Note that gmx_pme_check_restrictions checked for this already.
808      */
809     if (pme->bUseThreads && (pme->overlap[0].comm_data.size() > 1))
810     {
811         gmx_incons(
812                 "More than one communication pulse required for grid overlap communication along "
813                 "the major dimension while using threads");
814     }
815
816     snew(pme->bsp_mod[XX], pme->nkx);
817     snew(pme->bsp_mod[YY], pme->nky);
818     snew(pme->bsp_mod[ZZ], pme->nkz);
819
820     pme->gpu     = pmeGpu; /* Carrying over the single GPU structure */
821     pme->runMode = runMode;
822
823     /* The required size of the interpolation grid, including overlap.
824      * The allocated size (pmegrid_n?) might be slightly larger.
825      */
826     pme->pmegrid_nx = pme->overlap[0].s2g1[pme->nodeid_major] - pme->overlap[0].s2g0[pme->nodeid_major];
827     pme->pmegrid_ny = pme->overlap[1].s2g1[pme->nodeid_minor] - pme->overlap[1].s2g0[pme->nodeid_minor];
828     pme->pmegrid_nz_base = pme->nkz;
829     pme->pmegrid_nz      = pme->pmegrid_nz_base + pme->pme_order - 1;
830     set_grid_alignment(&pme->pmegrid_nz, pme->pme_order);
831     pme->pmegrid_start_ix = pme->overlap[0].s2g0[pme->nodeid_major];
832     pme->pmegrid_start_iy = pme->overlap[1].s2g0[pme->nodeid_minor];
833     pme->pmegrid_start_iz = 0;
834
835     make_gridindex_to_localindex(
836             pme->nkx, pme->pmegrid_start_ix, pme->pmegrid_nx - (pme->pme_order - 1), &pme->nnx, &pme->fshx);
837     make_gridindex_to_localindex(
838             pme->nky, pme->pmegrid_start_iy, pme->pmegrid_ny - (pme->pme_order - 1), &pme->nny, &pme->fshy);
839     make_gridindex_to_localindex(
840             pme->nkz, pme->pmegrid_start_iz, pme->pmegrid_nz_base, &pme->nnz, &pme->fshz);
841
842     pme->spline_work = make_pme_spline_work(pme->pme_order);
843
844     ndata[0] = pme->nkx;
845     ndata[1] = pme->nky;
846     ndata[2] = pme->nkz;
847     /* It doesn't matter if we allocate too many grids here,
848      * we only allocate and use the ones we need.
849      */
850     if (pme->doLJ)
851     {
852         pme->ngrids = ((ir->ljpme_combination_rule == LongRangeVdW::LB) ? DO_Q_AND_LJ_LB : DO_Q_AND_LJ);
853     }
854     else
855     {
856         pme->ngrids = DO_Q;
857     }
858     snew(pme->fftgrid, pme->ngrids);
859     snew(pme->cfftgrid, pme->ngrids);
860     snew(pme->pfft_setup, pme->ngrids);
861
862     for (i = 0; i < pme->ngrids; ++i)
863     {
864         if ((i < DO_Q && pme->doCoulomb && (i == 0 || bFreeEnergy_q))
865             || (i >= DO_Q && pme->doLJ
866                 && (i == 2 || bFreeEnergy_lj || ir->ljpme_combination_rule == LongRangeVdW::LB)))
867         {
868             pmegrids_init(&pme->pmegrid[i],
869                           pme->pmegrid_nx,
870                           pme->pmegrid_ny,
871                           pme->pmegrid_nz,
872                           pme->pmegrid_nz_base,
873                           pme->pme_order,
874                           pme->bUseThreads,
875                           pme->nthread,
876                           pme->overlap[0].s2g1[pme->nodeid_major]
877                                   - pme->overlap[0].s2g0[pme->nodeid_major + 1],
878                           pme->overlap[1].s2g1[pme->nodeid_minor]
879                                   - pme->overlap[1].s2g0[pme->nodeid_minor + 1]);
880             /* This routine will allocate the grid data to fit the FFTs */
881             const auto allocateRealGridForGpu = (pme->runMode == PmeRunMode::Mixed)
882                                                         ? gmx::PinningPolicy::PinnedIfSupported
883                                                         : gmx::PinningPolicy::CannotBePinned;
884             gmx_parallel_3dfft_init(&pme->pfft_setup[i],
885                                     ndata,
886                                     &pme->fftgrid[i],
887                                     &pme->cfftgrid[i],
888                                     pme->mpi_comm_d,
889                                     bReproducible,
890                                     pme->nthread,
891                                     allocateRealGridForGpu);
892         }
893     }
894
895     if (!pme->bP3M)
896     {
897         /* Use plain SPME B-spline interpolation */
898         make_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
899     }
900     else
901     {
902         /* Use the P3M grid-optimized influence function */
903         make_p3m_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
904     }
905
906     /* Use atc[0] for spreading */
907     const int firstDimIndex   = (numPmeDomains.x > 1 ? 0 : 1);
908     MPI_Comm  mpiCommFirstDim = (pme->nnodes > 1 ? pme->mpi_comm_d[firstDimIndex] : MPI_COMM_NULL);
909     bool      doSpread        = true;
910     pme->atc.emplace_back(mpiCommFirstDim, pme->nthread, pme->pme_order, firstDimIndex, doSpread);
911     if (pme->ndecompdim >= 2)
912     {
913         const int secondDimIndex = 1;
914         doSpread                 = false;
915         pme->atc.emplace_back(pme->mpi_comm_d[1], pme->nthread, pme->pme_order, secondDimIndex, doSpread);
916     }
917
918     // Initial check of validity of the input for running on the GPU
919     if (pme->runMode != PmeRunMode::CPU)
920     {
921         std::string errorString;
922         bool        canRunOnGpu = pme_gpu_check_restrictions(pme.get(), &errorString);
923         if (!canRunOnGpu)
924         {
925             GMX_THROW(gmx::NotImplementedError(errorString));
926         }
927         pme_gpu_reinit(pme.get(), deviceContext, deviceStream, pmeGpuProgram);
928     }
929     else
930     {
931         GMX_ASSERT(pme->gpu == nullptr, "Should not have PME GPU object when PME is on a CPU.");
932     }
933
934
935     pme_init_all_work(&pme->solve_work, pme->nthread, pme->nkx);
936
937     // no exception was thrown during the init, so we hand over the PME structure handle
938     return pme.release();
939 }
940
941 void gmx_pme_reinit(struct gmx_pme_t** pmedata,
942                     const t_commrec*   cr,
943                     struct gmx_pme_t*  pme_src,
944                     const t_inputrec*  ir,
945                     const ivec         grid_size,
946                     real               ewaldcoeff_q,
947                     real               ewaldcoeff_lj)
948 {
949     // Create a copy of t_inputrec fields that are used in gmx_pme_init().
950     // TODO: This would be better as just copying a sub-structure that contains
951     // all the PME parameters and nothing else.
952     t_inputrec irc;
953     irc.pbcType                = ir->pbcType;
954     irc.coulombtype            = ir->coulombtype;
955     irc.vdwtype                = ir->vdwtype;
956     irc.efep                   = ir->efep;
957     irc.pme_order              = ir->pme_order;
958     irc.epsilon_r              = ir->epsilon_r;
959     irc.ljpme_combination_rule = ir->ljpme_combination_rule;
960     irc.nkx                    = grid_size[XX];
961     irc.nky                    = grid_size[YY];
962     irc.nkz                    = grid_size[ZZ];
963
964     try
965     {
966         // This is reinit. Any logging should have been done at first init.
967         // Here we should avoid writing notes for settings the user did not
968         // set directly.
969         const gmx::MDLogger dummyLogger;
970         GMX_ASSERT(pmedata, "Invalid PME pointer");
971         NumPmeDomains numPmeDomains = { pme_src->nnodes_major, pme_src->nnodes_minor };
972         *pmedata                    = gmx_pme_init(cr,
973                                 numPmeDomains,
974                                 &irc,
975                                 pme_src->bFEP_q,
976                                 pme_src->bFEP_lj,
977                                 FALSE,
978                                 ewaldcoeff_q,
979                                 ewaldcoeff_lj,
980                                 pme_src->nthread,
981                                 pme_src->runMode,
982                                 pme_src->gpu,
983                                 nullptr,
984                                 nullptr,
985                                 nullptr,
986                                 dummyLogger);
987         /* When running PME on the CPU not using domain decomposition,
988          * the atom data is allocated once only in gmx_pme_(re)init().
989          */
990         if (!pme_src->gpu && pme_src->nnodes == 1)
991         {
992             gmx_pme_reinit_atoms(*pmedata, pme_src->atc[0].numAtoms(), nullptr, nullptr);
993         }
994         // TODO this is mostly passing around current values
995     }
996     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
997
998     /* We can easily reuse the allocated pme grids in pme_src */
999     reuse_pmegrids(&pme_src->pmegrid[PME_GRID_QA], &(*pmedata)->pmegrid[PME_GRID_QA]);
1000     /* We would like to reuse the fft grids, but that's harder */
1001 }
1002
1003 void gmx_pme_calc_energy(gmx_pme_t* pme, gmx::ArrayRef<const gmx::RVec> x, gmx::ArrayRef<const real> q, real* V)
1004 {
1005     pmegrids_t* grid;
1006
1007     if (pme->nnodes > 1)
1008     {
1009         gmx_incons("gmx_pme_calc_energy called in parallel");
1010     }
1011     if (pme->bFEP_q)
1012     {
1013         gmx_incons("gmx_pme_calc_energy with free energy");
1014     }
1015
1016     if (!pme->atc_energy)
1017     {
1018         pme->atc_energy = std::make_unique<PmeAtomComm>(MPI_COMM_NULL, 1, pme->pme_order, 0, true);
1019     }
1020     PmeAtomComm* atc = pme->atc_energy.get();
1021     atc->setNumAtoms(x.ssize());
1022     atc->x           = x;
1023     atc->coefficient = q;
1024
1025     /* We only use the A-charges grid */
1026     grid = &pme->pmegrid[PME_GRID_QA];
1027
1028     /* Only calculate the spline coefficients, don't actually spread */
1029     spread_on_grid(pme, atc, nullptr, TRUE, FALSE, pme->fftgrid[PME_GRID_QA], FALSE, PME_GRID_QA);
1030
1031     *V = gather_energy_bsplines(pme, grid->grid.grid, atc);
1032 }
1033
1034 /*! \brief Calculate initial Lorentz-Berthelot coefficients for LJ-PME */
1035 static void calc_initial_lb_coeffs(gmx::ArrayRef<real> coefficient, const real* local_c6, const real* local_sigma)
1036 {
1037     for (gmx::index i = 0; i < coefficient.ssize(); ++i)
1038     {
1039         real sigma4    = local_sigma[i];
1040         sigma4         = sigma4 * sigma4;
1041         sigma4         = sigma4 * sigma4;
1042         coefficient[i] = local_c6[i] / sigma4;
1043     }
1044 }
1045
1046 /*! \brief Calculate next Lorentz-Berthelot coefficients for LJ-PME */
1047 static void calc_next_lb_coeffs(gmx::ArrayRef<real> coefficient, const real* local_sigma)
1048 {
1049     for (gmx::index i = 0; i < coefficient.ssize(); ++i)
1050     {
1051         coefficient[i] *= local_sigma[i];
1052     }
1053 }
1054
1055 int gmx_pme_do(struct gmx_pme_t*              pme,
1056                gmx::ArrayRef<const gmx::RVec> coordinates,
1057                gmx::ArrayRef<gmx::RVec>       forces,
1058                real                           chargeA[],
1059                real                           chargeB[],
1060                real                           c6A[],
1061                real                           c6B[],
1062                real                           sigmaA[],
1063                real                           sigmaB[],
1064                const matrix                   box,
1065                const t_commrec*               cr,
1066                int                            maxshift_x,
1067                int                            maxshift_y,
1068                t_nrnb*                        nrnb,
1069                gmx_wallcycle*                 wcycle,
1070                matrix                         vir_q,
1071                matrix                         vir_lj,
1072                real*                          energy_q,
1073                real*                          energy_lj,
1074                real                           lambda_q,
1075                real                           lambda_lj,
1076                real*                          dvdlambda_q,
1077                real*                          dvdlambda_lj,
1078                const gmx::StepWorkload&       stepWork)
1079 {
1080     GMX_ASSERT(pme->runMode == PmeRunMode::CPU,
1081                "gmx_pme_do should not be called on the GPU PME run.");
1082
1083     int                  d, npme, grid_index, max_grid_index;
1084     PmeAtomComm&         atc         = pme->atc[0];
1085     pmegrids_t*          pmegrid     = nullptr;
1086     real*                grid        = nullptr;
1087     real*                coefficient = nullptr;
1088     PmeOutput            output[2]; // The second is used for the B state with FEP
1089     real                 scale, lambda;
1090     gmx_bool             bClearF;
1091     gmx_parallel_3dfft_t pfft_setup;
1092     real*                fftgrid;
1093     t_complex*           cfftgrid;
1094     int                  thread;
1095     gmx_bool             bFirst, bDoSplines;
1096     int                  fep_state;
1097     int                  fep_states_lj = pme->bFEP_lj ? 2 : 1;
1098     // There's no support for computing energy without virial, or vice versa
1099     const bool computeEnergyAndVirial = (stepWork.computeEnergy || stepWork.computeVirial);
1100
1101     /* We could be passing lambda!=0 while no q or LJ is actually perturbed */
1102     if (!pme->bFEP_q)
1103     {
1104         lambda_q = 0;
1105     }
1106     if (!pme->bFEP_lj)
1107     {
1108         lambda_lj = 0;
1109     }
1110
1111     assert(pme->nnodes > 0);
1112     assert(pme->nnodes == 1 || pme->ndecompdim > 0);
1113
1114     if (pme->nnodes > 1)
1115     {
1116         atc.pd.resize(coordinates.ssize());
1117         for (int d = pme->ndecompdim - 1; d >= 0; d--)
1118         {
1119             PmeAtomComm& atc = pme->atc[d];
1120             atc.maxshift     = (atc.dimind == 0 ? maxshift_x : maxshift_y);
1121         }
1122     }
1123     else
1124     {
1125         GMX_ASSERT(coordinates.ssize() == atc.numAtoms(), "We expect atc.numAtoms() coordinates");
1126         GMX_ASSERT(forces.ssize() >= atc.numAtoms(),
1127                    "We need a force buffer with at least atc.numAtoms() elements");
1128
1129         atc.x = coordinates;
1130         atc.f = forces;
1131     }
1132
1133     matrix scaledBox;
1134     pme->boxScaler->scaleBox(box, scaledBox);
1135
1136     gmx::invertBoxMatrix(scaledBox, pme->recipbox);
1137     bFirst = TRUE;
1138
1139     /* For simplicity, we construct the splines for all particles if
1140      * more than one PME calculations is needed. Some optimization
1141      * could be done by keeping track of which atoms have splines
1142      * constructed, and construct new splines on each pass for atoms
1143      * that don't yet have them.
1144      */
1145
1146     bDoSplines = pme->bFEP || (pme->doCoulomb && pme->doLJ);
1147
1148     /* We need a maximum of four separate PME calculations:
1149      * grid_index=0: Coulomb PME with charges from state A
1150      * grid_index=1: Coulomb PME with charges from state B
1151      * grid_index=2: LJ PME with C6 from state A
1152      * grid_index=3: LJ PME with C6 from state B
1153      * For Lorentz-Berthelot combination rules, a separate loop is used to
1154      * calculate all the terms
1155      */
1156
1157     /* If we are doing LJ-PME with LB, we only do Q here */
1158     max_grid_index = (pme->ljpme_combination_rule == LongRangeVdW::LB) ? DO_Q : DO_Q_AND_LJ;
1159
1160     for (grid_index = 0; grid_index < max_grid_index; ++grid_index)
1161     {
1162         /* Check if we should do calculations at this grid_index
1163          * If grid_index is odd we should be doing FEP
1164          * If grid_index < 2 we should be doing electrostatic PME
1165          * If grid_index >= 2 we should be doing LJ-PME
1166          */
1167         if ((grid_index < DO_Q && (!pme->doCoulomb || (grid_index == 1 && !pme->bFEP_q)))
1168             || (grid_index >= DO_Q && (!pme->doLJ || (grid_index == 3 && !pme->bFEP_lj))))
1169         {
1170             continue;
1171         }
1172         /* Unpack structure */
1173         pmegrid    = &pme->pmegrid[grid_index];
1174         fftgrid    = pme->fftgrid[grid_index];
1175         cfftgrid   = pme->cfftgrid[grid_index];
1176         pfft_setup = pme->pfft_setup[grid_index];
1177         switch (grid_index)
1178         {
1179             case 0: coefficient = chargeA; break;
1180             case 1: coefficient = chargeB; break;
1181             case 2: coefficient = c6A; break;
1182             case 3: coefficient = c6B; break;
1183         }
1184
1185         grid = pmegrid->grid.grid;
1186
1187         if (debug)
1188         {
1189             fprintf(debug, "PME: number of ranks = %d, rank = %d\n", cr->nnodes, cr->nodeid);
1190             fprintf(debug, "Grid = %p\n", static_cast<void*>(grid));
1191             if (grid == nullptr)
1192             {
1193                 gmx_fatal(FARGS, "No grid!");
1194             }
1195         }
1196
1197         if (pme->nnodes == 1)
1198         {
1199             atc.coefficient = gmx::arrayRefFromArray(coefficient, coordinates.size());
1200         }
1201         else
1202         {
1203             wallcycle_start(wcycle, ewcPME_REDISTXF);
1204             do_redist_pos_coeffs(pme, cr, bFirst, coordinates, coefficient);
1205
1206             wallcycle_stop(wcycle, ewcPME_REDISTXF);
1207         }
1208
1209         if (debug)
1210         {
1211             fprintf(debug, "Rank= %6d, pme local particles=%6d\n", cr->nodeid, atc.numAtoms());
1212         }
1213
1214         wallcycle_start(wcycle, ewcPME_SPREAD);
1215
1216         /* Spread the coefficients on a grid */
1217         spread_on_grid(pme, &atc, pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1218
1219         if (bFirst)
1220         {
1221             inc_nrnb(nrnb, eNR_WEIGHTS, DIM * atc.numAtoms());
1222         }
1223         inc_nrnb(nrnb, eNR_SPREADBSP, pme->pme_order * pme->pme_order * pme->pme_order * atc.numAtoms());
1224
1225         if (!pme->bUseThreads)
1226         {
1227             wrap_periodic_pmegrid(pme, grid);
1228
1229             /* sum contributions to local grid from other nodes */
1230             if (pme->nnodes > 1)
1231             {
1232                 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1233             }
1234
1235             copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1236         }
1237
1238         wallcycle_stop(wcycle, ewcPME_SPREAD);
1239
1240         /* TODO If the OpenMP and single-threaded implementations
1241            converge, then spread_on_grid() and
1242            copy_pmegrid_to_fftgrid() will perhaps live in the same
1243            source file.
1244         */
1245
1246         /* Here we start a large thread parallel region */
1247 #pragma omp parallel num_threads(pme->nthread) private(thread)
1248         {
1249             try
1250             {
1251                 thread = gmx_omp_get_thread_num();
1252                 int loop_count;
1253
1254                 /* do 3d-fft */
1255                 if (thread == 0)
1256                 {
1257                     wallcycle_start(wcycle, ewcPME_FFT);
1258                 }
1259                 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX, thread, wcycle);
1260                 if (thread == 0)
1261                 {
1262                     wallcycle_stop(wcycle, ewcPME_FFT);
1263                 }
1264
1265                 /* solve in k-space for our local cells */
1266                 if (thread == 0)
1267                 {
1268                     wallcycle_start(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
1269                 }
1270                 if (grid_index < DO_Q)
1271                 {
1272                     loop_count = solve_pme_yzx(pme,
1273                                                cfftgrid,
1274                                                scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ],
1275                                                computeEnergyAndVirial,
1276                                                pme->nthread,
1277                                                thread);
1278                 }
1279                 else
1280                 {
1281                     loop_count =
1282                             solve_pme_lj_yzx(pme,
1283                                              &cfftgrid,
1284                                              FALSE,
1285                                              scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ],
1286                                              computeEnergyAndVirial,
1287                                              pme->nthread,
1288                                              thread);
1289                 }
1290
1291                 if (thread == 0)
1292                 {
1293                     wallcycle_stop(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
1294                     inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1295                 }
1296
1297                 /* do 3d-invfft */
1298                 if (thread == 0)
1299                 {
1300                     wallcycle_start(wcycle, ewcPME_FFT);
1301                 }
1302                 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL, thread, wcycle);
1303                 if (thread == 0)
1304                 {
1305                     wallcycle_stop(wcycle, ewcPME_FFT);
1306
1307
1308                     if (pme->nodeid == 0)
1309                     {
1310                         real ntot = pme->nkx * pme->nky * pme->nkz;
1311                         npme      = static_cast<int>(ntot * std::log(ntot) / std::log(2.0));
1312                         inc_nrnb(nrnb, eNR_FFT, 2 * npme);
1313                     }
1314
1315                     /* Note: this wallcycle region is closed below
1316                        outside an OpenMP region, so take care if
1317                        refactoring code here. */
1318                     wallcycle_start(wcycle, ewcPME_GATHER);
1319                 }
1320
1321                 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1322             }
1323             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1324         }
1325         /* End of thread parallel section.
1326          * With MPI we have to synchronize here before gmx_sum_qgrid_dd.
1327          */
1328
1329         /* distribute local grid to all nodes */
1330         if (pme->nnodes > 1)
1331         {
1332             gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1333         }
1334
1335         unwrap_periodic_pmegrid(pme, grid);
1336
1337         if (stepWork.computeForces)
1338         {
1339             /* interpolate forces for our local atoms */
1340
1341
1342             /* If we are running without parallelization,
1343              * atc->f is the actual force array, not a buffer,
1344              * therefore we should not clear it.
1345              */
1346             lambda  = grid_index < DO_Q ? lambda_q : lambda_lj;
1347             bClearF = (bFirst && PAR(cr));
1348 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1349             for (thread = 0; thread < pme->nthread; thread++)
1350             {
1351                 try
1352                 {
1353                     gather_f_bsplines(pme,
1354                                       grid,
1355                                       bClearF,
1356                                       &atc,
1357                                       &atc.spline[thread],
1358                                       pme->bFEP ? (grid_index % 2 == 0 ? 1.0 - lambda : lambda) : 1.0);
1359                 }
1360                 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1361             }
1362
1363
1364             inc_nrnb(nrnb, eNR_GATHERFBSP, pme->pme_order * pme->pme_order * pme->pme_order * atc.numAtoms());
1365             /* Note: this wallcycle region is opened above inside an OpenMP
1366                region, so take care if refactoring code here. */
1367             wallcycle_stop(wcycle, ewcPME_GATHER);
1368         }
1369
1370         if (computeEnergyAndVirial)
1371         {
1372             /* This should only be called on the master thread
1373              * and after the threads have synchronized.
1374              */
1375             if (grid_index < 2)
1376             {
1377                 get_pme_ener_vir_q(pme->solve_work, pme->nthread, &output[grid_index % 2]);
1378             }
1379             else
1380             {
1381                 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output[grid_index % 2]);
1382             }
1383         }
1384         bFirst = FALSE;
1385     } /* of grid_index-loop */
1386
1387     /* For Lorentz-Berthelot combination rules in LJ-PME, we need to calculate
1388      * seven terms. */
1389
1390     if (pme->doLJ && pme->ljpme_combination_rule == LongRangeVdW::LB)
1391     {
1392         /* Loop over A- and B-state if we are doing FEP */
1393         for (fep_state = 0; fep_state < fep_states_lj; ++fep_state)
1394         {
1395             real *local_c6 = nullptr, *local_sigma = nullptr, *RedistC6 = nullptr, *RedistSigma = nullptr;
1396             gmx::ArrayRef<real> coefficientBuffer;
1397             if (pme->nnodes == 1)
1398             {
1399                 pme->lb_buf1.resize(atc.numAtoms());
1400                 coefficientBuffer = pme->lb_buf1;
1401                 switch (fep_state)
1402                 {
1403                     case 0:
1404                         local_c6    = c6A;
1405                         local_sigma = sigmaA;
1406                         break;
1407                     case 1:
1408                         local_c6    = c6B;
1409                         local_sigma = sigmaB;
1410                         break;
1411                     default: gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1412                 }
1413             }
1414             else
1415             {
1416                 coefficientBuffer = atc.coefficientBuffer;
1417                 switch (fep_state)
1418                 {
1419                     case 0:
1420                         RedistC6    = c6A;
1421                         RedistSigma = sigmaA;
1422                         break;
1423                     case 1:
1424                         RedistC6    = c6B;
1425                         RedistSigma = sigmaB;
1426                         break;
1427                     default: gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1428                 }
1429                 wallcycle_start(wcycle, ewcPME_REDISTXF);
1430
1431                 do_redist_pos_coeffs(pme, cr, bFirst, coordinates, RedistC6);
1432                 pme->lb_buf1.resize(atc.numAtoms());
1433                 pme->lb_buf2.resize(atc.numAtoms());
1434                 local_c6 = pme->lb_buf1.data();
1435                 for (int i = 0; i < atc.numAtoms(); ++i)
1436                 {
1437                     local_c6[i] = atc.coefficient[i];
1438                 }
1439
1440                 do_redist_pos_coeffs(pme, cr, FALSE, coordinates, RedistSigma);
1441                 local_sigma = pme->lb_buf2.data();
1442                 for (int i = 0; i < atc.numAtoms(); ++i)
1443                 {
1444                     local_sigma[i] = atc.coefficient[i];
1445                 }
1446
1447                 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1448             }
1449             atc.coefficient = coefficientBuffer;
1450             calc_initial_lb_coeffs(coefficientBuffer, local_c6, local_sigma);
1451
1452             /*Seven terms in LJ-PME with LB, grid_index < 2 reserved for electrostatics*/
1453             for (grid_index = 2; grid_index < 9; ++grid_index)
1454             {
1455                 /* Unpack structure */
1456                 pmegrid    = &pme->pmegrid[grid_index];
1457                 fftgrid    = pme->fftgrid[grid_index];
1458                 pfft_setup = pme->pfft_setup[grid_index];
1459                 calc_next_lb_coeffs(coefficientBuffer, local_sigma);
1460                 grid = pmegrid->grid.grid;
1461
1462                 wallcycle_start(wcycle, ewcPME_SPREAD);
1463                 /* Spread the c6 on a grid */
1464                 spread_on_grid(pme, &atc, pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1465
1466                 if (bFirst)
1467                 {
1468                     inc_nrnb(nrnb, eNR_WEIGHTS, DIM * atc.numAtoms());
1469                 }
1470
1471                 inc_nrnb(nrnb,
1472                          eNR_SPREADBSP,
1473                          pme->pme_order * pme->pme_order * pme->pme_order * atc.numAtoms());
1474                 if (pme->nthread == 1)
1475                 {
1476                     wrap_periodic_pmegrid(pme, grid);
1477                     /* sum contributions to local grid from other nodes */
1478                     if (pme->nnodes > 1)
1479                     {
1480                         gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1481                     }
1482                     copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1483                 }
1484                 wallcycle_stop(wcycle, ewcPME_SPREAD);
1485
1486                 /*Here we start a large thread parallel region*/
1487 #pragma omp parallel num_threads(pme->nthread) private(thread)
1488                 {
1489                     try
1490                     {
1491                         thread = gmx_omp_get_thread_num();
1492                         /* do 3d-fft */
1493                         if (thread == 0)
1494                         {
1495                             wallcycle_start(wcycle, ewcPME_FFT);
1496                         }
1497
1498                         gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX, thread, wcycle);
1499                         if (thread == 0)
1500                         {
1501                             wallcycle_stop(wcycle, ewcPME_FFT);
1502                         }
1503                     }
1504                     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1505                 }
1506                 bFirst = FALSE;
1507             }
1508             /* solve in k-space for our local cells */
1509 #pragma omp parallel num_threads(pme->nthread) private(thread)
1510             {
1511                 try
1512                 {
1513                     int loop_count;
1514                     thread = gmx_omp_get_thread_num();
1515                     if (thread == 0)
1516                     {
1517                         wallcycle_start(wcycle, ewcLJPME);
1518                     }
1519
1520                     loop_count =
1521                             solve_pme_lj_yzx(pme,
1522                                              &pme->cfftgrid[2],
1523                                              TRUE,
1524                                              scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ],
1525                                              computeEnergyAndVirial,
1526                                              pme->nthread,
1527                                              thread);
1528                     if (thread == 0)
1529                     {
1530                         wallcycle_stop(wcycle, ewcLJPME);
1531                         inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1532                     }
1533                 }
1534                 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1535             }
1536
1537             if (computeEnergyAndVirial)
1538             {
1539                 /* This should only be called on the master thread and
1540                  * after the threads have synchronized.
1541                  */
1542                 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output[fep_state]);
1543             }
1544
1545             bFirst = !pme->doCoulomb;
1546             calc_initial_lb_coeffs(coefficientBuffer, local_c6, local_sigma);
1547             for (grid_index = 8; grid_index >= 2; --grid_index)
1548             {
1549                 /* Unpack structure */
1550                 pmegrid    = &pme->pmegrid[grid_index];
1551                 fftgrid    = pme->fftgrid[grid_index];
1552                 pfft_setup = pme->pfft_setup[grid_index];
1553                 grid       = pmegrid->grid.grid;
1554                 calc_next_lb_coeffs(coefficientBuffer, local_sigma);
1555 #pragma omp parallel num_threads(pme->nthread) private(thread)
1556                 {
1557                     try
1558                     {
1559                         thread = gmx_omp_get_thread_num();
1560                         /* do 3d-invfft */
1561                         if (thread == 0)
1562                         {
1563                             wallcycle_start(wcycle, ewcPME_FFT);
1564                         }
1565
1566                         gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL, thread, wcycle);
1567                         if (thread == 0)
1568                         {
1569                             wallcycle_stop(wcycle, ewcPME_FFT);
1570
1571
1572                             if (pme->nodeid == 0)
1573                             {
1574                                 real ntot = pme->nkx * pme->nky * pme->nkz;
1575                                 npme      = static_cast<int>(ntot * std::log(ntot) / std::log(2.0));
1576                                 inc_nrnb(nrnb, eNR_FFT, 2 * npme);
1577                             }
1578                             wallcycle_start(wcycle, ewcPME_GATHER);
1579                         }
1580
1581                         copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1582                     }
1583                     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1584                 } /*#pragma omp parallel*/
1585
1586                 /* distribute local grid to all nodes */
1587                 if (pme->nnodes > 1)
1588                 {
1589                     gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1590                 }
1591
1592                 unwrap_periodic_pmegrid(pme, grid);
1593
1594                 if (stepWork.computeForces)
1595                 {
1596                     /* interpolate forces for our local atoms */
1597                     bClearF = (bFirst && PAR(cr));
1598                     scale   = pme->bFEP ? (fep_state < 1 ? 1.0 - lambda_lj : lambda_lj) : 1.0;
1599                     scale *= lb_scale_factor[grid_index - 2];
1600
1601 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1602                     for (thread = 0; thread < pme->nthread; thread++)
1603                     {
1604                         try
1605                         {
1606                             gather_f_bsplines(
1607                                     pme, grid, bClearF, &pme->atc[0], &pme->atc[0].spline[thread], scale);
1608                         }
1609                         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1610                     }
1611
1612
1613                     inc_nrnb(nrnb,
1614                              eNR_GATHERFBSP,
1615                              pme->pme_order * pme->pme_order * pme->pme_order * pme->atc[0].numAtoms());
1616                 }
1617                 wallcycle_stop(wcycle, ewcPME_GATHER);
1618
1619                 bFirst = FALSE;
1620             } /* for (grid_index = 8; grid_index >= 2; --grid_index) */
1621         }     /* for (fep_state = 0; fep_state < fep_states_lj; ++fep_state) */
1622     }         /* if (pme->doLJ && pme->ljpme_combination_rule == LongRangeVdW::LB) */
1623
1624     if (stepWork.computeForces && pme->nnodes > 1)
1625     {
1626         wallcycle_start(wcycle, ewcPME_REDISTXF);
1627         for (d = 0; d < pme->ndecompdim; d++)
1628         {
1629             gmx::ArrayRef<gmx::RVec> forcesRef;
1630             if (d == pme->ndecompdim - 1)
1631             {
1632                 const size_t numAtoms = coordinates.size();
1633                 GMX_ASSERT(forces.size() >= numAtoms, "Need at least numAtoms forces");
1634                 forcesRef = forces.subArray(0, numAtoms);
1635             }
1636             else
1637             {
1638                 forcesRef = pme->atc[d + 1].f;
1639             }
1640             if (DOMAINDECOMP(cr))
1641             {
1642                 dd_pmeredist_f(pme, &pme->atc[d], forcesRef, d == pme->ndecompdim - 1 && pme->bPPnode);
1643             }
1644         }
1645
1646         wallcycle_stop(wcycle, ewcPME_REDISTXF);
1647     }
1648
1649     if (computeEnergyAndVirial)
1650     {
1651         if (pme->doCoulomb)
1652         {
1653             if (!pme->bFEP_q)
1654             {
1655                 *energy_q = output[0].coulombEnergy_;
1656                 m_add(vir_q, output[0].coulombVirial_, vir_q);
1657             }
1658             else
1659             {
1660                 *energy_q = (1.0 - lambda_q) * output[0].coulombEnergy_ + lambda_q * output[1].coulombEnergy_;
1661                 *dvdlambda_q += output[1].coulombEnergy_ - output[0].coulombEnergy_;
1662                 for (int i = 0; i < DIM; i++)
1663                 {
1664                     for (int j = 0; j < DIM; j++)
1665                     {
1666                         vir_q[i][j] += (1.0 - lambda_q) * output[0].coulombVirial_[i][j]
1667                                        + lambda_q * output[1].coulombVirial_[i][j];
1668                     }
1669                 }
1670             }
1671             if (debug)
1672             {
1673                 fprintf(debug, "Electrostatic PME mesh energy: %g\n", *energy_q);
1674             }
1675         }
1676         else
1677         {
1678             *energy_q = 0;
1679         }
1680
1681         if (pme->doLJ)
1682         {
1683             if (!pme->bFEP_lj)
1684             {
1685                 *energy_lj = output[0].lennardJonesEnergy_;
1686                 m_add(vir_lj, output[0].lennardJonesVirial_, vir_lj);
1687             }
1688             else
1689             {
1690                 *energy_lj = (1.0 - lambda_lj) * output[0].lennardJonesEnergy_
1691                              + lambda_lj * output[1].lennardJonesEnergy_;
1692                 *dvdlambda_lj += output[1].lennardJonesEnergy_ - output[0].lennardJonesEnergy_;
1693                 for (int i = 0; i < DIM; i++)
1694                 {
1695                     for (int j = 0; j < DIM; j++)
1696                     {
1697                         vir_lj[i][j] += (1.0 - lambda_lj) * output[0].lennardJonesVirial_[i][j]
1698                                         + lambda_lj * output[1].lennardJonesVirial_[i][j];
1699                     }
1700                 }
1701             }
1702             if (debug)
1703             {
1704                 fprintf(debug, "Lennard-Jones PME mesh energy: %g\n", *energy_lj);
1705             }
1706         }
1707         else
1708         {
1709             *energy_lj = 0;
1710         }
1711     }
1712     return 0;
1713 }
1714
1715 void gmx_pme_destroy(gmx_pme_t* pme)
1716 {
1717     if (!pme)
1718     {
1719         return;
1720     }
1721
1722     delete pme->boxScaler;
1723
1724     sfree(pme->nnx);
1725     sfree(pme->nny);
1726     sfree(pme->nnz);
1727     sfree(pme->fshx);
1728     sfree(pme->fshy);
1729     sfree(pme->fshz);
1730
1731     for (int i = 0; i < pme->ngrids; ++i)
1732     {
1733         pmegrids_destroy(&pme->pmegrid[i]);
1734     }
1735     if (pme->pfft_setup)
1736     {
1737         for (int i = 0; i < pme->ngrids; ++i)
1738         {
1739             gmx_parallel_3dfft_destroy(pme->pfft_setup[i]);
1740         }
1741     }
1742     sfree(pme->fftgrid);
1743     sfree(pme->cfftgrid);
1744     sfree(pme->pfft_setup);
1745
1746     for (int i = 0; i < DIM; i++)
1747     {
1748         sfree(pme->bsp_mod[i]);
1749     }
1750
1751     sfree(pme->bufv);
1752     sfree(pme->bufr);
1753
1754     if (pme->solve_work)
1755     {
1756         pme_free_all_work(&pme->solve_work, pme->nthread);
1757     }
1758
1759     sfree(pme->sum_qgrid_tmp);
1760     sfree(pme->sum_qgrid_dd_tmp);
1761
1762     destroy_pme_spline_work(pme->spline_work);
1763
1764     if (pme->gpu != nullptr)
1765     {
1766         pme_gpu_destroy(pme->gpu);
1767     }
1768
1769     delete pme;
1770 }
1771
1772 void gmx_pme_reinit_atoms(gmx_pme_t* pme, const int numAtoms, const real* chargesA, const real* chargesB)
1773 {
1774     if (pme->gpu != nullptr)
1775     {
1776         GMX_ASSERT(!(pme->bFEP_q && chargesB == nullptr),
1777                    "B state charges must be specified if running Coulomb FEP on the GPU");
1778         pme_gpu_reinit_atoms(pme->gpu, numAtoms, chargesA, pme->bFEP_q ? chargesB : nullptr);
1779     }
1780     else
1781     {
1782         pme->atc[0].setNumAtoms(numAtoms);
1783         // TODO: set the charges here as well
1784     }
1785 }
1786
1787 bool gmx_pme_grid_matches(const gmx_pme_t& pme, const ivec grid_size)
1788 {
1789     return (pme.nkx == grid_size[XX] && pme.nky == grid_size[YY] && pme.nkz == grid_size[ZZ]);
1790 }
1791
1792 void gmx::SeparatePmeRanksPermitted::disablePmeRanks(const std::string& reason)
1793 {
1794     permitSeparatePmeRanks_ = false;
1795
1796     if (!reason.empty())
1797     {
1798         reasons_.push_back(reason);
1799     }
1800 }
1801
1802 bool gmx::SeparatePmeRanksPermitted::permitSeparatePmeRanks() const
1803 {
1804     return permitSeparatePmeRanks_;
1805 }
1806
1807 std::string gmx::SeparatePmeRanksPermitted::reasonsWhyDisabled() const
1808 {
1809     return joinStrings(reasons_, "; ");
1810 }