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 The GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 * \brief This file contains function definitions necessary for
41 * computing energies and forces for the PME long-ranged part (Coulomb
44 * \author Erik Lindahl <erik@kth.se>
45 * \author Berk Hess <hess@kth.se>
46 * \ingroup module_ewald
48 /* IMPORTANT FOR DEVELOPERS:
50 * Triclinic pme stuff isn't entirely trivial, and we've experienced
51 * some bugs during development (many of them due to me). To avoid
52 * this in the future, please check the following things if you make
53 * changes in this file:
55 * 1. You should obtain identical (at least to the PME precision)
56 * energies, forces, and virial for
57 * a rectangular box and a triclinic one where the z (or y) axis is
58 * tilted a whole box side. For instance you could use these boxes:
60 * rectangular triclinic
65 * 2. You should check the energy conservation in a triclinic box.
67 * It might seem an overkill, but better safe than sorry.
86 #include "gromacs/domdec/domdec.h"
87 #include "gromacs/ewald/ewald_utils.h"
88 #include "gromacs/fft/parallel_3dfft.h"
89 #include "gromacs/fileio/pdbio.h"
90 #include "gromacs/gmxlib/network.h"
91 #include "gromacs/gmxlib/nrnb.h"
92 #include "gromacs/hardware/hw_info.h"
93 #include "gromacs/math/gmxcomplex.h"
94 #include "gromacs/math/invertmatrix.h"
95 #include "gromacs/math/units.h"
96 #include "gromacs/math/vec.h"
97 #include "gromacs/math/vectypes.h"
98 #include "gromacs/mdtypes/commrec.h"
99 #include "gromacs/mdtypes/forcerec.h"
100 #include "gromacs/mdtypes/inputrec.h"
101 #include "gromacs/mdtypes/md_enums.h"
102 #include "gromacs/pbcutil/pbc.h"
103 #include "gromacs/timing/cyclecounter.h"
104 #include "gromacs/timing/wallcycle.h"
105 #include "gromacs/timing/walltime_accounting.h"
106 #include "gromacs/topology/topology.h"
107 #include "gromacs/utility/basedefinitions.h"
108 #include "gromacs/utility/cstringutil.h"
109 #include "gromacs/utility/exceptions.h"
110 #include "gromacs/utility/fatalerror.h"
111 #include "gromacs/utility/gmxmpi.h"
112 #include "gromacs/utility/gmxomp.h"
113 #include "gromacs/utility/logger.h"
114 #include "gromacs/utility/real.h"
115 #include "gromacs/utility/smalloc.h"
116 #include "gromacs/utility/stringutil.h"
117 #include "gromacs/utility/unique_cptr.h"
119 #include "calculate_spline_moduli.h"
120 #include "pme_gather.h"
121 #include "pme_gpu_internal.h"
122 #include "pme_grid.h"
123 #include "pme_internal.h"
124 #include "pme_redistribute.h"
125 #include "pme_solve.h"
126 #include "pme_spline_work.h"
127 #include "pme_spread.h"
129 /*! \brief Help build a descriptive message in \c error if there are
130 * \c errorReasons why PME on GPU is not supported.
132 * \returns Whether the lack of errorReasons indicate there is support. */
133 static bool addMessageIfNotSupported(const std::list<std::string>& errorReasons, std::string* error)
135 bool isSupported = errorReasons.empty();
136 if (!isSupported && error)
138 std::string regressionTestMarker = "PME GPU does not support";
139 // this prefix is tested for in the regression tests script gmxtest.pl
140 *error = regressionTestMarker;
141 if (errorReasons.size() == 1)
143 *error += " " + errorReasons.back();
147 *error += ": " + gmx::joinStrings(errorReasons, "; ");
154 bool pme_gpu_supports_build(std::string* error)
156 std::list<std::string> errorReasons;
159 errorReasons.emplace_back("a double-precision build");
161 if (GMX_GPU == GMX_GPU_NONE)
163 errorReasons.emplace_back("a non-GPU build");
165 return addMessageIfNotSupported(errorReasons, error);
168 bool pme_gpu_supports_hardware(const gmx_hw_info_t gmx_unused& hwinfo, std::string* error)
170 std::list<std::string> errorReasons;
172 if (GMX_GPU == GMX_GPU_OPENCL)
175 errorReasons.emplace_back("Apple OS X operating system");
178 return addMessageIfNotSupported(errorReasons, error);
181 bool pme_gpu_supports_input(const t_inputrec& ir, const gmx_mtop_t& mtop, std::string* error)
183 std::list<std::string> errorReasons;
184 if (!EEL_PME(ir.coulombtype))
186 errorReasons.emplace_back("systems that do not use PME for electrostatics");
188 if (ir.pme_order != 4)
190 errorReasons.emplace_back("interpolation orders other than 4");
192 if (ir.efep != efepNO)
194 if (gmx_mtop_has_perturbed_charges(mtop))
196 errorReasons.emplace_back(
197 "free energy calculations with perturbed charges (multiple grids)");
200 if (EVDW_PME(ir.vdwtype))
202 errorReasons.emplace_back("Lennard-Jones PME");
204 if (!EI_DYNAMICS(ir.eI))
206 errorReasons.emplace_back("not a dynamical integrator");
208 return addMessageIfNotSupported(errorReasons, error);
211 /*! \brief \libinternal
212 * Finds out if PME with given inputs is possible to run on GPU.
213 * This function is an internal final check, validating the whole PME structure on creation,
214 * but it still duplicates the preliminary checks from the above (externally exposed) pme_gpu_supports_input() - just in case.
216 * \param[in] pme The PME structure.
217 * \param[out] error The error message if the input is not supported on GPU.
218 * \returns True if this PME input is possible to run on GPU, false otherwise.
220 static bool pme_gpu_check_restrictions(const gmx_pme_t* pme, std::string* error)
222 std::list<std::string> errorReasons;
223 if (pme->nnodes != 1)
225 errorReasons.emplace_back("PME decomposition");
227 if (pme->pme_order != 4)
229 errorReasons.emplace_back("interpolation orders other than 4");
233 errorReasons.emplace_back("free energy calculations (multiple grids)");
237 errorReasons.emplace_back("Lennard-Jones PME");
241 errorReasons.emplace_back("double precision");
243 if (GMX_GPU == GMX_GPU_NONE)
245 errorReasons.emplace_back("non-GPU build of GROMACS");
248 return addMessageIfNotSupported(errorReasons, error);
251 PmeRunMode pme_run_mode(const gmx_pme_t* pme)
253 GMX_ASSERT(pme != nullptr, "Expecting valid PME data pointer");
257 gmx::PinningPolicy pme_get_pinning_policy()
259 return gmx::PinningPolicy::PinnedIfSupported;
262 /*! \brief Number of bytes in a cache line.
264 * Must also be a multiple of the SIMD and SIMD4 register size, to
265 * preserve alignment.
267 const int gmxCacheLineSize = 64;
269 //! Set up coordinate communication
270 static void setup_coordinate_communication(PmeAtomComm* atc)
278 for (i = 1; i <= nslab / 2; i++)
280 fw = (atc->nodeid + i) % nslab;
281 bw = (atc->nodeid - i + nslab) % nslab;
284 atc->slabCommSetup[n].node_dest = fw;
285 atc->slabCommSetup[n].node_src = bw;
290 atc->slabCommSetup[n].node_dest = bw;
291 atc->slabCommSetup[n].node_src = fw;
297 /*! \brief Round \p n up to the next multiple of \p f */
298 static int mult_up(int n, int f)
300 return ((n + f - 1) / f) * f;
303 /*! \brief Return estimate of the load imbalance from the PME grid not being a good match for the number of PME ranks */
304 static double estimate_pme_load_imbalance(struct gmx_pme_t* pme)
309 nma = pme->nnodes_major;
310 nmi = pme->nnodes_minor;
312 n1 = mult_up(pme->nkx, nma) * mult_up(pme->nky, nmi) * pme->nkz;
313 n2 = mult_up(pme->nkx, nma) * mult_up(pme->nkz, nmi) * pme->nky;
314 n3 = mult_up(pme->nky, nma) * mult_up(pme->nkz, nmi) * pme->nkx;
316 /* pme_solve is roughly double the cost of an fft */
318 return (n1 + n2 + 3 * n3) / static_cast<double>(6 * pme->nkx * pme->nky * pme->nkz);
323 PmeAtomComm::PmeAtomComm(MPI_Comm PmeMpiCommunicator,
324 const int numThreads,
327 const bool doSpread) :
334 if (PmeMpiCommunicator != MPI_COMM_NULL)
336 mpi_comm = PmeMpiCommunicator;
338 MPI_Comm_size(mpi_comm, &nslab);
339 MPI_Comm_rank(mpi_comm, &nodeid);
344 fprintf(debug, "For PME atom communication in dimind %d: nslab %d rank %d\n", dimind, nslab, nodeid);
349 slabCommSetup.resize(nslab);
350 setup_coordinate_communication(this);
352 count_thread.resize(nthread);
353 for (auto& countThread : count_thread)
355 countThread.resize(nslab);
361 threadMap.resize(nthread);
363 # pragma omp parallel for num_threads(nthread) schedule(static)
364 for (int thread = 0; thread < nthread; thread++)
368 /* Allocate buffer with padding to avoid cache polution */
369 threadMap[thread].nBuffer.resize(nthread + 2 * gmxCacheLineSize);
370 threadMap[thread].n = threadMap[thread].nBuffer.data() + gmxCacheLineSize;
372 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
379 /*! \brief Initialize data structure for communication */
380 static void init_overlap_comm(pme_overlap_t* ol, int norder, MPI_Comm comm, int nnodes, int nodeid, int ndata, int commplainsize)
388 /* Linear translation of the PME grid won't affect reciprocal space
389 * calculations, so to optimize we only interpolate "upwards",
390 * which also means we only have to consider overlap in one direction.
391 * I.e., particles on this node might also be spread to grid indices
392 * that belong to higher nodes (modulo nnodes)
395 ol->s2g0.resize(ol->nnodes + 1);
396 ol->s2g1.resize(ol->nnodes);
399 fprintf(debug, "PME slab boundaries:");
401 for (int i = 0; i < nnodes; i++)
403 /* s2g0 the local interpolation grid start.
404 * s2g1 the local interpolation grid end.
405 * Since in calc_pidx we divide particles, and not grid lines,
406 * spatially uniform along dimension x or y, we need to round
407 * s2g0 down and s2g1 up.
409 ol->s2g0[i] = (i * ndata + 0) / nnodes;
410 ol->s2g1[i] = ((i + 1) * ndata + nnodes - 1) / nnodes + norder - 1;
414 fprintf(debug, " %3d %3d", ol->s2g0[i], ol->s2g1[i]);
417 ol->s2g0[nnodes] = ndata;
420 fprintf(debug, "\n");
423 /* Determine with how many nodes we need to communicate the grid overlap */
424 int testRankCount = 0;
429 for (int i = 0; i < nnodes; i++)
431 if ((i + testRankCount < nnodes && ol->s2g1[i] > ol->s2g0[i + testRankCount])
432 || (i + testRankCount >= nnodes && ol->s2g1[i] > ol->s2g0[i + testRankCount - nnodes] + ndata))
437 } while (bCont && testRankCount < nnodes);
439 ol->comm_data.resize(testRankCount - 1);
442 for (size_t b = 0; b < ol->comm_data.size(); b++)
444 pme_grid_comm_t* pgc = &ol->comm_data[b];
447 pgc->send_id = (ol->nodeid + (b + 1)) % ol->nnodes;
448 int fft_start = ol->s2g0[pgc->send_id];
449 int fft_end = ol->s2g0[pgc->send_id + 1];
450 if (pgc->send_id < nodeid)
455 int send_index1 = ol->s2g1[nodeid];
456 send_index1 = std::min(send_index1, fft_end);
457 pgc->send_index0 = fft_start;
458 pgc->send_nindex = std::max(0, send_index1 - pgc->send_index0);
459 ol->send_size += pgc->send_nindex;
461 /* We always start receiving to the first index of our slab */
462 pgc->recv_id = (ol->nodeid - (b + 1) + ol->nnodes) % ol->nnodes;
463 fft_start = ol->s2g0[ol->nodeid];
464 fft_end = ol->s2g0[ol->nodeid + 1];
465 int recv_index1 = ol->s2g1[pgc->recv_id];
466 if (pgc->recv_id > nodeid)
468 recv_index1 -= ndata;
470 recv_index1 = std::min(recv_index1, fft_end);
471 pgc->recv_index0 = fft_start;
472 pgc->recv_nindex = std::max(0, recv_index1 - pgc->recv_index0);
476 /* Communicate the buffer sizes to receive */
478 for (size_t b = 0; b < ol->comm_data.size(); b++)
480 MPI_Sendrecv(&ol->send_size, 1, MPI_INT, ol->comm_data[b].send_id, b, &ol->comm_data[b].recv_size,
481 1, MPI_INT, ol->comm_data[b].recv_id, b, ol->mpi_comm, &stat);
485 /* For non-divisible grid we need pme_order iso pme_order-1 */
486 ol->sendbuf.resize(norder * commplainsize);
487 ol->recvbuf.resize(norder * commplainsize);
490 int minimalPmeGridSize(int pmeOrder)
492 /* The actual grid size limitations are:
493 * serial: >= pme_order
494 * DD, no OpenMP: >= 2*(pme_order - 1)
495 * DD, OpenMP: >= pme_order + 1
496 * But we use the maximum for simplicity since in practice there is not
497 * much performance difference between pme_order and 2*(pme_order -1).
499 int minimalSize = 2 * (pmeOrder - 1);
501 GMX_RELEASE_ASSERT(pmeOrder >= 3, "pmeOrder has to be >= 3");
502 GMX_RELEASE_ASSERT(minimalSize >= pmeOrder + 1, "The grid size should be >= pmeOrder + 1");
507 bool gmx_pme_check_restrictions(int pme_order, int nkx, int nky, int nkz, int numPmeDomainsAlongX, bool useThreads, bool errorsAreFatal)
509 if (pme_order > PME_ORDER_MAX)
516 std::string message = gmx::formatString(
517 "pme_order (%d) is larger than the maximum allowed value (%d). Modify and "
518 "recompile the code if you really need such a high order.",
519 pme_order, PME_ORDER_MAX);
520 GMX_THROW(gmx::InconsistentInputError(message));
523 const int minGridSize = minimalPmeGridSize(pme_order);
524 if (nkx < minGridSize || nky < minGridSize || nkz < minGridSize)
530 std::string message =
531 gmx::formatString("The PME grid sizes need to be >= 2*(pme_order-1) (%d)", minGridSize);
532 GMX_THROW(gmx::InconsistentInputError(message));
535 /* Check for a limitation of the (current) sum_fftgrid_dd code.
536 * We only allow multiple communication pulses in dim 1, not in dim 0.
539 && (nkx < numPmeDomainsAlongX * pme_order && nkx != numPmeDomainsAlongX * (pme_order - 1)))
546 "The number of PME grid lines per rank along x is %g. But when using OpenMP "
547 "threads, the number of grid lines per rank along x should be >= pme_order (%d) "
548 "or = pmeorder-1. To resolve this issue, use fewer ranks along x (and possibly "
549 "more along y and/or z) by specifying -dd manually.",
550 nkx / static_cast<double>(numPmeDomainsAlongX), pme_order);
556 /*! \brief Round \p enumerator */
557 static int div_round_up(int enumerator, int denominator)
559 return (enumerator + denominator - 1) / denominator;
562 gmx_pme_t* gmx_pme_init(const t_commrec* cr,
563 const NumPmeDomains& numPmeDomains,
564 const t_inputrec* ir,
565 gmx_bool bFreeEnergy_q,
566 gmx_bool bFreeEnergy_lj,
567 gmx_bool bReproducible,
573 const gmx_device_info_t* gpuInfo,
574 PmeGpuProgramHandle pmeGpuProgram,
575 const gmx::MDLogger& mdlog)
577 int use_threads, sum_use_threads, i;
582 fprintf(debug, "Creating PME data structures.\n");
585 gmx::unique_cptr<gmx_pme_t, gmx_pme_destroy> pme(new gmx_pme_t());
587 pme->sum_qgrid_tmp = nullptr;
588 pme->sum_qgrid_dd_tmp = nullptr;
595 pme->nnodes_major = numPmeDomains.x;
596 pme->nnodes_minor = numPmeDomains.y;
598 if (numPmeDomains.x * numPmeDomains.y > 1)
600 pme->mpi_comm = cr->mpi_comm_mygroup;
603 MPI_Comm_rank(pme->mpi_comm, &pme->nodeid);
604 MPI_Comm_size(pme->mpi_comm, &pme->nnodes);
606 if (pme->nnodes != numPmeDomains.x * numPmeDomains.y)
608 gmx_incons("PME rank count mismatch");
613 pme->mpi_comm = MPI_COMM_NULL;
616 if (pme->nnodes == 1)
618 pme->mpi_comm_d[0] = MPI_COMM_NULL;
619 pme->mpi_comm_d[1] = MPI_COMM_NULL;
621 pme->nodeid_major = 0;
622 pme->nodeid_minor = 0;
626 if (numPmeDomains.y == 1)
628 pme->mpi_comm_d[0] = pme->mpi_comm;
629 pme->mpi_comm_d[1] = MPI_COMM_NULL;
631 pme->nodeid_major = pme->nodeid;
632 pme->nodeid_minor = 0;
634 else if (numPmeDomains.x == 1)
636 pme->mpi_comm_d[0] = MPI_COMM_NULL;
637 pme->mpi_comm_d[1] = pme->mpi_comm;
639 pme->nodeid_major = 0;
640 pme->nodeid_minor = pme->nodeid;
644 if (pme->nnodes % numPmeDomains.x != 0)
647 "For 2D PME decomposition, #PME ranks must be divisible by the number of "
653 MPI_Comm_split(pme->mpi_comm, pme->nodeid % numPmeDomains.y, pme->nodeid,
654 &pme->mpi_comm_d[0]); /* My communicator along major dimension */
655 MPI_Comm_split(pme->mpi_comm, pme->nodeid / numPmeDomains.y, pme->nodeid,
656 &pme->mpi_comm_d[1]); /* My communicator along minor dimension */
658 MPI_Comm_rank(pme->mpi_comm_d[0], &pme->nodeid_major);
659 MPI_Comm_size(pme->mpi_comm_d[0], &pme->nnodes_major);
660 MPI_Comm_rank(pme->mpi_comm_d[1], &pme->nodeid_minor);
661 MPI_Comm_size(pme->mpi_comm_d[1], &pme->nnodes_minor);
665 // cr is always initialized if there is a a PP rank, so we can safely assume
666 // that when it is not, like in ewald tests, we not on a PP rank.
667 pme->bPPnode = ((cr != nullptr && cr->duty != 0) && thisRankHasDuty(cr, DUTY_PP));
669 pme->nthread = nthread;
671 /* Check if any of the PME MPI ranks uses threads */
672 use_threads = (pme->nthread > 1 ? 1 : 0);
676 MPI_Allreduce(&use_threads, &sum_use_threads, 1, MPI_INT, MPI_SUM, pme->mpi_comm);
681 sum_use_threads = use_threads;
683 pme->bUseThreads = (sum_use_threads > 0);
685 if (ir->ePBC == epbcSCREW)
687 gmx_fatal(FARGS, "pme does not (yet) work with pbc = screw");
691 * It is likely that the current gmx_pme_do() routine supports calculating
692 * only Coulomb or LJ while gmx_pme_init() configures for both,
693 * but that has never been tested.
694 * It is likely that the current gmx_pme_do() routine supports calculating,
695 * not calculating free-energy for Coulomb and/or LJ while gmx_pme_init()
696 * configures with free-energy, but that has never been tested.
698 pme->doCoulomb = EEL_PME(ir->coulombtype);
699 pme->doLJ = EVDW_PME(ir->vdwtype);
700 pme->bFEP_q = ((ir->efep != efepNO) && bFreeEnergy_q);
701 pme->bFEP_lj = ((ir->efep != efepNO) && bFreeEnergy_lj);
702 pme->bFEP = (pme->bFEP_q || pme->bFEP_lj);
706 pme->bP3M = (ir->coulombtype == eelP3M_AD || getenv("GMX_PME_P3M") != nullptr);
707 pme->pme_order = ir->pme_order;
708 pme->ewaldcoeff_q = ewaldcoeff_q;
709 pme->ewaldcoeff_lj = ewaldcoeff_lj;
711 /* Always constant electrostatics coefficients */
712 pme->epsilon_r = ir->epsilon_r;
714 /* Always constant LJ coefficients */
715 pme->ljpme_combination_rule = ir->ljpme_combination_rule;
717 // The box requires scaling with nwalls = 2, we store that condition as well
718 // as the scaling factor
719 delete pme->boxScaler;
720 pme->boxScaler = new EwaldBoxZScaler(*ir);
722 /* If we violate restrictions, generate a fatal error here */
723 gmx_pme_check_restrictions(pme->pme_order, pme->nkx, pme->nky, pme->nkz, pme->nnodes_major,
724 pme->bUseThreads, true);
731 MPI_Type_contiguous(DIM, GMX_MPI_REAL, &(pme->rvec_mpi));
732 MPI_Type_commit(&(pme->rvec_mpi));
735 /* Note that the coefficient spreading and force gathering, which usually
736 * takes about the same amount of time as FFT+solve_pme,
737 * is always fully load balanced
738 * (unless the coefficient distribution is inhomogeneous).
741 imbal = estimate_pme_load_imbalance(pme.get());
742 if (imbal >= 1.2 && pme->nodeid_major == 0 && pme->nodeid_minor == 0)
744 GMX_LOG(mdlog.warning)
746 .appendTextFormatted(
747 "NOTE: The load imbalance in PME FFT and solve is %d%%.\n"
748 " For optimal PME load balancing\n"
749 " PME grid_x (%d) and grid_y (%d) should be divisible by "
752 " and PME grid_y (%d) and grid_z (%d) should be divisible by "
755 gmx::roundToInt((imbal - 1) * 100), pme->nkx, pme->nky,
756 pme->nnodes_major, pme->nky, pme->nkz, pme->nnodes_minor);
760 /* For non-divisible grid we need pme_order iso pme_order-1 */
761 /* In sum_qgrid_dd x overlap is copied in place: take padding into account.
762 * y is always copied through a buffer: we don't need padding in z,
763 * but we do need the overlap in x because of the communication order.
765 init_overlap_comm(&pme->overlap[0], pme->pme_order, pme->mpi_comm_d[0], pme->nnodes_major,
766 pme->nodeid_major, pme->nkx,
767 (div_round_up(pme->nky, pme->nnodes_minor) + pme->pme_order)
768 * (pme->nkz + pme->pme_order - 1));
770 /* Along overlap dim 1 we can send in multiple pulses in sum_fftgrid_dd.
771 * We do this with an offset buffer of equal size, so we need to allocate
772 * extra for the offset. That's what the (+1)*pme->nkz is for.
774 init_overlap_comm(&pme->overlap[1], pme->pme_order, pme->mpi_comm_d[1], pme->nnodes_minor,
775 pme->nodeid_minor, pme->nky,
776 (div_round_up(pme->nkx, pme->nnodes_major) + pme->pme_order + 1) * pme->nkz);
778 /* Double-check for a limitation of the (current) sum_fftgrid_dd code.
779 * Note that gmx_pme_check_restrictions checked for this already.
781 if (pme->bUseThreads && (pme->overlap[0].comm_data.size() > 1))
784 "More than one communication pulse required for grid overlap communication along "
785 "the major dimension while using threads");
788 snew(pme->bsp_mod[XX], pme->nkx);
789 snew(pme->bsp_mod[YY], pme->nky);
790 snew(pme->bsp_mod[ZZ], pme->nkz);
792 pme->gpu = pmeGpu; /* Carrying over the single GPU structure */
793 pme->runMode = runMode;
795 /* The required size of the interpolation grid, including overlap.
796 * The allocated size (pmegrid_n?) might be slightly larger.
798 pme->pmegrid_nx = pme->overlap[0].s2g1[pme->nodeid_major] - pme->overlap[0].s2g0[pme->nodeid_major];
799 pme->pmegrid_ny = pme->overlap[1].s2g1[pme->nodeid_minor] - pme->overlap[1].s2g0[pme->nodeid_minor];
800 pme->pmegrid_nz_base = pme->nkz;
801 pme->pmegrid_nz = pme->pmegrid_nz_base + pme->pme_order - 1;
802 set_grid_alignment(&pme->pmegrid_nz, pme->pme_order);
803 pme->pmegrid_start_ix = pme->overlap[0].s2g0[pme->nodeid_major];
804 pme->pmegrid_start_iy = pme->overlap[1].s2g0[pme->nodeid_minor];
805 pme->pmegrid_start_iz = 0;
807 make_gridindex_to_localindex(pme->nkx, pme->pmegrid_start_ix,
808 pme->pmegrid_nx - (pme->pme_order - 1), &pme->nnx, &pme->fshx);
809 make_gridindex_to_localindex(pme->nky, pme->pmegrid_start_iy,
810 pme->pmegrid_ny - (pme->pme_order - 1), &pme->nny, &pme->fshy);
811 make_gridindex_to_localindex(pme->nkz, pme->pmegrid_start_iz, pme->pmegrid_nz_base, &pme->nnz,
814 pme->spline_work = make_pme_spline_work(pme->pme_order);
819 /* It doesn't matter if we allocate too many grids here,
820 * we only allocate and use the ones we need.
824 pme->ngrids = ((ir->ljpme_combination_rule == eljpmeLB) ? DO_Q_AND_LJ_LB : DO_Q_AND_LJ);
830 snew(pme->fftgrid, pme->ngrids);
831 snew(pme->cfftgrid, pme->ngrids);
832 snew(pme->pfft_setup, pme->ngrids);
834 for (i = 0; i < pme->ngrids; ++i)
836 if ((i < DO_Q && pme->doCoulomb && (i == 0 || bFreeEnergy_q))
837 || (i >= DO_Q && pme->doLJ
838 && (i == 2 || bFreeEnergy_lj || ir->ljpme_combination_rule == eljpmeLB)))
840 pmegrids_init(&pme->pmegrid[i], pme->pmegrid_nx, pme->pmegrid_ny, pme->pmegrid_nz,
841 pme->pmegrid_nz_base, pme->pme_order, pme->bUseThreads, pme->nthread,
842 pme->overlap[0].s2g1[pme->nodeid_major]
843 - pme->overlap[0].s2g0[pme->nodeid_major + 1],
844 pme->overlap[1].s2g1[pme->nodeid_minor]
845 - pme->overlap[1].s2g0[pme->nodeid_minor + 1]);
846 /* This routine will allocate the grid data to fit the FFTs */
847 const auto allocateRealGridForGpu = (pme->runMode == PmeRunMode::Mixed)
848 ? gmx::PinningPolicy::PinnedIfSupported
849 : gmx::PinningPolicy::CannotBePinned;
850 gmx_parallel_3dfft_init(&pme->pfft_setup[i], ndata, &pme->fftgrid[i], &pme->cfftgrid[i],
851 pme->mpi_comm_d, bReproducible, pme->nthread, allocateRealGridForGpu);
857 /* Use plain SPME B-spline interpolation */
858 make_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
862 /* Use the P3M grid-optimized influence function */
863 make_p3m_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
866 /* Use atc[0] for spreading */
867 const int firstDimIndex = (numPmeDomains.x > 1 ? 0 : 1);
868 MPI_Comm mpiCommFirstDim = (pme->nnodes > 1 ? pme->mpi_comm_d[firstDimIndex] : MPI_COMM_NULL);
869 bool doSpread = true;
870 pme->atc.emplace_back(mpiCommFirstDim, pme->nthread, pme->pme_order, firstDimIndex, doSpread);
871 if (pme->ndecompdim >= 2)
873 const int secondDimIndex = 1;
875 pme->atc.emplace_back(pme->mpi_comm_d[1], pme->nthread, pme->pme_order, secondDimIndex, doSpread);
878 if (pme_gpu_active(pme.get()))
882 // Initial check of validity of the data
883 std::string errorString;
884 bool canRunOnGpu = pme_gpu_check_restrictions(pme.get(), &errorString);
887 GMX_THROW(gmx::NotImplementedError(errorString));
891 pme_gpu_reinit(pme.get(), gpuInfo, pmeGpuProgram);
894 pme_init_all_work(&pme->solve_work, pme->nthread, pme->nkx);
896 // no exception was thrown during the init, so we hand over the PME structure handle
897 return pme.release();
900 void gmx_pme_reinit(struct gmx_pme_t** pmedata,
902 struct gmx_pme_t* pme_src,
903 const t_inputrec* ir,
904 const ivec grid_size,
908 // Create a copy of t_inputrec fields that are used in gmx_pme_init().
909 // TODO: This would be better as just copying a sub-structure that contains
910 // all the PME parameters and nothing else.
913 irc.coulombtype = ir->coulombtype;
914 irc.vdwtype = ir->vdwtype;
916 irc.pme_order = ir->pme_order;
917 irc.epsilon_r = ir->epsilon_r;
918 irc.ljpme_combination_rule = ir->ljpme_combination_rule;
919 irc.nkx = grid_size[XX];
920 irc.nky = grid_size[YY];
921 irc.nkz = grid_size[ZZ];
925 // This is reinit. Any logging should have been done at first init.
926 // Here we should avoid writing notes for settings the user did not
928 const gmx::MDLogger dummyLogger;
929 GMX_ASSERT(pmedata, "Invalid PME pointer");
930 NumPmeDomains numPmeDomains = { pme_src->nnodes_major, pme_src->nnodes_minor };
931 *pmedata = gmx_pme_init(cr, numPmeDomains, &irc, pme_src->bFEP_q, pme_src->bFEP_lj, FALSE,
932 ewaldcoeff_q, ewaldcoeff_lj, pme_src->nthread, pme_src->runMode,
933 pme_src->gpu, nullptr, nullptr, dummyLogger);
934 /* When running PME on the CPU not using domain decomposition,
935 * the atom data is allocated once only in gmx_pme_(re)init().
937 if (!pme_src->gpu && pme_src->nnodes == 1)
939 gmx_pme_reinit_atoms(*pmedata, pme_src->atc[0].numAtoms(), nullptr);
941 // TODO this is mostly passing around current values
943 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
945 /* We can easily reuse the allocated pme grids in pme_src */
946 reuse_pmegrids(&pme_src->pmegrid[PME_GRID_QA], &(*pmedata)->pmegrid[PME_GRID_QA]);
947 /* We would like to reuse the fft grids, but that's harder */
950 void gmx_pme_calc_energy(gmx_pme_t* pme, gmx::ArrayRef<const gmx::RVec> x, gmx::ArrayRef<const real> q, real* V)
956 gmx_incons("gmx_pme_calc_energy called in parallel");
960 gmx_incons("gmx_pme_calc_energy with free energy");
963 if (!pme->atc_energy)
965 pme->atc_energy = std::make_unique<PmeAtomComm>(MPI_COMM_NULL, 1, pme->pme_order, 0, true);
967 PmeAtomComm* atc = pme->atc_energy.get();
968 atc->setNumAtoms(x.ssize());
970 atc->coefficient = q;
972 /* We only use the A-charges grid */
973 grid = &pme->pmegrid[PME_GRID_QA];
975 /* Only calculate the spline coefficients, don't actually spread */
976 spread_on_grid(pme, atc, nullptr, TRUE, FALSE, pme->fftgrid[PME_GRID_QA], FALSE, PME_GRID_QA);
978 *V = gather_energy_bsplines(pme, grid->grid.grid, atc);
981 /*! \brief Calculate initial Lorentz-Berthelot coefficients for LJ-PME */
982 static void calc_initial_lb_coeffs(gmx::ArrayRef<real> coefficient, const real* local_c6, const real* local_sigma)
984 for (gmx::index i = 0; i < coefficient.ssize(); ++i)
986 real sigma4 = local_sigma[i];
987 sigma4 = sigma4 * sigma4;
988 sigma4 = sigma4 * sigma4;
989 coefficient[i] = local_c6[i] / sigma4;
993 /*! \brief Calculate next Lorentz-Berthelot coefficients for LJ-PME */
994 static void calc_next_lb_coeffs(gmx::ArrayRef<real> coefficient, const real* local_sigma)
996 for (gmx::index i = 0; i < coefficient.ssize(); ++i)
998 coefficient[i] *= local_sigma[i];
1002 int gmx_pme_do(struct gmx_pme_t* pme,
1003 gmx::ArrayRef<const gmx::RVec> coordinates,
1004 gmx::ArrayRef<gmx::RVec> forces,
1012 const t_commrec* cr,
1016 gmx_wallcycle* wcycle,
1027 GMX_ASSERT(pme->runMode == PmeRunMode::CPU,
1028 "gmx_pme_do should not be called on the GPU PME run.");
1030 int d, npme, grid_index, max_grid_index;
1031 PmeAtomComm& atc = pme->atc[0];
1032 pmegrids_t* pmegrid = nullptr;
1033 real* grid = nullptr;
1034 real* coefficient = nullptr;
1035 PmeOutput output[2]; // The second is used for the B state with FEP
1038 gmx_parallel_3dfft_t pfft_setup;
1040 t_complex* cfftgrid;
1042 gmx_bool bFirst, bDoSplines;
1044 int fep_states_lj = pme->bFEP_lj ? 2 : 1;
1045 const gmx_bool bCalcEnerVir = (flags & GMX_PME_CALC_ENER_VIR) != 0;
1046 const gmx_bool bBackFFT = (flags & (GMX_PME_CALC_F | GMX_PME_CALC_POT)) != 0;
1047 const gmx_bool bCalcF = (flags & GMX_PME_CALC_F) != 0;
1049 /* We could be passing lambda!=0 while no q or LJ is actually perturbed */
1059 assert(pme->nnodes > 0);
1060 assert(pme->nnodes == 1 || pme->ndecompdim > 0);
1062 if (pme->nnodes > 1)
1064 atc.pd.resize(coordinates.ssize());
1065 for (int d = pme->ndecompdim - 1; d >= 0; d--)
1067 PmeAtomComm& atc = pme->atc[d];
1068 atc.maxshift = (atc.dimind == 0 ? maxshift_x : maxshift_y);
1073 GMX_ASSERT(coordinates.ssize() == atc.numAtoms(), "We expect atc.numAtoms() coordinates");
1074 GMX_ASSERT(forces.ssize() >= atc.numAtoms(),
1075 "We need a force buffer with at least atc.numAtoms() elements");
1077 atc.x = coordinates;
1082 pme->boxScaler->scaleBox(box, scaledBox);
1084 gmx::invertBoxMatrix(scaledBox, pme->recipbox);
1087 /* For simplicity, we construct the splines for all particles if
1088 * more than one PME calculations is needed. Some optimization
1089 * could be done by keeping track of which atoms have splines
1090 * constructed, and construct new splines on each pass for atoms
1091 * that don't yet have them.
1094 bDoSplines = pme->bFEP || (pme->doCoulomb && pme->doLJ);
1096 /* We need a maximum of four separate PME calculations:
1097 * grid_index=0: Coulomb PME with charges from state A
1098 * grid_index=1: Coulomb PME with charges from state B
1099 * grid_index=2: LJ PME with C6 from state A
1100 * grid_index=3: LJ PME with C6 from state B
1101 * For Lorentz-Berthelot combination rules, a separate loop is used to
1102 * calculate all the terms
1105 /* If we are doing LJ-PME with LB, we only do Q here */
1106 max_grid_index = (pme->ljpme_combination_rule == eljpmeLB) ? DO_Q : DO_Q_AND_LJ;
1108 for (grid_index = 0; grid_index < max_grid_index; ++grid_index)
1110 /* Check if we should do calculations at this grid_index
1111 * If grid_index is odd we should be doing FEP
1112 * If grid_index < 2 we should be doing electrostatic PME
1113 * If grid_index >= 2 we should be doing LJ-PME
1115 if ((grid_index < DO_Q && (!pme->doCoulomb || (grid_index == 1 && !pme->bFEP_q)))
1116 || (grid_index >= DO_Q && (!pme->doLJ || (grid_index == 3 && !pme->bFEP_lj))))
1120 /* Unpack structure */
1121 pmegrid = &pme->pmegrid[grid_index];
1122 fftgrid = pme->fftgrid[grid_index];
1123 cfftgrid = pme->cfftgrid[grid_index];
1124 pfft_setup = pme->pfft_setup[grid_index];
1127 case 0: coefficient = chargeA; break;
1128 case 1: coefficient = chargeB; break;
1129 case 2: coefficient = c6A; break;
1130 case 3: coefficient = c6B; break;
1133 grid = pmegrid->grid.grid;
1137 fprintf(debug, "PME: number of ranks = %d, rank = %d\n", cr->nnodes, cr->nodeid);
1138 fprintf(debug, "Grid = %p\n", static_cast<void*>(grid));
1139 if (grid == nullptr)
1141 gmx_fatal(FARGS, "No grid!");
1145 if (pme->nnodes == 1)
1147 atc.coefficient = gmx::arrayRefFromArray(coefficient, coordinates.size());
1151 wallcycle_start(wcycle, ewcPME_REDISTXF);
1152 do_redist_pos_coeffs(pme, cr, bFirst, coordinates, coefficient);
1154 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1159 fprintf(debug, "Rank= %6d, pme local particles=%6d\n", cr->nodeid, atc.numAtoms());
1162 if (flags & GMX_PME_SPREAD)
1164 wallcycle_start(wcycle, ewcPME_SPREAD);
1166 /* Spread the coefficients on a grid */
1167 spread_on_grid(pme, &atc, pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1171 inc_nrnb(nrnb, eNR_WEIGHTS, DIM * atc.numAtoms());
1173 inc_nrnb(nrnb, eNR_SPREADBSP,
1174 pme->pme_order * pme->pme_order * pme->pme_order * atc.numAtoms());
1176 if (!pme->bUseThreads)
1178 wrap_periodic_pmegrid(pme, grid);
1180 /* sum contributions to local grid from other nodes */
1181 if (pme->nnodes > 1)
1183 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1186 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1189 wallcycle_stop(wcycle, ewcPME_SPREAD);
1191 /* TODO If the OpenMP and single-threaded implementations
1192 converge, then spread_on_grid() and
1193 copy_pmegrid_to_fftgrid() will perhaps live in the same
1198 /* Here we start a large thread parallel region */
1199 #pragma omp parallel num_threads(pme->nthread) private(thread)
1203 thread = gmx_omp_get_thread_num();
1204 if (flags & GMX_PME_SOLVE)
1211 wallcycle_start(wcycle, ewcPME_FFT);
1213 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX, thread, wcycle);
1216 wallcycle_stop(wcycle, ewcPME_FFT);
1219 /* solve in k-space for our local cells */
1222 wallcycle_start(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
1224 if (grid_index < DO_Q)
1226 loop_count = solve_pme_yzx(
1227 pme, cfftgrid, scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ],
1228 bCalcEnerVir, pme->nthread, thread);
1232 loop_count = solve_pme_lj_yzx(
1233 pme, &cfftgrid, FALSE,
1234 scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ],
1235 bCalcEnerVir, pme->nthread, thread);
1240 wallcycle_stop(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
1241 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1250 wallcycle_start(wcycle, ewcPME_FFT);
1252 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL, thread, wcycle);
1255 wallcycle_stop(wcycle, ewcPME_FFT);
1258 if (pme->nodeid == 0)
1260 real ntot = pme->nkx * pme->nky * pme->nkz;
1261 npme = static_cast<int>(ntot * std::log(ntot) / std::log(2.0));
1262 inc_nrnb(nrnb, eNR_FFT, 2 * npme);
1265 /* Note: this wallcycle region is closed below
1266 outside an OpenMP region, so take care if
1267 refactoring code here. */
1268 wallcycle_start(wcycle, ewcPME_GATHER);
1271 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1274 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1276 /* End of thread parallel section.
1277 * With MPI we have to synchronize here before gmx_sum_qgrid_dd.
1282 /* distribute local grid to all nodes */
1283 if (pme->nnodes > 1)
1285 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1288 unwrap_periodic_pmegrid(pme, grid);
1293 /* interpolate forces for our local atoms */
1296 /* If we are running without parallelization,
1297 * atc->f is the actual force array, not a buffer,
1298 * therefore we should not clear it.
1300 lambda = grid_index < DO_Q ? lambda_q : lambda_lj;
1301 bClearF = (bFirst && PAR(cr));
1302 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1303 for (thread = 0; thread < pme->nthread; thread++)
1307 gather_f_bsplines(pme, grid, bClearF, &atc, &atc.spline[thread],
1308 pme->bFEP ? (grid_index % 2 == 0 ? 1.0 - lambda : lambda) : 1.0);
1310 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1314 inc_nrnb(nrnb, eNR_GATHERFBSP,
1315 pme->pme_order * pme->pme_order * pme->pme_order * atc.numAtoms());
1316 /* Note: this wallcycle region is opened above inside an OpenMP
1317 region, so take care if refactoring code here. */
1318 wallcycle_stop(wcycle, ewcPME_GATHER);
1323 /* This should only be called on the master thread
1324 * and after the threads have synchronized.
1328 get_pme_ener_vir_q(pme->solve_work, pme->nthread, &output[grid_index % 2]);
1332 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output[grid_index % 2]);
1336 } /* of grid_index-loop */
1338 /* For Lorentz-Berthelot combination rules in LJ-PME, we need to calculate
1341 if (pme->doLJ && pme->ljpme_combination_rule == eljpmeLB)
1343 /* Loop over A- and B-state if we are doing FEP */
1344 for (fep_state = 0; fep_state < fep_states_lj; ++fep_state)
1346 real *local_c6 = nullptr, *local_sigma = nullptr, *RedistC6 = nullptr, *RedistSigma = nullptr;
1347 gmx::ArrayRef<real> coefficientBuffer;
1348 if (pme->nnodes == 1)
1350 pme->lb_buf1.resize(atc.numAtoms());
1351 coefficientBuffer = pme->lb_buf1;
1356 local_sigma = sigmaA;
1360 local_sigma = sigmaB;
1362 default: gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1367 coefficientBuffer = atc.coefficientBuffer;
1372 RedistSigma = sigmaA;
1376 RedistSigma = sigmaB;
1378 default: gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1380 wallcycle_start(wcycle, ewcPME_REDISTXF);
1382 do_redist_pos_coeffs(pme, cr, bFirst, coordinates, RedistC6);
1383 pme->lb_buf1.resize(atc.numAtoms());
1384 pme->lb_buf2.resize(atc.numAtoms());
1385 local_c6 = pme->lb_buf1.data();
1386 for (int i = 0; i < atc.numAtoms(); ++i)
1388 local_c6[i] = atc.coefficient[i];
1391 do_redist_pos_coeffs(pme, cr, FALSE, coordinates, RedistSigma);
1392 local_sigma = pme->lb_buf2.data();
1393 for (int i = 0; i < atc.numAtoms(); ++i)
1395 local_sigma[i] = atc.coefficient[i];
1398 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1400 atc.coefficient = coefficientBuffer;
1401 calc_initial_lb_coeffs(coefficientBuffer, local_c6, local_sigma);
1403 /*Seven terms in LJ-PME with LB, grid_index < 2 reserved for electrostatics*/
1404 for (grid_index = 2; grid_index < 9; ++grid_index)
1406 /* Unpack structure */
1407 pmegrid = &pme->pmegrid[grid_index];
1408 fftgrid = pme->fftgrid[grid_index];
1409 pfft_setup = pme->pfft_setup[grid_index];
1410 calc_next_lb_coeffs(coefficientBuffer, local_sigma);
1411 grid = pmegrid->grid.grid;
1413 if (flags & GMX_PME_SPREAD)
1415 wallcycle_start(wcycle, ewcPME_SPREAD);
1416 /* Spread the c6 on a grid */
1417 spread_on_grid(pme, &atc, pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1421 inc_nrnb(nrnb, eNR_WEIGHTS, DIM * atc.numAtoms());
1424 inc_nrnb(nrnb, eNR_SPREADBSP,
1425 pme->pme_order * pme->pme_order * pme->pme_order * atc.numAtoms());
1426 if (pme->nthread == 1)
1428 wrap_periodic_pmegrid(pme, grid);
1429 /* sum contributions to local grid from other nodes */
1430 if (pme->nnodes > 1)
1432 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1434 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1436 wallcycle_stop(wcycle, ewcPME_SPREAD);
1438 /*Here we start a large thread parallel region*/
1439 #pragma omp parallel num_threads(pme->nthread) private(thread)
1443 thread = gmx_omp_get_thread_num();
1444 if (flags & GMX_PME_SOLVE)
1449 wallcycle_start(wcycle, ewcPME_FFT);
1452 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX, thread, wcycle);
1455 wallcycle_stop(wcycle, ewcPME_FFT);
1459 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1463 if (flags & GMX_PME_SOLVE)
1465 /* solve in k-space for our local cells */
1466 #pragma omp parallel num_threads(pme->nthread) private(thread)
1471 thread = gmx_omp_get_thread_num();
1474 wallcycle_start(wcycle, ewcLJPME);
1477 loop_count = solve_pme_lj_yzx(
1478 pme, &pme->cfftgrid[2], TRUE,
1479 scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ],
1480 bCalcEnerVir, pme->nthread, thread);
1483 wallcycle_stop(wcycle, ewcLJPME);
1484 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1487 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1493 /* This should only be called on the master thread and
1494 * after the threads have synchronized.
1496 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output[fep_state]);
1501 bFirst = !pme->doCoulomb;
1502 calc_initial_lb_coeffs(coefficientBuffer, local_c6, local_sigma);
1503 for (grid_index = 8; grid_index >= 2; --grid_index)
1505 /* Unpack structure */
1506 pmegrid = &pme->pmegrid[grid_index];
1507 fftgrid = pme->fftgrid[grid_index];
1508 pfft_setup = pme->pfft_setup[grid_index];
1509 grid = pmegrid->grid.grid;
1510 calc_next_lb_coeffs(coefficientBuffer, local_sigma);
1511 #pragma omp parallel num_threads(pme->nthread) private(thread)
1515 thread = gmx_omp_get_thread_num();
1519 wallcycle_start(wcycle, ewcPME_FFT);
1522 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL, thread, wcycle);
1525 wallcycle_stop(wcycle, ewcPME_FFT);
1528 if (pme->nodeid == 0)
1530 real ntot = pme->nkx * pme->nky * pme->nkz;
1531 npme = static_cast<int>(ntot * std::log(ntot) / std::log(2.0));
1532 inc_nrnb(nrnb, eNR_FFT, 2 * npme);
1534 wallcycle_start(wcycle, ewcPME_GATHER);
1537 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1539 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1540 } /*#pragma omp parallel*/
1542 /* distribute local grid to all nodes */
1543 if (pme->nnodes > 1)
1545 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1548 unwrap_periodic_pmegrid(pme, grid);
1552 /* interpolate forces for our local atoms */
1553 bClearF = (bFirst && PAR(cr));
1554 scale = pme->bFEP ? (fep_state < 1 ? 1.0 - lambda_lj : lambda_lj) : 1.0;
1555 scale *= lb_scale_factor[grid_index - 2];
1557 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1558 for (thread = 0; thread < pme->nthread; thread++)
1562 gather_f_bsplines(pme, grid, bClearF, &pme->atc[0],
1563 &pme->atc[0].spline[thread], scale);
1565 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1569 inc_nrnb(nrnb, eNR_GATHERFBSP,
1570 pme->pme_order * pme->pme_order * pme->pme_order * pme->atc[0].numAtoms());
1572 wallcycle_stop(wcycle, ewcPME_GATHER);
1575 } /* for (grid_index = 8; grid_index >= 2; --grid_index) */
1577 } /* for (fep_state = 0; fep_state < fep_states_lj; ++fep_state) */
1578 } /* if ((flags & GMX_PME_DO_LJ) && pme->ljpme_combination_rule == eljpmeLB) */
1580 if (bCalcF && pme->nnodes > 1)
1582 wallcycle_start(wcycle, ewcPME_REDISTXF);
1583 for (d = 0; d < pme->ndecompdim; d++)
1585 gmx::ArrayRef<gmx::RVec> forcesRef;
1586 if (d == pme->ndecompdim - 1)
1588 const size_t numAtoms = coordinates.size();
1589 GMX_ASSERT(forces.size() >= numAtoms, "Need at least numAtoms forces");
1590 forcesRef = forces.subArray(0, numAtoms);
1594 forcesRef = pme->atc[d + 1].f;
1596 if (DOMAINDECOMP(cr))
1598 dd_pmeredist_f(pme, &pme->atc[d], forcesRef, d == pme->ndecompdim - 1 && pme->bPPnode);
1602 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1611 *energy_q = output[0].coulombEnergy_;
1612 m_add(vir_q, output[0].coulombVirial_, vir_q);
1616 *energy_q = (1.0 - lambda_q) * output[0].coulombEnergy_ + lambda_q * output[1].coulombEnergy_;
1617 *dvdlambda_q += output[1].coulombEnergy_ - output[0].coulombEnergy_;
1618 for (int i = 0; i < DIM; i++)
1620 for (int j = 0; j < DIM; j++)
1622 vir_q[i][j] += (1.0 - lambda_q) * output[0].coulombVirial_[i][j]
1623 + lambda_q * output[1].coulombVirial_[i][j];
1629 fprintf(debug, "Electrostatic PME mesh energy: %g\n", *energy_q);
1641 *energy_lj = output[0].lennardJonesEnergy_;
1642 m_add(vir_lj, output[0].lennardJonesVirial_, vir_lj);
1646 *energy_lj = (1.0 - lambda_lj) * output[0].lennardJonesEnergy_
1647 + lambda_lj * output[1].lennardJonesEnergy_;
1648 *dvdlambda_lj += output[1].lennardJonesEnergy_ - output[0].lennardJonesEnergy_;
1649 for (int i = 0; i < DIM; i++)
1651 for (int j = 0; j < DIM; j++)
1653 vir_lj[i][j] += (1.0 - lambda_lj) * output[0].lennardJonesVirial_[i][j]
1654 + lambda_lj * output[1].lennardJonesVirial_[i][j];
1660 fprintf(debug, "Lennard-Jones PME mesh energy: %g\n", *energy_lj);
1671 void gmx_pme_destroy(gmx_pme_t* pme)
1678 delete pme->boxScaler;
1687 for (int i = 0; i < pme->ngrids; ++i)
1689 pmegrids_destroy(&pme->pmegrid[i]);
1691 if (pme->pfft_setup)
1693 for (int i = 0; i < pme->ngrids; ++i)
1695 gmx_parallel_3dfft_destroy(pme->pfft_setup[i]);
1698 sfree(pme->fftgrid);
1699 sfree(pme->cfftgrid);
1700 sfree(pme->pfft_setup);
1702 for (int i = 0; i < DIM; i++)
1704 sfree(pme->bsp_mod[i]);
1710 if (pme->solve_work)
1712 pme_free_all_work(&pme->solve_work, pme->nthread);
1715 sfree(pme->sum_qgrid_tmp);
1716 sfree(pme->sum_qgrid_dd_tmp);
1718 destroy_pme_spline_work(pme->spline_work);
1720 if (pme_gpu_active(pme) && pme->gpu)
1722 pme_gpu_destroy(pme->gpu);
1728 void gmx_pme_reinit_atoms(gmx_pme_t* pme, const int numAtoms, const real* charges)
1730 if (pme_gpu_active(pme))
1732 pme_gpu_reinit_atoms(pme->gpu, numAtoms, charges);
1736 pme->atc[0].setNumAtoms(numAtoms);
1737 // TODO: set the charges here as well