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