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,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief This file contains function definitions necessary for
40 * computing energies and forces for the PME long-ranged part (Coulomb
43 * \author Erik Lindahl <erik@kth.se>
44 * \author Berk Hess <hess@kth.se>
45 * \ingroup module_ewald
47 /* IMPORTANT FOR DEVELOPERS:
49 * Triclinic pme stuff isn't entirely trivial, and we've experienced
50 * some bugs during development (many of them due to me). To avoid
51 * this in the future, please check the following things if you make
52 * changes in this file:
54 * 1. You should obtain identical (at least to the PME precision)
55 * energies, forces, and virial for
56 * a rectangular box and a triclinic one where the z (or y) axis is
57 * tilted a whole box side. For instance you could use these boxes:
59 * rectangular triclinic
64 * 2. You should check the energy conservation in a triclinic box.
66 * It might seem an overkill, but better safe than sorry.
85 #include "gromacs/domdec/domdec.h"
86 #include "gromacs/ewald/ewald-utils.h"
87 #include "gromacs/fft/parallel_3dfft.h"
88 #include "gromacs/fileio/pdbio.h"
89 #include "gromacs/gmxlib/network.h"
90 #include "gromacs/gmxlib/nrnb.h"
91 #include "gromacs/hardware/hw_info.h"
92 #include "gromacs/math/gmxcomplex.h"
93 #include "gromacs/math/invertmatrix.h"
94 #include "gromacs/math/units.h"
95 #include "gromacs/math/vec.h"
96 #include "gromacs/math/vectypes.h"
97 #include "gromacs/mdtypes/commrec.h"
98 #include "gromacs/mdtypes/forcerec.h"
99 #include "gromacs/mdtypes/inputrec.h"
100 #include "gromacs/mdtypes/md_enums.h"
101 #include "gromacs/pbcutil/pbc.h"
102 #include "gromacs/timing/cyclecounter.h"
103 #include "gromacs/timing/wallcycle.h"
104 #include "gromacs/timing/walltime_accounting.h"
105 #include "gromacs/topology/topology.h"
106 #include "gromacs/utility/basedefinitions.h"
107 #include "gromacs/utility/cstringutil.h"
108 #include "gromacs/utility/exceptions.h"
109 #include "gromacs/utility/fatalerror.h"
110 #include "gromacs/utility/gmxmpi.h"
111 #include "gromacs/utility/gmxomp.h"
112 #include "gromacs/utility/logger.h"
113 #include "gromacs/utility/real.h"
114 #include "gromacs/utility/smalloc.h"
115 #include "gromacs/utility/stringutil.h"
116 #include "gromacs/utility/unique_cptr.h"
118 #include "calculate-spline-moduli.h"
119 #include "pme-gather.h"
120 #include "pme-gpu-internal.h"
121 #include "pme-grid.h"
122 #include "pme-internal.h"
123 #include "pme-redistribute.h"
124 #include "pme-solve.h"
125 #include "pme-spline-work.h"
126 #include "pme-spread.h"
128 /*! \brief Help build a descriptive message in \c error if there are
129 * \c errorReasons why PME on GPU is not supported.
131 * \returns Whether the lack of errorReasons indicate there is support. */
133 addMessageIfNotSupported(const std::list<std::string> &errorReasons,
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");
162 if (GMX_GPU == GMX_GPU_NONE)
164 errorReasons.emplace_back("a non-GPU build");
166 return addMessageIfNotSupported(errorReasons, error);
169 bool pme_gpu_supports_hardware(const gmx_hw_info_t &hwinfo,
172 std::list<std::string> errorReasons;
173 if (GMX_GPU == GMX_GPU_OPENCL)
175 if (!areAllGpuDevicesFromAmd(hwinfo.gpu_info))
177 errorReasons.emplace_back("non-AMD devices");
180 return addMessageIfNotSupported(errorReasons, error);
183 bool pme_gpu_supports_input(const t_inputrec &ir, const gmx_mtop_t &mtop, std::string *error)
185 std::list<std::string> errorReasons;
186 if (!EEL_PME(ir.coulombtype))
188 errorReasons.emplace_back("systems that do not use PME for electrostatics");
190 if (ir.pme_order != 4)
192 errorReasons.emplace_back("interpolation orders other than 4");
194 if (ir.efep != efepNO)
196 if (gmx_mtop_has_perturbed_charges(mtop))
198 errorReasons.emplace_back("free energy calculations with perturbed charges (multiple grids)");
201 if (EVDW_PME(ir.vdwtype))
203 errorReasons.emplace_back("Lennard-Jones PME");
205 if (ir.cutoff_scheme == ecutsGROUP)
207 errorReasons.emplace_back("group cutoff scheme");
209 if (!EI_DYNAMICS(ir.eI))
211 errorReasons.emplace_back("not a dynamical integrator");
213 return addMessageIfNotSupported(errorReasons, error);
216 /*! \brief \libinternal
217 * Finds out if PME with given inputs is possible to run on GPU.
218 * This function is an internal final check, validating the whole PME structure on creation,
219 * but it still duplicates the preliminary checks from the above (externally exposed) pme_gpu_supports_input() - just in case.
221 * \param[in] pme The PME structure.
222 * \param[out] error The error message if the input is not supported on GPU.
223 * \returns True if this PME input is possible to run on GPU, false otherwise.
225 static bool pme_gpu_check_restrictions(const gmx_pme_t *pme, std::string *error)
227 std::list<std::string> errorReasons;
228 if (pme->nnodes != 1)
230 errorReasons.emplace_back("PME decomposition");
232 if (pme->pme_order != 4)
234 errorReasons.emplace_back("interpolation orders other than 4");
238 errorReasons.emplace_back("free energy calculations (multiple grids)");
242 errorReasons.emplace_back("Lennard-Jones PME");
246 errorReasons.emplace_back("double precision");
248 if (GMX_GPU == GMX_GPU_NONE)
250 errorReasons.emplace_back("non-GPU build of GROMACS");
253 return addMessageIfNotSupported(errorReasons, error);
256 PmeRunMode pme_run_mode(const gmx_pme_t *pme)
258 GMX_ASSERT(pme != nullptr, "Expecting valid PME data pointer");
262 gmx::PinningPolicy pme_get_pinning_policy()
264 return gmx::PinningPolicy::PinnedIfSupported;
267 /*! \brief Number of bytes in a cache line.
269 * Must also be a multiple of the SIMD and SIMD4 register size, to
270 * preserve alignment.
272 const int gmxCacheLineSize = 64;
274 //! Set up coordinate communication
275 static void setup_coordinate_communication(pme_atomcomm_t *atc)
283 for (i = 1; i <= nslab/2; i++)
285 fw = (atc->nodeid + i) % nslab;
286 bw = (atc->nodeid - i + nslab) % nslab;
289 atc->node_dest[n] = fw;
290 atc->node_src[n] = bw;
295 atc->node_dest[n] = bw;
296 atc->node_src[n] = fw;
302 /*! \brief Round \p n up to the next multiple of \p f */
303 static int mult_up(int n, int f)
305 return ((n + f - 1)/f)*f;
308 /*! \brief Return estimate of the load imbalance from the PME grid not being a good match for the number of PME ranks */
309 static double estimate_pme_load_imbalance(struct gmx_pme_t *pme)
314 nma = pme->nnodes_major;
315 nmi = pme->nnodes_minor;
317 n1 = mult_up(pme->nkx, nma)*mult_up(pme->nky, nmi)*pme->nkz;
318 n2 = mult_up(pme->nkx, nma)*mult_up(pme->nkz, nmi)*pme->nky;
319 n3 = mult_up(pme->nky, nma)*mult_up(pme->nkz, nmi)*pme->nkx;
321 /* pme_solve is roughly double the cost of an fft */
323 return (n1 + n2 + 3*n3)/static_cast<double>(6*pme->nkx*pme->nky*pme->nkz);
326 /*! \brief Initialize atom communication data structure */
327 static void init_atomcomm(struct gmx_pme_t *pme, pme_atomcomm_t *atc,
328 int dimind, gmx_bool bSpread)
332 atc->dimind = dimind;
339 atc->mpi_comm = pme->mpi_comm_d[dimind];
340 MPI_Comm_size(atc->mpi_comm, &atc->nslab);
341 MPI_Comm_rank(atc->mpi_comm, &atc->nodeid);
345 fprintf(debug, "For PME atom communication in dimind %d: nslab %d rank %d\n", atc->dimind, atc->nslab, atc->nodeid);
349 atc->bSpread = bSpread;
350 atc->pme_order = pme->pme_order;
354 snew(atc->node_dest, atc->nslab);
355 snew(atc->node_src, atc->nslab);
356 setup_coordinate_communication(atc);
358 snew(atc->count_thread, pme->nthread);
359 for (thread = 0; thread < pme->nthread; thread++)
361 snew(atc->count_thread[thread], atc->nslab);
363 atc->count = atc->count_thread[0];
364 snew(atc->rcount, atc->nslab);
365 snew(atc->buf_index, atc->nslab);
368 atc->nthread = pme->nthread;
369 if (atc->nthread > 1)
371 snew(atc->thread_plist, atc->nthread);
373 snew(atc->spline, atc->nthread);
374 for (thread = 0; thread < atc->nthread; thread++)
376 if (atc->nthread > 1)
378 snew(atc->thread_plist[thread].n, atc->nthread+2*gmxCacheLineSize);
379 atc->thread_plist[thread].n += gmxCacheLineSize;
384 /*! \brief Destroy an atom communication data structure and its child structs */
385 static void destroy_atomcomm(pme_atomcomm_t *atc)
390 sfree(atc->node_dest);
391 sfree(atc->node_src);
392 for (int i = 0; i < atc->nthread; i++)
394 sfree(atc->count_thread[i]);
396 sfree(atc->count_thread);
398 sfree(atc->buf_index);
401 sfree(atc->coefficient);
407 sfree(atc->thread_idx);
408 for (int i = 0; i < atc->nthread; i++)
410 if (atc->nthread > 1)
412 int *n_ptr = atc->thread_plist[i].n - gmxCacheLineSize;
414 sfree(atc->thread_plist[i].i);
416 sfree(atc->spline[i].ind);
417 for (int d = 0; d < ZZ; d++)
419 sfree(atc->spline[i].theta[d]);
420 sfree(atc->spline[i].dtheta[d]);
422 sfree_aligned(atc->spline[i].ptr_dtheta_z);
423 sfree_aligned(atc->spline[i].ptr_theta_z);
425 if (atc->nthread > 1)
427 sfree(atc->thread_plist);
432 /*! \brief Initialize data structure for communication */
434 init_overlap_comm(pme_overlap_t * ol,
454 /* Linear translation of the PME grid won't affect reciprocal space
455 * calculations, so to optimize we only interpolate "upwards",
456 * which also means we only have to consider overlap in one direction.
457 * I.e., particles on this node might also be spread to grid indices
458 * that belong to higher nodes (modulo nnodes)
461 ol->s2g0.resize(ol->nnodes + 1);
462 ol->s2g1.resize(ol->nnodes);
465 fprintf(debug, "PME slab boundaries:");
467 for (int i = 0; i < nnodes; i++)
469 /* s2g0 the local interpolation grid start.
470 * s2g1 the local interpolation grid end.
471 * Since in calc_pidx we divide particles, and not grid lines,
472 * spatially uniform along dimension x or y, we need to round
473 * s2g0 down and s2g1 up.
475 ol->s2g0[i] = (i * ndata + 0) / nnodes;
476 ol->s2g1[i] = ((i + 1) * ndata + nnodes - 1) / nnodes + norder - 1;
480 fprintf(debug, " %3d %3d", ol->s2g0[i], ol->s2g1[i]);
483 ol->s2g0[nnodes] = ndata;
486 fprintf(debug, "\n");
489 /* Determine with how many nodes we need to communicate the grid overlap */
490 int testRankCount = 0;
495 for (int i = 0; i < nnodes; i++)
497 if ((i + testRankCount < nnodes && ol->s2g1[i] > ol->s2g0[i + testRankCount]) ||
498 (i + testRankCount >= nnodes && ol->s2g1[i] > ol->s2g0[i + testRankCount - nnodes] + ndata))
504 while (bCont && testRankCount < nnodes);
506 ol->comm_data.resize(testRankCount - 1);
509 for (size_t b = 0; b < ol->comm_data.size(); b++)
511 pme_grid_comm_t *pgc = &ol->comm_data[b];
514 pgc->send_id = (ol->nodeid + (b + 1)) % ol->nnodes;
515 int fft_start = ol->s2g0[pgc->send_id];
516 int fft_end = ol->s2g0[pgc->send_id + 1];
517 if (pgc->send_id < nodeid)
522 int send_index1 = ol->s2g1[nodeid];
523 send_index1 = std::min(send_index1, fft_end);
524 pgc->send_index0 = fft_start;
525 pgc->send_nindex = std::max(0, send_index1 - pgc->send_index0);
526 ol->send_size += pgc->send_nindex;
528 /* We always start receiving to the first index of our slab */
529 pgc->recv_id = (ol->nodeid - (b + 1) + ol->nnodes) % ol->nnodes;
530 fft_start = ol->s2g0[ol->nodeid];
531 fft_end = ol->s2g0[ol->nodeid + 1];
532 int recv_index1 = ol->s2g1[pgc->recv_id];
533 if (pgc->recv_id > nodeid)
535 recv_index1 -= ndata;
537 recv_index1 = std::min(recv_index1, fft_end);
538 pgc->recv_index0 = fft_start;
539 pgc->recv_nindex = std::max(0, recv_index1 - pgc->recv_index0);
543 /* Communicate the buffer sizes to receive */
544 for (size_t b = 0; b < ol->comm_data.size(); b++)
546 MPI_Sendrecv(&ol->send_size, 1, MPI_INT, ol->comm_data[b].send_id, b,
547 &ol->comm_data[b].recv_size, 1, MPI_INT, ol->comm_data[b].recv_id, b,
548 ol->mpi_comm, &stat);
552 /* For non-divisible grid we need pme_order iso pme_order-1 */
553 ol->sendbuf.resize(norder * commplainsize);
554 ol->recvbuf.resize(norder * commplainsize);
557 int minimalPmeGridSize(int pmeOrder)
559 /* The actual grid size limitations are:
560 * serial: >= pme_order
561 * DD, no OpenMP: >= 2*(pme_order - 1)
562 * DD, OpenMP: >= pme_order + 1
563 * But we use the maximum for simplicity since in practice there is not
564 * much performance difference between pme_order and 2*(pme_order -1).
566 int minimalSize = 2*(pmeOrder - 1);
568 GMX_RELEASE_ASSERT(pmeOrder >= 3, "pmeOrder has to be >= 3");
569 GMX_RELEASE_ASSERT(minimalSize >= pmeOrder + 1, "The grid size should be >= pmeOrder + 1");
574 bool gmx_pme_check_restrictions(int pme_order,
575 int nkx, int nky, int nkz,
576 int numPmeDomainsAlongX,
580 if (pme_order > PME_ORDER_MAX)
587 std::string message = gmx::formatString(
588 "pme_order (%d) is larger than the maximum allowed value (%d). Modify and recompile the code if you really need such a high order.",
589 pme_order, PME_ORDER_MAX);
590 GMX_THROW(gmx::InconsistentInputError(message));
593 const int minGridSize = minimalPmeGridSize(pme_order);
594 if (nkx < minGridSize ||
602 std::string message = gmx::formatString(
603 "The PME grid sizes need to be >= 2*(pme_order-1) (%d)",
605 GMX_THROW(gmx::InconsistentInputError(message));
608 /* Check for a limitation of the (current) sum_fftgrid_dd code.
609 * We only allow multiple communication pulses in dim 1, not in dim 0.
611 if (useThreads && (nkx < numPmeDomainsAlongX*pme_order &&
612 nkx != numPmeDomainsAlongX*(pme_order - 1)))
618 gmx_fatal(FARGS, "The number of PME grid lines per rank along x is %g. But when using OpenMP threads, the number of grid lines per rank along x should be >= pme_order (%d) or = pmeorder-1. To resolve this issue, use fewer ranks along x (and possibly more along y and/or z) by specifying -dd manually.",
619 nkx/static_cast<double>(numPmeDomainsAlongX), pme_order);
625 /*! \brief Round \p enumerator */
626 static int div_round_up(int enumerator, int denominator)
628 return (enumerator + denominator - 1)/denominator;
631 gmx_pme_t *gmx_pme_init(const t_commrec *cr,
632 const NumPmeDomains &numPmeDomains,
633 const t_inputrec *ir,
635 gmx_bool bFreeEnergy_q,
636 gmx_bool bFreeEnergy_lj,
637 gmx_bool bReproducible,
643 const gmx_device_info_t *gpuInfo,
644 PmeGpuProgramHandle pmeGpuProgram,
645 const gmx::MDLogger & /*mdlog*/)
647 int use_threads, sum_use_threads, i;
652 fprintf(debug, "Creating PME data structures.\n");
655 gmx::unique_cptr<gmx_pme_t, gmx_pme_destroy> pme(new gmx_pme_t());
657 pme->sum_qgrid_tmp = nullptr;
658 pme->sum_qgrid_dd_tmp = nullptr;
665 pme->nnodes_major = numPmeDomains.x;
666 pme->nnodes_minor = numPmeDomains.y;
669 if (numPmeDomains.x*numPmeDomains.y > 1)
671 pme->mpi_comm = cr->mpi_comm_mygroup;
673 MPI_Comm_rank(pme->mpi_comm, &pme->nodeid);
674 MPI_Comm_size(pme->mpi_comm, &pme->nnodes);
675 if (pme->nnodes != numPmeDomains.x*numPmeDomains.y)
677 gmx_incons("PME rank count mismatch");
682 pme->mpi_comm = MPI_COMM_NULL;
686 if (pme->nnodes == 1)
689 pme->mpi_comm_d[0] = MPI_COMM_NULL;
690 pme->mpi_comm_d[1] = MPI_COMM_NULL;
693 pme->nodeid_major = 0;
694 pme->nodeid_minor = 0;
696 pme->mpi_comm_d[0] = pme->mpi_comm_d[1] = MPI_COMM_NULL;
701 if (numPmeDomains.y == 1)
704 pme->mpi_comm_d[0] = pme->mpi_comm;
705 pme->mpi_comm_d[1] = MPI_COMM_NULL;
708 pme->nodeid_major = pme->nodeid;
709 pme->nodeid_minor = 0;
712 else if (numPmeDomains.x == 1)
715 pme->mpi_comm_d[0] = MPI_COMM_NULL;
716 pme->mpi_comm_d[1] = pme->mpi_comm;
719 pme->nodeid_major = 0;
720 pme->nodeid_minor = pme->nodeid;
724 if (pme->nnodes % numPmeDomains.x != 0)
726 gmx_incons("For 2D PME decomposition, #PME ranks must be divisible by the number of domains along x");
731 MPI_Comm_split(pme->mpi_comm, pme->nodeid % numPmeDomains.y,
732 pme->nodeid, &pme->mpi_comm_d[0]); /* My communicator along major dimension */
733 MPI_Comm_split(pme->mpi_comm, pme->nodeid/numPmeDomains.y,
734 pme->nodeid, &pme->mpi_comm_d[1]); /* My communicator along minor dimension */
736 MPI_Comm_rank(pme->mpi_comm_d[0], &pme->nodeid_major);
737 MPI_Comm_size(pme->mpi_comm_d[0], &pme->nnodes_major);
738 MPI_Comm_rank(pme->mpi_comm_d[1], &pme->nodeid_minor);
739 MPI_Comm_size(pme->mpi_comm_d[1], &pme->nnodes_minor);
742 pme->bPPnode = thisRankHasDuty(cr, DUTY_PP);
745 pme->nthread = nthread;
747 /* Check if any of the PME MPI ranks uses threads */
748 use_threads = (pme->nthread > 1 ? 1 : 0);
752 MPI_Allreduce(&use_threads, &sum_use_threads, 1, MPI_INT,
753 MPI_SUM, pme->mpi_comm);
758 sum_use_threads = use_threads;
760 pme->bUseThreads = (sum_use_threads > 0);
762 if (ir->ePBC == epbcSCREW)
764 gmx_fatal(FARGS, "pme does not (yet) work with pbc = screw");
768 * It is likely that the current gmx_pme_do() routine supports calculating
769 * only Coulomb or LJ while gmx_pme_init() configures for both,
770 * but that has never been tested.
771 * It is likely that the current gmx_pme_do() routine supports calculating,
772 * not calculating free-energy for Coulomb and/or LJ while gmx_pme_init()
773 * configures with free-energy, but that has never been tested.
775 pme->doCoulomb = EEL_PME(ir->coulombtype);
776 pme->doLJ = EVDW_PME(ir->vdwtype);
777 pme->bFEP_q = ((ir->efep != efepNO) && bFreeEnergy_q);
778 pme->bFEP_lj = ((ir->efep != efepNO) && bFreeEnergy_lj);
779 pme->bFEP = (pme->bFEP_q || pme->bFEP_lj);
783 pme->bP3M = (ir->coulombtype == eelP3M_AD || getenv("GMX_PME_P3M") != nullptr);
784 pme->pme_order = ir->pme_order;
785 pme->ewaldcoeff_q = ewaldcoeff_q;
786 pme->ewaldcoeff_lj = ewaldcoeff_lj;
788 /* Always constant electrostatics coefficients */
789 pme->epsilon_r = ir->epsilon_r;
791 /* Always constant LJ coefficients */
792 pme->ljpme_combination_rule = ir->ljpme_combination_rule;
794 // The box requires scaling with nwalls = 2, we store that condition as well
795 // as the scaling factor
796 delete pme->boxScaler;
797 pme->boxScaler = new EwaldBoxZScaler(*ir);
799 /* If we violate restrictions, generate a fatal error here */
800 gmx_pme_check_restrictions(pme->pme_order,
801 pme->nkx, pme->nky, pme->nkz,
811 MPI_Type_contiguous(DIM, GMX_MPI_REAL, &(pme->rvec_mpi));
812 MPI_Type_commit(&(pme->rvec_mpi));
815 /* Note that the coefficient spreading and force gathering, which usually
816 * takes about the same amount of time as FFT+solve_pme,
817 * is always fully load balanced
818 * (unless the coefficient distribution is inhomogeneous).
821 imbal = estimate_pme_load_imbalance(pme.get());
822 if (imbal >= 1.2 && pme->nodeid_major == 0 && pme->nodeid_minor == 0)
826 "NOTE: The load imbalance in PME FFT and solve is %d%%.\n"
827 " For optimal PME load balancing\n"
828 " PME grid_x (%d) and grid_y (%d) should be divisible by #PME_ranks_x (%d)\n"
829 " and PME grid_y (%d) and grid_z (%d) should be divisible by #PME_ranks_y (%d)\n"
831 gmx::roundToInt((imbal-1)*100),
832 pme->nkx, pme->nky, pme->nnodes_major,
833 pme->nky, pme->nkz, pme->nnodes_minor);
837 /* For non-divisible grid we need pme_order iso pme_order-1 */
838 /* In sum_qgrid_dd x overlap is copied in place: take padding into account.
839 * y is always copied through a buffer: we don't need padding in z,
840 * but we do need the overlap in x because of the communication order.
842 init_overlap_comm(&pme->overlap[0], pme->pme_order,
846 pme->nnodes_major, pme->nodeid_major,
848 (div_round_up(pme->nky, pme->nnodes_minor)+pme->pme_order)*(pme->nkz+pme->pme_order-1));
850 /* Along overlap dim 1 we can send in multiple pulses in sum_fftgrid_dd.
851 * We do this with an offset buffer of equal size, so we need to allocate
852 * extra for the offset. That's what the (+1)*pme->nkz is for.
854 init_overlap_comm(&pme->overlap[1], pme->pme_order,
858 pme->nnodes_minor, pme->nodeid_minor,
860 (div_round_up(pme->nkx, pme->nnodes_major)+pme->pme_order+1)*pme->nkz);
862 /* Double-check for a limitation of the (current) sum_fftgrid_dd code.
863 * Note that gmx_pme_check_restrictions checked for this already.
865 if (pme->bUseThreads && (pme->overlap[0].comm_data.size() > 1))
867 gmx_incons("More than one communication pulse required for grid overlap communication along the major dimension while using threads");
870 snew(pme->bsp_mod[XX], pme->nkx);
871 snew(pme->bsp_mod[YY], pme->nky);
872 snew(pme->bsp_mod[ZZ], pme->nkz);
874 pme->gpu = pmeGpu; /* Carrying over the single GPU structure */
875 pme->runMode = runMode;
877 /* The required size of the interpolation grid, including overlap.
878 * The allocated size (pmegrid_n?) might be slightly larger.
880 pme->pmegrid_nx = pme->overlap[0].s2g1[pme->nodeid_major] -
881 pme->overlap[0].s2g0[pme->nodeid_major];
882 pme->pmegrid_ny = pme->overlap[1].s2g1[pme->nodeid_minor] -
883 pme->overlap[1].s2g0[pme->nodeid_minor];
884 pme->pmegrid_nz_base = pme->nkz;
885 pme->pmegrid_nz = pme->pmegrid_nz_base + pme->pme_order - 1;
886 set_grid_alignment(&pme->pmegrid_nz, pme->pme_order);
887 pme->pmegrid_start_ix = pme->overlap[0].s2g0[pme->nodeid_major];
888 pme->pmegrid_start_iy = pme->overlap[1].s2g0[pme->nodeid_minor];
889 pme->pmegrid_start_iz = 0;
891 make_gridindex_to_localindex(pme->nkx,
892 pme->pmegrid_start_ix,
893 pme->pmegrid_nx - (pme->pme_order-1),
894 &pme->nnx, &pme->fshx);
895 make_gridindex_to_localindex(pme->nky,
896 pme->pmegrid_start_iy,
897 pme->pmegrid_ny - (pme->pme_order-1),
898 &pme->nny, &pme->fshy);
899 make_gridindex_to_localindex(pme->nkz,
900 pme->pmegrid_start_iz,
901 pme->pmegrid_nz_base,
902 &pme->nnz, &pme->fshz);
904 pme->spline_work = make_pme_spline_work(pme->pme_order);
909 /* It doesn't matter if we allocate too many grids here,
910 * we only allocate and use the ones we need.
914 pme->ngrids = ((ir->ljpme_combination_rule == eljpmeLB) ? DO_Q_AND_LJ_LB : DO_Q_AND_LJ);
920 snew(pme->fftgrid, pme->ngrids);
921 snew(pme->cfftgrid, pme->ngrids);
922 snew(pme->pfft_setup, pme->ngrids);
924 for (i = 0; i < pme->ngrids; ++i)
926 if ((i < DO_Q && pme->doCoulomb && (i == 0 ||
928 (i >= DO_Q && pme->doLJ && (i == 2 ||
930 ir->ljpme_combination_rule == eljpmeLB)))
932 pmegrids_init(&pme->pmegrid[i],
933 pme->pmegrid_nx, pme->pmegrid_ny, pme->pmegrid_nz,
934 pme->pmegrid_nz_base,
938 pme->overlap[0].s2g1[pme->nodeid_major]-pme->overlap[0].s2g0[pme->nodeid_major+1],
939 pme->overlap[1].s2g1[pme->nodeid_minor]-pme->overlap[1].s2g0[pme->nodeid_minor+1]);
940 /* This routine will allocate the grid data to fit the FFTs */
941 const auto allocateRealGridForGpu = (pme->runMode == PmeRunMode::Mixed) ? gmx::PinningPolicy::PinnedIfSupported : gmx::PinningPolicy::CannotBePinned;
942 gmx_parallel_3dfft_init(&pme->pfft_setup[i], ndata,
943 &pme->fftgrid[i], &pme->cfftgrid[i],
945 bReproducible, pme->nthread, allocateRealGridForGpu);
952 /* Use plain SPME B-spline interpolation */
953 make_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
957 /* Use the P3M grid-optimized influence function */
958 make_p3m_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
961 /* Use atc[0] for spreading */
962 init_atomcomm(pme.get(), &pme->atc[0], numPmeDomains.x > 1 ? 0 : 1, TRUE);
963 if (pme->ndecompdim >= 2)
965 init_atomcomm(pme.get(), &pme->atc[1], 1, FALSE);
968 if (pme->nnodes == 1)
970 pme->atc[0].n = homenr;
971 pme_realloc_atomcomm_things(&pme->atc[0]);
974 pme->lb_buf1 = nullptr;
975 pme->lb_buf2 = nullptr;
976 pme->lb_buf_nalloc = 0;
978 if (pme_gpu_active(pme.get()))
982 // Initial check of validity of the data
983 std::string errorString;
984 bool canRunOnGpu = pme_gpu_check_restrictions(pme.get(), &errorString);
987 GMX_THROW(gmx::NotImplementedError(errorString));
991 pme_gpu_reinit(pme.get(), gpuInfo, pmeGpuProgram);
994 pme_init_all_work(&pme->solve_work, pme->nthread, pme->nkx);
996 // no exception was thrown during the init, so we hand over the PME structure handle
997 return pme.release();
1000 void gmx_pme_reinit(struct gmx_pme_t **pmedata,
1001 const t_commrec *cr,
1002 struct gmx_pme_t * pme_src,
1003 const t_inputrec * ir,
1004 const ivec grid_size,
1010 // Create a copy of t_inputrec fields that are used in gmx_pme_init().
1011 // TODO: This would be better as just copying a sub-structure that contains
1012 // all the PME parameters and nothing else.
1014 irc.ePBC = ir->ePBC;
1015 irc.coulombtype = ir->coulombtype;
1016 irc.vdwtype = ir->vdwtype;
1017 irc.efep = ir->efep;
1018 irc.pme_order = ir->pme_order;
1019 irc.epsilon_r = ir->epsilon_r;
1020 irc.ljpme_combination_rule = ir->ljpme_combination_rule;
1021 irc.nkx = grid_size[XX];
1022 irc.nky = grid_size[YY];
1023 irc.nkz = grid_size[ZZ];
1025 if (pme_src->nnodes == 1)
1027 homenr = pme_src->atc[0].n;
1036 const gmx::MDLogger dummyLogger;
1037 // This is reinit which is currently only changing grid size/coefficients,
1038 // so we don't expect the actual logging.
1039 // TODO: when PME is an object, it should take reference to mdlog on construction and save it.
1040 GMX_ASSERT(pmedata, "Invalid PME pointer");
1041 NumPmeDomains numPmeDomains = { pme_src->nnodes_major, pme_src->nnodes_minor };
1042 *pmedata = gmx_pme_init(cr, numPmeDomains,
1043 &irc, homenr, pme_src->bFEP_q, pme_src->bFEP_lj, FALSE, ewaldcoeff_q, ewaldcoeff_lj,
1044 pme_src->nthread, pme_src->runMode, pme_src->gpu, nullptr, nullptr, dummyLogger);
1045 //TODO this is mostly passing around current values
1047 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1049 /* We can easily reuse the allocated pme grids in pme_src */
1050 reuse_pmegrids(&pme_src->pmegrid[PME_GRID_QA], &(*pmedata)->pmegrid[PME_GRID_QA]);
1051 /* We would like to reuse the fft grids, but that's harder */
1054 void gmx_pme_calc_energy(struct gmx_pme_t *pme, int n, rvec *x, real *q, real *V)
1056 pme_atomcomm_t *atc;
1059 if (pme->nnodes > 1)
1061 gmx_incons("gmx_pme_calc_energy called in parallel");
1065 gmx_incons("gmx_pme_calc_energy with free energy");
1068 atc = &pme->atc_energy;
1070 if (atc->spline == nullptr)
1072 snew(atc->spline, atc->nthread);
1075 atc->bSpread = TRUE;
1076 atc->pme_order = pme->pme_order;
1078 pme_realloc_atomcomm_things(atc);
1080 atc->coefficient = q;
1082 /* We only use the A-charges grid */
1083 grid = &pme->pmegrid[PME_GRID_QA];
1085 /* Only calculate the spline coefficients, don't actually spread */
1086 spread_on_grid(pme, atc, nullptr, TRUE, FALSE, pme->fftgrid[PME_GRID_QA], FALSE, PME_GRID_QA);
1088 *V = gather_energy_bsplines(pme, grid->grid.grid, atc);
1091 /*! \brief Calculate initial Lorentz-Berthelot coefficients for LJ-PME */
1093 calc_initial_lb_coeffs(struct gmx_pme_t *pme, const real *local_c6, const real *local_sigma)
1096 for (i = 0; i < pme->atc[0].n; ++i)
1099 sigma4 = local_sigma[i];
1100 sigma4 = sigma4*sigma4;
1101 sigma4 = sigma4*sigma4;
1102 pme->atc[0].coefficient[i] = local_c6[i] / sigma4;
1106 /*! \brief Calculate next Lorentz-Berthelot coefficients for LJ-PME */
1108 calc_next_lb_coeffs(struct gmx_pme_t *pme, const real *local_sigma)
1112 for (i = 0; i < pme->atc[0].n; ++i)
1114 pme->atc[0].coefficient[i] *= local_sigma[i];
1118 int gmx_pme_do(struct gmx_pme_t *pme,
1119 int start, int homenr,
1121 real chargeA[], real chargeB[],
1122 real c6A[], real c6B[],
1123 real sigmaA[], real sigmaB[],
1124 matrix box, const t_commrec *cr,
1125 int maxshift_x, int maxshift_y,
1126 t_nrnb *nrnb, gmx_wallcycle *wcycle,
1127 matrix vir_q, matrix vir_lj,
1128 real *energy_q, real *energy_lj,
1129 real lambda_q, real lambda_lj,
1130 real *dvdlambda_q, real *dvdlambda_lj,
1133 GMX_ASSERT(pme->runMode == PmeRunMode::CPU, "gmx_pme_do should not be called on the GPU PME run.");
1135 int d, i, j, npme, grid_index, max_grid_index;
1137 pme_atomcomm_t *atc = nullptr;
1138 pmegrids_t *pmegrid = nullptr;
1139 real *grid = nullptr;
1141 real *coefficient = nullptr;
1142 PmeOutput output[2]; // The second is used for the B state with FEP
1145 gmx_parallel_3dfft_t pfft_setup;
1147 t_complex * cfftgrid;
1149 gmx_bool bFirst, bDoSplines;
1151 int fep_states_lj = pme->bFEP_lj ? 2 : 1;
1152 const gmx_bool bCalcEnerVir = (flags & GMX_PME_CALC_ENER_VIR) != 0;
1153 const gmx_bool bBackFFT = (flags & (GMX_PME_CALC_F | GMX_PME_CALC_POT)) != 0;
1154 const gmx_bool bCalcF = (flags & GMX_PME_CALC_F) != 0;
1156 /* We could be passing lambda!=1 while no q or LJ is actually perturbed */
1166 assert(pme->nnodes > 0);
1167 assert(pme->nnodes == 1 || pme->ndecompdim > 0);
1169 if (pme->nnodes > 1)
1173 if (atc->npd > atc->pd_nalloc)
1175 atc->pd_nalloc = over_alloc_dd(atc->npd);
1176 srenew(atc->pd, atc->pd_nalloc);
1178 for (d = pme->ndecompdim-1; d >= 0; d--)
1181 atc->maxshift = (atc->dimind == 0 ? maxshift_x : maxshift_y);
1187 /* This could be necessary for TPI */
1188 pme->atc[0].n = homenr;
1189 if (DOMAINDECOMP(cr))
1191 pme_realloc_atomcomm_things(atc);
1198 pme->boxScaler->scaleBox(box, scaledBox);
1200 gmx::invertBoxMatrix(scaledBox, pme->recipbox);
1203 /* For simplicity, we construct the splines for all particles if
1204 * more than one PME calculations is needed. Some optimization
1205 * could be done by keeping track of which atoms have splines
1206 * constructed, and construct new splines on each pass for atoms
1207 * that don't yet have them.
1210 bDoSplines = pme->bFEP || (pme->doCoulomb && pme->doLJ);
1212 /* We need a maximum of four separate PME calculations:
1213 * grid_index=0: Coulomb PME with charges from state A
1214 * grid_index=1: Coulomb PME with charges from state B
1215 * grid_index=2: LJ PME with C6 from state A
1216 * grid_index=3: LJ PME with C6 from state B
1217 * For Lorentz-Berthelot combination rules, a separate loop is used to
1218 * calculate all the terms
1221 /* If we are doing LJ-PME with LB, we only do Q here */
1222 max_grid_index = (pme->ljpme_combination_rule == eljpmeLB) ? DO_Q : DO_Q_AND_LJ;
1224 for (grid_index = 0; grid_index < max_grid_index; ++grid_index)
1226 /* Check if we should do calculations at this grid_index
1227 * If grid_index is odd we should be doing FEP
1228 * If grid_index < 2 we should be doing electrostatic PME
1229 * If grid_index >= 2 we should be doing LJ-PME
1231 if ((grid_index < DO_Q && (!pme->doCoulomb ||
1232 (grid_index == 1 && !pme->bFEP_q))) ||
1233 (grid_index >= DO_Q && (!pme->doLJ ||
1234 (grid_index == 3 && !pme->bFEP_lj))))
1238 /* Unpack structure */
1239 pmegrid = &pme->pmegrid[grid_index];
1240 fftgrid = pme->fftgrid[grid_index];
1241 cfftgrid = pme->cfftgrid[grid_index];
1242 pfft_setup = pme->pfft_setup[grid_index];
1245 case 0: coefficient = chargeA + start; break;
1246 case 1: coefficient = chargeB + start; break;
1247 case 2: coefficient = c6A + start; break;
1248 case 3: coefficient = c6B + start; break;
1251 grid = pmegrid->grid.grid;
1255 fprintf(debug, "PME: number of ranks = %d, rank = %d\n",
1256 cr->nnodes, cr->nodeid);
1257 fprintf(debug, "Grid = %p\n", static_cast<void*>(grid));
1258 if (grid == nullptr)
1260 gmx_fatal(FARGS, "No grid!");
1264 if (pme->nnodes == 1)
1266 atc->coefficient = coefficient;
1270 wallcycle_start(wcycle, ewcPME_REDISTXF);
1271 do_redist_pos_coeffs(pme, cr, start, homenr, bFirst, x, coefficient);
1273 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1278 fprintf(debug, "Rank= %6d, pme local particles=%6d\n",
1279 cr->nodeid, atc->n);
1282 if (flags & GMX_PME_SPREAD)
1284 wallcycle_start(wcycle, ewcPME_SPREAD);
1286 /* Spread the coefficients on a grid */
1287 spread_on_grid(pme, &pme->atc[0], pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1291 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
1293 inc_nrnb(nrnb, eNR_SPREADBSP,
1294 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
1296 if (!pme->bUseThreads)
1298 wrap_periodic_pmegrid(pme, grid);
1300 /* sum contributions to local grid from other nodes */
1302 if (pme->nnodes > 1)
1304 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1308 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1311 wallcycle_stop(wcycle, ewcPME_SPREAD);
1313 /* TODO If the OpenMP and single-threaded implementations
1314 converge, then spread_on_grid() and
1315 copy_pmegrid_to_fftgrid() will perhaps live in the same
1320 /* Here we start a large thread parallel region */
1321 #pragma omp parallel num_threads(pme->nthread) private(thread)
1325 thread = gmx_omp_get_thread_num();
1326 if (flags & GMX_PME_SOLVE)
1333 wallcycle_start(wcycle, ewcPME_FFT);
1335 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
1339 wallcycle_stop(wcycle, ewcPME_FFT);
1342 /* solve in k-space for our local cells */
1345 wallcycle_start(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
1347 if (grid_index < DO_Q)
1350 solve_pme_yzx(pme, cfftgrid,
1351 scaledBox[XX][XX]*scaledBox[YY][YY]*scaledBox[ZZ][ZZ],
1353 pme->nthread, thread);
1358 solve_pme_lj_yzx(pme, &cfftgrid, FALSE,
1359 scaledBox[XX][XX]*scaledBox[YY][YY]*scaledBox[ZZ][ZZ],
1361 pme->nthread, thread);
1366 wallcycle_stop(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
1367 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1376 wallcycle_start(wcycle, ewcPME_FFT);
1378 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
1382 wallcycle_stop(wcycle, ewcPME_FFT);
1385 if (pme->nodeid == 0)
1387 real ntot = pme->nkx*pme->nky*pme->nkz;
1388 npme = static_cast<int>(ntot*std::log(ntot)/std::log(2.0));
1389 inc_nrnb(nrnb, eNR_FFT, 2*npme);
1392 /* Note: this wallcycle region is closed below
1393 outside an OpenMP region, so take care if
1394 refactoring code here. */
1395 wallcycle_start(wcycle, ewcPME_GATHER);
1398 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1400 } GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1402 /* End of thread parallel section.
1403 * With MPI we have to synchronize here before gmx_sum_qgrid_dd.
1408 /* distribute local grid to all nodes */
1410 if (pme->nnodes > 1)
1412 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1416 unwrap_periodic_pmegrid(pme, grid);
1421 /* interpolate forces for our local atoms */
1424 /* If we are running without parallelization,
1425 * atc->f is the actual force array, not a buffer,
1426 * therefore we should not clear it.
1428 lambda = grid_index < DO_Q ? lambda_q : lambda_lj;
1429 bClearF = (bFirst && PAR(cr));
1430 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1431 for (thread = 0; thread < pme->nthread; thread++)
1435 gather_f_bsplines(pme, grid, bClearF, atc,
1436 &atc->spline[thread],
1437 pme->bFEP ? (grid_index % 2 == 0 ? 1.0-lambda : lambda) : 1.0);
1439 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1443 inc_nrnb(nrnb, eNR_GATHERFBSP,
1444 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
1445 /* Note: this wallcycle region is opened above inside an OpenMP
1446 region, so take care if refactoring code here. */
1447 wallcycle_stop(wcycle, ewcPME_GATHER);
1452 /* This should only be called on the master thread
1453 * and after the threads have synchronized.
1457 get_pme_ener_vir_q(pme->solve_work, pme->nthread, &output[grid_index % 2]);
1461 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output[grid_index % 2]);
1465 } /* of grid_index-loop */
1467 /* For Lorentz-Berthelot combination rules in LJ-PME, we need to calculate
1470 if (pme->doLJ && pme->ljpme_combination_rule == eljpmeLB)
1472 /* Loop over A- and B-state if we are doing FEP */
1473 for (fep_state = 0; fep_state < fep_states_lj; ++fep_state)
1475 real *local_c6 = nullptr, *local_sigma = nullptr, *RedistC6 = nullptr, *RedistSigma = nullptr;
1476 if (pme->nnodes == 1)
1478 if (pme->lb_buf1 == nullptr)
1480 pme->lb_buf_nalloc = pme->atc[0].n;
1481 snew(pme->lb_buf1, pme->lb_buf_nalloc);
1483 pme->atc[0].coefficient = pme->lb_buf1;
1488 local_sigma = sigmaA;
1492 local_sigma = sigmaB;
1495 gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1505 RedistSigma = sigmaA;
1509 RedistSigma = sigmaB;
1512 gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1514 wallcycle_start(wcycle, ewcPME_REDISTXF);
1516 do_redist_pos_coeffs(pme, cr, start, homenr, bFirst, x, RedistC6);
1517 if (pme->lb_buf_nalloc < atc->n)
1519 pme->lb_buf_nalloc = atc->nalloc;
1520 srenew(pme->lb_buf1, pme->lb_buf_nalloc);
1521 srenew(pme->lb_buf2, pme->lb_buf_nalloc);
1523 local_c6 = pme->lb_buf1;
1524 for (i = 0; i < atc->n; ++i)
1526 local_c6[i] = atc->coefficient[i];
1529 do_redist_pos_coeffs(pme, cr, start, homenr, FALSE, x, RedistSigma);
1530 local_sigma = pme->lb_buf2;
1531 for (i = 0; i < atc->n; ++i)
1533 local_sigma[i] = atc->coefficient[i];
1536 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1538 calc_initial_lb_coeffs(pme, local_c6, local_sigma);
1540 /*Seven terms in LJ-PME with LB, grid_index < 2 reserved for electrostatics*/
1541 for (grid_index = 2; grid_index < 9; ++grid_index)
1543 /* Unpack structure */
1544 pmegrid = &pme->pmegrid[grid_index];
1545 fftgrid = pme->fftgrid[grid_index];
1546 pfft_setup = pme->pfft_setup[grid_index];
1547 calc_next_lb_coeffs(pme, local_sigma);
1548 grid = pmegrid->grid.grid;
1550 if (flags & GMX_PME_SPREAD)
1552 wallcycle_start(wcycle, ewcPME_SPREAD);
1553 /* Spread the c6 on a grid */
1554 spread_on_grid(pme, &pme->atc[0], pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1558 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
1561 inc_nrnb(nrnb, eNR_SPREADBSP,
1562 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
1563 if (pme->nthread == 1)
1565 wrap_periodic_pmegrid(pme, grid);
1566 /* sum contributions to local grid from other nodes */
1568 if (pme->nnodes > 1)
1570 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1573 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1575 wallcycle_stop(wcycle, ewcPME_SPREAD);
1577 /*Here we start a large thread parallel region*/
1578 #pragma omp parallel num_threads(pme->nthread) private(thread)
1582 thread = gmx_omp_get_thread_num();
1583 if (flags & GMX_PME_SOLVE)
1588 wallcycle_start(wcycle, ewcPME_FFT);
1591 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
1595 wallcycle_stop(wcycle, ewcPME_FFT);
1599 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1603 if (flags & GMX_PME_SOLVE)
1605 /* solve in k-space for our local cells */
1606 #pragma omp parallel num_threads(pme->nthread) private(thread)
1611 thread = gmx_omp_get_thread_num();
1614 wallcycle_start(wcycle, ewcLJPME);
1618 solve_pme_lj_yzx(pme, &pme->cfftgrid[2], TRUE,
1619 scaledBox[XX][XX]*scaledBox[YY][YY]*scaledBox[ZZ][ZZ],
1621 pme->nthread, thread);
1624 wallcycle_stop(wcycle, ewcLJPME);
1625 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1628 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1634 /* This should only be called on the master thread and
1635 * after the threads have synchronized.
1637 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output[fep_state]);
1642 bFirst = !pme->doCoulomb;
1643 calc_initial_lb_coeffs(pme, local_c6, local_sigma);
1644 for (grid_index = 8; grid_index >= 2; --grid_index)
1646 /* Unpack structure */
1647 pmegrid = &pme->pmegrid[grid_index];
1648 fftgrid = pme->fftgrid[grid_index];
1649 pfft_setup = pme->pfft_setup[grid_index];
1650 grid = pmegrid->grid.grid;
1651 calc_next_lb_coeffs(pme, local_sigma);
1652 #pragma omp parallel num_threads(pme->nthread) private(thread)
1656 thread = gmx_omp_get_thread_num();
1660 wallcycle_start(wcycle, ewcPME_FFT);
1663 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
1667 wallcycle_stop(wcycle, ewcPME_FFT);
1670 if (pme->nodeid == 0)
1672 real ntot = pme->nkx*pme->nky*pme->nkz;
1673 npme = static_cast<int>(ntot*std::log(ntot)/std::log(2.0));
1674 inc_nrnb(nrnb, eNR_FFT, 2*npme);
1676 wallcycle_start(wcycle, ewcPME_GATHER);
1679 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1681 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1682 } /*#pragma omp parallel*/
1684 /* distribute local grid to all nodes */
1686 if (pme->nnodes > 1)
1688 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1692 unwrap_periodic_pmegrid(pme, grid);
1696 /* interpolate forces for our local atoms */
1697 bClearF = (bFirst && PAR(cr));
1698 scale = pme->bFEP ? (fep_state < 1 ? 1.0-lambda_lj : lambda_lj) : 1.0;
1699 scale *= lb_scale_factor[grid_index-2];
1701 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1702 for (thread = 0; thread < pme->nthread; thread++)
1706 gather_f_bsplines(pme, grid, bClearF, &pme->atc[0],
1707 &pme->atc[0].spline[thread],
1710 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1714 inc_nrnb(nrnb, eNR_GATHERFBSP,
1715 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
1717 wallcycle_stop(wcycle, ewcPME_GATHER);
1720 } /* for (grid_index = 8; grid_index >= 2; --grid_index) */
1722 } /* for (fep_state = 0; fep_state < fep_states_lj; ++fep_state) */
1723 } /* if ((flags & GMX_PME_DO_LJ) && pme->ljpme_combination_rule == eljpmeLB) */
1725 if (bCalcF && pme->nnodes > 1)
1727 wallcycle_start(wcycle, ewcPME_REDISTXF);
1728 for (d = 0; d < pme->ndecompdim; d++)
1731 if (d == pme->ndecompdim - 1)
1738 n_d = pme->atc[d+1].n;
1739 f_d = pme->atc[d+1].f;
1741 if (DOMAINDECOMP(cr))
1743 dd_pmeredist_f(pme, atc, n_d, f_d,
1744 d == pme->ndecompdim-1 && pme->bPPnode);
1748 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1757 *energy_q = output[0].coulombEnergy_;
1758 m_add(vir_q, output[0].coulombVirial_, vir_q);
1762 *energy_q = (1.0-lambda_q)*output[0].coulombEnergy_ + lambda_q*output[1].coulombEnergy_;
1763 *dvdlambda_q += output[1].coulombEnergy_ - output[0].coulombEnergy_;
1764 for (i = 0; i < DIM; i++)
1766 for (j = 0; j < DIM; j++)
1768 vir_q[i][j] += (1.0-lambda_q)*output[0].coulombVirial_[i][j] +
1769 lambda_q*output[1].coulombVirial_[i][j];
1775 fprintf(debug, "Electrostatic PME mesh energy: %g\n", *energy_q);
1787 *energy_lj = output[0].lennardJonesEnergy_;
1788 m_add(vir_lj, output[0].lennardJonesVirial_, vir_lj);
1792 *energy_lj = (1.0-lambda_lj)*output[0].lennardJonesEnergy_ + lambda_lj*output[1].lennardJonesEnergy_;
1793 *dvdlambda_lj += output[1].lennardJonesEnergy_ - output[0].lennardJonesEnergy_;
1794 for (i = 0; i < DIM; i++)
1796 for (j = 0; j < DIM; j++)
1798 vir_lj[i][j] += (1.0-lambda_lj)*output[0].lennardJonesVirial_[i][j] + lambda_lj*output[1].lennardJonesVirial_[i][j];
1804 fprintf(debug, "Lennard-Jones PME mesh energy: %g\n", *energy_lj);
1815 void gmx_pme_destroy(gmx_pme_t *pme)
1822 delete pme->boxScaler;
1831 for (int i = 0; i < pme->ngrids; ++i)
1833 pmegrids_destroy(&pme->pmegrid[i]);
1835 if (pme->pfft_setup)
1837 for (int i = 0; i < pme->ngrids; ++i)
1839 gmx_parallel_3dfft_destroy(pme->pfft_setup[i]);
1842 sfree(pme->fftgrid);
1843 sfree(pme->cfftgrid);
1844 sfree(pme->pfft_setup);
1846 for (int i = 0; i < std::max(1, pme->ndecompdim); i++) //pme->atc[0] is always allocated
1848 destroy_atomcomm(&pme->atc[i]);
1851 for (int i = 0; i < DIM; i++)
1853 sfree(pme->bsp_mod[i]);
1856 sfree(pme->lb_buf1);
1857 sfree(pme->lb_buf2);
1862 if (pme->solve_work)
1864 pme_free_all_work(&pme->solve_work, pme->nthread);
1867 sfree(pme->sum_qgrid_tmp);
1868 sfree(pme->sum_qgrid_dd_tmp);
1870 destroy_pme_spline_work(pme->spline_work);
1872 if (pme_gpu_active(pme) && pme->gpu)
1874 pme_gpu_destroy(pme->gpu);
1880 void gmx_pme_reinit_atoms(const gmx_pme_t *pme, const int nAtoms, const real *charges)
1882 if (pme_gpu_active(pme))
1884 pme_gpu_reinit_atoms(pme->gpu, nAtoms, charges);
1886 // TODO: handle the CPU case here; handle the whole t_mdatoms