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/stringutil.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 /*! \brief Help build a descriptive message in \c error if there are
131 * \c errorReasons why PME on GPU is not supported.
133 * \returns Whether the lack of errorReasons indicate there is support. */
134 static bool addMessageIfNotSupported(const std::list<std::string>& errorReasons, std::string* error)
136 bool isSupported = errorReasons.empty();
137 if (!isSupported && error)
139 std::string regressionTestMarker = "PME GPU does not support";
140 // this prefix is tested for in the regression tests script gmxtest.pl
141 *error = regressionTestMarker;
142 if (errorReasons.size() == 1)
144 *error += " " + errorReasons.back();
148 *error += ": " + gmx::joinStrings(errorReasons, "; ");
155 bool pme_gpu_supports_build(std::string* error)
157 std::list<std::string> errorReasons;
160 errorReasons.emplace_back("a double-precision build");
164 errorReasons.emplace_back("a non-GPU build");
168 errorReasons.emplace_back("SYCL build"); // SYCL-TODO
170 return addMessageIfNotSupported(errorReasons, error);
173 bool pme_gpu_supports_hardware(const gmx_hw_info_t gmx_unused& hwinfo, std::string* error)
175 std::list<std::string> errorReasons;
180 errorReasons.emplace_back("Apple OS X operating system");
183 return addMessageIfNotSupported(errorReasons, error);
186 bool pme_gpu_supports_input(const t_inputrec& ir, std::string* error)
188 std::list<std::string> errorReasons;
189 if (!EEL_PME(ir.coulombtype))
191 errorReasons.emplace_back("systems that do not use PME for electrostatics");
193 if (ir.pme_order != 4)
195 errorReasons.emplace_back("interpolation orders other than 4");
197 if (EVDW_PME(ir.vdwtype))
199 errorReasons.emplace_back("Lennard-Jones PME");
201 if (!EI_DYNAMICS(ir.eI))
203 errorReasons.emplace_back(
204 "Cannot compute PME interactions on a GPU, because PME GPU requires a dynamical "
205 "integrator (md, sd, etc).");
207 return addMessageIfNotSupported(errorReasons, error);
210 /*! \brief \libinternal
211 * Finds out if PME with given inputs is possible to run on GPU.
212 * This function is an internal final check, validating the whole PME structure on creation,
213 * but it still duplicates the preliminary checks from the above (externally exposed) pme_gpu_supports_input() - just in case.
215 * \param[in] pme The PME structure.
216 * \param[out] error The error message if the input is not supported on GPU.
217 * \returns True if this PME input is possible to run on GPU, false otherwise.
219 static bool pme_gpu_check_restrictions(const gmx_pme_t* pme, std::string* error)
221 std::list<std::string> errorReasons;
222 if (pme->nnodes != 1)
224 errorReasons.emplace_back("PME decomposition");
226 if (pme->pme_order != 4)
228 errorReasons.emplace_back("interpolation orders other than 4");
232 errorReasons.emplace_back("Lennard-Jones PME");
236 errorReasons.emplace_back("double precision");
240 errorReasons.emplace_back("non-GPU build of GROMACS");
244 errorReasons.emplace_back("SYCL build of GROMACS"); // SYCL-TODO
246 return addMessageIfNotSupported(errorReasons, error);
249 PmeRunMode pme_run_mode(const gmx_pme_t* pme)
251 GMX_ASSERT(pme != nullptr, "Expecting valid PME data pointer");
255 gmx::PinningPolicy pme_get_pinning_policy()
257 return gmx::PinningPolicy::PinnedIfSupported;
260 /*! \brief Number of bytes in a cache line.
262 * Must also be a multiple of the SIMD and SIMD4 register size, to
263 * preserve alignment.
265 const int gmxCacheLineSize = 64;
267 //! Set up coordinate communication
268 static void setup_coordinate_communication(PmeAtomComm* atc)
276 for (i = 1; i <= nslab / 2; i++)
278 fw = (atc->nodeid + i) % nslab;
279 bw = (atc->nodeid - i + nslab) % nslab;
282 atc->slabCommSetup[n].node_dest = fw;
283 atc->slabCommSetup[n].node_src = bw;
288 atc->slabCommSetup[n].node_dest = bw;
289 atc->slabCommSetup[n].node_src = fw;
295 /*! \brief Round \p n up to the next multiple of \p f */
296 static int mult_up(int n, int f)
298 return ((n + f - 1) / f) * f;
301 /*! \brief Return estimate of the load imbalance from the PME grid not being a good match for the number of PME ranks */
302 static double estimate_pme_load_imbalance(struct gmx_pme_t* pme)
307 nma = pme->nnodes_major;
308 nmi = pme->nnodes_minor;
310 n1 = mult_up(pme->nkx, nma) * mult_up(pme->nky, nmi) * pme->nkz;
311 n2 = mult_up(pme->nkx, nma) * mult_up(pme->nkz, nmi) * pme->nky;
312 n3 = mult_up(pme->nky, nma) * mult_up(pme->nkz, nmi) * pme->nkx;
314 /* pme_solve is roughly double the cost of an fft */
316 return (n1 + n2 + 3 * n3) / static_cast<double>(6 * pme->nkx * pme->nky * pme->nkz);
321 PmeAtomComm::PmeAtomComm(MPI_Comm PmeMpiCommunicator,
322 const int numThreads,
325 const bool doSpread) :
332 if (PmeMpiCommunicator != MPI_COMM_NULL)
334 mpi_comm = PmeMpiCommunicator;
336 MPI_Comm_size(mpi_comm, &nslab);
337 MPI_Comm_rank(mpi_comm, &nodeid);
342 fprintf(debug, "For PME atom communication in dimind %d: nslab %d rank %d\n", dimind, nslab, nodeid);
347 slabCommSetup.resize(nslab);
348 setup_coordinate_communication(this);
350 count_thread.resize(nthread);
351 for (auto& countThread : count_thread)
353 countThread.resize(nslab);
359 threadMap.resize(nthread);
361 # pragma omp parallel for num_threads(nthread) schedule(static)
362 for (int thread = 0; thread < nthread; thread++)
366 /* Allocate buffer with padding to avoid cache polution */
367 threadMap[thread].nBuffer.resize(nthread + 2 * gmxCacheLineSize);
368 threadMap[thread].n = threadMap[thread].nBuffer.data() + gmxCacheLineSize;
370 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
377 /*! \brief Initialize data structure for communication */
378 static void init_overlap_comm(pme_overlap_t* ol, int norder, MPI_Comm comm, int nnodes, int nodeid, int ndata, int commplainsize)
386 /* Linear translation of the PME grid won't affect reciprocal space
387 * calculations, so to optimize we only interpolate "upwards",
388 * which also means we only have to consider overlap in one direction.
389 * I.e., particles on this node might also be spread to grid indices
390 * that belong to higher nodes (modulo nnodes)
393 ol->s2g0.resize(ol->nnodes + 1);
394 ol->s2g1.resize(ol->nnodes);
397 fprintf(debug, "PME slab boundaries:");
399 for (int i = 0; i < nnodes; i++)
401 /* s2g0 the local interpolation grid start.
402 * s2g1 the local interpolation grid end.
403 * Since in calc_pidx we divide particles, and not grid lines,
404 * spatially uniform along dimension x or y, we need to round
405 * s2g0 down and s2g1 up.
407 ol->s2g0[i] = (i * ndata + 0) / nnodes;
408 ol->s2g1[i] = ((i + 1) * ndata + nnodes - 1) / nnodes + norder - 1;
412 fprintf(debug, " %3d %3d", ol->s2g0[i], ol->s2g1[i]);
415 ol->s2g0[nnodes] = ndata;
418 fprintf(debug, "\n");
421 /* Determine with how many nodes we need to communicate the grid overlap */
422 int testRankCount = 0;
427 for (int i = 0; i < nnodes; i++)
429 if ((i + testRankCount < nnodes && ol->s2g1[i] > ol->s2g0[i + testRankCount])
430 || (i + testRankCount >= nnodes && ol->s2g1[i] > ol->s2g0[i + testRankCount - nnodes] + ndata))
435 } while (bCont && testRankCount < nnodes);
437 ol->comm_data.resize(testRankCount - 1);
440 for (size_t b = 0; b < ol->comm_data.size(); b++)
442 pme_grid_comm_t* pgc = &ol->comm_data[b];
445 pgc->send_id = (ol->nodeid + (b + 1)) % ol->nnodes;
446 int fft_start = ol->s2g0[pgc->send_id];
447 int fft_end = ol->s2g0[pgc->send_id + 1];
448 if (pgc->send_id < nodeid)
453 int send_index1 = ol->s2g1[nodeid];
454 send_index1 = std::min(send_index1, fft_end);
455 pgc->send_index0 = fft_start;
456 pgc->send_nindex = std::max(0, send_index1 - pgc->send_index0);
457 ol->send_size += pgc->send_nindex;
459 /* We always start receiving to the first index of our slab */
460 pgc->recv_id = (ol->nodeid - (b + 1) + ol->nnodes) % ol->nnodes;
461 fft_start = ol->s2g0[ol->nodeid];
462 fft_end = ol->s2g0[ol->nodeid + 1];
463 int recv_index1 = ol->s2g1[pgc->recv_id];
464 if (pgc->recv_id > nodeid)
466 recv_index1 -= ndata;
468 recv_index1 = std::min(recv_index1, fft_end);
469 pgc->recv_index0 = fft_start;
470 pgc->recv_nindex = std::max(0, recv_index1 - pgc->recv_index0);
474 /* Communicate the buffer sizes to receive */
476 for (size_t b = 0; b < ol->comm_data.size(); b++)
478 MPI_Sendrecv(&ol->send_size,
481 ol->comm_data[b].send_id,
483 &ol->comm_data[b].recv_size,
486 ol->comm_data[b].recv_id,
493 /* For non-divisible grid we need pme_order iso pme_order-1 */
494 ol->sendbuf.resize(norder * commplainsize);
495 ol->recvbuf.resize(norder * commplainsize);
498 int minimalPmeGridSize(int pmeOrder)
500 /* The actual grid size limitations are:
501 * serial: >= pme_order
502 * DD, no OpenMP: >= 2*(pme_order - 1)
503 * DD, OpenMP: >= pme_order + 1
504 * But we use the maximum for simplicity since in practice there is not
505 * much performance difference between pme_order and 2*(pme_order -1).
507 int minimalSize = 2 * (pmeOrder - 1);
509 GMX_RELEASE_ASSERT(pmeOrder >= 3, "pmeOrder has to be >= 3");
510 GMX_RELEASE_ASSERT(minimalSize >= pmeOrder + 1, "The grid size should be >= pmeOrder + 1");
515 bool gmx_pme_check_restrictions(int pme_order, int nkx, int nky, int nkz, int numPmeDomainsAlongX, bool useThreads, bool errorsAreFatal)
517 if (pme_order > PME_ORDER_MAX)
524 std::string message = gmx::formatString(
525 "pme_order (%d) is larger than the maximum allowed value (%d). Modify and "
526 "recompile the code if you really need such a high order.",
529 GMX_THROW(gmx::InconsistentInputError(message));
532 const int minGridSize = minimalPmeGridSize(pme_order);
533 if (nkx < minGridSize || nky < minGridSize || nkz < minGridSize)
539 std::string message =
540 gmx::formatString("The PME grid sizes need to be >= 2*(pme_order-1) (%d)", minGridSize);
541 GMX_THROW(gmx::InconsistentInputError(message));
544 /* Check for a limitation of the (current) sum_fftgrid_dd code.
545 * We only allow multiple communication pulses in dim 1, not in dim 0.
548 && (nkx < numPmeDomainsAlongX * pme_order && nkx != numPmeDomainsAlongX * (pme_order - 1)))
555 "The number of PME grid lines per rank along x is %g. But when using OpenMP "
556 "threads, the number of grid lines per rank along x should be >= pme_order (%d) "
557 "or = pmeorder-1. To resolve this issue, use fewer ranks along x (and possibly "
558 "more along y and/or z) by specifying -dd manually.",
559 nkx / static_cast<double>(numPmeDomainsAlongX),
566 /*! \brief Round \p enumerator */
567 static int div_round_up(int enumerator, int denominator)
569 return (enumerator + denominator - 1) / denominator;
572 gmx_pme_t* gmx_pme_init(const t_commrec* cr,
573 const NumPmeDomains& numPmeDomains,
574 const t_inputrec* ir,
575 gmx_bool bFreeEnergy_q,
576 gmx_bool bFreeEnergy_lj,
577 gmx_bool bReproducible,
583 const DeviceContext* deviceContext,
584 const DeviceStream* deviceStream,
585 const PmeGpuProgram* pmeGpuProgram,
586 const gmx::MDLogger& mdlog)
588 int use_threads, sum_use_threads, i;
593 fprintf(debug, "Creating PME data structures.\n");
596 gmx::unique_cptr<gmx_pme_t, gmx_pme_destroy> pme(new gmx_pme_t());
598 pme->sum_qgrid_tmp = nullptr;
599 pme->sum_qgrid_dd_tmp = nullptr;
606 pme->nnodes_major = numPmeDomains.x;
607 pme->nnodes_minor = numPmeDomains.y;
609 if (numPmeDomains.x * numPmeDomains.y > 1)
611 pme->mpi_comm = cr->mpi_comm_mygroup;
614 MPI_Comm_rank(pme->mpi_comm, &pme->nodeid);
615 MPI_Comm_size(pme->mpi_comm, &pme->nnodes);
617 if (pme->nnodes != numPmeDomains.x * numPmeDomains.y)
619 gmx_incons("PME rank count mismatch");
624 pme->mpi_comm = MPI_COMM_NULL;
627 if (pme->nnodes == 1)
629 pme->mpi_comm_d[0] = MPI_COMM_NULL;
630 pme->mpi_comm_d[1] = MPI_COMM_NULL;
632 pme->nodeid_major = 0;
633 pme->nodeid_minor = 0;
637 if (numPmeDomains.y == 1)
639 pme->mpi_comm_d[0] = pme->mpi_comm;
640 pme->mpi_comm_d[1] = MPI_COMM_NULL;
642 pme->nodeid_major = pme->nodeid;
643 pme->nodeid_minor = 0;
645 else if (numPmeDomains.x == 1)
647 pme->mpi_comm_d[0] = MPI_COMM_NULL;
648 pme->mpi_comm_d[1] = pme->mpi_comm;
650 pme->nodeid_major = 0;
651 pme->nodeid_minor = pme->nodeid;
655 if (pme->nnodes % numPmeDomains.x != 0)
658 "For 2D PME decomposition, #PME ranks must be divisible by the number of "
664 MPI_Comm_split(pme->mpi_comm,
665 pme->nodeid % numPmeDomains.y,
667 &pme->mpi_comm_d[0]); /* My communicator along major dimension */
668 MPI_Comm_split(pme->mpi_comm,
669 pme->nodeid / numPmeDomains.y,
671 &pme->mpi_comm_d[1]); /* My communicator along minor dimension */
673 MPI_Comm_rank(pme->mpi_comm_d[0], &pme->nodeid_major);
674 MPI_Comm_size(pme->mpi_comm_d[0], &pme->nnodes_major);
675 MPI_Comm_rank(pme->mpi_comm_d[1], &pme->nodeid_minor);
676 MPI_Comm_size(pme->mpi_comm_d[1], &pme->nnodes_minor);
680 // cr is always initialized if there is a a PP rank, so we can safely assume
681 // that when it is not, like in ewald tests, we not on a PP rank.
682 pme->bPPnode = ((cr != nullptr && cr->duty != 0) && thisRankHasDuty(cr, DUTY_PP));
684 pme->nthread = nthread;
686 /* Check if any of the PME MPI ranks uses threads */
687 use_threads = (pme->nthread > 1 ? 1 : 0);
691 MPI_Allreduce(&use_threads, &sum_use_threads, 1, MPI_INT, MPI_SUM, pme->mpi_comm);
696 sum_use_threads = use_threads;
698 pme->bUseThreads = (sum_use_threads > 0);
700 if (ir->pbcType == PbcType::Screw)
702 gmx_fatal(FARGS, "pme does not (yet) work with pbc = screw");
706 * It is likely that the current gmx_pme_do() routine supports calculating
707 * only Coulomb or LJ while gmx_pme_init() configures for both,
708 * but that has never been tested.
709 * It is likely that the current gmx_pme_do() routine supports calculating,
710 * not calculating free-energy for Coulomb and/or LJ while gmx_pme_init()
711 * configures with free-energy, but that has never been tested.
713 pme->doCoulomb = EEL_PME(ir->coulombtype);
714 pme->doLJ = EVDW_PME(ir->vdwtype);
715 pme->bFEP_q = ((ir->efep != FreeEnergyPerturbationType::No) && bFreeEnergy_q);
716 pme->bFEP_lj = ((ir->efep != FreeEnergyPerturbationType::No) && bFreeEnergy_lj);
717 pme->bFEP = (pme->bFEP_q || pme->bFEP_lj);
721 pme->bP3M = (ir->coulombtype == CoulombInteractionType::P3mAD || getenv("GMX_PME_P3M") != nullptr);
722 pme->pme_order = ir->pme_order;
723 pme->ewaldcoeff_q = ewaldcoeff_q;
724 pme->ewaldcoeff_lj = ewaldcoeff_lj;
726 /* Always constant electrostatics coefficients */
727 pme->epsilon_r = ir->epsilon_r;
729 /* Always constant LJ coefficients */
730 pme->ljpme_combination_rule = ir->ljpme_combination_rule;
732 // The box requires scaling with nwalls = 2, we store that condition as well
733 // as the scaling factor
734 delete pme->boxScaler;
735 pme->boxScaler = new EwaldBoxZScaler(*ir);
737 /* If we violate restrictions, generate a fatal error here */
738 gmx_pme_check_restrictions(
739 pme->pme_order, pme->nkx, pme->nky, pme->nkz, pme->nnodes_major, pme->bUseThreads, true);
746 MPI_Type_contiguous(DIM, GMX_MPI_REAL, &(pme->rvec_mpi));
747 MPI_Type_commit(&(pme->rvec_mpi));
750 /* Note that the coefficient spreading and force gathering, which usually
751 * takes about the same amount of time as FFT+solve_pme,
752 * is always fully load balanced
753 * (unless the coefficient distribution is inhomogeneous).
756 imbal = estimate_pme_load_imbalance(pme.get());
757 if (imbal >= 1.2 && pme->nodeid_major == 0 && pme->nodeid_minor == 0)
759 GMX_LOG(mdlog.warning)
761 .appendTextFormatted(
762 "NOTE: The load imbalance in PME FFT and solve is %d%%.\n"
763 " For optimal PME load balancing\n"
764 " PME grid_x (%d) and grid_y (%d) should be divisible by "
767 " and PME grid_y (%d) and grid_z (%d) should be divisible by "
770 gmx::roundToInt((imbal - 1) * 100),
780 /* For non-divisible grid we need pme_order iso pme_order-1 */
781 /* In sum_qgrid_dd x overlap is copied in place: take padding into account.
782 * y is always copied through a buffer: we don't need padding in z,
783 * but we do need the overlap in x because of the communication order.
785 init_overlap_comm(&pme->overlap[0],
791 (div_round_up(pme->nky, pme->nnodes_minor) + pme->pme_order)
792 * (pme->nkz + pme->pme_order - 1));
794 /* Along overlap dim 1 we can send in multiple pulses in sum_fftgrid_dd.
795 * We do this with an offset buffer of equal size, so we need to allocate
796 * extra for the offset. That's what the (+1)*pme->nkz is for.
798 init_overlap_comm(&pme->overlap[1],
804 (div_round_up(pme->nkx, pme->nnodes_major) + pme->pme_order + 1) * pme->nkz);
806 /* Double-check for a limitation of the (current) sum_fftgrid_dd code.
807 * Note that gmx_pme_check_restrictions checked for this already.
809 if (pme->bUseThreads && (pme->overlap[0].comm_data.size() > 1))
812 "More than one communication pulse required for grid overlap communication along "
813 "the major dimension while using threads");
816 snew(pme->bsp_mod[XX], pme->nkx);
817 snew(pme->bsp_mod[YY], pme->nky);
818 snew(pme->bsp_mod[ZZ], pme->nkz);
820 pme->gpu = pmeGpu; /* Carrying over the single GPU structure */
821 pme->runMode = runMode;
823 /* The required size of the interpolation grid, including overlap.
824 * The allocated size (pmegrid_n?) might be slightly larger.
826 pme->pmegrid_nx = pme->overlap[0].s2g1[pme->nodeid_major] - pme->overlap[0].s2g0[pme->nodeid_major];
827 pme->pmegrid_ny = pme->overlap[1].s2g1[pme->nodeid_minor] - pme->overlap[1].s2g0[pme->nodeid_minor];
828 pme->pmegrid_nz_base = pme->nkz;
829 pme->pmegrid_nz = pme->pmegrid_nz_base + pme->pme_order - 1;
830 set_grid_alignment(&pme->pmegrid_nz, pme->pme_order);
831 pme->pmegrid_start_ix = pme->overlap[0].s2g0[pme->nodeid_major];
832 pme->pmegrid_start_iy = pme->overlap[1].s2g0[pme->nodeid_minor];
833 pme->pmegrid_start_iz = 0;
835 make_gridindex_to_localindex(
836 pme->nkx, pme->pmegrid_start_ix, pme->pmegrid_nx - (pme->pme_order - 1), &pme->nnx, &pme->fshx);
837 make_gridindex_to_localindex(
838 pme->nky, pme->pmegrid_start_iy, pme->pmegrid_ny - (pme->pme_order - 1), &pme->nny, &pme->fshy);
839 make_gridindex_to_localindex(
840 pme->nkz, pme->pmegrid_start_iz, pme->pmegrid_nz_base, &pme->nnz, &pme->fshz);
842 pme->spline_work = make_pme_spline_work(pme->pme_order);
847 /* It doesn't matter if we allocate too many grids here,
848 * we only allocate and use the ones we need.
852 pme->ngrids = ((ir->ljpme_combination_rule == LongRangeVdW::LB) ? DO_Q_AND_LJ_LB : DO_Q_AND_LJ);
858 snew(pme->fftgrid, pme->ngrids);
859 snew(pme->cfftgrid, pme->ngrids);
860 snew(pme->pfft_setup, pme->ngrids);
862 for (i = 0; i < pme->ngrids; ++i)
864 if ((i < DO_Q && pme->doCoulomb && (i == 0 || bFreeEnergy_q))
865 || (i >= DO_Q && pme->doLJ
866 && (i == 2 || bFreeEnergy_lj || ir->ljpme_combination_rule == LongRangeVdW::LB)))
868 pmegrids_init(&pme->pmegrid[i],
872 pme->pmegrid_nz_base,
876 pme->overlap[0].s2g1[pme->nodeid_major]
877 - pme->overlap[0].s2g0[pme->nodeid_major + 1],
878 pme->overlap[1].s2g1[pme->nodeid_minor]
879 - pme->overlap[1].s2g0[pme->nodeid_minor + 1]);
880 /* This routine will allocate the grid data to fit the FFTs */
881 const auto allocateRealGridForGpu = (pme->runMode == PmeRunMode::Mixed)
882 ? gmx::PinningPolicy::PinnedIfSupported
883 : gmx::PinningPolicy::CannotBePinned;
884 gmx_parallel_3dfft_init(&pme->pfft_setup[i],
891 allocateRealGridForGpu);
897 /* Use plain SPME B-spline interpolation */
898 make_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
902 /* Use the P3M grid-optimized influence function */
903 make_p3m_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
906 /* Use atc[0] for spreading */
907 const int firstDimIndex = (numPmeDomains.x > 1 ? 0 : 1);
908 MPI_Comm mpiCommFirstDim = (pme->nnodes > 1 ? pme->mpi_comm_d[firstDimIndex] : MPI_COMM_NULL);
909 bool doSpread = true;
910 pme->atc.emplace_back(mpiCommFirstDim, pme->nthread, pme->pme_order, firstDimIndex, doSpread);
911 if (pme->ndecompdim >= 2)
913 const int secondDimIndex = 1;
915 pme->atc.emplace_back(pme->mpi_comm_d[1], pme->nthread, pme->pme_order, secondDimIndex, doSpread);
918 // Initial check of validity of the input for running on the GPU
919 if (pme->runMode != PmeRunMode::CPU)
921 std::string errorString;
922 bool canRunOnGpu = pme_gpu_check_restrictions(pme.get(), &errorString);
925 GMX_THROW(gmx::NotImplementedError(errorString));
927 pme_gpu_reinit(pme.get(), deviceContext, deviceStream, pmeGpuProgram);
931 GMX_ASSERT(pme->gpu == nullptr, "Should not have PME GPU object when PME is on a CPU.");
935 pme_init_all_work(&pme->solve_work, pme->nthread, pme->nkx);
937 // no exception was thrown during the init, so we hand over the PME structure handle
938 return pme.release();
941 void gmx_pme_reinit(struct gmx_pme_t** pmedata,
943 struct gmx_pme_t* pme_src,
944 const t_inputrec* ir,
945 const ivec grid_size,
949 // Create a copy of t_inputrec fields that are used in gmx_pme_init().
950 // TODO: This would be better as just copying a sub-structure that contains
951 // all the PME parameters and nothing else.
953 irc.pbcType = ir->pbcType;
954 irc.coulombtype = ir->coulombtype;
955 irc.vdwtype = ir->vdwtype;
957 irc.pme_order = ir->pme_order;
958 irc.epsilon_r = ir->epsilon_r;
959 irc.ljpme_combination_rule = ir->ljpme_combination_rule;
960 irc.nkx = grid_size[XX];
961 irc.nky = grid_size[YY];
962 irc.nkz = grid_size[ZZ];
966 // This is reinit. Any logging should have been done at first init.
967 // Here we should avoid writing notes for settings the user did not
969 const gmx::MDLogger dummyLogger;
970 GMX_ASSERT(pmedata, "Invalid PME pointer");
971 NumPmeDomains numPmeDomains = { pme_src->nnodes_major, pme_src->nnodes_minor };
972 *pmedata = gmx_pme_init(cr,
987 /* When running PME on the CPU not using domain decomposition,
988 * the atom data is allocated once only in gmx_pme_(re)init().
990 if (!pme_src->gpu && pme_src->nnodes == 1)
992 gmx_pme_reinit_atoms(*pmedata, pme_src->atc[0].numAtoms(), nullptr, nullptr);
994 // TODO this is mostly passing around current values
996 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
998 /* We can easily reuse the allocated pme grids in pme_src */
999 reuse_pmegrids(&pme_src->pmegrid[PME_GRID_QA], &(*pmedata)->pmegrid[PME_GRID_QA]);
1000 /* We would like to reuse the fft grids, but that's harder */
1003 void gmx_pme_calc_energy(gmx_pme_t* pme, gmx::ArrayRef<const gmx::RVec> x, gmx::ArrayRef<const real> q, real* V)
1007 if (pme->nnodes > 1)
1009 gmx_incons("gmx_pme_calc_energy called in parallel");
1013 gmx_incons("gmx_pme_calc_energy with free energy");
1016 if (!pme->atc_energy)
1018 pme->atc_energy = std::make_unique<PmeAtomComm>(MPI_COMM_NULL, 1, pme->pme_order, 0, true);
1020 PmeAtomComm* atc = pme->atc_energy.get();
1021 atc->setNumAtoms(x.ssize());
1023 atc->coefficient = q;
1025 /* We only use the A-charges grid */
1026 grid = &pme->pmegrid[PME_GRID_QA];
1028 /* Only calculate the spline coefficients, don't actually spread */
1029 spread_on_grid(pme, atc, nullptr, TRUE, FALSE, pme->fftgrid[PME_GRID_QA], FALSE, PME_GRID_QA);
1031 *V = gather_energy_bsplines(pme, grid->grid.grid, atc);
1034 /*! \brief Calculate initial Lorentz-Berthelot coefficients for LJ-PME */
1035 static void calc_initial_lb_coeffs(gmx::ArrayRef<real> coefficient,
1036 gmx::ArrayRef<const real> local_c6,
1037 gmx::ArrayRef<const real> local_sigma)
1039 for (gmx::index i = 0; i < coefficient.ssize(); ++i)
1041 real sigma4 = local_sigma[i];
1042 sigma4 = sigma4 * sigma4;
1043 sigma4 = sigma4 * sigma4;
1044 coefficient[i] = local_c6[i] / sigma4;
1048 /*! \brief Calculate next Lorentz-Berthelot coefficients for LJ-PME */
1049 static void calc_next_lb_coeffs(gmx::ArrayRef<real> coefficient, gmx::ArrayRef<const real> local_sigma)
1051 for (gmx::index i = 0; i < coefficient.ssize(); ++i)
1053 coefficient[i] *= local_sigma[i];
1057 int gmx_pme_do(struct gmx_pme_t* pme,
1058 gmx::ArrayRef<const gmx::RVec> coordinates,
1059 gmx::ArrayRef<gmx::RVec> forces,
1060 gmx::ArrayRef<const real> chargeA,
1061 gmx::ArrayRef<const real> chargeB,
1062 gmx::ArrayRef<const real> c6A,
1063 gmx::ArrayRef<const real> c6B,
1064 gmx::ArrayRef<const real> sigmaA,
1065 gmx::ArrayRef<const real> sigmaB,
1067 const t_commrec* cr,
1071 gmx_wallcycle* wcycle,
1080 const gmx::StepWorkload& stepWork)
1082 GMX_ASSERT(pme->runMode == PmeRunMode::CPU,
1083 "gmx_pme_do should not be called on the GPU PME run.");
1085 int d, npme, grid_index, max_grid_index;
1086 PmeAtomComm& atc = pme->atc[0];
1087 pmegrids_t* pmegrid = nullptr;
1088 real* grid = nullptr;
1089 gmx::ArrayRef<const real> coefficient;
1090 PmeOutput output[2]; // The second is used for the B state with FEP
1093 gmx_parallel_3dfft_t pfft_setup;
1095 t_complex* cfftgrid;
1097 gmx_bool bFirst, bDoSplines;
1099 int fep_states_lj = pme->bFEP_lj ? 2 : 1;
1100 // There's no support for computing energy without virial, or vice versa
1101 const bool computeEnergyAndVirial = (stepWork.computeEnergy || stepWork.computeVirial);
1103 /* We could be passing lambda!=0 while no q or LJ is actually perturbed */
1113 assert(pme->nnodes > 0);
1114 assert(pme->nnodes == 1 || pme->ndecompdim > 0);
1116 if (pme->nnodes > 1)
1118 atc.pd.resize(coordinates.ssize());
1119 for (int d = pme->ndecompdim - 1; d >= 0; d--)
1121 PmeAtomComm& atc = pme->atc[d];
1122 atc.maxshift = (atc.dimind == 0 ? maxshift_x : maxshift_y);
1127 GMX_ASSERT(coordinates.ssize() == atc.numAtoms(), "We expect atc.numAtoms() coordinates");
1128 GMX_ASSERT(forces.ssize() >= atc.numAtoms(),
1129 "We need a force buffer with at least atc.numAtoms() elements");
1131 atc.x = coordinates;
1136 pme->boxScaler->scaleBox(box, scaledBox);
1138 gmx::invertBoxMatrix(scaledBox, pme->recipbox);
1141 /* For simplicity, we construct the splines for all particles if
1142 * more than one PME calculations is needed. Some optimization
1143 * could be done by keeping track of which atoms have splines
1144 * constructed, and construct new splines on each pass for atoms
1145 * that don't yet have them.
1148 bDoSplines = pme->bFEP || (pme->doCoulomb && pme->doLJ);
1150 /* We need a maximum of four separate PME calculations:
1151 * grid_index=0: Coulomb PME with charges from state A
1152 * grid_index=1: Coulomb PME with charges from state B
1153 * grid_index=2: LJ PME with C6 from state A
1154 * grid_index=3: LJ PME with C6 from state B
1155 * For Lorentz-Berthelot combination rules, a separate loop is used to
1156 * calculate all the terms
1159 /* If we are doing LJ-PME with LB, we only do Q here */
1160 max_grid_index = (pme->ljpme_combination_rule == LongRangeVdW::LB) ? DO_Q : DO_Q_AND_LJ;
1162 for (grid_index = 0; grid_index < max_grid_index; ++grid_index)
1164 /* Check if we should do calculations at this grid_index
1165 * If grid_index is odd we should be doing FEP
1166 * If grid_index < 2 we should be doing electrostatic PME
1167 * If grid_index >= 2 we should be doing LJ-PME
1169 if ((grid_index < DO_Q && (!pme->doCoulomb || (grid_index == 1 && !pme->bFEP_q)))
1170 || (grid_index >= DO_Q && (!pme->doLJ || (grid_index == 3 && !pme->bFEP_lj))))
1174 /* Unpack structure */
1175 pmegrid = &pme->pmegrid[grid_index];
1176 fftgrid = pme->fftgrid[grid_index];
1177 cfftgrid = pme->cfftgrid[grid_index];
1178 pfft_setup = pme->pfft_setup[grid_index];
1181 case 0: coefficient = chargeA; break;
1182 case 1: coefficient = chargeB; break;
1183 case 2: coefficient = c6A; break;
1184 case 3: coefficient = c6B; break;
1187 grid = pmegrid->grid.grid;
1191 fprintf(debug, "PME: number of ranks = %d, rank = %d\n", cr->nnodes, cr->nodeid);
1192 fprintf(debug, "Grid = %p\n", static_cast<void*>(grid));
1193 if (grid == nullptr)
1195 gmx_fatal(FARGS, "No grid!");
1199 if (pme->nnodes == 1)
1201 atc.coefficient = coefficient;
1205 wallcycle_start(wcycle, ewcPME_REDISTXF);
1206 do_redist_pos_coeffs(pme, cr, bFirst, coordinates, coefficient);
1208 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1213 fprintf(debug, "Rank= %6d, pme local particles=%6d\n", cr->nodeid, atc.numAtoms());
1216 wallcycle_start(wcycle, ewcPME_SPREAD);
1218 /* Spread the coefficients on a grid */
1219 spread_on_grid(pme, &atc, pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1223 inc_nrnb(nrnb, eNR_WEIGHTS, DIM * atc.numAtoms());
1225 inc_nrnb(nrnb, eNR_SPREADBSP, pme->pme_order * pme->pme_order * pme->pme_order * atc.numAtoms());
1227 if (!pme->bUseThreads)
1229 wrap_periodic_pmegrid(pme, grid);
1231 /* sum contributions to local grid from other nodes */
1232 if (pme->nnodes > 1)
1234 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1237 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1240 wallcycle_stop(wcycle, ewcPME_SPREAD);
1242 /* TODO If the OpenMP and single-threaded implementations
1243 converge, then spread_on_grid() and
1244 copy_pmegrid_to_fftgrid() will perhaps live in the same
1248 /* Here we start a large thread parallel region */
1249 #pragma omp parallel num_threads(pme->nthread) private(thread)
1253 thread = gmx_omp_get_thread_num();
1259 wallcycle_start(wcycle, ewcPME_FFT);
1261 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX, thread, wcycle);
1264 wallcycle_stop(wcycle, ewcPME_FFT);
1267 /* solve in k-space for our local cells */
1270 wallcycle_start(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
1272 if (grid_index < DO_Q)
1274 loop_count = solve_pme_yzx(pme,
1276 scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ],
1277 computeEnergyAndVirial,
1284 solve_pme_lj_yzx(pme,
1287 scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ],
1288 computeEnergyAndVirial,
1295 wallcycle_stop(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
1296 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1302 wallcycle_start(wcycle, ewcPME_FFT);
1304 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL, thread, wcycle);
1307 wallcycle_stop(wcycle, ewcPME_FFT);
1310 if (pme->nodeid == 0)
1312 real ntot = pme->nkx * pme->nky * pme->nkz;
1313 npme = static_cast<int>(ntot * std::log(ntot) / std::log(2.0));
1314 inc_nrnb(nrnb, eNR_FFT, 2 * npme);
1317 /* Note: this wallcycle region is closed below
1318 outside an OpenMP region, so take care if
1319 refactoring code here. */
1320 wallcycle_start(wcycle, ewcPME_GATHER);
1323 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1325 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1327 /* End of thread parallel section.
1328 * With MPI we have to synchronize here before gmx_sum_qgrid_dd.
1331 /* distribute local grid to all nodes */
1332 if (pme->nnodes > 1)
1334 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1337 unwrap_periodic_pmegrid(pme, grid);
1339 if (stepWork.computeForces)
1341 /* interpolate forces for our local atoms */
1344 /* If we are running without parallelization,
1345 * atc->f is the actual force array, not a buffer,
1346 * therefore we should not clear it.
1348 lambda = grid_index < DO_Q ? lambda_q : lambda_lj;
1349 bClearF = (bFirst && PAR(cr));
1350 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1351 for (thread = 0; thread < pme->nthread; thread++)
1355 gather_f_bsplines(pme,
1359 &atc.spline[thread],
1360 pme->bFEP ? (grid_index % 2 == 0 ? 1.0 - lambda : lambda) : 1.0);
1362 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1366 inc_nrnb(nrnb, eNR_GATHERFBSP, pme->pme_order * pme->pme_order * pme->pme_order * atc.numAtoms());
1367 /* Note: this wallcycle region is opened above inside an OpenMP
1368 region, so take care if refactoring code here. */
1369 wallcycle_stop(wcycle, ewcPME_GATHER);
1372 if (computeEnergyAndVirial)
1374 /* This should only be called on the master thread
1375 * and after the threads have synchronized.
1379 get_pme_ener_vir_q(pme->solve_work, pme->nthread, &output[grid_index % 2]);
1383 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output[grid_index % 2]);
1387 } /* of grid_index-loop */
1389 /* For Lorentz-Berthelot combination rules in LJ-PME, we need to calculate
1392 if (pme->doLJ && pme->ljpme_combination_rule == LongRangeVdW::LB)
1394 /* Loop over A- and B-state if we are doing FEP */
1395 for (fep_state = 0; fep_state < fep_states_lj; ++fep_state)
1397 std::vector<real> local_c6;
1398 std::vector<real> local_sigma;
1399 gmx::ArrayRef<const real> RedistC6;
1400 gmx::ArrayRef<const real> RedistSigma;
1401 gmx::ArrayRef<real> coefficientBuffer;
1402 if (pme->nnodes == 1)
1404 pme->lb_buf1.resize(atc.numAtoms());
1405 coefficientBuffer = pme->lb_buf1;
1409 local_c6.assign(c6A.begin(), c6A.end());
1410 local_sigma.assign(sigmaA.begin(), sigmaA.end());
1413 local_c6.assign(c6B.begin(), c6B.end());
1414 local_sigma.assign(sigmaB.begin(), sigmaB.end());
1416 default: gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1421 coefficientBuffer = atc.coefficientBuffer;
1426 RedistSigma = sigmaA;
1430 RedistSigma = sigmaB;
1432 default: gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1434 wallcycle_start(wcycle, ewcPME_REDISTXF);
1436 do_redist_pos_coeffs(pme, cr, bFirst, coordinates, RedistC6);
1437 pme->lb_buf1.resize(atc.numAtoms());
1438 pme->lb_buf2.resize(atc.numAtoms());
1439 local_c6.assign(pme->lb_buf1.begin(), pme->lb_buf1.end());
1440 for (int i = 0; i < atc.numAtoms(); ++i)
1442 local_c6[i] = atc.coefficient[i];
1445 do_redist_pos_coeffs(pme, cr, FALSE, coordinates, RedistSigma);
1446 local_sigma.assign(pme->lb_buf2.begin(), pme->lb_buf2.end());
1447 for (int i = 0; i < atc.numAtoms(); ++i)
1449 local_sigma[i] = atc.coefficient[i];
1452 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1454 atc.coefficient = coefficientBuffer;
1455 calc_initial_lb_coeffs(coefficientBuffer, local_c6, local_sigma);
1457 /*Seven terms in LJ-PME with LB, grid_index < 2 reserved for electrostatics*/
1458 for (grid_index = 2; grid_index < 9; ++grid_index)
1460 /* Unpack structure */
1461 pmegrid = &pme->pmegrid[grid_index];
1462 fftgrid = pme->fftgrid[grid_index];
1463 pfft_setup = pme->pfft_setup[grid_index];
1464 calc_next_lb_coeffs(coefficientBuffer, local_sigma);
1465 grid = pmegrid->grid.grid;
1467 wallcycle_start(wcycle, ewcPME_SPREAD);
1468 /* Spread the c6 on a grid */
1469 spread_on_grid(pme, &atc, pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1473 inc_nrnb(nrnb, eNR_WEIGHTS, DIM * atc.numAtoms());
1478 pme->pme_order * pme->pme_order * pme->pme_order * atc.numAtoms());
1479 if (pme->nthread == 1)
1481 wrap_periodic_pmegrid(pme, grid);
1482 /* sum contributions to local grid from other nodes */
1483 if (pme->nnodes > 1)
1485 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1487 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1489 wallcycle_stop(wcycle, ewcPME_SPREAD);
1491 /*Here we start a large thread parallel region*/
1492 #pragma omp parallel num_threads(pme->nthread) private(thread)
1496 thread = gmx_omp_get_thread_num();
1500 wallcycle_start(wcycle, ewcPME_FFT);
1503 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX, thread, wcycle);
1506 wallcycle_stop(wcycle, ewcPME_FFT);
1509 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1513 /* solve in k-space for our local cells */
1514 #pragma omp parallel num_threads(pme->nthread) private(thread)
1519 thread = gmx_omp_get_thread_num();
1522 wallcycle_start(wcycle, ewcLJPME);
1526 solve_pme_lj_yzx(pme,
1529 scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ],
1530 computeEnergyAndVirial,
1535 wallcycle_stop(wcycle, ewcLJPME);
1536 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1539 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1542 if (computeEnergyAndVirial)
1544 /* This should only be called on the master thread and
1545 * after the threads have synchronized.
1547 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output[fep_state]);
1550 bFirst = !pme->doCoulomb;
1551 calc_initial_lb_coeffs(coefficientBuffer, local_c6, local_sigma);
1552 for (grid_index = 8; grid_index >= 2; --grid_index)
1554 /* Unpack structure */
1555 pmegrid = &pme->pmegrid[grid_index];
1556 fftgrid = pme->fftgrid[grid_index];
1557 pfft_setup = pme->pfft_setup[grid_index];
1558 grid = pmegrid->grid.grid;
1559 calc_next_lb_coeffs(coefficientBuffer, local_sigma);
1560 #pragma omp parallel num_threads(pme->nthread) private(thread)
1564 thread = gmx_omp_get_thread_num();
1568 wallcycle_start(wcycle, ewcPME_FFT);
1571 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL, thread, wcycle);
1574 wallcycle_stop(wcycle, ewcPME_FFT);
1577 if (pme->nodeid == 0)
1579 real ntot = pme->nkx * pme->nky * pme->nkz;
1580 npme = static_cast<int>(ntot * std::log(ntot) / std::log(2.0));
1581 inc_nrnb(nrnb, eNR_FFT, 2 * npme);
1583 wallcycle_start(wcycle, ewcPME_GATHER);
1586 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1588 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1589 } /*#pragma omp parallel*/
1591 /* distribute local grid to all nodes */
1592 if (pme->nnodes > 1)
1594 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1597 unwrap_periodic_pmegrid(pme, grid);
1599 if (stepWork.computeForces)
1601 /* interpolate forces for our local atoms */
1602 bClearF = (bFirst && PAR(cr));
1603 scale = pme->bFEP ? (fep_state < 1 ? 1.0 - lambda_lj : lambda_lj) : 1.0;
1604 scale *= lb_scale_factor[grid_index - 2];
1606 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1607 for (thread = 0; thread < pme->nthread; thread++)
1612 pme, grid, bClearF, &pme->atc[0], &pme->atc[0].spline[thread], scale);
1614 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1620 pme->pme_order * pme->pme_order * pme->pme_order * pme->atc[0].numAtoms());
1622 wallcycle_stop(wcycle, ewcPME_GATHER);
1625 } /* for (grid_index = 8; grid_index >= 2; --grid_index) */
1626 } /* for (fep_state = 0; fep_state < fep_states_lj; ++fep_state) */
1627 } /* if (pme->doLJ && pme->ljpme_combination_rule == LongRangeVdW::LB) */
1629 if (stepWork.computeForces && pme->nnodes > 1)
1631 wallcycle_start(wcycle, ewcPME_REDISTXF);
1632 for (d = 0; d < pme->ndecompdim; d++)
1634 gmx::ArrayRef<gmx::RVec> forcesRef;
1635 if (d == pme->ndecompdim - 1)
1637 const size_t numAtoms = coordinates.size();
1638 GMX_ASSERT(forces.size() >= numAtoms, "Need at least numAtoms forces");
1639 forcesRef = forces.subArray(0, numAtoms);
1643 forcesRef = pme->atc[d + 1].f;
1645 if (DOMAINDECOMP(cr))
1647 dd_pmeredist_f(pme, &pme->atc[d], forcesRef, d == pme->ndecompdim - 1 && pme->bPPnode);
1651 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1654 if (computeEnergyAndVirial)
1660 *energy_q = output[0].coulombEnergy_;
1661 m_add(vir_q, output[0].coulombVirial_, vir_q);
1665 *energy_q = (1.0 - lambda_q) * output[0].coulombEnergy_ + lambda_q * output[1].coulombEnergy_;
1666 *dvdlambda_q += output[1].coulombEnergy_ - output[0].coulombEnergy_;
1667 for (int i = 0; i < DIM; i++)
1669 for (int j = 0; j < DIM; j++)
1671 vir_q[i][j] += (1.0 - lambda_q) * output[0].coulombVirial_[i][j]
1672 + lambda_q * output[1].coulombVirial_[i][j];
1678 fprintf(debug, "Electrostatic PME mesh energy: %g\n", *energy_q);
1690 *energy_lj = output[0].lennardJonesEnergy_;
1691 m_add(vir_lj, output[0].lennardJonesVirial_, vir_lj);
1695 *energy_lj = (1.0 - lambda_lj) * output[0].lennardJonesEnergy_
1696 + lambda_lj * output[1].lennardJonesEnergy_;
1697 *dvdlambda_lj += output[1].lennardJonesEnergy_ - output[0].lennardJonesEnergy_;
1698 for (int i = 0; i < DIM; i++)
1700 for (int j = 0; j < DIM; j++)
1702 vir_lj[i][j] += (1.0 - lambda_lj) * output[0].lennardJonesVirial_[i][j]
1703 + lambda_lj * output[1].lennardJonesVirial_[i][j];
1709 fprintf(debug, "Lennard-Jones PME mesh energy: %g\n", *energy_lj);
1720 void gmx_pme_destroy(gmx_pme_t* pme)
1727 delete pme->boxScaler;
1736 for (int i = 0; i < pme->ngrids; ++i)
1738 pmegrids_destroy(&pme->pmegrid[i]);
1740 if (pme->pfft_setup)
1742 for (int i = 0; i < pme->ngrids; ++i)
1744 gmx_parallel_3dfft_destroy(pme->pfft_setup[i]);
1747 sfree(pme->fftgrid);
1748 sfree(pme->cfftgrid);
1749 sfree(pme->pfft_setup);
1751 for (int i = 0; i < DIM; i++)
1753 sfree(pme->bsp_mod[i]);
1759 if (pme->solve_work)
1761 pme_free_all_work(&pme->solve_work, pme->nthread);
1764 sfree(pme->sum_qgrid_tmp);
1765 sfree(pme->sum_qgrid_dd_tmp);
1767 destroy_pme_spline_work(pme->spline_work);
1769 if (pme->gpu != nullptr)
1771 pme_gpu_destroy(pme->gpu);
1777 void gmx_pme_reinit_atoms(gmx_pme_t* pme, const int numAtoms, const real* chargesA, const real* chargesB)
1779 if (pme->gpu != nullptr)
1781 GMX_ASSERT(!(pme->bFEP_q && chargesB == nullptr),
1782 "B state charges must be specified if running Coulomb FEP on the GPU");
1783 pme_gpu_reinit_atoms(pme->gpu, numAtoms, chargesA, pme->bFEP_q ? chargesB : nullptr);
1787 pme->atc[0].setNumAtoms(numAtoms);
1788 // TODO: set the charges here as well
1792 bool gmx_pme_grid_matches(const gmx_pme_t& pme, const ivec grid_size)
1794 return (pme.nkx == grid_size[XX] && pme.nky == grid_size[YY] && pme.nkz == grid_size[ZZ]);
1797 void gmx::SeparatePmeRanksPermitted::disablePmeRanks(const std::string& reason)
1799 permitSeparatePmeRanks_ = false;
1801 if (!reason.empty())
1803 reasons_.push_back(reason);
1807 bool gmx::SeparatePmeRanksPermitted::permitSeparatePmeRanks() const
1809 return permitSeparatePmeRanks_;
1812 std::string gmx::SeparatePmeRanksPermitted::reasonsWhyDisabled() const
1814 return joinStrings(reasons_, "; ");