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, 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/math/gmxcomplex.h"
92 #include "gromacs/math/invertmatrix.h"
93 #include "gromacs/math/units.h"
94 #include "gromacs/math/vec.h"
95 #include "gromacs/math/vectypes.h"
96 #include "gromacs/mdtypes/commrec.h"
97 #include "gromacs/mdtypes/forcerec.h"
98 #include "gromacs/mdtypes/inputrec.h"
99 #include "gromacs/mdtypes/md_enums.h"
100 #include "gromacs/pbcutil/pbc.h"
101 #include "gromacs/timing/cyclecounter.h"
102 #include "gromacs/timing/wallcycle.h"
103 #include "gromacs/timing/walltime_accounting.h"
104 #include "gromacs/topology/topology.h"
105 #include "gromacs/utility/basedefinitions.h"
106 #include "gromacs/utility/exceptions.h"
107 #include "gromacs/utility/fatalerror.h"
108 #include "gromacs/utility/gmxmpi.h"
109 #include "gromacs/utility/gmxomp.h"
110 #include "gromacs/utility/logger.h"
111 #include "gromacs/utility/real.h"
112 #include "gromacs/utility/smalloc.h"
113 #include "gromacs/utility/stringutil.h"
114 #include "gromacs/utility/unique_cptr.h"
116 #include "calculate-spline-moduli.h"
117 #include "pme-gather.h"
118 #include "pme-gpu-internal.h"
119 #include "pme-grid.h"
120 #include "pme-internal.h"
121 #include "pme-redistribute.h"
122 #include "pme-solve.h"
123 #include "pme-spline-work.h"
124 #include "pme-spread.h"
126 /*! \brief Help build a descriptive message in \c error if there are
127 * \c errorReasons why PME on GPU is not supported.
129 * \returns Whether the lack of errorReasons indicate there is support. */
131 addMessageIfNotSupported(const std::list<std::string> &errorReasons,
134 bool foundErrorReasons = errorReasons.empty();
135 if (!foundErrorReasons && error)
137 std::string regressionTestMarker = "PME GPU does not support";
138 // this prefix is tested for in the regression tests script gmxtest.pl
139 *error = regressionTestMarker + ": " + gmx::joinStrings(errorReasons, "; ") + ".";
141 return foundErrorReasons;
144 bool pme_gpu_supports_build(std::string *error)
146 std::list<std::string> errorReasons;
149 errorReasons.emplace_back("double precision");
151 if (GMX_GPU != GMX_GPU_CUDA)
153 errorReasons.emplace_back("non-CUDA build of GROMACS");
155 return addMessageIfNotSupported(errorReasons, error);
158 bool pme_gpu_supports_input(const t_inputrec &ir, const gmx_mtop_t &mtop, std::string *error)
160 std::list<std::string> errorReasons;
161 if (!EEL_PME(ir.coulombtype))
163 errorReasons.emplace_back("systems that do not use PME for electrostatics");
165 if (ir.pme_order != 4)
167 errorReasons.emplace_back("interpolation orders other than 4");
169 if (ir.efep != efepNO)
171 if (gmx_mtop_has_perturbed_charges(mtop))
173 errorReasons.emplace_back("free energy calculations with perturbed charges (multiple grids)");
176 if (EVDW_PME(ir.vdwtype))
178 errorReasons.emplace_back("Lennard-Jones PME");
180 if (ir.cutoff_scheme == ecutsGROUP)
182 errorReasons.emplace_back("group cutoff scheme");
184 if (!EI_DYNAMICS(ir.eI))
186 errorReasons.emplace_back("not a dynamical integrator");
188 return addMessageIfNotSupported(errorReasons, error);
191 /*! \brief \libinternal
192 * Finds out if PME with given inputs is possible to run on GPU.
193 * This function is an internal final check, validating the whole PME structure on creation,
194 * but it still duplicates the preliminary checks from the above (externally exposed) pme_gpu_supports_input() - just in case.
196 * \param[in] pme The PME structure.
197 * \param[out] error The error message if the input is not supported on GPU.
198 * \returns True if this PME input is possible to run on GPU, false otherwise.
200 static bool pme_gpu_check_restrictions(const gmx_pme_t *pme, std::string *error)
202 std::list<std::string> errorReasons;
203 if (pme->nnodes != 1)
205 errorReasons.emplace_back("PME decomposition");
207 if (pme->pme_order != 4)
209 errorReasons.emplace_back("interpolation orders other than 4");
213 errorReasons.emplace_back("free energy calculations (multiple grids)");
217 errorReasons.emplace_back("Lennard-Jones PME");
221 errorReasons.emplace_back("double precision");
223 if (GMX_GPU != GMX_GPU_CUDA)
225 errorReasons.emplace_back("non-CUDA build of GROMACS");
228 return addMessageIfNotSupported(errorReasons, error);
231 PmeRunMode pme_run_mode(const gmx_pme_t *pme)
233 GMX_ASSERT(pme != nullptr, "Expecting valid PME data pointer");
237 gmx::PinningPolicy pme_get_pinning_policy()
239 return gmx::PinningPolicy::PinnedIfSupported;
242 /*! \brief Number of bytes in a cache line.
244 * Must also be a multiple of the SIMD and SIMD4 register size, to
245 * preserve alignment.
247 const int gmxCacheLineSize = 64;
249 //! Set up coordinate communication
250 static void setup_coordinate_communication(pme_atomcomm_t *atc)
258 for (i = 1; i <= nslab/2; i++)
260 fw = (atc->nodeid + i) % nslab;
261 bw = (atc->nodeid - i + nslab) % nslab;
264 atc->node_dest[n] = fw;
265 atc->node_src[n] = bw;
270 atc->node_dest[n] = bw;
271 atc->node_src[n] = fw;
277 /*! \brief Round \p n up to the next multiple of \p f */
278 static int mult_up(int n, int f)
280 return ((n + f - 1)/f)*f;
283 /*! \brief Return estimate of the load imbalance from the PME grid not being a good match for the number of PME ranks */
284 static double estimate_pme_load_imbalance(struct gmx_pme_t *pme)
289 nma = pme->nnodes_major;
290 nmi = pme->nnodes_minor;
292 n1 = mult_up(pme->nkx, nma)*mult_up(pme->nky, nmi)*pme->nkz;
293 n2 = mult_up(pme->nkx, nma)*mult_up(pme->nkz, nmi)*pme->nky;
294 n3 = mult_up(pme->nky, nma)*mult_up(pme->nkz, nmi)*pme->nkx;
296 /* pme_solve is roughly double the cost of an fft */
298 return (n1 + n2 + 3*n3)/static_cast<double>(6*pme->nkx*pme->nky*pme->nkz);
301 /*! \brief Initialize atom communication data structure */
302 static void init_atomcomm(struct gmx_pme_t *pme, pme_atomcomm_t *atc,
303 int dimind, gmx_bool bSpread)
307 atc->dimind = dimind;
314 atc->mpi_comm = pme->mpi_comm_d[dimind];
315 MPI_Comm_size(atc->mpi_comm, &atc->nslab);
316 MPI_Comm_rank(atc->mpi_comm, &atc->nodeid);
320 fprintf(debug, "For PME atom communication in dimind %d: nslab %d rank %d\n", atc->dimind, atc->nslab, atc->nodeid);
324 atc->bSpread = bSpread;
325 atc->pme_order = pme->pme_order;
329 snew(atc->node_dest, atc->nslab);
330 snew(atc->node_src, atc->nslab);
331 setup_coordinate_communication(atc);
333 snew(atc->count_thread, pme->nthread);
334 for (thread = 0; thread < pme->nthread; thread++)
336 snew(atc->count_thread[thread], atc->nslab);
338 atc->count = atc->count_thread[0];
339 snew(atc->rcount, atc->nslab);
340 snew(atc->buf_index, atc->nslab);
343 atc->nthread = pme->nthread;
344 if (atc->nthread > 1)
346 snew(atc->thread_plist, atc->nthread);
348 snew(atc->spline, atc->nthread);
349 for (thread = 0; thread < atc->nthread; thread++)
351 if (atc->nthread > 1)
353 snew(atc->thread_plist[thread].n, atc->nthread+2*gmxCacheLineSize);
354 atc->thread_plist[thread].n += gmxCacheLineSize;
359 /*! \brief Destroy an atom communication data structure and its child structs */
360 static void destroy_atomcomm(pme_atomcomm_t *atc)
365 sfree(atc->node_dest);
366 sfree(atc->node_src);
367 for (int i = 0; i < atc->nthread; i++)
369 sfree(atc->count_thread[i]);
371 sfree(atc->count_thread);
373 sfree(atc->buf_index);
376 sfree(atc->coefficient);
382 sfree(atc->thread_idx);
383 for (int i = 0; i < atc->nthread; i++)
385 if (atc->nthread > 1)
387 int *n_ptr = atc->thread_plist[i].n - gmxCacheLineSize;
389 sfree(atc->thread_plist[i].i);
391 sfree(atc->spline[i].ind);
392 for (int d = 0; d < ZZ; d++)
394 sfree(atc->spline[i].theta[d]);
395 sfree(atc->spline[i].dtheta[d]);
397 sfree_aligned(atc->spline[i].ptr_dtheta_z);
398 sfree_aligned(atc->spline[i].ptr_theta_z);
400 if (atc->nthread > 1)
402 sfree(atc->thread_plist);
407 /*! \brief Initialize data structure for communication */
409 init_overlap_comm(pme_overlap_t * ol,
429 /* Linear translation of the PME grid won't affect reciprocal space
430 * calculations, so to optimize we only interpolate "upwards",
431 * which also means we only have to consider overlap in one direction.
432 * I.e., particles on this node might also be spread to grid indices
433 * that belong to higher nodes (modulo nnodes)
436 ol->s2g0.resize(ol->nnodes + 1);
437 ol->s2g1.resize(ol->nnodes);
440 fprintf(debug, "PME slab boundaries:");
442 for (int i = 0; i < nnodes; i++)
444 /* s2g0 the local interpolation grid start.
445 * s2g1 the local interpolation grid end.
446 * Since in calc_pidx we divide particles, and not grid lines,
447 * spatially uniform along dimension x or y, we need to round
448 * s2g0 down and s2g1 up.
450 ol->s2g0[i] = (i * ndata + 0) / nnodes;
451 ol->s2g1[i] = ((i + 1) * ndata + nnodes - 1) / nnodes + norder - 1;
455 fprintf(debug, " %3d %3d", ol->s2g0[i], ol->s2g1[i]);
458 ol->s2g0[nnodes] = ndata;
461 fprintf(debug, "\n");
464 /* Determine with how many nodes we need to communicate the grid overlap */
465 int testRankCount = 0;
470 for (int i = 0; i < nnodes; i++)
472 if ((i + testRankCount < nnodes && ol->s2g1[i] > ol->s2g0[i + testRankCount]) ||
473 (i + testRankCount >= nnodes && ol->s2g1[i] > ol->s2g0[i + testRankCount - nnodes] + ndata))
479 while (bCont && testRankCount < nnodes);
481 ol->comm_data.resize(testRankCount - 1);
484 for (size_t b = 0; b < ol->comm_data.size(); b++)
486 pme_grid_comm_t *pgc = &ol->comm_data[b];
489 pgc->send_id = (ol->nodeid + (b + 1)) % ol->nnodes;
490 int fft_start = ol->s2g0[pgc->send_id];
491 int fft_end = ol->s2g0[pgc->send_id + 1];
492 if (pgc->send_id < nodeid)
497 int send_index1 = ol->s2g1[nodeid];
498 send_index1 = std::min(send_index1, fft_end);
499 pgc->send_index0 = fft_start;
500 pgc->send_nindex = std::max(0, send_index1 - pgc->send_index0);
501 ol->send_size += pgc->send_nindex;
503 /* We always start receiving to the first index of our slab */
504 pgc->recv_id = (ol->nodeid - (b + 1) + ol->nnodes) % ol->nnodes;
505 fft_start = ol->s2g0[ol->nodeid];
506 fft_end = ol->s2g0[ol->nodeid + 1];
507 int recv_index1 = ol->s2g1[pgc->recv_id];
508 if (pgc->recv_id > nodeid)
510 recv_index1 -= ndata;
512 recv_index1 = std::min(recv_index1, fft_end);
513 pgc->recv_index0 = fft_start;
514 pgc->recv_nindex = std::max(0, recv_index1 - pgc->recv_index0);
518 /* Communicate the buffer sizes to receive */
519 for (size_t b = 0; b < ol->comm_data.size(); b++)
521 MPI_Sendrecv(&ol->send_size, 1, MPI_INT, ol->comm_data[b].send_id, b,
522 &ol->comm_data[b].recv_size, 1, MPI_INT, ol->comm_data[b].recv_id, b,
523 ol->mpi_comm, &stat);
527 /* For non-divisible grid we need pme_order iso pme_order-1 */
528 ol->sendbuf.resize(norder * commplainsize);
529 ol->recvbuf.resize(norder * commplainsize);
532 int minimalPmeGridSize(int pmeOrder)
534 /* The actual grid size limitations are:
535 * serial: >= pme_order
536 * DD, no OpenMP: >= 2*(pme_order - 1)
537 * DD, OpenMP: >= pme_order + 1
538 * But we use the maximum for simplicity since in practice there is not
539 * much performance difference between pme_order and 2*(pme_order -1).
541 int minimalSize = 2*(pmeOrder - 1);
543 GMX_RELEASE_ASSERT(pmeOrder >= 3, "pmeOrder has to be >= 3");
544 GMX_RELEASE_ASSERT(minimalSize >= pmeOrder + 1, "The grid size should be >= pmeOrder + 1");
549 bool gmx_pme_check_restrictions(int pme_order,
550 int nkx, int nky, int nkz,
551 int numPmeDomainsAlongX,
555 if (pme_order > PME_ORDER_MAX)
562 std::string message = gmx::formatString(
563 "pme_order (%d) is larger than the maximum allowed value (%d). Modify and recompile the code if you really need such a high order.",
564 pme_order, PME_ORDER_MAX);
565 GMX_THROW(gmx::InconsistentInputError(message));
568 const int minGridSize = minimalPmeGridSize(pme_order);
569 if (nkx < minGridSize ||
577 std::string message = gmx::formatString(
578 "The PME grid sizes need to be >= 2*(pme_order-1) (%d)",
580 GMX_THROW(gmx::InconsistentInputError(message));
583 /* Check for a limitation of the (current) sum_fftgrid_dd code.
584 * We only allow multiple communication pulses in dim 1, not in dim 0.
586 if (useThreads && (nkx < numPmeDomainsAlongX*pme_order &&
587 nkx != numPmeDomainsAlongX*(pme_order - 1)))
593 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.",
594 nkx/static_cast<double>(numPmeDomainsAlongX), pme_order);
600 /*! \brief Round \p enumerator */
601 static int div_round_up(int enumerator, int denominator)
603 return (enumerator + denominator - 1)/denominator;
606 gmx_pme_t *gmx_pme_init(const t_commrec *cr,
607 const NumPmeDomains &numPmeDomains,
608 const t_inputrec *ir,
610 gmx_bool bFreeEnergy_q,
611 gmx_bool bFreeEnergy_lj,
612 gmx_bool bReproducible,
618 gmx_device_info_t *gpuInfo,
619 PmeGpuProgramHandle pmeGpuProgram,
620 const gmx::MDLogger & /*mdlog*/)
622 int use_threads, sum_use_threads, i;
627 fprintf(debug, "Creating PME data structures.\n");
630 gmx::unique_cptr<gmx_pme_t, gmx_pme_destroy> pme(new gmx_pme_t());
632 pme->sum_qgrid_tmp = nullptr;
633 pme->sum_qgrid_dd_tmp = nullptr;
640 pme->nnodes_major = numPmeDomains.x;
641 pme->nnodes_minor = numPmeDomains.y;
644 if (numPmeDomains.x*numPmeDomains.y > 1)
646 pme->mpi_comm = cr->mpi_comm_mygroup;
648 MPI_Comm_rank(pme->mpi_comm, &pme->nodeid);
649 MPI_Comm_size(pme->mpi_comm, &pme->nnodes);
650 if (pme->nnodes != numPmeDomains.x*numPmeDomains.y)
652 gmx_incons("PME rank count mismatch");
657 pme->mpi_comm = MPI_COMM_NULL;
661 if (pme->nnodes == 1)
664 pme->mpi_comm_d[0] = MPI_COMM_NULL;
665 pme->mpi_comm_d[1] = MPI_COMM_NULL;
668 pme->nodeid_major = 0;
669 pme->nodeid_minor = 0;
671 pme->mpi_comm_d[0] = pme->mpi_comm_d[1] = MPI_COMM_NULL;
676 if (numPmeDomains.y == 1)
679 pme->mpi_comm_d[0] = pme->mpi_comm;
680 pme->mpi_comm_d[1] = MPI_COMM_NULL;
683 pme->nodeid_major = pme->nodeid;
684 pme->nodeid_minor = 0;
687 else if (numPmeDomains.x == 1)
690 pme->mpi_comm_d[0] = MPI_COMM_NULL;
691 pme->mpi_comm_d[1] = pme->mpi_comm;
694 pme->nodeid_major = 0;
695 pme->nodeid_minor = pme->nodeid;
699 if (pme->nnodes % numPmeDomains.x != 0)
701 gmx_incons("For 2D PME decomposition, #PME ranks must be divisible by the number of domains along x");
706 MPI_Comm_split(pme->mpi_comm, pme->nodeid % numPmeDomains.y,
707 pme->nodeid, &pme->mpi_comm_d[0]); /* My communicator along major dimension */
708 MPI_Comm_split(pme->mpi_comm, pme->nodeid/numPmeDomains.y,
709 pme->nodeid, &pme->mpi_comm_d[1]); /* My communicator along minor dimension */
711 MPI_Comm_rank(pme->mpi_comm_d[0], &pme->nodeid_major);
712 MPI_Comm_size(pme->mpi_comm_d[0], &pme->nnodes_major);
713 MPI_Comm_rank(pme->mpi_comm_d[1], &pme->nodeid_minor);
714 MPI_Comm_size(pme->mpi_comm_d[1], &pme->nnodes_minor);
717 pme->bPPnode = thisRankHasDuty(cr, DUTY_PP);
720 pme->nthread = nthread;
722 /* Check if any of the PME MPI ranks uses threads */
723 use_threads = (pme->nthread > 1 ? 1 : 0);
727 MPI_Allreduce(&use_threads, &sum_use_threads, 1, MPI_INT,
728 MPI_SUM, pme->mpi_comm);
733 sum_use_threads = use_threads;
735 pme->bUseThreads = (sum_use_threads > 0);
737 if (ir->ePBC == epbcSCREW)
739 gmx_fatal(FARGS, "pme does not (yet) work with pbc = screw");
743 * It is likely that the current gmx_pme_do() routine supports calculating
744 * only Coulomb or LJ while gmx_pme_init() configures for both,
745 * but that has never been tested.
746 * It is likely that the current gmx_pme_do() routine supports calculating,
747 * not calculating free-energy for Coulomb and/or LJ while gmx_pme_init()
748 * configures with free-energy, but that has never been tested.
750 pme->doCoulomb = EEL_PME(ir->coulombtype);
751 pme->doLJ = EVDW_PME(ir->vdwtype);
752 pme->bFEP_q = ((ir->efep != efepNO) && bFreeEnergy_q);
753 pme->bFEP_lj = ((ir->efep != efepNO) && bFreeEnergy_lj);
754 pme->bFEP = (pme->bFEP_q || pme->bFEP_lj);
758 pme->bP3M = (ir->coulombtype == eelP3M_AD || getenv("GMX_PME_P3M") != nullptr);
759 pme->pme_order = ir->pme_order;
760 pme->ewaldcoeff_q = ewaldcoeff_q;
761 pme->ewaldcoeff_lj = ewaldcoeff_lj;
763 /* Always constant electrostatics coefficients */
764 pme->epsilon_r = ir->epsilon_r;
766 /* Always constant LJ coefficients */
767 pme->ljpme_combination_rule = ir->ljpme_combination_rule;
769 // The box requires scaling with nwalls = 2, we store that condition as well
770 // as the scaling factor
771 delete pme->boxScaler;
772 pme->boxScaler = new EwaldBoxZScaler(*ir);
774 /* If we violate restrictions, generate a fatal error here */
775 gmx_pme_check_restrictions(pme->pme_order,
776 pme->nkx, pme->nky, pme->nkz,
786 MPI_Type_contiguous(DIM, GMX_MPI_REAL, &(pme->rvec_mpi));
787 MPI_Type_commit(&(pme->rvec_mpi));
790 /* Note that the coefficient spreading and force gathering, which usually
791 * takes about the same amount of time as FFT+solve_pme,
792 * is always fully load balanced
793 * (unless the coefficient distribution is inhomogeneous).
796 imbal = estimate_pme_load_imbalance(pme.get());
797 if (imbal >= 1.2 && pme->nodeid_major == 0 && pme->nodeid_minor == 0)
801 "NOTE: The load imbalance in PME FFT and solve is %d%%.\n"
802 " For optimal PME load balancing\n"
803 " PME grid_x (%d) and grid_y (%d) should be divisible by #PME_ranks_x (%d)\n"
804 " and PME grid_y (%d) and grid_z (%d) should be divisible by #PME_ranks_y (%d)\n"
806 static_cast<int>((imbal-1)*100 + 0.5),
807 pme->nkx, pme->nky, pme->nnodes_major,
808 pme->nky, pme->nkz, pme->nnodes_minor);
812 /* For non-divisible grid we need pme_order iso pme_order-1 */
813 /* In sum_qgrid_dd x overlap is copied in place: take padding into account.
814 * y is always copied through a buffer: we don't need padding in z,
815 * but we do need the overlap in x because of the communication order.
817 init_overlap_comm(&pme->overlap[0], pme->pme_order,
821 pme->nnodes_major, pme->nodeid_major,
823 (div_round_up(pme->nky, pme->nnodes_minor)+pme->pme_order)*(pme->nkz+pme->pme_order-1));
825 /* Along overlap dim 1 we can send in multiple pulses in sum_fftgrid_dd.
826 * We do this with an offset buffer of equal size, so we need to allocate
827 * extra for the offset. That's what the (+1)*pme->nkz is for.
829 init_overlap_comm(&pme->overlap[1], pme->pme_order,
833 pme->nnodes_minor, pme->nodeid_minor,
835 (div_round_up(pme->nkx, pme->nnodes_major)+pme->pme_order+1)*pme->nkz);
837 /* Double-check for a limitation of the (current) sum_fftgrid_dd code.
838 * Note that gmx_pme_check_restrictions checked for this already.
840 if (pme->bUseThreads && (pme->overlap[0].comm_data.size() > 1))
842 gmx_incons("More than one communication pulse required for grid overlap communication along the major dimension while using threads");
845 snew(pme->bsp_mod[XX], pme->nkx);
846 snew(pme->bsp_mod[YY], pme->nky);
847 snew(pme->bsp_mod[ZZ], pme->nkz);
849 pme->gpu = pmeGpu; /* Carrying over the single GPU structure */
850 pme->runMode = runMode;
852 /* The required size of the interpolation grid, including overlap.
853 * The allocated size (pmegrid_n?) might be slightly larger.
855 pme->pmegrid_nx = pme->overlap[0].s2g1[pme->nodeid_major] -
856 pme->overlap[0].s2g0[pme->nodeid_major];
857 pme->pmegrid_ny = pme->overlap[1].s2g1[pme->nodeid_minor] -
858 pme->overlap[1].s2g0[pme->nodeid_minor];
859 pme->pmegrid_nz_base = pme->nkz;
860 pme->pmegrid_nz = pme->pmegrid_nz_base + pme->pme_order - 1;
861 set_grid_alignment(&pme->pmegrid_nz, pme->pme_order);
862 pme->pmegrid_start_ix = pme->overlap[0].s2g0[pme->nodeid_major];
863 pme->pmegrid_start_iy = pme->overlap[1].s2g0[pme->nodeid_minor];
864 pme->pmegrid_start_iz = 0;
866 make_gridindex_to_localindex(pme->nkx,
867 pme->pmegrid_start_ix,
868 pme->pmegrid_nx - (pme->pme_order-1),
869 &pme->nnx, &pme->fshx);
870 make_gridindex_to_localindex(pme->nky,
871 pme->pmegrid_start_iy,
872 pme->pmegrid_ny - (pme->pme_order-1),
873 &pme->nny, &pme->fshy);
874 make_gridindex_to_localindex(pme->nkz,
875 pme->pmegrid_start_iz,
876 pme->pmegrid_nz_base,
877 &pme->nnz, &pme->fshz);
879 pme->spline_work = make_pme_spline_work(pme->pme_order);
884 /* It doesn't matter if we allocate too many grids here,
885 * we only allocate and use the ones we need.
889 pme->ngrids = ((ir->ljpme_combination_rule == eljpmeLB) ? DO_Q_AND_LJ_LB : DO_Q_AND_LJ);
895 snew(pme->fftgrid, pme->ngrids);
896 snew(pme->cfftgrid, pme->ngrids);
897 snew(pme->pfft_setup, pme->ngrids);
899 for (i = 0; i < pme->ngrids; ++i)
901 if ((i < DO_Q && pme->doCoulomb && (i == 0 ||
903 (i >= DO_Q && pme->doLJ && (i == 2 ||
905 ir->ljpme_combination_rule == eljpmeLB)))
907 pmegrids_init(&pme->pmegrid[i],
908 pme->pmegrid_nx, pme->pmegrid_ny, pme->pmegrid_nz,
909 pme->pmegrid_nz_base,
913 pme->overlap[0].s2g1[pme->nodeid_major]-pme->overlap[0].s2g0[pme->nodeid_major+1],
914 pme->overlap[1].s2g1[pme->nodeid_minor]-pme->overlap[1].s2g0[pme->nodeid_minor+1]);
915 /* This routine will allocate the grid data to fit the FFTs */
916 const auto allocateRealGridForGpu = (pme->runMode == PmeRunMode::Mixed) ? gmx::PinningPolicy::PinnedIfSupported : gmx::PinningPolicy::CannotBePinned;
917 gmx_parallel_3dfft_init(&pme->pfft_setup[i], ndata,
918 &pme->fftgrid[i], &pme->cfftgrid[i],
920 bReproducible, pme->nthread, allocateRealGridForGpu);
927 /* Use plain SPME B-spline interpolation */
928 make_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
932 /* Use the P3M grid-optimized influence function */
933 make_p3m_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
936 /* Use atc[0] for spreading */
937 init_atomcomm(pme.get(), &pme->atc[0], numPmeDomains.x > 1 ? 0 : 1, TRUE);
938 if (pme->ndecompdim >= 2)
940 init_atomcomm(pme.get(), &pme->atc[1], 1, FALSE);
943 if (pme->nnodes == 1)
945 pme->atc[0].n = homenr;
946 pme_realloc_atomcomm_things(&pme->atc[0]);
949 pme->lb_buf1 = nullptr;
950 pme->lb_buf2 = nullptr;
951 pme->lb_buf_nalloc = 0;
953 if (pme_gpu_active(pme.get()))
957 // Initial check of validity of the data
958 std::string errorString;
959 bool canRunOnGpu = pme_gpu_check_restrictions(pme.get(), &errorString);
962 GMX_THROW(gmx::NotImplementedError(errorString));
966 pme_gpu_reinit(pme.get(), gpuInfo, pmeGpuProgram);
969 pme_init_all_work(&pme->solve_work, pme->nthread, pme->nkx);
971 // no exception was thrown during the init, so we hand over the PME structure handle
972 return pme.release();
975 void gmx_pme_reinit(struct gmx_pme_t **pmedata,
977 struct gmx_pme_t * pme_src,
978 const t_inputrec * ir,
979 const ivec grid_size,
985 // Create a copy of t_inputrec fields that are used in gmx_pme_init().
986 // TODO: This would be better as just copying a sub-structure that contains
987 // all the PME parameters and nothing else.
990 irc.coulombtype = ir->coulombtype;
991 irc.vdwtype = ir->vdwtype;
993 irc.pme_order = ir->pme_order;
994 irc.epsilon_r = ir->epsilon_r;
995 irc.ljpme_combination_rule = ir->ljpme_combination_rule;
996 irc.nkx = grid_size[XX];
997 irc.nky = grid_size[YY];
998 irc.nkz = grid_size[ZZ];
1000 if (pme_src->nnodes == 1)
1002 homenr = pme_src->atc[0].n;
1011 const gmx::MDLogger dummyLogger;
1012 // This is reinit which is currently only changing grid size/coefficients,
1013 // so we don't expect the actual logging.
1014 // TODO: when PME is an object, it should take reference to mdlog on construction and save it.
1015 GMX_ASSERT(pmedata, "Invalid PME pointer");
1016 NumPmeDomains numPmeDomains = { pme_src->nnodes_major, pme_src->nnodes_minor };
1017 *pmedata = gmx_pme_init(cr, numPmeDomains,
1018 &irc, homenr, pme_src->bFEP_q, pme_src->bFEP_lj, FALSE, ewaldcoeff_q, ewaldcoeff_lj,
1019 pme_src->nthread, pme_src->runMode, pme_src->gpu, nullptr, nullptr, dummyLogger);
1020 //TODO this is mostly passing around current values
1022 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1024 /* We can easily reuse the allocated pme grids in pme_src */
1025 reuse_pmegrids(&pme_src->pmegrid[PME_GRID_QA], &(*pmedata)->pmegrid[PME_GRID_QA]);
1026 /* We would like to reuse the fft grids, but that's harder */
1029 void gmx_pme_calc_energy(struct gmx_pme_t *pme, int n, rvec *x, real *q, real *V)
1031 pme_atomcomm_t *atc;
1034 if (pme->nnodes > 1)
1036 gmx_incons("gmx_pme_calc_energy called in parallel");
1040 gmx_incons("gmx_pme_calc_energy with free energy");
1043 atc = &pme->atc_energy;
1045 if (atc->spline == nullptr)
1047 snew(atc->spline, atc->nthread);
1050 atc->bSpread = TRUE;
1051 atc->pme_order = pme->pme_order;
1053 pme_realloc_atomcomm_things(atc);
1055 atc->coefficient = q;
1057 /* We only use the A-charges grid */
1058 grid = &pme->pmegrid[PME_GRID_QA];
1060 /* Only calculate the spline coefficients, don't actually spread */
1061 spread_on_grid(pme, atc, nullptr, TRUE, FALSE, pme->fftgrid[PME_GRID_QA], FALSE, PME_GRID_QA);
1063 *V = gather_energy_bsplines(pme, grid->grid.grid, atc);
1066 /*! \brief Calculate initial Lorentz-Berthelot coefficients for LJ-PME */
1068 calc_initial_lb_coeffs(struct gmx_pme_t *pme, const real *local_c6, const real *local_sigma)
1071 for (i = 0; i < pme->atc[0].n; ++i)
1074 sigma4 = local_sigma[i];
1075 sigma4 = sigma4*sigma4;
1076 sigma4 = sigma4*sigma4;
1077 pme->atc[0].coefficient[i] = local_c6[i] / sigma4;
1081 /*! \brief Calculate next Lorentz-Berthelot coefficients for LJ-PME */
1083 calc_next_lb_coeffs(struct gmx_pme_t *pme, const real *local_sigma)
1087 for (i = 0; i < pme->atc[0].n; ++i)
1089 pme->atc[0].coefficient[i] *= local_sigma[i];
1093 int gmx_pme_do(struct gmx_pme_t *pme,
1094 int start, int homenr,
1096 real chargeA[], real chargeB[],
1097 real c6A[], real c6B[],
1098 real sigmaA[], real sigmaB[],
1099 matrix box, const t_commrec *cr,
1100 int maxshift_x, int maxshift_y,
1101 t_nrnb *nrnb, gmx_wallcycle *wcycle,
1102 matrix vir_q, matrix vir_lj,
1103 real *energy_q, real *energy_lj,
1104 real lambda_q, real lambda_lj,
1105 real *dvdlambda_q, real *dvdlambda_lj,
1108 GMX_ASSERT(pme->runMode == PmeRunMode::CPU, "gmx_pme_do should not be called on the GPU PME run.");
1110 int d, i, j, npme, grid_index, max_grid_index;
1112 pme_atomcomm_t *atc = nullptr;
1113 pmegrids_t *pmegrid = nullptr;
1114 real *grid = nullptr;
1116 real *coefficient = nullptr;
1121 gmx_parallel_3dfft_t pfft_setup;
1123 t_complex * cfftgrid;
1125 gmx_bool bFirst, bDoSplines;
1127 int fep_states_lj = pme->bFEP_lj ? 2 : 1;
1128 const gmx_bool bCalcEnerVir = (flags & GMX_PME_CALC_ENER_VIR) != 0;
1129 const gmx_bool bBackFFT = (flags & (GMX_PME_CALC_F | GMX_PME_CALC_POT)) != 0;
1130 const gmx_bool bCalcF = (flags & GMX_PME_CALC_F) != 0;
1132 assert(pme->nnodes > 0);
1133 assert(pme->nnodes == 1 || pme->ndecompdim > 0);
1135 if (pme->nnodes > 1)
1139 if (atc->npd > atc->pd_nalloc)
1141 atc->pd_nalloc = over_alloc_dd(atc->npd);
1142 srenew(atc->pd, atc->pd_nalloc);
1144 for (d = pme->ndecompdim-1; d >= 0; d--)
1147 atc->maxshift = (atc->dimind == 0 ? maxshift_x : maxshift_y);
1153 /* This could be necessary for TPI */
1154 pme->atc[0].n = homenr;
1155 if (DOMAINDECOMP(cr))
1157 pme_realloc_atomcomm_things(atc);
1164 pme->boxScaler->scaleBox(box, scaledBox);
1166 gmx::invertBoxMatrix(scaledBox, pme->recipbox);
1169 /* For simplicity, we construct the splines for all particles if
1170 * more than one PME calculations is needed. Some optimization
1171 * could be done by keeping track of which atoms have splines
1172 * constructed, and construct new splines on each pass for atoms
1173 * that don't yet have them.
1176 bDoSplines = pme->bFEP || (pme->doCoulomb && pme->doLJ);
1178 /* We need a maximum of four separate PME calculations:
1179 * grid_index=0: Coulomb PME with charges from state A
1180 * grid_index=1: Coulomb PME with charges from state B
1181 * grid_index=2: LJ PME with C6 from state A
1182 * grid_index=3: LJ PME with C6 from state B
1183 * For Lorentz-Berthelot combination rules, a separate loop is used to
1184 * calculate all the terms
1187 /* If we are doing LJ-PME with LB, we only do Q here */
1188 max_grid_index = (pme->ljpme_combination_rule == eljpmeLB) ? DO_Q : DO_Q_AND_LJ;
1190 for (grid_index = 0; grid_index < max_grid_index; ++grid_index)
1192 /* Check if we should do calculations at this grid_index
1193 * If grid_index is odd we should be doing FEP
1194 * If grid_index < 2 we should be doing electrostatic PME
1195 * If grid_index >= 2 we should be doing LJ-PME
1197 if ((grid_index < DO_Q && (!pme->doCoulomb ||
1198 (grid_index == 1 && !pme->bFEP_q))) ||
1199 (grid_index >= DO_Q && (!pme->doLJ ||
1200 (grid_index == 3 && !pme->bFEP_lj))))
1204 /* Unpack structure */
1205 pmegrid = &pme->pmegrid[grid_index];
1206 fftgrid = pme->fftgrid[grid_index];
1207 cfftgrid = pme->cfftgrid[grid_index];
1208 pfft_setup = pme->pfft_setup[grid_index];
1211 case 0: coefficient = chargeA + start; break;
1212 case 1: coefficient = chargeB + start; break;
1213 case 2: coefficient = c6A + start; break;
1214 case 3: coefficient = c6B + start; break;
1217 grid = pmegrid->grid.grid;
1221 fprintf(debug, "PME: number of ranks = %d, rank = %d\n",
1222 cr->nnodes, cr->nodeid);
1223 fprintf(debug, "Grid = %p\n", static_cast<void*>(grid));
1224 if (grid == nullptr)
1226 gmx_fatal(FARGS, "No grid!");
1230 if (pme->nnodes == 1)
1232 atc->coefficient = coefficient;
1236 wallcycle_start(wcycle, ewcPME_REDISTXF);
1237 do_redist_pos_coeffs(pme, cr, start, homenr, bFirst, x, coefficient);
1239 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1244 fprintf(debug, "Rank= %6d, pme local particles=%6d\n",
1245 cr->nodeid, atc->n);
1248 if (flags & GMX_PME_SPREAD)
1250 wallcycle_start(wcycle, ewcPME_SPREAD);
1252 /* Spread the coefficients on a grid */
1253 spread_on_grid(pme, &pme->atc[0], pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1257 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
1259 inc_nrnb(nrnb, eNR_SPREADBSP,
1260 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
1262 if (!pme->bUseThreads)
1264 wrap_periodic_pmegrid(pme, grid);
1266 /* sum contributions to local grid from other nodes */
1268 if (pme->nnodes > 1)
1270 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1274 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1277 wallcycle_stop(wcycle, ewcPME_SPREAD);
1279 /* TODO If the OpenMP and single-threaded implementations
1280 converge, then spread_on_grid() and
1281 copy_pmegrid_to_fftgrid() will perhaps live in the same
1286 /* Here we start a large thread parallel region */
1287 #pragma omp parallel num_threads(pme->nthread) private(thread)
1291 thread = gmx_omp_get_thread_num();
1292 if (flags & GMX_PME_SOLVE)
1299 wallcycle_start(wcycle, ewcPME_FFT);
1301 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
1305 wallcycle_stop(wcycle, ewcPME_FFT);
1308 /* solve in k-space for our local cells */
1311 wallcycle_start(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
1313 if (grid_index < DO_Q)
1316 solve_pme_yzx(pme, cfftgrid,
1317 scaledBox[XX][XX]*scaledBox[YY][YY]*scaledBox[ZZ][ZZ],
1319 pme->nthread, thread);
1324 solve_pme_lj_yzx(pme, &cfftgrid, FALSE,
1325 scaledBox[XX][XX]*scaledBox[YY][YY]*scaledBox[ZZ][ZZ],
1327 pme->nthread, thread);
1332 wallcycle_stop(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
1333 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1342 wallcycle_start(wcycle, ewcPME_FFT);
1344 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
1348 wallcycle_stop(wcycle, ewcPME_FFT);
1351 if (pme->nodeid == 0)
1353 real ntot = pme->nkx*pme->nky*pme->nkz;
1354 npme = static_cast<int>(ntot*std::log(ntot)/std::log(2.0));
1355 inc_nrnb(nrnb, eNR_FFT, 2*npme);
1358 /* Note: this wallcycle region is closed below
1359 outside an OpenMP region, so take care if
1360 refactoring code here. */
1361 wallcycle_start(wcycle, ewcPME_GATHER);
1364 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1366 } GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1368 /* End of thread parallel section.
1369 * With MPI we have to synchronize here before gmx_sum_qgrid_dd.
1374 /* distribute local grid to all nodes */
1376 if (pme->nnodes > 1)
1378 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1382 unwrap_periodic_pmegrid(pme, grid);
1387 /* interpolate forces for our local atoms */
1390 /* If we are running without parallelization,
1391 * atc->f is the actual force array, not a buffer,
1392 * therefore we should not clear it.
1394 lambda = grid_index < DO_Q ? lambda_q : lambda_lj;
1395 bClearF = (bFirst && PAR(cr));
1396 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1397 for (thread = 0; thread < pme->nthread; thread++)
1401 gather_f_bsplines(pme, grid, bClearF, atc,
1402 &atc->spline[thread],
1403 pme->bFEP ? (grid_index % 2 == 0 ? 1.0-lambda : lambda) : 1.0);
1405 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1409 inc_nrnb(nrnb, eNR_GATHERFBSP,
1410 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
1411 /* Note: this wallcycle region is opened above inside an OpenMP
1412 region, so take care if refactoring code here. */
1413 wallcycle_stop(wcycle, ewcPME_GATHER);
1418 /* This should only be called on the master thread
1419 * and after the threads have synchronized.
1423 get_pme_ener_vir_q(pme->solve_work, pme->nthread, &energy_AB[grid_index], vir_AB[grid_index]);
1427 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &energy_AB[grid_index], vir_AB[grid_index]);
1431 } /* of grid_index-loop */
1433 /* For Lorentz-Berthelot combination rules in LJ-PME, we need to calculate
1436 if (pme->doLJ && pme->ljpme_combination_rule == eljpmeLB)
1438 /* Loop over A- and B-state if we are doing FEP */
1439 for (fep_state = 0; fep_state < fep_states_lj; ++fep_state)
1441 real *local_c6 = nullptr, *local_sigma = nullptr, *RedistC6 = nullptr, *RedistSigma = nullptr;
1442 if (pme->nnodes == 1)
1444 if (pme->lb_buf1 == nullptr)
1446 pme->lb_buf_nalloc = pme->atc[0].n;
1447 snew(pme->lb_buf1, pme->lb_buf_nalloc);
1449 pme->atc[0].coefficient = pme->lb_buf1;
1454 local_sigma = sigmaA;
1458 local_sigma = sigmaB;
1461 gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1471 RedistSigma = sigmaA;
1475 RedistSigma = sigmaB;
1478 gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
1480 wallcycle_start(wcycle, ewcPME_REDISTXF);
1482 do_redist_pos_coeffs(pme, cr, start, homenr, bFirst, x, RedistC6);
1483 if (pme->lb_buf_nalloc < atc->n)
1485 pme->lb_buf_nalloc = atc->nalloc;
1486 srenew(pme->lb_buf1, pme->lb_buf_nalloc);
1487 srenew(pme->lb_buf2, pme->lb_buf_nalloc);
1489 local_c6 = pme->lb_buf1;
1490 for (i = 0; i < atc->n; ++i)
1492 local_c6[i] = atc->coefficient[i];
1495 do_redist_pos_coeffs(pme, cr, start, homenr, FALSE, x, RedistSigma);
1496 local_sigma = pme->lb_buf2;
1497 for (i = 0; i < atc->n; ++i)
1499 local_sigma[i] = atc->coefficient[i];
1502 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1504 calc_initial_lb_coeffs(pme, local_c6, local_sigma);
1506 /*Seven terms in LJ-PME with LB, grid_index < 2 reserved for electrostatics*/
1507 for (grid_index = 2; grid_index < 9; ++grid_index)
1509 /* Unpack structure */
1510 pmegrid = &pme->pmegrid[grid_index];
1511 fftgrid = pme->fftgrid[grid_index];
1512 pfft_setup = pme->pfft_setup[grid_index];
1513 calc_next_lb_coeffs(pme, local_sigma);
1514 grid = pmegrid->grid.grid;
1516 if (flags & GMX_PME_SPREAD)
1518 wallcycle_start(wcycle, ewcPME_SPREAD);
1519 /* Spread the c6 on a grid */
1520 spread_on_grid(pme, &pme->atc[0], pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
1524 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
1527 inc_nrnb(nrnb, eNR_SPREADBSP,
1528 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
1529 if (pme->nthread == 1)
1531 wrap_periodic_pmegrid(pme, grid);
1532 /* sum contributions to local grid from other nodes */
1534 if (pme->nnodes > 1)
1536 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
1539 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
1541 wallcycle_stop(wcycle, ewcPME_SPREAD);
1543 /*Here we start a large thread parallel region*/
1544 #pragma omp parallel num_threads(pme->nthread) private(thread)
1548 thread = gmx_omp_get_thread_num();
1549 if (flags & GMX_PME_SOLVE)
1554 wallcycle_start(wcycle, ewcPME_FFT);
1557 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
1561 wallcycle_stop(wcycle, ewcPME_FFT);
1565 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1569 if (flags & GMX_PME_SOLVE)
1571 /* solve in k-space for our local cells */
1572 #pragma omp parallel num_threads(pme->nthread) private(thread)
1577 thread = gmx_omp_get_thread_num();
1580 wallcycle_start(wcycle, ewcLJPME);
1584 solve_pme_lj_yzx(pme, &pme->cfftgrid[2], TRUE,
1585 scaledBox[XX][XX]*scaledBox[YY][YY]*scaledBox[ZZ][ZZ],
1587 pme->nthread, thread);
1590 wallcycle_stop(wcycle, ewcLJPME);
1591 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
1594 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1600 /* This should only be called on the master thread and
1601 * after the threads have synchronized.
1603 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &energy_AB[2+fep_state], vir_AB[2+fep_state]);
1608 bFirst = !pme->doCoulomb;
1609 calc_initial_lb_coeffs(pme, local_c6, local_sigma);
1610 for (grid_index = 8; grid_index >= 2; --grid_index)
1612 /* Unpack structure */
1613 pmegrid = &pme->pmegrid[grid_index];
1614 fftgrid = pme->fftgrid[grid_index];
1615 pfft_setup = pme->pfft_setup[grid_index];
1616 grid = pmegrid->grid.grid;
1617 calc_next_lb_coeffs(pme, local_sigma);
1618 #pragma omp parallel num_threads(pme->nthread) private(thread)
1622 thread = gmx_omp_get_thread_num();
1626 wallcycle_start(wcycle, ewcPME_FFT);
1629 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
1633 wallcycle_stop(wcycle, ewcPME_FFT);
1636 if (pme->nodeid == 0)
1638 real ntot = pme->nkx*pme->nky*pme->nkz;
1639 npme = static_cast<int>(ntot*std::log(ntot)/std::log(2.0));
1640 inc_nrnb(nrnb, eNR_FFT, 2*npme);
1642 wallcycle_start(wcycle, ewcPME_GATHER);
1645 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
1647 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1648 } /*#pragma omp parallel*/
1650 /* distribute local grid to all nodes */
1652 if (pme->nnodes > 1)
1654 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
1658 unwrap_periodic_pmegrid(pme, grid);
1662 /* interpolate forces for our local atoms */
1663 bClearF = (bFirst && PAR(cr));
1664 scale = pme->bFEP ? (fep_state < 1 ? 1.0-lambda_lj : lambda_lj) : 1.0;
1665 scale *= lb_scale_factor[grid_index-2];
1667 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1668 for (thread = 0; thread < pme->nthread; thread++)
1672 gather_f_bsplines(pme, grid, bClearF, &pme->atc[0],
1673 &pme->atc[0].spline[thread],
1676 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1680 inc_nrnb(nrnb, eNR_GATHERFBSP,
1681 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
1683 wallcycle_stop(wcycle, ewcPME_GATHER);
1686 } /* for (grid_index = 8; grid_index >= 2; --grid_index) */
1688 } /* for (fep_state = 0; fep_state < fep_states_lj; ++fep_state) */
1689 } /* if ((flags & GMX_PME_DO_LJ) && pme->ljpme_combination_rule == eljpmeLB) */
1691 if (bCalcF && pme->nnodes > 1)
1693 wallcycle_start(wcycle, ewcPME_REDISTXF);
1694 for (d = 0; d < pme->ndecompdim; d++)
1697 if (d == pme->ndecompdim - 1)
1704 n_d = pme->atc[d+1].n;
1705 f_d = pme->atc[d+1].f;
1707 if (DOMAINDECOMP(cr))
1709 dd_pmeredist_f(pme, atc, n_d, f_d,
1710 d == pme->ndecompdim-1 && pme->bPPnode);
1714 wallcycle_stop(wcycle, ewcPME_REDISTXF);
1723 *energy_q = energy_AB[0];
1724 m_add(vir_q, vir_AB[0], vir_q);
1728 *energy_q = (1.0-lambda_q)*energy_AB[0] + lambda_q*energy_AB[1];
1729 *dvdlambda_q += energy_AB[1] - energy_AB[0];
1730 for (i = 0; i < DIM; i++)
1732 for (j = 0; j < DIM; j++)
1734 vir_q[i][j] += (1.0-lambda_q)*vir_AB[0][i][j] +
1735 lambda_q*vir_AB[1][i][j];
1741 fprintf(debug, "Electrostatic PME mesh energy: %g\n", *energy_q);
1753 *energy_lj = energy_AB[2];
1754 m_add(vir_lj, vir_AB[2], vir_lj);
1758 *energy_lj = (1.0-lambda_lj)*energy_AB[2] + lambda_lj*energy_AB[3];
1759 *dvdlambda_lj += energy_AB[3] - energy_AB[2];
1760 for (i = 0; i < DIM; i++)
1762 for (j = 0; j < DIM; j++)
1764 vir_lj[i][j] += (1.0-lambda_lj)*vir_AB[2][i][j] + lambda_lj*vir_AB[3][i][j];
1770 fprintf(debug, "Lennard-Jones PME mesh energy: %g\n", *energy_lj);
1781 void gmx_pme_destroy(gmx_pme_t *pme)
1788 delete pme->boxScaler;
1797 for (int i = 0; i < pme->ngrids; ++i)
1799 pmegrids_destroy(&pme->pmegrid[i]);
1801 if (pme->pfft_setup)
1803 for (int i = 0; i < pme->ngrids; ++i)
1805 gmx_parallel_3dfft_destroy(pme->pfft_setup[i]);
1808 sfree(pme->fftgrid);
1809 sfree(pme->cfftgrid);
1810 sfree(pme->pfft_setup);
1812 for (int i = 0; i < std::max(1, pme->ndecompdim); i++) //pme->atc[0] is always allocated
1814 destroy_atomcomm(&pme->atc[i]);
1817 for (int i = 0; i < DIM; i++)
1819 sfree(pme->bsp_mod[i]);
1822 sfree(pme->lb_buf1);
1823 sfree(pme->lb_buf2);
1828 if (pme->solve_work)
1830 pme_free_all_work(&pme->solve_work, pme->nthread);
1833 sfree(pme->sum_qgrid_tmp);
1834 sfree(pme->sum_qgrid_dd_tmp);
1836 destroy_pme_spline_work(pme->spline_work);
1838 if (pme_gpu_active(pme) && pme->gpu)
1840 pme_gpu_destroy(pme->gpu);
1846 void gmx_pme_reinit_atoms(const gmx_pme_t *pme, const int nAtoms, const real *charges)
1848 if (pme_gpu_active(pme))
1850 pme_gpu_reinit_atoms(pme->gpu, nAtoms, charges);
1852 // TODO: handle the CPU case here; handle the whole t_mdatoms