2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017, by the GROMACS development team, 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.
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.
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.
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.
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.
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.
39 * \brief This file contains function definitions necessary for
40 * computing energies and forces for the PME long-ranged part (Coulomb
43 * \author Erik Lindahl <erik@kth.se>
44 * \author Berk Hess <hess@kth.se>
45 * \ingroup module_ewald
47 /* IMPORTANT FOR DEVELOPERS:
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:
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:
59 * rectangular triclinic
64 * 2. You should check the energy conservation in a triclinic box.
66 * It might seem an overkill, but better safe than sorry.
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"
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"
124 /*! \brief Number of bytes in a cache line.
126 * Must also be a multiple of the SIMD and SIMD4 register size, to
127 * preserve alignment.
129 const int gmxCacheLineSize = 64;
131 //! Set up coordinate communication
132 static void setup_coordinate_communication(pme_atomcomm_t *atc)
140 for (i = 1; i <= nslab/2; i++)
142 fw = (atc->nodeid + i) % nslab;
143 bw = (atc->nodeid - i + nslab) % nslab;
146 atc->node_dest[n] = fw;
147 atc->node_src[n] = bw;
152 atc->node_dest[n] = bw;
153 atc->node_src[n] = fw;
159 /*! \brief Round \p n up to the next multiple of \p f */
160 static int mult_up(int n, int f)
162 return ((n + f - 1)/f)*f;
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)
171 nma = pme->nnodes_major;
172 nmi = pme->nnodes_minor;
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;
178 /* pme_solve is roughly double the cost of an fft */
180 return (n1 + n2 + 3*n3)/(double)(6*pme->nkx*pme->nky*pme->nkz);
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)
189 atc->dimind = dimind;
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);
202 fprintf(debug, "For PME atom communication in dimind %d: nslab %d rank %d\n", atc->dimind, atc->nslab, atc->nodeid);
206 atc->bSpread = bSpread;
207 atc->pme_order = pme->pme_order;
211 snew(atc->node_dest, atc->nslab);
212 snew(atc->node_src, atc->nslab);
213 setup_coordinate_communication(atc);
215 snew(atc->count_thread, pme->nthread);
216 for (thread = 0; thread < pme->nthread; thread++)
218 snew(atc->count_thread[thread], atc->nslab);
220 atc->count = atc->count_thread[0];
221 snew(atc->rcount, atc->nslab);
222 snew(atc->buf_index, atc->nslab);
225 atc->nthread = pme->nthread;
226 if (atc->nthread > 1)
228 snew(atc->thread_plist, atc->nthread);
230 snew(atc->spline, atc->nthread);
231 for (thread = 0; thread < atc->nthread; thread++)
233 if (atc->nthread > 1)
235 snew(atc->thread_plist[thread].n, atc->nthread+2*gmxCacheLineSize);
236 atc->thread_plist[thread].n += gmxCacheLineSize;
241 /*! \brief Destroy an atom communication data structure and its child structs */
242 static void destroy_atomcomm(pme_atomcomm_t *atc)
247 sfree(atc->node_dest);
248 sfree(atc->node_src);
249 for (int i = 0; i < atc->nthread; i++)
251 sfree(atc->count_thread[i]);
253 sfree(atc->count_thread);
255 sfree(atc->buf_index);
258 sfree(atc->coefficient);
264 sfree(atc->thread_idx);
265 for (int i = 0; i < atc->nthread; i++)
267 if (atc->nthread > 1)
269 int *n_ptr = atc->thread_plist[i].n - gmxCacheLineSize;
271 sfree(atc->thread_plist[i].i);
273 sfree(atc->spline[i].ind);
274 for (int d = 0; d < ZZ; d++)
276 sfree(atc->spline[i].theta[d]);
277 sfree(atc->spline[i].dtheta[d]);
279 sfree_aligned(atc->spline[i].ptr_dtheta_z);
280 sfree_aligned(atc->spline[i].ptr_theta_z);
282 if (atc->nthread > 1)
284 sfree(atc->thread_plist);
289 /*! \brief Initialize data structure for communication */
291 init_overlap_comm(pme_overlap_t * ol,
302 pme_grid_comm_t *pgc;
304 int fft_start, fft_end, send_index1, recv_index1;
314 /* Linear translation of the PME grid won't affect reciprocal space
315 * calculations, so to optimize we only interpolate "upwards",
316 * which also means we only have to consider overlap in one direction.
317 * I.e., particles on this node might also be spread to grid indices
318 * that belong to higher nodes (modulo nnodes)
321 snew(ol->s2g0, ol->nnodes+1);
322 snew(ol->s2g1, ol->nnodes);
325 fprintf(debug, "PME slab boundaries:");
327 for (i = 0; i < nnodes; i++)
329 /* s2g0 the local interpolation grid start.
330 * s2g1 the local interpolation grid end.
331 * Since in calc_pidx we divide particles, and not grid lines,
332 * spatially uniform along dimension x or y, we need to round
333 * s2g0 down and s2g1 up.
335 ol->s2g0[i] = ( i *ndata + 0 )/nnodes;
336 ol->s2g1[i] = ((i+1)*ndata + nnodes-1)/nnodes + norder - 1;
340 fprintf(debug, " %3d %3d", ol->s2g0[i], ol->s2g1[i]);
343 ol->s2g0[nnodes] = ndata;
346 fprintf(debug, "\n");
349 /* Determine with how many nodes we need to communicate the grid overlap */
355 for (i = 0; i < nnodes; i++)
357 if ((i+b < nnodes && ol->s2g1[i] > ol->s2g0[i+b]) ||
358 (i+b >= nnodes && ol->s2g1[i] > ol->s2g0[i+b-nnodes] + ndata))
364 while (bCont && b < nnodes);
365 ol->noverlap_nodes = b - 1;
367 snew(ol->send_id, ol->noverlap_nodes);
368 snew(ol->recv_id, ol->noverlap_nodes);
369 for (b = 0; b < ol->noverlap_nodes; b++)
371 ol->send_id[b] = (ol->nodeid + (b + 1)) % ol->nnodes;
372 ol->recv_id[b] = (ol->nodeid - (b + 1) + ol->nnodes) % ol->nnodes;
374 snew(ol->comm_data, ol->noverlap_nodes);
377 for (b = 0; b < ol->noverlap_nodes; b++)
379 pgc = &ol->comm_data[b];
381 fft_start = ol->s2g0[ol->send_id[b]];
382 fft_end = ol->s2g0[ol->send_id[b]+1];
383 if (ol->send_id[b] < nodeid)
388 send_index1 = ol->s2g1[nodeid];
389 send_index1 = std::min(send_index1, fft_end);
390 pgc->send_index0 = fft_start;
391 pgc->send_nindex = std::max(0, send_index1 - pgc->send_index0);
392 ol->send_size += pgc->send_nindex;
394 /* We always start receiving to the first index of our slab */
395 fft_start = ol->s2g0[ol->nodeid];
396 fft_end = ol->s2g0[ol->nodeid+1];
397 recv_index1 = ol->s2g1[ol->recv_id[b]];
398 if (ol->recv_id[b] > nodeid)
400 recv_index1 -= ndata;
402 recv_index1 = std::min(recv_index1, fft_end);
403 pgc->recv_index0 = fft_start;
404 pgc->recv_nindex = std::max(0, recv_index1 - pgc->recv_index0);
408 /* Communicate the buffer sizes to receive */
409 for (b = 0; b < ol->noverlap_nodes; b++)
411 MPI_Sendrecv(&ol->send_size, 1, MPI_INT, ol->send_id[b], b,
412 &ol->comm_data[b].recv_size, 1, MPI_INT, ol->recv_id[b], b,
413 ol->mpi_comm, &stat);
417 /* For non-divisible grid we need pme_order iso pme_order-1 */
418 snew(ol->sendbuf, norder*commplainsize);
419 snew(ol->recvbuf, norder*commplainsize);
422 /*! \brief Destroy data structure for communication */
424 destroy_overlap_comm(const pme_overlap_t *ol)
430 sfree(ol->comm_data);
435 int minimalPmeGridSize(int pmeOrder)
437 /* The actual grid size limitations are:
438 * serial: >= pme_order
439 * DD, no OpenMP: >= 2*(pme_order - 1)
440 * DD, OpenMP: >= pme_order + 1
441 * But we use the maximum for simplicity since in practice there is not
442 * much performance difference between pme_order and 2*(pme_order -1).
444 int minimalSize = 2*(pmeOrder - 1);
446 GMX_RELEASE_ASSERT(pmeOrder >= 3, "pmeOrder has to be >= 3");
447 GMX_RELEASE_ASSERT(minimalSize >= pmeOrder + 1, "The grid size should be >= pmeOrder + 1");
452 bool gmx_pme_check_restrictions(int pme_order,
453 int nkx, int nky, int nkz,
458 if (pme_order > PME_ORDER_MAX)
465 std::string message = gmx::formatString(
466 "pme_order (%d) is larger than the maximum allowed value (%d). Modify and recompile the code if you really need such a high order.",
467 pme_order, PME_ORDER_MAX);
468 GMX_THROW(InconsistentInputError(message));
471 const int minGridSize = minimalPmeGridSize(pme_order);
472 if (nkx < minGridSize ||
480 std::string message = gmx::formatString(
481 "The PME grid sizes need to be >= 2*(pme_order-1) (%d)",
483 GMX_THROW(InconsistentInputError(message));
486 /* Check for a limitation of the (current) sum_fftgrid_dd code.
487 * We only allow multiple communication pulses in dim 1, not in dim 0.
489 if (useThreads && (nkx < nnodes_major*pme_order &&
490 nkx != nnodes_major*(pme_order - 1)))
496 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.",
497 nkx/(double)nnodes_major, pme_order);
503 /*! \brief Round \p enumerator */
504 static int div_round_up(int enumerator, int denominator)
506 return (enumerator + denominator - 1)/denominator;
509 gmx_pme_t *gmx_pme_init(const t_commrec *cr,
512 const t_inputrec *ir,
514 gmx_bool bFreeEnergy_q,
515 gmx_bool bFreeEnergy_lj,
516 gmx_bool bReproducible,
522 gmx_device_info_t *gpuInfo,
523 const gmx::MDLogger & /*mdlog*/)
525 int use_threads, sum_use_threads, i;
530 fprintf(debug, "Creating PME data structures.\n");
533 unique_cptr<gmx_pme_t, gmx_pme_destroy> pme(new gmx_pme_t());
535 pme->sum_qgrid_tmp = nullptr;
536 pme->sum_qgrid_dd_tmp = nullptr;
543 pme->nnodes_major = nnodes_major;
544 pme->nnodes_minor = nnodes_minor;
547 if (nnodes_major*nnodes_minor > 1)
549 pme->mpi_comm = cr->mpi_comm_mygroup;
551 MPI_Comm_rank(pme->mpi_comm, &pme->nodeid);
552 MPI_Comm_size(pme->mpi_comm, &pme->nnodes);
553 if (pme->nnodes != nnodes_major*nnodes_minor)
555 gmx_incons("PME rank count mismatch");
560 pme->mpi_comm = MPI_COMM_NULL;
564 if (pme->nnodes == 1)
567 pme->mpi_comm_d[0] = MPI_COMM_NULL;
568 pme->mpi_comm_d[1] = MPI_COMM_NULL;
571 pme->nodeid_major = 0;
572 pme->nodeid_minor = 0;
574 pme->mpi_comm_d[0] = pme->mpi_comm_d[1] = MPI_COMM_NULL;
579 if (nnodes_minor == 1)
582 pme->mpi_comm_d[0] = pme->mpi_comm;
583 pme->mpi_comm_d[1] = MPI_COMM_NULL;
586 pme->nodeid_major = pme->nodeid;
587 pme->nodeid_minor = 0;
590 else if (nnodes_major == 1)
593 pme->mpi_comm_d[0] = MPI_COMM_NULL;
594 pme->mpi_comm_d[1] = pme->mpi_comm;
597 pme->nodeid_major = 0;
598 pme->nodeid_minor = pme->nodeid;
602 if (pme->nnodes % nnodes_major != 0)
604 gmx_incons("For 2D PME decomposition, #PME ranks must be divisible by the number of ranks in the major dimension");
609 MPI_Comm_split(pme->mpi_comm, pme->nodeid % nnodes_minor,
610 pme->nodeid, &pme->mpi_comm_d[0]); /* My communicator along major dimension */
611 MPI_Comm_split(pme->mpi_comm, pme->nodeid/nnodes_minor,
612 pme->nodeid, &pme->mpi_comm_d[1]); /* My communicator along minor dimension */
614 MPI_Comm_rank(pme->mpi_comm_d[0], &pme->nodeid_major);
615 MPI_Comm_size(pme->mpi_comm_d[0], &pme->nnodes_major);
616 MPI_Comm_rank(pme->mpi_comm_d[1], &pme->nodeid_minor);
617 MPI_Comm_size(pme->mpi_comm_d[1], &pme->nnodes_minor);
620 pme->bPPnode = thisRankHasDuty(cr, DUTY_PP);
623 pme->nthread = nthread;
625 /* Check if any of the PME MPI ranks uses threads */
626 use_threads = (pme->nthread > 1 ? 1 : 0);
630 MPI_Allreduce(&use_threads, &sum_use_threads, 1, MPI_INT,
631 MPI_SUM, pme->mpi_comm);
636 sum_use_threads = use_threads;
638 pme->bUseThreads = (sum_use_threads > 0);
640 if (ir->ePBC == epbcSCREW)
642 gmx_fatal(FARGS, "pme does not (yet) work with pbc = screw");
646 * It is likely that the current gmx_pme_do() routine supports calculating
647 * only Coulomb or LJ while gmx_pme_init() configures for both,
648 * but that has never been tested.
649 * It is likely that the current gmx_pme_do() routine supports calculating,
650 * not calculating free-energy for Coulomb and/or LJ while gmx_pme_init()
651 * configures with free-energy, but that has never been tested.
653 pme->doCoulomb = EEL_PME(ir->coulombtype);
654 pme->doLJ = EVDW_PME(ir->vdwtype);
655 pme->bFEP_q = ((ir->efep != efepNO) && bFreeEnergy_q);
656 pme->bFEP_lj = ((ir->efep != efepNO) && bFreeEnergy_lj);
657 pme->bFEP = (pme->bFEP_q || pme->bFEP_lj);
661 pme->bP3M = (ir->coulombtype == eelP3M_AD || getenv("GMX_PME_P3M") != nullptr);
662 pme->pme_order = ir->pme_order;
663 pme->ewaldcoeff_q = ewaldcoeff_q;
664 pme->ewaldcoeff_lj = ewaldcoeff_lj;
666 /* Always constant electrostatics coefficients */
667 pme->epsilon_r = ir->epsilon_r;
669 /* Always constant LJ coefficients */
670 pme->ljpme_combination_rule = ir->ljpme_combination_rule;
672 // The box requires scaling with nwalls = 2, we store that condition as well
673 // as the scaling factor
674 delete pme->boxScaler;
675 pme->boxScaler = new EwaldBoxZScaler(*ir);
677 /* If we violate restrictions, generate a fatal error here */
678 gmx_pme_check_restrictions(pme->pme_order,
679 pme->nkx, pme->nky, pme->nkz,
689 MPI_Type_contiguous(DIM, GMX_MPI_REAL, &(pme->rvec_mpi));
690 MPI_Type_commit(&(pme->rvec_mpi));
693 /* Note that the coefficient spreading and force gathering, which usually
694 * takes about the same amount of time as FFT+solve_pme,
695 * is always fully load balanced
696 * (unless the coefficient distribution is inhomogeneous).
699 imbal = estimate_pme_load_imbalance(pme.get());
700 if (imbal >= 1.2 && pme->nodeid_major == 0 && pme->nodeid_minor == 0)
704 "NOTE: The load imbalance in PME FFT and solve is %d%%.\n"
705 " For optimal PME load balancing\n"
706 " PME grid_x (%d) and grid_y (%d) should be divisible by #PME_ranks_x (%d)\n"
707 " and PME grid_y (%d) and grid_z (%d) should be divisible by #PME_ranks_y (%d)\n"
709 (int)((imbal-1)*100 + 0.5),
710 pme->nkx, pme->nky, pme->nnodes_major,
711 pme->nky, pme->nkz, pme->nnodes_minor);
715 /* For non-divisible grid we need pme_order iso pme_order-1 */
716 /* In sum_qgrid_dd x overlap is copied in place: take padding into account.
717 * y is always copied through a buffer: we don't need padding in z,
718 * but we do need the overlap in x because of the communication order.
720 init_overlap_comm(&pme->overlap[0], pme->pme_order,
724 pme->nnodes_major, pme->nodeid_major,
726 (div_round_up(pme->nky, pme->nnodes_minor)+pme->pme_order)*(pme->nkz+pme->pme_order-1));
728 /* Along overlap dim 1 we can send in multiple pulses in sum_fftgrid_dd.
729 * We do this with an offset buffer of equal size, so we need to allocate
730 * extra for the offset. That's what the (+1)*pme->nkz is for.
732 init_overlap_comm(&pme->overlap[1], pme->pme_order,
736 pme->nnodes_minor, pme->nodeid_minor,
738 (div_round_up(pme->nkx, pme->nnodes_major)+pme->pme_order+1)*pme->nkz);
740 /* Double-check for a limitation of the (current) sum_fftgrid_dd code.
741 * Note that gmx_pme_check_restrictions checked for this already.
743 if (pme->bUseThreads && pme->overlap[0].noverlap_nodes > 1)
745 gmx_incons("More than one communication pulse required for grid overlap communication along the major dimension while using threads");
748 snew(pme->bsp_mod[XX], pme->nkx);
749 snew(pme->bsp_mod[YY], pme->nky);
750 snew(pme->bsp_mod[ZZ], pme->nkz);
752 pme->gpu = pmeGPU; /* Carrying over the single GPU structure */
753 pme->runMode = runMode;
755 /* The required size of the interpolation grid, including overlap.
756 * The allocated size (pmegrid_n?) might be slightly larger.
758 pme->pmegrid_nx = pme->overlap[0].s2g1[pme->nodeid_major] -
759 pme->overlap[0].s2g0[pme->nodeid_major];
760 pme->pmegrid_ny = pme->overlap[1].s2g1[pme->nodeid_minor] -
761 pme->overlap[1].s2g0[pme->nodeid_minor];
762 pme->pmegrid_nz_base = pme->nkz;
763 pme->pmegrid_nz = pme->pmegrid_nz_base + pme->pme_order - 1;
764 set_grid_alignment(&pme->pmegrid_nz, pme->pme_order);
765 pme->pmegrid_start_ix = pme->overlap[0].s2g0[pme->nodeid_major];
766 pme->pmegrid_start_iy = pme->overlap[1].s2g0[pme->nodeid_minor];
767 pme->pmegrid_start_iz = 0;
769 make_gridindex_to_localindex(pme->nkx,
770 pme->pmegrid_start_ix,
771 pme->pmegrid_nx - (pme->pme_order-1),
772 &pme->nnx, &pme->fshx);
773 make_gridindex_to_localindex(pme->nky,
774 pme->pmegrid_start_iy,
775 pme->pmegrid_ny - (pme->pme_order-1),
776 &pme->nny, &pme->fshy);
777 make_gridindex_to_localindex(pme->nkz,
778 pme->pmegrid_start_iz,
779 pme->pmegrid_nz_base,
780 &pme->nnz, &pme->fshz);
782 pme->spline_work = make_pme_spline_work(pme->pme_order);
787 /* It doesn't matter if we allocate too many grids here,
788 * we only allocate and use the ones we need.
792 pme->ngrids = ((ir->ljpme_combination_rule == eljpmeLB) ? DO_Q_AND_LJ_LB : DO_Q_AND_LJ);
798 snew(pme->fftgrid, pme->ngrids);
799 snew(pme->cfftgrid, pme->ngrids);
800 snew(pme->pfft_setup, pme->ngrids);
802 for (i = 0; i < pme->ngrids; ++i)
804 if ((i < DO_Q && pme->doCoulomb && (i == 0 ||
806 (i >= DO_Q && pme->doLJ && (i == 2 ||
808 ir->ljpme_combination_rule == eljpmeLB)))
810 pmegrids_init(&pme->pmegrid[i],
811 pme->pmegrid_nx, pme->pmegrid_ny, pme->pmegrid_nz,
812 pme->pmegrid_nz_base,
816 pme->overlap[0].s2g1[pme->nodeid_major]-pme->overlap[0].s2g0[pme->nodeid_major+1],
817 pme->overlap[1].s2g1[pme->nodeid_minor]-pme->overlap[1].s2g0[pme->nodeid_minor+1]);
818 /* This routine will allocate the grid data to fit the FFTs */
819 const auto allocateRealGridForGpu = (pme->runMode == PmeRunMode::Hybrid) ? gmx::PinningPolicy::CanBePinned : gmx::PinningPolicy::CannotBePinned;
820 gmx_parallel_3dfft_init(&pme->pfft_setup[i], ndata,
821 &pme->fftgrid[i], &pme->cfftgrid[i],
823 bReproducible, pme->nthread, allocateRealGridForGpu);
830 /* Use plain SPME B-spline interpolation */
831 make_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
835 /* Use the P3M grid-optimized influence function */
836 make_p3m_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
839 /* Use atc[0] for spreading */
840 init_atomcomm(pme.get(), &pme->atc[0], nnodes_major > 1 ? 0 : 1, TRUE);
841 if (pme->ndecompdim >= 2)
843 init_atomcomm(pme.get(), &pme->atc[1], 1, FALSE);
846 if (pme->nnodes == 1)
848 pme->atc[0].n = homenr;
849 pme_realloc_atomcomm_things(&pme->atc[0]);
852 pme->lb_buf1 = nullptr;
853 pme->lb_buf2 = nullptr;
854 pme->lb_buf_nalloc = 0;
856 pme_gpu_reinit(pme.get(), gpuInfo);
858 pme_init_all_work(&pme->solve_work, pme->nthread, pme->nkx);
860 // no exception was thrown during the init, so we hand over the PME structure handle
861 return pme.release();
864 void gmx_pme_reinit(struct gmx_pme_t **pmedata,
866 struct gmx_pme_t * pme_src,
867 const t_inputrec * ir,
868 const ivec grid_size,
874 // Create a copy of t_inputrec fields that are used in gmx_pme_init().
875 // TODO: This would be better as just copying a sub-structure that contains
876 // all the PME parameters and nothing else.
879 irc.coulombtype = ir->coulombtype;
880 irc.vdwtype = ir->vdwtype;
882 irc.pme_order = ir->pme_order;
883 irc.epsilon_r = ir->epsilon_r;
884 irc.ljpme_combination_rule = ir->ljpme_combination_rule;
885 irc.nkx = grid_size[XX];
886 irc.nky = grid_size[YY];
887 irc.nkz = grid_size[ZZ];
889 if (pme_src->nnodes == 1)
891 homenr = pme_src->atc[0].n;
900 const gmx::MDLogger dummyLogger;
901 // This is reinit which is currently only changing grid size/coefficients,
902 // so we don't expect the actual logging.
903 // TODO: when PME is an object, it should take reference to mdlog on construction and save it.
904 GMX_ASSERT(pmedata, "Invalid PME pointer");
905 *pmedata = gmx_pme_init(cr, pme_src->nnodes_major, pme_src->nnodes_minor,
906 &irc, homenr, pme_src->bFEP_q, pme_src->bFEP_lj, FALSE, ewaldcoeff_q, ewaldcoeff_lj,
907 pme_src->nthread, pme_src->runMode, pme_src->gpu, nullptr, dummyLogger);
908 //TODO this is mostly passing around current values
910 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
912 /* We can easily reuse the allocated pme grids in pme_src */
913 reuse_pmegrids(&pme_src->pmegrid[PME_GRID_QA], &(*pmedata)->pmegrid[PME_GRID_QA]);
914 /* We would like to reuse the fft grids, but that's harder */
917 void gmx_pme_calc_energy(struct gmx_pme_t *pme, int n, rvec *x, real *q, real *V)
924 gmx_incons("gmx_pme_calc_energy called in parallel");
928 gmx_incons("gmx_pme_calc_energy with free energy");
931 atc = &pme->atc_energy;
933 if (atc->spline == nullptr)
935 snew(atc->spline, atc->nthread);
939 atc->pme_order = pme->pme_order;
941 pme_realloc_atomcomm_things(atc);
943 atc->coefficient = q;
945 /* We only use the A-charges grid */
946 grid = &pme->pmegrid[PME_GRID_QA];
948 /* Only calculate the spline coefficients, don't actually spread */
949 spread_on_grid(pme, atc, nullptr, TRUE, FALSE, pme->fftgrid[PME_GRID_QA], FALSE, PME_GRID_QA);
951 *V = gather_energy_bsplines(pme, grid->grid.grid, atc);
954 /*! \brief Calculate initial Lorentz-Berthelot coefficients for LJ-PME */
956 calc_initial_lb_coeffs(struct gmx_pme_t *pme, real *local_c6, real *local_sigma)
959 for (i = 0; i < pme->atc[0].n; ++i)
962 sigma4 = local_sigma[i];
963 sigma4 = sigma4*sigma4;
964 sigma4 = sigma4*sigma4;
965 pme->atc[0].coefficient[i] = local_c6[i] / sigma4;
969 /*! \brief Calculate next Lorentz-Berthelot coefficients for LJ-PME */
971 calc_next_lb_coeffs(struct gmx_pme_t *pme, real *local_sigma)
975 for (i = 0; i < pme->atc[0].n; ++i)
977 pme->atc[0].coefficient[i] *= local_sigma[i];
981 int gmx_pme_do(struct gmx_pme_t *pme,
982 int start, int homenr,
984 real chargeA[], real chargeB[],
985 real c6A[], real c6B[],
986 real sigmaA[], real sigmaB[],
987 matrix box, t_commrec *cr,
988 int maxshift_x, int maxshift_y,
989 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
990 matrix vir_q, matrix vir_lj,
991 real *energy_q, real *energy_lj,
992 real lambda_q, real lambda_lj,
993 real *dvdlambda_q, real *dvdlambda_lj,
996 GMX_ASSERT(pme->runMode == PmeRunMode::CPU, "gmx_pme_do should not be called on the GPU PME run.");
998 int d, i, j, npme, grid_index, max_grid_index;
1000 pme_atomcomm_t *atc = nullptr;
1001 pmegrids_t *pmegrid = nullptr;
1002 real *grid = nullptr;
1004 real *coefficient = nullptr;
1009 gmx_parallel_3dfft_t pfft_setup;
1011 t_complex * cfftgrid;
1013 gmx_bool bFirst, bDoSplines;
1015 int fep_states_lj = pme->bFEP_lj ? 2 : 1;
1016 const gmx_bool bCalcEnerVir = flags & GMX_PME_CALC_ENER_VIR;
1017 const gmx_bool bBackFFT = flags & (GMX_PME_CALC_F | GMX_PME_CALC_POT);
1018 const gmx_bool bCalcF = flags & GMX_PME_CALC_F;
1020 assert(pme->nnodes > 0);
1021 assert(pme->nnodes == 1 || pme->ndecompdim > 0);
1023 if (pme->nnodes > 1)
1027 if (atc->npd > atc->pd_nalloc)
1029 atc->pd_nalloc = over_alloc_dd(atc->npd);
1030 srenew(atc->pd, atc->pd_nalloc);
1032 for (d = pme->ndecompdim-1; d >= 0; d--)
1035 atc->maxshift = (atc->dimind == 0 ? maxshift_x : maxshift_y);
1041 /* This could be necessary for TPI */
1042 pme->atc[0].n = homenr;
1043 if (DOMAINDECOMP(cr))
1045 pme_realloc_atomcomm_things(atc);
1052 pme->boxScaler->scaleBox(box, scaledBox);
1054 gmx::invertBoxMatrix(scaledBox, pme->recipbox);
1057 /* For simplicity, we construct the splines for all particles if
1058 * more than one PME calculations is needed. Some optimization
1059 * could be done by keeping track of which atoms have splines
1060 * constructed, and construct new splines on each pass for atoms
1061 * that don't yet have them.
1064 bDoSplines = pme->bFEP || (pme->doCoulomb && pme->doLJ);
1066 /* We need a maximum of four separate PME calculations:
1067 * grid_index=0: Coulomb PME with charges from state A
1068 * grid_index=1: Coulomb PME with charges from state B
1069 * grid_index=2: LJ PME with C6 from state A
1070 * grid_index=3: LJ PME with C6 from state B
1071 * For Lorentz-Berthelot combination rules, a separate loop is used to
1072 * calculate all the terms
1075 /* If we are doing LJ-PME with LB, we only do Q here */
1076 max_grid_index = (pme->ljpme_combination_rule == eljpmeLB) ? DO_Q : DO_Q_AND_LJ;
1078 for (grid_index = 0; grid_index < max_grid_index; ++grid_index)
1080 /* Check if we should do calculations at this grid_index
1081 * If grid_index is odd we should be doing FEP
1082 * If grid_index < 2 we should be doing electrostatic PME
1083 * If grid_index >= 2 we should be doing LJ-PME
1085 if ((grid_index < DO_Q && (!pme->doCoulomb ||
1086 (grid_index == 1 && !pme->bFEP_q))) ||
1087 (grid_index >= DO_Q && (!pme->doLJ ||
1088 (grid_index == 3 && !pme->bFEP_lj))))
1092 /* Unpack structure */
1093 pmegrid = &pme->pmegrid[grid_index];
1094 fftgrid = pme->fftgrid[grid_index];
1095 cfftgrid = pme->cfftgrid[grid_index];
1096 pfft_setup = pme->pfft_setup[grid_index];
1099 case 0: coefficient = chargeA + start; break;
1100 case 1: coefficient = chargeB + start; break;
1101 case 2: coefficient = c6A + start; break;
1102 case 3: coefficient = c6B + start; break;
1105 grid = pmegrid->grid.grid;
1109 fprintf(debug, "PME: number of ranks = %d, rank = %d\n",
1110 cr->nnodes, cr->nodeid);
1111 fprintf(debug, "Grid = %p\n", (void*)grid);
1112 if (grid == nullptr)
1114 gmx_fatal(FARGS, "No grid!");
1119 if (pme->nnodes == 1)
1121 atc->coefficient = coefficient;
1125 wallcycle_start(wcycle, ewcPME_REDISTXF);
1126 do_redist_pos_coeffs(pme, cr, start, homenr, bFirst, x, coefficient);
1129 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1134 fprintf(debug, "Rank= %6d, pme local particles=%6d\n",
1135 cr->nodeid, atc->n);
1138 if (flags & GMX_PME_SPREAD)
1140 wallcycle_start(wcycle, ewcPME_SPREAD);
1142 /* Spread the coefficients on a grid */
1143 spread_on_grid(pme, &pme->atc[0], pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1147 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
1149 inc_nrnb(nrnb, eNR_SPREADBSP,
1150 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
1152 if (!pme->bUseThreads)
1154 wrap_periodic_pmegrid(pme, grid);
1156 /* sum contributions to local grid from other nodes */
1158 if (pme->nnodes > 1)
1160 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1165 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1168 wallcycle_stop(wcycle, ewcPME_SPREAD);
1170 /* TODO If the OpenMP and single-threaded implementations
1171 converge, then spread_on_grid() and
1172 copy_pmegrid_to_fftgrid() will perhaps live in the same
1177 /* Here we start a large thread parallel region */
1178 #pragma omp parallel num_threads(pme->nthread) private(thread)
1182 thread = gmx_omp_get_thread_num();
1183 if (flags & GMX_PME_SOLVE)
1190 wallcycle_start(wcycle, ewcPME_FFT);
1192 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
1196 wallcycle_stop(wcycle, ewcPME_FFT);
1200 /* solve in k-space for our local cells */
1203 wallcycle_start(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
1205 if (grid_index < DO_Q)
1208 solve_pme_yzx(pme, cfftgrid,
1209 scaledBox[XX][XX]*scaledBox[YY][YY]*scaledBox[ZZ][ZZ],
1211 pme->nthread, thread);
1216 solve_pme_lj_yzx(pme, &cfftgrid, FALSE,
1217 scaledBox[XX][XX]*scaledBox[YY][YY]*scaledBox[ZZ][ZZ],
1219 pme->nthread, thread);
1224 wallcycle_stop(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
1226 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1236 wallcycle_start(wcycle, ewcPME_FFT);
1238 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
1242 wallcycle_stop(wcycle, ewcPME_FFT);
1246 if (pme->nodeid == 0)
1248 real ntot = pme->nkx*pme->nky*pme->nkz;
1249 npme = static_cast<int>(ntot*std::log(ntot)/std::log(2.0));
1250 inc_nrnb(nrnb, eNR_FFT, 2*npme);
1253 /* Note: this wallcycle region is closed below
1254 outside an OpenMP region, so take care if
1255 refactoring code here. */
1256 wallcycle_start(wcycle, ewcPME_GATHER);
1259 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1261 } GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1263 /* End of thread parallel section.
1264 * With MPI we have to synchronize here before gmx_sum_qgrid_dd.
1269 /* distribute local grid to all nodes */
1271 if (pme->nnodes > 1)
1273 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1278 unwrap_periodic_pmegrid(pme, grid);
1283 /* interpolate forces for our local atoms */
1287 /* If we are running without parallelization,
1288 * atc->f is the actual force array, not a buffer,
1289 * therefore we should not clear it.
1291 lambda = grid_index < DO_Q ? lambda_q : lambda_lj;
1292 bClearF = (bFirst && PAR(cr));
1293 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1294 for (thread = 0; thread < pme->nthread; thread++)
1298 gather_f_bsplines(pme, grid, bClearF, atc,
1299 &atc->spline[thread],
1300 pme->bFEP ? (grid_index % 2 == 0 ? 1.0-lambda : lambda) : 1.0);
1302 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1307 inc_nrnb(nrnb, eNR_GATHERFBSP,
1308 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
1309 /* Note: this wallcycle region is opened above inside an OpenMP
1310 region, so take care if refactoring code here. */
1311 wallcycle_stop(wcycle, ewcPME_GATHER);
1316 /* This should only be called on the master thread
1317 * and after the threads have synchronized.
1321 get_pme_ener_vir_q(pme->solve_work, pme->nthread, &energy_AB[grid_index], vir_AB[grid_index]);
1325 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &energy_AB[grid_index], vir_AB[grid_index]);
1329 } /* of grid_index-loop */
1331 /* For Lorentz-Berthelot combination rules in LJ-PME, we need to calculate
1334 if (pme->doLJ && pme->ljpme_combination_rule == eljpmeLB)
1336 /* Loop over A- and B-state if we are doing FEP */
1337 for (fep_state = 0; fep_state < fep_states_lj; ++fep_state)
1339 real *local_c6 = nullptr, *local_sigma = nullptr, *RedistC6 = nullptr, *RedistSigma = nullptr;
1340 if (pme->nnodes == 1)
1342 if (pme->lb_buf1 == nullptr)
1344 pme->lb_buf_nalloc = pme->atc[0].n;
1345 snew(pme->lb_buf1, pme->lb_buf_nalloc);
1347 pme->atc[0].coefficient = pme->lb_buf1;
1352 local_sigma = sigmaA;
1356 local_sigma = sigmaB;
1359 gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1369 RedistSigma = sigmaA;
1373 RedistSigma = sigmaB;
1376 gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1378 wallcycle_start(wcycle, ewcPME_REDISTXF);
1380 do_redist_pos_coeffs(pme, cr, start, homenr, bFirst, x, RedistC6);
1381 if (pme->lb_buf_nalloc < atc->n)
1383 pme->lb_buf_nalloc = atc->nalloc;
1384 srenew(pme->lb_buf1, pme->lb_buf_nalloc);
1385 srenew(pme->lb_buf2, pme->lb_buf_nalloc);
1387 local_c6 = pme->lb_buf1;
1388 for (i = 0; i < atc->n; ++i)
1390 local_c6[i] = atc->coefficient[i];
1394 do_redist_pos_coeffs(pme, cr, start, homenr, FALSE, x, RedistSigma);
1395 local_sigma = pme->lb_buf2;
1396 for (i = 0; i < atc->n; ++i)
1398 local_sigma[i] = atc->coefficient[i];
1402 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1404 calc_initial_lb_coeffs(pme, local_c6, local_sigma);
1406 /*Seven terms in LJ-PME with LB, grid_index < 2 reserved for electrostatics*/
1407 for (grid_index = 2; grid_index < 9; ++grid_index)
1409 /* Unpack structure */
1410 pmegrid = &pme->pmegrid[grid_index];
1411 fftgrid = pme->fftgrid[grid_index];
1412 pfft_setup = pme->pfft_setup[grid_index];
1413 calc_next_lb_coeffs(pme, local_sigma);
1414 grid = pmegrid->grid.grid;
1417 if (flags & GMX_PME_SPREAD)
1419 wallcycle_start(wcycle, ewcPME_SPREAD);
1420 /* Spread the c6 on a grid */
1421 spread_on_grid(pme, &pme->atc[0], pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1425 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
1428 inc_nrnb(nrnb, eNR_SPREADBSP,
1429 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
1430 if (pme->nthread == 1)
1432 wrap_periodic_pmegrid(pme, grid);
1433 /* sum contributions to local grid from other nodes */
1435 if (pme->nnodes > 1)
1437 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1441 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1443 wallcycle_stop(wcycle, ewcPME_SPREAD);
1445 /*Here we start a large thread parallel region*/
1446 #pragma omp parallel num_threads(pme->nthread) private(thread)
1450 thread = gmx_omp_get_thread_num();
1451 if (flags & GMX_PME_SOLVE)
1456 wallcycle_start(wcycle, ewcPME_FFT);
1459 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
1463 wallcycle_stop(wcycle, ewcPME_FFT);
1468 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1472 if (flags & GMX_PME_SOLVE)
1474 /* solve in k-space for our local cells */
1475 #pragma omp parallel num_threads(pme->nthread) private(thread)
1480 thread = gmx_omp_get_thread_num();
1483 wallcycle_start(wcycle, ewcLJPME);
1487 solve_pme_lj_yzx(pme, &pme->cfftgrid[2], TRUE,
1488 scaledBox[XX][XX]*scaledBox[YY][YY]*scaledBox[ZZ][ZZ],
1490 pme->nthread, thread);
1493 wallcycle_stop(wcycle, ewcLJPME);
1495 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1498 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1504 /* This should only be called on the master thread and
1505 * after the threads have synchronized.
1507 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &energy_AB[2+fep_state], vir_AB[2+fep_state]);
1512 bFirst = !pme->doCoulomb;
1513 calc_initial_lb_coeffs(pme, local_c6, local_sigma);
1514 for (grid_index = 8; grid_index >= 2; --grid_index)
1516 /* Unpack structure */
1517 pmegrid = &pme->pmegrid[grid_index];
1518 fftgrid = pme->fftgrid[grid_index];
1519 pfft_setup = pme->pfft_setup[grid_index];
1520 grid = pmegrid->grid.grid;
1521 calc_next_lb_coeffs(pme, local_sigma);
1523 #pragma omp parallel num_threads(pme->nthread) private(thread)
1527 thread = gmx_omp_get_thread_num();
1532 wallcycle_start(wcycle, ewcPME_FFT);
1535 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
1539 wallcycle_stop(wcycle, ewcPME_FFT);
1543 if (pme->nodeid == 0)
1545 real ntot = pme->nkx*pme->nky*pme->nkz;
1546 npme = static_cast<int>(ntot*std::log(ntot)/std::log(2.0));
1547 inc_nrnb(nrnb, eNR_FFT, 2*npme);
1549 wallcycle_start(wcycle, ewcPME_GATHER);
1552 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1554 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1555 } /*#pragma omp parallel*/
1557 /* distribute local grid to all nodes */
1559 if (pme->nnodes > 1)
1561 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1566 unwrap_periodic_pmegrid(pme, grid);
1570 /* interpolate forces for our local atoms */
1572 bClearF = (bFirst && PAR(cr));
1573 scale = pme->bFEP ? (fep_state < 1 ? 1.0-lambda_lj : lambda_lj) : 1.0;
1574 scale *= lb_scale_factor[grid_index-2];
1576 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1577 for (thread = 0; thread < pme->nthread; thread++)
1581 gather_f_bsplines(pme, grid, bClearF, &pme->atc[0],
1582 &pme->atc[0].spline[thread],
1585 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1590 inc_nrnb(nrnb, eNR_GATHERFBSP,
1591 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
1593 wallcycle_stop(wcycle, ewcPME_GATHER);
1596 } /* for (grid_index = 8; grid_index >= 2; --grid_index) */
1598 } /* for (fep_state = 0; fep_state < fep_states_lj; ++fep_state) */
1599 } /* if ((flags & GMX_PME_DO_LJ) && pme->ljpme_combination_rule == eljpmeLB) */
1601 if (bCalcF && pme->nnodes > 1)
1603 wallcycle_start(wcycle, ewcPME_REDISTXF);
1604 for (d = 0; d < pme->ndecompdim; d++)
1607 if (d == pme->ndecompdim - 1)
1614 n_d = pme->atc[d+1].n;
1615 f_d = pme->atc[d+1].f;
1617 if (DOMAINDECOMP(cr))
1619 dd_pmeredist_f(pme, atc, n_d, f_d,
1620 d == pme->ndecompdim-1 && pme->bPPnode);
1624 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1634 *energy_q = energy_AB[0];
1635 m_add(vir_q, vir_AB[0], vir_q);
1639 *energy_q = (1.0-lambda_q)*energy_AB[0] + lambda_q*energy_AB[1];
1640 *dvdlambda_q += energy_AB[1] - energy_AB[0];
1641 for (i = 0; i < DIM; i++)
1643 for (j = 0; j < DIM; j++)
1645 vir_q[i][j] += (1.0-lambda_q)*vir_AB[0][i][j] +
1646 lambda_q*vir_AB[1][i][j];
1652 fprintf(debug, "Electrostatic PME mesh energy: %g\n", *energy_q);
1664 *energy_lj = energy_AB[2];
1665 m_add(vir_lj, vir_AB[2], vir_lj);
1669 *energy_lj = (1.0-lambda_lj)*energy_AB[2] + lambda_lj*energy_AB[3];
1670 *dvdlambda_lj += energy_AB[3] - energy_AB[2];
1671 for (i = 0; i < DIM; i++)
1673 for (j = 0; j < DIM; j++)
1675 vir_lj[i][j] += (1.0-lambda_lj)*vir_AB[2][i][j] + lambda_lj*vir_AB[3][i][j];
1681 fprintf(debug, "Lennard-Jones PME mesh energy: %g\n", *energy_lj);
1692 void gmx_pme_destroy(gmx_pme_t *pme)
1699 delete pme->boxScaler;
1708 for (int i = 0; i < pme->ngrids; ++i)
1710 pmegrids_destroy(&pme->pmegrid[i]);
1712 if (pme->pfft_setup)
1714 for (int i = 0; i < pme->ngrids; ++i)
1716 gmx_parallel_3dfft_destroy(pme->pfft_setup[i]);
1719 sfree(pme->fftgrid);
1720 sfree(pme->cfftgrid);
1721 sfree(pme->pfft_setup);
1723 for (int i = 0; i < std::max(1, pme->ndecompdim); i++) //pme->atc[0] is always allocated
1725 destroy_atomcomm(&pme->atc[i]);
1728 for (int i = 0; i < DIM; i++)
1730 sfree(pme->bsp_mod[i]);
1733 destroy_overlap_comm(&pme->overlap[0]);
1734 destroy_overlap_comm(&pme->overlap[1]);
1736 sfree(pme->lb_buf1);
1737 sfree(pme->lb_buf2);
1742 if (pme->solve_work)
1744 pme_free_all_work(&pme->solve_work, pme->nthread);
1747 sfree(pme->sum_qgrid_tmp);
1748 sfree(pme->sum_qgrid_dd_tmp);
1750 if (pme_gpu_active(pme) && pme->gpu)
1752 pme_gpu_destroy(pme->gpu);
1758 void gmx_pme_reinit_atoms(const gmx_pme_t *pme, const int nAtoms, const real *charges)
1760 if (pme_gpu_active(pme))
1762 pme_gpu_reinit_atoms(pme->gpu, nAtoms, charges);
1764 // TODO: handle the CPU case here; handle the whole t_mdatoms