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