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