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,2021, 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/mdtypes/simulation_workload.h"
103 #include "gromacs/pbcutil/pbc.h"
104 #include "gromacs/timing/cyclecounter.h"
105 #include "gromacs/timing/wallcycle.h"
106 #include "gromacs/timing/walltime_accounting.h"
107 #include "gromacs/topology/topology.h"
108 #include "gromacs/utility/basedefinitions.h"
109 #include "gromacs/utility/cstringutil.h"
110 #include "gromacs/utility/exceptions.h"
111 #include "gromacs/utility/fatalerror.h"
112 #include "gromacs/utility/gmxmpi.h"
113 #include "gromacs/utility/gmxomp.h"
114 #include "gromacs/utility/logger.h"
115 #include "gromacs/utility/real.h"
116 #include "gromacs/utility/smalloc.h"
117 #include "gromacs/utility/message_string_collector.h"
118 #include "gromacs/utility/unique_cptr.h"
120 #include "calculate_spline_moduli.h"
121 #include "pme_gather.h"
122 #include "pme_gpu_internal.h"
123 #include "pme_grid.h"
124 #include "pme_internal.h"
125 #include "pme_redistribute.h"
126 #include "pme_solve.h"
127 #include "pme_spline_work.h"
128 #include "pme_spread.h"
130 //NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
131 bool g_allowPmeWithSyclForTesting = false;
133 bool pme_gpu_supports_build(std::string* error)
135 gmx::MessageStringCollector errorReasons;
136 // Before changing the prefix string, make sure that it is not searched for in regression tests.
137 errorReasons.startContext("PME GPU does not support:");
138 errorReasons.appendIf(GMX_DOUBLE, "Double-precision build of GROMACS.");
139 errorReasons.appendIf(!GMX_GPU, "Non-GPU build of GROMACS.");
140 errorReasons.appendIf(GMX_GPU_SYCL && !g_allowPmeWithSyclForTesting, "SYCL build."); // SYCL-TODO
141 errorReasons.finishContext();
142 if (error != nullptr)
144 *error = errorReasons.toString();
146 return errorReasons.isEmpty();
149 bool pme_gpu_supports_hardware(const gmx_hw_info_t gmx_unused& hwinfo, std::string* error)
151 gmx::MessageStringCollector errorReasons;
152 // Before changing the prefix string, make sure that it is not searched for in regression tests.
153 errorReasons.startContext("PME GPU does not support:");
155 errorReasons.appendIf(GMX_GPU_OPENCL, "Apple OS X operating system");
157 errorReasons.finishContext();
158 if (error != nullptr)
160 *error = errorReasons.toString();
162 return errorReasons.isEmpty();
165 bool pme_gpu_supports_input(const t_inputrec& ir, std::string* error)
167 gmx::MessageStringCollector errorReasons;
168 // Before changing the prefix string, make sure that it is not searched for in regression tests.
169 errorReasons.startContext("PME GPU does not support:");
170 errorReasons.appendIf(!EEL_PME(ir.coulombtype),
171 "Systems that do not use PME for electrostatics.");
172 errorReasons.appendIf((ir.pme_order != 4), "Interpolation orders other than 4.");
173 errorReasons.appendIf(EVDW_PME(ir.vdwtype), "Lennard-Jones PME.");
174 errorReasons.appendIf(!EI_DYNAMICS(ir.eI), "Non-dynamical integrator (use md, sd, etc).");
175 errorReasons.finishContext();
176 if (error != nullptr)
178 *error = errorReasons.toString();
180 return errorReasons.isEmpty();
183 /*! \brief \libinternal
184 * Finds out if PME with given inputs is possible to run on GPU.
185 * This function is an internal final check, validating the whole PME structure on creation,
186 * but it still duplicates the preliminary checks from the above (externally exposed) pme_gpu_supports_input() - just in case.
188 * \param[in] pme The PME structure.
189 * \param[out] error The error message if the input is not supported on GPU.
190 * \returns True if this PME input is possible to run on GPU, false otherwise.
192 static bool pme_gpu_check_restrictions(const gmx_pme_t* pme, std::string* error)
194 gmx::MessageStringCollector errorReasons;
195 // Before changing the prefix string, make sure that it is not searched for in regression tests.
196 errorReasons.startContext("PME GPU does not support:");
197 errorReasons.appendIf((pme->nnodes != 1), "PME decomposition.");
198 errorReasons.appendIf((pme->pme_order != 4), "interpolation orders other than 4.");
199 errorReasons.appendIf(pme->doLJ, "Lennard-Jones PME.");
200 errorReasons.appendIf(GMX_DOUBLE, "Double precision build of GROMACS.");
201 errorReasons.appendIf(!GMX_GPU, "Non-GPU build of GROMACS.");
202 errorReasons.appendIf(GMX_GPU_SYCL && !g_allowPmeWithSyclForTesting, "SYCL build of GROMACS."); // SYCL-TODO
203 errorReasons.finishContext();
204 if (error != nullptr)
206 *error = errorReasons.toString();
208 return errorReasons.isEmpty();
211 PmeRunMode pme_run_mode(const gmx_pme_t* pme)
213 GMX_ASSERT(pme != nullptr, "Expecting valid PME data pointer");
217 gmx::PinningPolicy pme_get_pinning_policy()
219 return gmx::PinningPolicy::PinnedIfSupported;
222 /*! \brief Number of bytes in a cache line.
224 * Must also be a multiple of the SIMD and SIMD4 register size, to
225 * preserve alignment.
227 const int gmxCacheLineSize = 64;
229 //! Set up coordinate communication
230 static void setup_coordinate_communication(PmeAtomComm* atc)
238 for (i = 1; i <= nslab / 2; i++)
240 fw = (atc->nodeid + i) % nslab;
241 bw = (atc->nodeid - i + nslab) % nslab;
244 atc->slabCommSetup[n].node_dest = fw;
245 atc->slabCommSetup[n].node_src = bw;
250 atc->slabCommSetup[n].node_dest = bw;
251 atc->slabCommSetup[n].node_src = fw;
257 /*! \brief Round \p n up to the next multiple of \p f */
258 static int mult_up(int n, int f)
260 return ((n + f - 1) / f) * f;
263 /*! \brief Return estimate of the load imbalance from the PME grid not being a good match for the number of PME ranks */
264 static double estimate_pme_load_imbalance(struct gmx_pme_t* pme)
269 nma = pme->nnodes_major;
270 nmi = pme->nnodes_minor;
272 n1 = mult_up(pme->nkx, nma) * mult_up(pme->nky, nmi) * pme->nkz;
273 n2 = mult_up(pme->nkx, nma) * mult_up(pme->nkz, nmi) * pme->nky;
274 n3 = mult_up(pme->nky, nma) * mult_up(pme->nkz, nmi) * pme->nkx;
276 /* pme_solve is roughly double the cost of an fft */
278 return (n1 + n2 + 3 * n3) / static_cast<double>(6 * pme->nkx * pme->nky * pme->nkz);
283 PmeAtomComm::PmeAtomComm(MPI_Comm PmeMpiCommunicator,
284 const int numThreads,
287 const bool doSpread) :
288 dimind(dimIndex), bSpread(doSpread), pme_order(pmeOrder), nthread(numThreads), spline(nthread)
290 if (PmeMpiCommunicator != MPI_COMM_NULL)
292 mpi_comm = PmeMpiCommunicator;
294 MPI_Comm_size(mpi_comm, &nslab);
295 MPI_Comm_rank(mpi_comm, &nodeid);
300 fprintf(debug, "For PME atom communication in dimind %d: nslab %d rank %d\n", dimind, nslab, nodeid);
305 slabCommSetup.resize(nslab);
306 setup_coordinate_communication(this);
308 count_thread.resize(nthread);
309 for (auto& countThread : count_thread)
311 countThread.resize(nslab);
317 threadMap.resize(nthread);
319 # pragma omp parallel for num_threads(nthread) schedule(static)
320 for (int thread = 0; thread < nthread; thread++)
324 /* Allocate buffer with padding to avoid cache polution */
325 threadMap[thread].nBuffer.resize(nthread + 2 * gmxCacheLineSize);
326 threadMap[thread].n = threadMap[thread].nBuffer.data() + gmxCacheLineSize;
328 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
335 /*! \brief Initialize data structure for communication */
336 static void init_overlap_comm(pme_overlap_t* ol, int norder, MPI_Comm comm, int nnodes, int nodeid, int ndata, int commplainsize)
344 /* Linear translation of the PME grid won't affect reciprocal space
345 * calculations, so to optimize we only interpolate "upwards",
346 * which also means we only have to consider overlap in one direction.
347 * I.e., particles on this node might also be spread to grid indices
348 * that belong to higher nodes (modulo nnodes)
351 ol->s2g0.resize(ol->nnodes + 1);
352 ol->s2g1.resize(ol->nnodes);
355 fprintf(debug, "PME slab boundaries:");
357 for (int i = 0; i < nnodes; i++)
359 /* s2g0 the local interpolation grid start.
360 * s2g1 the local interpolation grid end.
361 * Since in calc_pidx we divide particles, and not grid lines,
362 * spatially uniform along dimension x or y, we need to round
363 * s2g0 down and s2g1 up.
365 ol->s2g0[i] = (i * ndata + 0) / nnodes;
366 ol->s2g1[i] = ((i + 1) * ndata + nnodes - 1) / nnodes + norder - 1;
370 fprintf(debug, " %3d %3d", ol->s2g0[i], ol->s2g1[i]);
373 ol->s2g0[nnodes] = ndata;
376 fprintf(debug, "\n");
379 /* Determine with how many nodes we need to communicate the grid overlap */
380 int testRankCount = 0;
385 for (int i = 0; i < nnodes; i++)
387 if ((i + testRankCount < nnodes && ol->s2g1[i] > ol->s2g0[i + testRankCount])
388 || (i + testRankCount >= nnodes && ol->s2g1[i] > ol->s2g0[i + testRankCount - nnodes] + ndata))
393 } while (bCont && testRankCount < nnodes);
395 ol->comm_data.resize(testRankCount - 1);
398 for (size_t b = 0; b < ol->comm_data.size(); b++)
400 pme_grid_comm_t* pgc = &ol->comm_data[b];
403 pgc->send_id = (ol->nodeid + (b + 1)) % ol->nnodes;
404 int fft_start = ol->s2g0[pgc->send_id];
405 int fft_end = ol->s2g0[pgc->send_id + 1];
406 if (pgc->send_id < nodeid)
411 int send_index1 = ol->s2g1[nodeid];
412 send_index1 = std::min(send_index1, fft_end);
413 pgc->send_index0 = fft_start;
414 pgc->send_nindex = std::max(0, send_index1 - pgc->send_index0);
415 ol->send_size += pgc->send_nindex;
417 /* We always start receiving to the first index of our slab */
418 pgc->recv_id = (ol->nodeid - (b + 1) + ol->nnodes) % ol->nnodes;
419 fft_start = ol->s2g0[ol->nodeid];
420 fft_end = ol->s2g0[ol->nodeid + 1];
421 int recv_index1 = ol->s2g1[pgc->recv_id];
422 if (pgc->recv_id > nodeid)
424 recv_index1 -= ndata;
426 recv_index1 = std::min(recv_index1, fft_end);
427 pgc->recv_index0 = fft_start;
428 pgc->recv_nindex = std::max(0, recv_index1 - pgc->recv_index0);
432 /* Communicate the buffer sizes to receive */
434 for (size_t b = 0; b < ol->comm_data.size(); b++)
436 MPI_Sendrecv(&ol->send_size,
439 ol->comm_data[b].send_id,
441 &ol->comm_data[b].recv_size,
444 ol->comm_data[b].recv_id,
451 /* For non-divisible grid we need pme_order iso pme_order-1 */
452 ol->sendbuf.resize(norder * commplainsize);
453 ol->recvbuf.resize(norder * commplainsize);
456 int minimalPmeGridSize(int pmeOrder)
458 /* The actual grid size limitations are:
459 * serial: >= pme_order
460 * DD, no OpenMP: >= 2*(pme_order - 1)
461 * DD, OpenMP: >= pme_order + 1
462 * But we use the maximum for simplicity since in practice there is not
463 * much performance difference between pme_order and 2*(pme_order -1).
465 int minimalSize = 2 * (pmeOrder - 1);
467 GMX_RELEASE_ASSERT(pmeOrder >= 3, "pmeOrder has to be >= 3");
468 GMX_RELEASE_ASSERT(minimalSize >= pmeOrder + 1, "The grid size should be >= pmeOrder + 1");
473 bool gmx_pme_check_restrictions(int pme_order, int nkx, int nky, int nkz, int numPmeDomainsAlongX, bool useThreads, bool errorsAreFatal)
475 if (pme_order > PME_ORDER_MAX)
482 std::string message = gmx::formatString(
483 "pme_order (%d) is larger than the maximum allowed value (%d). Modify and "
484 "recompile the code if you really need such a high order.",
487 GMX_THROW(gmx::InconsistentInputError(message));
490 const int minGridSize = minimalPmeGridSize(pme_order);
491 if (nkx < minGridSize || nky < minGridSize || nkz < minGridSize)
497 std::string message =
498 gmx::formatString("The PME grid sizes need to be >= 2*(pme_order-1) (%d)", minGridSize);
499 GMX_THROW(gmx::InconsistentInputError(message));
502 /* Check for a limitation of the (current) sum_fftgrid_dd code.
503 * We only allow multiple communication pulses in dim 1, not in dim 0.
506 && (nkx < numPmeDomainsAlongX * pme_order && nkx != numPmeDomainsAlongX * (pme_order - 1)))
513 "The number of PME grid lines per rank along x is %g. But when using OpenMP "
514 "threads, the number of grid lines per rank along x should be >= pme_order (%d) "
515 "or = pmeorder-1. To resolve this issue, use fewer ranks along x (and possibly "
516 "more along y and/or z) by specifying -dd manually.",
517 nkx / static_cast<double>(numPmeDomainsAlongX),
524 /*! \brief Round \p enumerator */
525 static int div_round_up(int enumerator, int denominator)
527 return (enumerator + denominator - 1) / denominator;
530 gmx_pme_t* gmx_pme_init(const t_commrec* cr,
531 const NumPmeDomains& numPmeDomains,
532 const t_inputrec* ir,
533 gmx_bool bFreeEnergy_q,
534 gmx_bool bFreeEnergy_lj,
535 gmx_bool bReproducible,
541 const DeviceContext* deviceContext,
542 const DeviceStream* deviceStream,
543 const PmeGpuProgram* pmeGpuProgram,
544 const gmx::MDLogger& mdlog)
546 int use_threads, sum_use_threads, i;
551 fprintf(debug, "Creating PME data structures.\n");
554 gmx::unique_cptr<gmx_pme_t, gmx_pme_destroy> pme(new gmx_pme_t());
561 pme->nnodes_major = numPmeDomains.x;
562 pme->nnodes_minor = numPmeDomains.y;
564 if (numPmeDomains.x * numPmeDomains.y > 1)
566 pme->mpi_comm = cr->mpi_comm_mygroup;
569 MPI_Comm_rank(pme->mpi_comm, &pme->nodeid);
570 MPI_Comm_size(pme->mpi_comm, &pme->nnodes);
572 if (pme->nnodes != numPmeDomains.x * numPmeDomains.y)
574 gmx_incons("PME rank count mismatch");
579 pme->mpi_comm = MPI_COMM_NULL;
582 if (pme->nnodes == 1)
584 pme->mpi_comm_d[0] = MPI_COMM_NULL;
585 pme->mpi_comm_d[1] = MPI_COMM_NULL;
587 pme->nodeid_major = 0;
588 pme->nodeid_minor = 0;
592 if (numPmeDomains.y == 1)
594 pme->mpi_comm_d[0] = pme->mpi_comm;
595 pme->mpi_comm_d[1] = MPI_COMM_NULL;
597 pme->nodeid_major = pme->nodeid;
598 pme->nodeid_minor = 0;
600 else if (numPmeDomains.x == 1)
602 pme->mpi_comm_d[0] = MPI_COMM_NULL;
603 pme->mpi_comm_d[1] = pme->mpi_comm;
605 pme->nodeid_major = 0;
606 pme->nodeid_minor = pme->nodeid;
610 if (pme->nnodes % numPmeDomains.x != 0)
613 "For 2D PME decomposition, #PME ranks must be divisible by the number of "
619 MPI_Comm_split(pme->mpi_comm,
620 pme->nodeid % numPmeDomains.y,
622 &pme->mpi_comm_d[0]); /* My communicator along major dimension */
623 MPI_Comm_split(pme->mpi_comm,
624 pme->nodeid / numPmeDomains.y,
626 &pme->mpi_comm_d[1]); /* My communicator along minor dimension */
628 MPI_Comm_rank(pme->mpi_comm_d[0], &pme->nodeid_major);
629 MPI_Comm_size(pme->mpi_comm_d[0], &pme->nnodes_major);
630 MPI_Comm_rank(pme->mpi_comm_d[1], &pme->nodeid_minor);
631 MPI_Comm_size(pme->mpi_comm_d[1], &pme->nnodes_minor);
635 // cr is always initialized if there is a a PP rank, so we can safely assume
636 // that when it is not, like in ewald tests, we not on a PP rank.
637 pme->bPPnode = ((cr != nullptr && cr->duty != 0) && thisRankHasDuty(cr, DUTY_PP));
639 pme->nthread = nthread;
641 /* Check if any of the PME MPI ranks uses threads */
642 use_threads = (pme->nthread > 1 ? 1 : 0);
646 MPI_Allreduce(&use_threads, &sum_use_threads, 1, MPI_INT, MPI_SUM, pme->mpi_comm);
651 sum_use_threads = use_threads;
653 pme->bUseThreads = (sum_use_threads > 0);
655 if (ir->pbcType == PbcType::Screw)
657 gmx_fatal(FARGS, "pme does not (yet) work with pbc = screw");
661 * It is likely that the current gmx_pme_do() routine supports calculating
662 * only Coulomb or LJ while gmx_pme_init() configures for both,
663 * but that has never been tested.
664 * It is likely that the current gmx_pme_do() routine supports calculating,
665 * not calculating free-energy for Coulomb and/or LJ while gmx_pme_init()
666 * configures with free-energy, but that has never been tested.
668 pme->doCoulomb = EEL_PME(ir->coulombtype);
669 pme->doLJ = EVDW_PME(ir->vdwtype);
670 pme->bFEP_q = ((ir->efep != FreeEnergyPerturbationType::No) && bFreeEnergy_q);
671 pme->bFEP_lj = ((ir->efep != FreeEnergyPerturbationType::No) && bFreeEnergy_lj);
672 pme->bFEP = (pme->bFEP_q || pme->bFEP_lj);
676 pme->bP3M = (ir->coulombtype == CoulombInteractionType::P3mAD || getenv("GMX_PME_P3M") != nullptr);
677 pme->pme_order = ir->pme_order;
678 pme->ewaldcoeff_q = ewaldcoeff_q;
679 pme->ewaldcoeff_lj = ewaldcoeff_lj;
681 /* Always constant electrostatics coefficients */
682 pme->epsilon_r = ir->epsilon_r;
684 /* Always constant LJ coefficients */
685 pme->ljpme_combination_rule = ir->ljpme_combination_rule;
687 // The box requires scaling with nwalls = 2, we store that condition as well
688 // as the scaling factor
689 pme->boxScaler = std::make_unique<EwaldBoxZScaler>(
690 EwaldBoxZScaler(inputrecPbcXY2Walls(ir), ir->wall_ewald_zfac));
692 /* If we violate restrictions, generate a fatal error here */
693 gmx_pme_check_restrictions(
694 pme->pme_order, pme->nkx, pme->nky, pme->nkz, pme->nnodes_major, pme->bUseThreads, true);
701 MPI_Type_contiguous(DIM, GMX_MPI_REAL, &(pme->rvec_mpi));
702 MPI_Type_commit(&(pme->rvec_mpi));
705 /* Note that the coefficient spreading and force gathering, which usually
706 * takes about the same amount of time as FFT+solve_pme,
707 * is always fully load balanced
708 * (unless the coefficient distribution is inhomogeneous).
711 imbal = estimate_pme_load_imbalance(pme.get());
712 if (imbal >= 1.2 && pme->nodeid_major == 0 && pme->nodeid_minor == 0)
714 GMX_LOG(mdlog.warning)
716 .appendTextFormatted(
717 "NOTE: The load imbalance in PME FFT and solve is %d%%.\n"
718 " For optimal PME load balancing\n"
719 " PME grid_x (%d) and grid_y (%d) should be divisible by "
722 " and PME grid_y (%d) and grid_z (%d) should be divisible by "
725 gmx::roundToInt((imbal - 1) * 100),
735 /* For non-divisible grid we need pme_order iso pme_order-1 */
736 /* In sum_qgrid_dd x overlap is copied in place: take padding into account.
737 * y is always copied through a buffer: we don't need padding in z,
738 * but we do need the overlap in x because of the communication order.
740 init_overlap_comm(&pme->overlap[0],
746 (div_round_up(pme->nky, pme->nnodes_minor) + pme->pme_order)
747 * (pme->nkz + pme->pme_order - 1));
749 /* Along overlap dim 1 we can send in multiple pulses in sum_fftgrid_dd.
750 * We do this with an offset buffer of equal size, so we need to allocate
751 * extra for the offset. That's what the (+1)*pme->nkz is for.
753 init_overlap_comm(&pme->overlap[1],
759 (div_round_up(pme->nkx, pme->nnodes_major) + pme->pme_order + 1) * pme->nkz);
761 /* Double-check for a limitation of the (current) sum_fftgrid_dd code.
762 * Note that gmx_pme_check_restrictions checked for this already.
764 if (pme->bUseThreads && (pme->overlap[0].comm_data.size() > 1))
767 "More than one communication pulse required for grid overlap communication along "
768 "the major dimension while using threads");
771 snew(pme->bsp_mod[XX], pme->nkx);
772 snew(pme->bsp_mod[YY], pme->nky);
773 snew(pme->bsp_mod[ZZ], pme->nkz);
775 pme->gpu = pmeGpu; /* Carrying over the single GPU structure */
776 pme->runMode = runMode;
778 /* The required size of the interpolation grid, including overlap.
779 * The allocated size (pmegrid_n?) might be slightly larger.
781 pme->pmegrid_nx = pme->overlap[0].s2g1[pme->nodeid_major] - pme->overlap[0].s2g0[pme->nodeid_major];
782 pme->pmegrid_ny = pme->overlap[1].s2g1[pme->nodeid_minor] - pme->overlap[1].s2g0[pme->nodeid_minor];
783 pme->pmegrid_nz_base = pme->nkz;
784 pme->pmegrid_nz = pme->pmegrid_nz_base + pme->pme_order - 1;
785 set_grid_alignment(&pme->pmegrid_nz, pme->pme_order);
786 pme->pmegrid_start_ix = pme->overlap[0].s2g0[pme->nodeid_major];
787 pme->pmegrid_start_iy = pme->overlap[1].s2g0[pme->nodeid_minor];
788 pme->pmegrid_start_iz = 0;
790 make_gridindex_to_localindex(
791 pme->nkx, pme->pmegrid_start_ix, pme->pmegrid_nx - (pme->pme_order - 1), &pme->nnx, &pme->fshx);
792 make_gridindex_to_localindex(
793 pme->nky, pme->pmegrid_start_iy, pme->pmegrid_ny - (pme->pme_order - 1), &pme->nny, &pme->fshy);
794 make_gridindex_to_localindex(
795 pme->nkz, pme->pmegrid_start_iz, pme->pmegrid_nz_base, &pme->nnz, &pme->fshz);
797 pme->spline_work = make_pme_spline_work(pme->pme_order);
802 /* It doesn't matter if we allocate too many grids here,
803 * we only allocate and use the ones we need.
807 pme->ngrids = ((ir->ljpme_combination_rule == LongRangeVdW::LB) ? DO_Q_AND_LJ_LB : DO_Q_AND_LJ);
813 snew(pme->fftgrid, pme->ngrids);
814 snew(pme->cfftgrid, pme->ngrids);
815 snew(pme->pfft_setup, pme->ngrids);
817 for (i = 0; i < pme->ngrids; ++i)
819 if ((i < DO_Q && pme->doCoulomb && (i == 0 || bFreeEnergy_q))
820 || (i >= DO_Q && pme->doLJ
821 && (i == 2 || bFreeEnergy_lj || ir->ljpme_combination_rule == LongRangeVdW::LB)))
823 pmegrids_init(&pme->pmegrid[i],
827 pme->pmegrid_nz_base,
831 pme->overlap[0].s2g1[pme->nodeid_major]
832 - pme->overlap[0].s2g0[pme->nodeid_major + 1],
833 pme->overlap[1].s2g1[pme->nodeid_minor]
834 - pme->overlap[1].s2g0[pme->nodeid_minor + 1]);
835 /* This routine will allocate the grid data to fit the FFTs */
836 const auto allocateRealGridForGpu = (pme->runMode == PmeRunMode::Mixed)
837 ? gmx::PinningPolicy::PinnedIfSupported
838 : gmx::PinningPolicy::CannotBePinned;
839 gmx_parallel_3dfft_init(&pme->pfft_setup[i],
846 allocateRealGridForGpu);
852 /* Use plain SPME B-spline interpolation */
853 make_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
857 /* Use the P3M grid-optimized influence function */
858 make_p3m_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
861 /* Use atc[0] for spreading */
862 const int firstDimIndex = (numPmeDomains.x > 1 ? 0 : 1);
863 MPI_Comm mpiCommFirstDim = (pme->nnodes > 1 ? pme->mpi_comm_d[firstDimIndex] : MPI_COMM_NULL);
864 bool doSpread = true;
865 pme->atc.emplace_back(mpiCommFirstDim, pme->nthread, pme->pme_order, firstDimIndex, doSpread);
866 if (pme->ndecompdim >= 2)
868 const int secondDimIndex = 1;
870 pme->atc.emplace_back(pme->mpi_comm_d[1], pme->nthread, pme->pme_order, secondDimIndex, doSpread);
873 // Initial check of validity of the input for running on the GPU
874 if (pme->runMode != PmeRunMode::CPU)
876 std::string errorString;
877 bool canRunOnGpu = pme_gpu_check_restrictions(pme.get(), &errorString);
880 GMX_THROW(gmx::NotImplementedError(errorString));
882 pme_gpu_reinit(pme.get(), deviceContext, deviceStream, pmeGpuProgram);
886 GMX_ASSERT(pme->gpu == nullptr, "Should not have PME GPU object when PME is on a CPU.");
890 pme_init_all_work(&pme->solve_work, pme->nthread, pme->nkx);
892 // no exception was thrown during the init, so we hand over the PME structure handle
893 return pme.release();
896 void gmx_pme_reinit(struct gmx_pme_t** pmedata,
898 struct gmx_pme_t* pme_src,
899 const t_inputrec* ir,
900 const ivec grid_size,
904 // Create a copy of t_inputrec fields that are used in gmx_pme_init().
905 // TODO: This would be better as just copying a sub-structure that contains
906 // all the PME parameters and nothing else.
908 irc.pbcType = ir->pbcType;
909 irc.coulombtype = ir->coulombtype;
910 irc.vdwtype = ir->vdwtype;
912 irc.pme_order = ir->pme_order;
913 irc.epsilon_r = ir->epsilon_r;
914 irc.ljpme_combination_rule = ir->ljpme_combination_rule;
915 irc.nkx = grid_size[XX];
916 irc.nky = grid_size[YY];
917 irc.nkz = grid_size[ZZ];
921 // This is reinit. Any logging should have been done at first init.
922 // Here we should avoid writing notes for settings the user did not
924 const gmx::MDLogger dummyLogger;
925 GMX_ASSERT(pmedata, "Invalid PME pointer");
926 NumPmeDomains numPmeDomains = { pme_src->nnodes_major, pme_src->nnodes_minor };
927 *pmedata = gmx_pme_init(cr,
942 /* When running PME on the CPU not using domain decomposition,
943 * the atom data is allocated once only in gmx_pme_(re)init().
945 if (!pme_src->gpu && pme_src->nnodes == 1)
947 gmx_pme_reinit_atoms(*pmedata, pme_src->atc[0].numAtoms(), {}, {});
949 // TODO this is mostly passing around current values
951 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
953 /* We can easily reuse the allocated pme grids in pme_src */
954 reuse_pmegrids(&pme_src->pmegrid[PME_GRID_QA], &(*pmedata)->pmegrid[PME_GRID_QA]);
955 /* We would like to reuse the fft grids, but that's harder */
958 real gmx_pme_calc_energy(gmx_pme_t* pme, gmx::ArrayRef<const gmx::RVec> x, gmx::ArrayRef<const real> q)
964 gmx_incons("gmx_pme_calc_energy called in parallel");
968 gmx_incons("gmx_pme_calc_energy with free energy");
971 if (!pme->atc_energy)
973 pme->atc_energy = std::make_unique<PmeAtomComm>(MPI_COMM_NULL, 1, pme->pme_order, 0, true);
975 PmeAtomComm* atc = pme->atc_energy.get();
976 atc->setNumAtoms(x.ssize());
978 atc->coefficient = q;
980 /* We only use the A-charges grid */
981 grid = &pme->pmegrid[PME_GRID_QA];
983 /* Only calculate the spline coefficients, don't actually spread */
984 spread_on_grid(pme, atc, nullptr, TRUE, FALSE, pme->fftgrid[PME_GRID_QA], FALSE, PME_GRID_QA);
986 return gather_energy_bsplines(pme, grid->grid.grid, atc);
989 /*! \brief Calculate initial Lorentz-Berthelot coefficients for LJ-PME */
990 static void calc_initial_lb_coeffs(gmx::ArrayRef<real> coefficient,
991 gmx::ArrayRef<const real> local_c6,
992 gmx::ArrayRef<const real> local_sigma)
994 for (gmx::index i = 0; i < coefficient.ssize(); ++i)
996 real sigma4 = local_sigma[i];
997 sigma4 = sigma4 * sigma4;
998 sigma4 = sigma4 * sigma4;
999 coefficient[i] = local_c6[i] / sigma4;
1003 /*! \brief Calculate next Lorentz-Berthelot coefficients for LJ-PME */
1004 static void calc_next_lb_coeffs(gmx::ArrayRef<real> coefficient, gmx::ArrayRef<const real> local_sigma)
1006 for (gmx::index i = 0; i < coefficient.ssize(); ++i)
1008 coefficient[i] *= local_sigma[i];
1012 int gmx_pme_do(struct gmx_pme_t* pme,
1013 gmx::ArrayRef<const gmx::RVec> coordinates,
1014 gmx::ArrayRef<gmx::RVec> forces,
1015 gmx::ArrayRef<const real> chargeA,
1016 gmx::ArrayRef<const real> chargeB,
1017 gmx::ArrayRef<const real> c6A,
1018 gmx::ArrayRef<const real> c6B,
1019 gmx::ArrayRef<const real> sigmaA,
1020 gmx::ArrayRef<const real> sigmaB,
1022 const t_commrec* cr,
1026 gmx_wallcycle* wcycle,
1035 const gmx::StepWorkload& stepWork)
1037 GMX_ASSERT(pme->runMode == PmeRunMode::CPU,
1038 "gmx_pme_do should not be called on the GPU PME run.");
1040 PmeAtomComm& atc = pme->atc[0];
1041 pmegrids_t* pmegrid = nullptr;
1042 real* grid = nullptr;
1043 gmx::ArrayRef<const real> coefficient;
1044 std::array<PmeOutput, 2> output; // The second is used for the B state with FEP
1045 gmx_parallel_3dfft_t pfft_setup;
1047 t_complex* cfftgrid;
1049 const int fep_states_lj = pme->bFEP_lj ? 2 : 1;
1050 // There's no support for computing energy without virial, or vice versa
1051 const bool computeEnergyAndVirial = (stepWork.computeEnergy || stepWork.computeVirial);
1053 /* We could be passing lambda!=0 while no q or LJ is actually perturbed */
1063 assert(pme->nnodes > 0);
1064 assert(pme->nnodes == 1 || pme->ndecompdim > 0);
1066 if (pme->nnodes > 1)
1068 atc.pd.resize(coordinates.ssize());
1069 for (int d = pme->ndecompdim - 1; d >= 0; d--)
1071 PmeAtomComm& atc = pme->atc[d];
1072 atc.maxshift = (atc.dimind == 0 ? maxshift_x : maxshift_y);
1077 GMX_ASSERT(coordinates.ssize() == atc.numAtoms(), "We expect atc.numAtoms() coordinates");
1078 GMX_ASSERT(forces.ssize() >= atc.numAtoms(),
1079 "We need a force buffer with at least atc.numAtoms() elements");
1081 atc.x = coordinates;
1086 pme->boxScaler->scaleBox(box, scaledBox);
1088 gmx::invertBoxMatrix(scaledBox, pme->recipbox);
1091 /* For simplicity, we construct the splines for all particles if
1092 * more than one PME calculations is needed. Some optimization
1093 * could be done by keeping track of which atoms have splines
1094 * constructed, and construct new splines on each pass for atoms
1095 * that don't yet have them.
1098 bool bDoSplines = pme->bFEP || (pme->doCoulomb && pme->doLJ);
1100 /* We need a maximum of four separate PME calculations:
1101 * grid_index=0: Coulomb PME with charges from state A
1102 * grid_index=1: Coulomb PME with charges from state B
1103 * grid_index=2: LJ PME with C6 from state A
1104 * grid_index=3: LJ PME with C6 from state B
1105 * For Lorentz-Berthelot combination rules, a separate loop is used to
1106 * calculate all the terms
1109 /* If we are doing LJ-PME with LB, we only do Q here */
1110 const int max_grid_index = (pme->ljpme_combination_rule == LongRangeVdW::LB) ? DO_Q : DO_Q_AND_LJ;
1112 for (int grid_index = 0; grid_index < max_grid_index; ++grid_index)
1114 /* Check if we should do calculations at this grid_index
1115 * If grid_index is odd we should be doing FEP
1116 * If grid_index < 2 we should be doing electrostatic PME
1117 * If grid_index >= 2 we should be doing LJ-PME
1119 if ((grid_index < DO_Q && (!pme->doCoulomb || (grid_index == 1 && !pme->bFEP_q)))
1120 || (grid_index >= DO_Q && (!pme->doLJ || (grid_index == 3 && !pme->bFEP_lj))))
1124 /* Unpack structure */
1125 pmegrid = &pme->pmegrid[grid_index];
1126 fftgrid = pme->fftgrid[grid_index];
1127 cfftgrid = pme->cfftgrid[grid_index];
1128 pfft_setup = pme->pfft_setup[grid_index];
1131 case 0: coefficient = chargeA; break;
1132 case 1: coefficient = chargeB; break;
1133 case 2: coefficient = c6A; break;
1134 case 3: coefficient = c6B; break;
1137 grid = pmegrid->grid.grid;
1139 if (pme->nnodes == 1)
1141 atc.coefficient = coefficient;
1145 wallcycle_start(wcycle, WallCycleCounter::PmeRedistXF);
1146 do_redist_pos_coeffs(pme, cr, bFirst, coordinates, coefficient);
1148 wallcycle_stop(wcycle, WallCycleCounter::PmeRedistXF);
1151 wallcycle_start(wcycle, WallCycleCounter::PmeSpread);
1153 /* Spread the coefficients on a grid */
1154 spread_on_grid(pme, &atc, pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1158 inc_nrnb(nrnb, eNR_WEIGHTS, DIM * atc.numAtoms());
1160 inc_nrnb(nrnb, eNR_SPREADBSP, pme->pme_order * pme->pme_order * pme->pme_order * atc.numAtoms());
1162 if (!pme->bUseThreads)
1164 wrap_periodic_pmegrid(pme, grid);
1166 /* sum contributions to local grid from other nodes */
1167 if (pme->nnodes > 1)
1169 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1172 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1175 wallcycle_stop(wcycle, WallCycleCounter::PmeSpread);
1177 /* TODO If the OpenMP and single-threaded implementations
1178 converge, then spread_on_grid() and
1179 copy_pmegrid_to_fftgrid() will perhaps live in the same
1183 /* Here we start a large thread parallel region */
1184 #pragma omp parallel num_threads(pme->nthread) private(thread)
1188 thread = gmx_omp_get_thread_num();
1194 wallcycle_start(wcycle, WallCycleCounter::PmeFft);
1196 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX, thread, wcycle);
1199 wallcycle_stop(wcycle, WallCycleCounter::PmeFft);
1202 /* solve in k-space for our local cells */
1207 (grid_index < DO_Q ? WallCycleCounter::PmeSolve : WallCycleCounter::LJPme));
1209 if (grid_index < DO_Q)
1211 loop_count = solve_pme_yzx(pme,
1213 scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ],
1214 computeEnergyAndVirial,
1221 solve_pme_lj_yzx(pme,
1224 scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ],
1225 computeEnergyAndVirial,
1234 (grid_index < DO_Q ? WallCycleCounter::PmeSolve : WallCycleCounter::LJPme));
1235 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1241 wallcycle_start(wcycle, WallCycleCounter::PmeFft);
1243 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL, thread, wcycle);
1246 wallcycle_stop(wcycle, WallCycleCounter::PmeFft);
1249 if (pme->nodeid == 0)
1251 real ntot = pme->nkx * pme->nky * pme->nkz;
1252 const int npme = static_cast<int>(ntot * std::log(ntot) / std::log(2.0));
1253 inc_nrnb(nrnb, eNR_FFT, 2 * npme);
1256 /* Note: this wallcycle region is closed below
1257 outside an OpenMP region, so take care if
1258 refactoring code here. */
1259 wallcycle_start(wcycle, WallCycleCounter::PmeGather);
1262 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1264 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1266 /* End of thread parallel section.
1267 * With MPI we have to synchronize here before gmx_sum_qgrid_dd.
1270 /* distribute local grid to all nodes */
1271 if (pme->nnodes > 1)
1273 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1276 unwrap_periodic_pmegrid(pme, grid);
1278 if (stepWork.computeForces)
1280 /* interpolate forces for our local atoms */
1283 /* If we are running without parallelization,
1284 * atc->f is the actual force array, not a buffer,
1285 * therefore we should not clear it.
1287 real lambda = grid_index < DO_Q ? lambda_q : lambda_lj;
1288 bClearF = (bFirst && PAR(cr));
1289 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1290 for (thread = 0; thread < pme->nthread; thread++)
1294 gather_f_bsplines(pme,
1298 &atc.spline[thread],
1299 pme->bFEP ? (grid_index % 2 == 0 ? 1.0 - lambda : lambda) : 1.0);
1301 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1305 inc_nrnb(nrnb, eNR_GATHERFBSP, pme->pme_order * pme->pme_order * pme->pme_order * atc.numAtoms());
1306 /* Note: this wallcycle region is opened above inside an OpenMP
1307 region, so take care if refactoring code here. */
1308 wallcycle_stop(wcycle, WallCycleCounter::PmeGather);
1311 if (computeEnergyAndVirial)
1313 /* This should only be called on the master thread
1314 * and after the threads have synchronized.
1318 get_pme_ener_vir_q(pme->solve_work, pme->nthread, &output[grid_index % 2]);
1322 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output[grid_index % 2]);
1326 } /* of grid_index-loop */
1328 /* For Lorentz-Berthelot combination rules in LJ-PME, we need to calculate
1331 if (pme->doLJ && pme->ljpme_combination_rule == LongRangeVdW::LB)
1333 /* Loop over A- and B-state if we are doing FEP */
1334 for (int fep_state = 0; fep_state < fep_states_lj; ++fep_state)
1336 std::vector<real> local_c6;
1337 std::vector<real> local_sigma;
1338 gmx::ArrayRef<const real> RedistC6;
1339 gmx::ArrayRef<const real> RedistSigma;
1340 gmx::ArrayRef<real> coefficientBuffer;
1341 if (pme->nnodes == 1)
1343 pme->lb_buf1.resize(atc.numAtoms());
1344 coefficientBuffer = pme->lb_buf1;
1348 local_c6.assign(c6A.begin(), c6A.end());
1349 local_sigma.assign(sigmaA.begin(), sigmaA.end());
1352 local_c6.assign(c6B.begin(), c6B.end());
1353 local_sigma.assign(sigmaB.begin(), sigmaB.end());
1355 default: gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1360 coefficientBuffer = atc.coefficientBuffer;
1365 RedistSigma = sigmaA;
1369 RedistSigma = sigmaB;
1371 default: gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1373 wallcycle_start(wcycle, WallCycleCounter::PmeRedistXF);
1375 do_redist_pos_coeffs(pme, cr, bFirst, coordinates, RedistC6);
1376 pme->lb_buf1.resize(atc.numAtoms());
1377 pme->lb_buf2.resize(atc.numAtoms());
1378 local_c6.assign(pme->lb_buf1.begin(), pme->lb_buf1.end());
1379 for (int i = 0; i < atc.numAtoms(); ++i)
1381 local_c6[i] = atc.coefficient[i];
1384 do_redist_pos_coeffs(pme, cr, FALSE, coordinates, RedistSigma);
1385 local_sigma.assign(pme->lb_buf2.begin(), pme->lb_buf2.end());
1386 for (int i = 0; i < atc.numAtoms(); ++i)
1388 local_sigma[i] = atc.coefficient[i];
1391 wallcycle_stop(wcycle, WallCycleCounter::PmeRedistXF);
1393 atc.coefficient = coefficientBuffer;
1394 calc_initial_lb_coeffs(coefficientBuffer, local_c6, local_sigma);
1396 /*Seven terms in LJ-PME with LB, grid_index < 2 reserved for electrostatics*/
1397 for (int grid_index = 2; grid_index < 9; ++grid_index)
1399 /* Unpack structure */
1400 pmegrid = &pme->pmegrid[grid_index];
1401 fftgrid = pme->fftgrid[grid_index];
1402 pfft_setup = pme->pfft_setup[grid_index];
1403 calc_next_lb_coeffs(coefficientBuffer, local_sigma);
1404 grid = pmegrid->grid.grid;
1406 wallcycle_start(wcycle, WallCycleCounter::PmeSpread);
1407 /* Spread the c6 on a grid */
1408 spread_on_grid(pme, &atc, pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1412 inc_nrnb(nrnb, eNR_WEIGHTS, DIM * atc.numAtoms());
1417 pme->pme_order * pme->pme_order * pme->pme_order * atc.numAtoms());
1418 if (pme->nthread == 1)
1420 wrap_periodic_pmegrid(pme, grid);
1421 /* sum contributions to local grid from other nodes */
1422 if (pme->nnodes > 1)
1424 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1426 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1428 wallcycle_stop(wcycle, WallCycleCounter::PmeSpread);
1430 /*Here we start a large thread parallel region*/
1431 #pragma omp parallel num_threads(pme->nthread) private(thread)
1435 thread = gmx_omp_get_thread_num();
1439 wallcycle_start(wcycle, WallCycleCounter::PmeFft);
1442 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX, thread, wcycle);
1445 wallcycle_stop(wcycle, WallCycleCounter::PmeFft);
1448 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1452 /* solve in k-space for our local cells */
1453 #pragma omp parallel num_threads(pme->nthread) private(thread)
1458 thread = gmx_omp_get_thread_num();
1461 wallcycle_start(wcycle, WallCycleCounter::LJPme);
1465 solve_pme_lj_yzx(pme,
1468 scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ],
1469 computeEnergyAndVirial,
1474 wallcycle_stop(wcycle, WallCycleCounter::LJPme);
1475 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1478 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1481 if (computeEnergyAndVirial)
1483 /* This should only be called on the master thread and
1484 * after the threads have synchronized.
1486 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output[fep_state]);
1489 bFirst = !pme->doCoulomb;
1490 calc_initial_lb_coeffs(coefficientBuffer, local_c6, local_sigma);
1491 for (int grid_index = 8; grid_index >= 2; --grid_index)
1493 /* Unpack structure */
1494 pmegrid = &pme->pmegrid[grid_index];
1495 fftgrid = pme->fftgrid[grid_index];
1496 pfft_setup = pme->pfft_setup[grid_index];
1497 grid = pmegrid->grid.grid;
1498 calc_next_lb_coeffs(coefficientBuffer, local_sigma);
1499 #pragma omp parallel num_threads(pme->nthread) private(thread)
1503 thread = gmx_omp_get_thread_num();
1507 wallcycle_start(wcycle, WallCycleCounter::PmeFft);
1510 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL, thread, wcycle);
1513 wallcycle_stop(wcycle, WallCycleCounter::PmeFft);
1516 if (pme->nodeid == 0)
1518 real ntot = pme->nkx * pme->nky * pme->nkz;
1519 const int npme = static_cast<int>(ntot * std::log(ntot) / std::log(2.0));
1520 inc_nrnb(nrnb, eNR_FFT, 2 * npme);
1522 wallcycle_start(wcycle, WallCycleCounter::PmeGather);
1525 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1527 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1528 } /*#pragma omp parallel*/
1530 /* distribute local grid to all nodes */
1531 if (pme->nnodes > 1)
1533 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1536 unwrap_periodic_pmegrid(pme, grid);
1538 if (stepWork.computeForces)
1540 /* interpolate forces for our local atoms */
1541 bClearF = (bFirst && PAR(cr));
1542 real scale = pme->bFEP ? (fep_state < 1 ? 1.0 - lambda_lj : lambda_lj) : 1.0;
1543 scale *= lb_scale_factor[grid_index - 2];
1545 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1546 for (thread = 0; thread < pme->nthread; thread++)
1551 pme, grid, bClearF, &pme->atc[0], &pme->atc[0].spline[thread], scale);
1553 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1559 pme->pme_order * pme->pme_order * pme->pme_order * pme->atc[0].numAtoms());
1561 wallcycle_stop(wcycle, WallCycleCounter::PmeGather);
1564 } /* for (grid_index = 8; grid_index >= 2; --grid_index) */
1565 } /* for (fep_state = 0; fep_state < fep_states_lj; ++fep_state) */
1566 } /* if (pme->doLJ && pme->ljpme_combination_rule == LongRangeVdW::LB) */
1568 if (stepWork.computeForces && pme->nnodes > 1)
1570 wallcycle_start(wcycle, WallCycleCounter::PmeRedistXF);
1571 for (int d = 0; d < pme->ndecompdim; d++)
1573 gmx::ArrayRef<gmx::RVec> forcesRef;
1574 if (d == pme->ndecompdim - 1)
1576 const size_t numAtoms = coordinates.size();
1577 GMX_ASSERT(forces.size() >= numAtoms, "Need at least numAtoms forces");
1578 forcesRef = forces.subArray(0, numAtoms);
1582 forcesRef = pme->atc[d + 1].f;
1584 if (haveDDAtomOrdering(*cr))
1586 dd_pmeredist_f(pme, &pme->atc[d], forcesRef, d == pme->ndecompdim - 1 && pme->bPPnode);
1590 wallcycle_stop(wcycle, WallCycleCounter::PmeRedistXF);
1593 if (computeEnergyAndVirial)
1599 *energy_q = output[0].coulombEnergy_;
1600 m_add(vir_q, output[0].coulombVirial_, vir_q);
1604 *energy_q = (1.0 - lambda_q) * output[0].coulombEnergy_ + lambda_q * output[1].coulombEnergy_;
1605 *dvdlambda_q += output[1].coulombEnergy_ - output[0].coulombEnergy_;
1606 for (int i = 0; i < DIM; i++)
1608 for (int j = 0; j < DIM; j++)
1610 vir_q[i][j] += (1.0 - lambda_q) * output[0].coulombVirial_[i][j]
1611 + lambda_q * output[1].coulombVirial_[i][j];
1625 *energy_lj = output[0].lennardJonesEnergy_;
1626 m_add(vir_lj, output[0].lennardJonesVirial_, vir_lj);
1630 *energy_lj = (1.0 - lambda_lj) * output[0].lennardJonesEnergy_
1631 + lambda_lj * output[1].lennardJonesEnergy_;
1632 *dvdlambda_lj += output[1].lennardJonesEnergy_ - output[0].lennardJonesEnergy_;
1633 for (int i = 0; i < DIM; i++)
1635 for (int j = 0; j < DIM; j++)
1637 vir_lj[i][j] += (1.0 - lambda_lj) * output[0].lennardJonesVirial_[i][j]
1638 + lambda_lj * output[1].lennardJonesVirial_[i][j];
1651 void gmx_pme_destroy(gmx_pme_t* pme)
1665 for (int i = 0; i < pme->ngrids; ++i)
1667 pmegrids_destroy(&pme->pmegrid[i]);
1669 if (pme->pfft_setup)
1671 for (int i = 0; i < pme->ngrids; ++i)
1673 gmx_parallel_3dfft_destroy(pme->pfft_setup[i]);
1676 sfree(pme->fftgrid);
1677 sfree(pme->cfftgrid);
1678 sfree(pme->pfft_setup);
1680 for (int i = 0; i < DIM; i++)
1682 sfree(pme->bsp_mod[i]);
1688 if (pme->solve_work)
1690 pme_free_all_work(&pme->solve_work, pme->nthread);
1693 destroy_pme_spline_work(pme->spline_work);
1695 if (pme->gpu != nullptr)
1697 pme_gpu_destroy(pme->gpu);
1703 void gmx_pme_reinit_atoms(gmx_pme_t* pme,
1705 gmx::ArrayRef<const real> chargesA,
1706 gmx::ArrayRef<const real> chargesB)
1708 if (pme->gpu != nullptr)
1710 GMX_ASSERT(!(pme->bFEP_q && chargesB.empty()),
1711 "B state charges must be specified if running Coulomb FEP on the GPU");
1712 pme_gpu_reinit_atoms(pme->gpu, numAtoms, chargesA.data(), pme->bFEP_q ? chargesB.data() : nullptr);
1716 pme->atc[0].setNumAtoms(numAtoms);
1717 // TODO: set the charges here as well
1721 bool gmx_pme_grid_matches(const gmx_pme_t& pme, const ivec grid_size)
1723 return (pme.nkx == grid_size[XX] && pme.nky == grid_size[YY] && pme.nkz == grid_size[ZZ]);
1726 void gmx::SeparatePmeRanksPermitted::disablePmeRanks(const std::string& reason)
1728 permitSeparatePmeRanks_ = false;
1730 if (!reason.empty())
1732 reasons_.push_back(reason);
1736 bool gmx::SeparatePmeRanksPermitted::permitSeparatePmeRanks() const
1738 return permitSeparatePmeRanks_;
1741 std::string gmx::SeparatePmeRanksPermitted::reasonsWhyDisabled() const
1743 return joinStrings(reasons_, "; ");