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, 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.
37 /* IMPORTANT FOR DEVELOPERS:
39 * Triclinic pme stuff isn't entirely trivial, and we've experienced
40 * some bugs during development (many of them due to me). To avoid
41 * this in the future, please check the following things if you make
42 * changes in this file:
44 * 1. You should obtain identical (at least to the PME precision)
45 * energies, forces, and virial for
46 * a rectangular box and a triclinic one where the z (or y) axis is
47 * tilted a whole box side. For instance you could use these boxes:
49 * rectangular triclinic
54 * 2. You should check the energy conservation in a triclinic box.
56 * It might seem an overkill, but better safe than sorry.
73 #include "gmx_fatal.h"
80 #include "gromacs/fft/parallel_3dfft.h"
81 #include "gromacs/fileio/futil.h"
82 #include "gromacs/fileio/pdbio.h"
83 #include "gromacs/math/gmxcomplex.h"
84 #include "gromacs/timing/cyclecounter.h"
85 #include "gromacs/timing/wallcycle.h"
86 #include "gromacs/utility/gmxmpi.h"
87 #include "gromacs/utility/gmxomp.h"
89 /* Include the SIMD macro file and then check for support */
90 #include "gromacs/simd/simd.h"
91 #include "gromacs/simd/simd_math.h"
92 #ifdef GMX_SIMD_HAVE_REAL
93 /* Turn on arbitrary width SIMD intrinsics for PME solve */
94 # define PME_SIMD_SOLVE
97 #define PME_GRID_QA 0 /* Gridindex for A-state for Q */
98 #define PME_GRID_C6A 2 /* Gridindex for A-state for LJ */
99 #define DO_Q 2 /* Electrostatic grids have index q<2 */
100 #define DO_Q_AND_LJ 4 /* non-LB LJ grids have index 2 <= q < 4 */
101 #define DO_Q_AND_LJ_LB 9 /* With LB rules we need a total of 2+7 grids */
103 /* Pascal triangle coefficients scaled with (1/2)^6 for LJ-PME with LB-rules */
104 const real lb_scale_factor[] = {
105 1.0/64, 6.0/64, 15.0/64, 20.0/64,
106 15.0/64, 6.0/64, 1.0/64
109 /* Pascal triangle coefficients used in solve_pme_lj_yzx, only need to do 4 calculations due to symmetry */
110 const real lb_scale_factor_symm[] = { 2.0/64, 12.0/64, 30.0/64, 20.0/64 };
112 /* Check if we have 4-wide SIMD macro support */
113 #if (defined GMX_SIMD4_HAVE_REAL)
114 /* Do PME spread and gather with 4-wide SIMD.
115 * NOTE: SIMD is only used with PME order 4 and 5 (which are the most common).
117 # define PME_SIMD4_SPREAD_GATHER
119 # if (defined GMX_SIMD_HAVE_LOADU) && (defined GMX_SIMD_HAVE_STOREU)
120 /* With PME-order=4 on x86, unaligned load+store is slightly faster
121 * than doubling all SIMD operations when using aligned load+store.
123 # define PME_SIMD4_UNALIGNED
128 /* #define PRT_FORCE */
129 /* conditions for on the fly time-measurement */
130 /* #define TAKETIME (step > 1 && timesteps < 10) */
131 #define TAKETIME FALSE
133 /* #define PME_TIME_THREADS */
136 #define mpi_type MPI_DOUBLE
138 #define mpi_type MPI_FLOAT
141 #ifdef PME_SIMD4_SPREAD_GATHER
142 # define SIMD4_ALIGNMENT (GMX_SIMD4_WIDTH*sizeof(real))
144 /* We can use any alignment, apart from 0, so we use 4 reals */
145 # define SIMD4_ALIGNMENT (4*sizeof(real))
148 /* GMX_CACHE_SEP should be a multiple of the SIMD and SIMD4 register size
149 * to preserve alignment.
151 #define GMX_CACHE_SEP 64
153 /* We only define a maximum to be able to use local arrays without allocation.
154 * An order larger than 12 should never be needed, even for test cases.
155 * If needed it can be changed here.
157 #define PME_ORDER_MAX 12
159 /* Internal datastructures */
165 int recv_size; /* Receive buffer width, used with OpenMP */
176 int *send_id, *recv_id;
177 int send_size; /* Send buffer width, used with OpenMP */
178 pme_grid_comm_t *comm_data;
184 int *n; /* Cumulative counts of the number of particles per thread */
185 int nalloc; /* Allocation size of i */
186 int *i; /* Particle indices ordered on thread index (n) */
200 int dimind; /* The index of the dimension, 0=x, 1=y */
207 int *node_dest; /* The nodes to send x and q to with DD */
208 int *node_src; /* The nodes to receive x and q from with DD */
209 int *buf_index; /* Index for commnode into the buffers */
216 int *count; /* The number of atoms to send to each node */
218 int *rcount; /* The number of atoms to receive */
225 gmx_bool bSpread; /* These coordinates are used for spreading */
228 rvec *fractx; /* Fractional coordinate relative to
229 * the lower cell boundary
232 int *thread_idx; /* Which thread should spread which charge */
233 thread_plist_t *thread_plist;
234 splinedata_t *spline;
241 ivec ci; /* The spatial location of this grid */
242 ivec n; /* The used size of *grid, including order-1 */
243 ivec offset; /* The grid offset from the full node grid */
244 int order; /* PME spreading order */
245 ivec s; /* The allocated size of *grid, s >= n */
246 real *grid; /* The grid local thread, size n */
250 pmegrid_t grid; /* The full node grid (non thread-local) */
251 int nthread; /* The number of threads operating on this grid */
252 ivec nc; /* The local spatial decomposition over the threads */
253 pmegrid_t *grid_th; /* Array of grids for each thread */
254 real *grid_all; /* Allocated array for the grids in *grid_th */
255 int **g2t; /* The grid to thread index */
256 ivec nthread_comm; /* The number of threads to communicate with */
260 #ifdef PME_SIMD4_SPREAD_GATHER
261 /* Masks for 4-wide SIMD aligned spreading and gathering */
262 gmx_simd4_bool_t mask_S0[6], mask_S1[6];
264 int dummy; /* C89 requires that struct has at least one member */
269 /* work data for solve_pme */
288 typedef struct gmx_pme {
289 int ndecompdim; /* The number of decomposition dimensions */
290 int nodeid; /* Our nodeid in mpi->mpi_comm */
293 int nnodes; /* The number of nodes doing PME */
298 MPI_Comm mpi_comm_d[2]; /* Indexed on dimension, 0=x, 1=y */
300 MPI_Datatype rvec_mpi; /* the pme vector's MPI type */
303 gmx_bool bUseThreads; /* Does any of the PME ranks have nthread>1 ? */
304 int nthread; /* The number of threads doing PME on our rank */
306 gmx_bool bPPnode; /* Node also does particle-particle forces */
308 gmx_bool bFEP; /* Compute Free energy contribution */
312 int nkx, nky, nkz; /* Grid dimensions */
314 gmx_bool bP3M; /* Do P3M: optimize the influence function */
318 int ljpme_combination_rule; /* Type of combination rule in LJ-PME */
320 int ngrids; /* number of grids we maintain for pmegrid, (c)fftgrid and pfft_setups*/
322 pmegrids_t pmegrid[DO_Q_AND_LJ_LB]; /* Grids on which we do spreading/interpolation,
323 * includes overlap Grid indices are ordered as
325 * 0: Coloumb PME, state A
326 * 1: Coloumb PME, state B
328 * This can probably be done in a better way
329 * but this simple hack works for now
331 /* The PME charge spreading grid sizes/strides, includes pme_order-1 */
332 int pmegrid_nx, pmegrid_ny, pmegrid_nz;
333 /* pmegrid_nz might be larger than strictly necessary to ensure
334 * memory alignment, pmegrid_nz_base gives the real base size.
337 /* The local PME grid starting indices */
338 int pmegrid_start_ix, pmegrid_start_iy, pmegrid_start_iz;
340 /* Work data for spreading and gathering */
341 pme_spline_work_t *spline_work;
343 real **fftgrid; /* Grids for FFT. With 1D FFT decomposition this can be a pointer */
344 /* inside the interpolation grid, but separate for 2D PME decomp. */
345 int fftgrid_nx, fftgrid_ny, fftgrid_nz;
347 t_complex **cfftgrid; /* Grids for complex FFT data */
349 int cfftgrid_nx, cfftgrid_ny, cfftgrid_nz;
351 gmx_parallel_3dfft_t *pfft_setup;
353 int *nnx, *nny, *nnz;
354 real *fshx, *fshy, *fshz;
356 pme_atomcomm_t atc[2]; /* Indexed on decomposition index */
359 /* Buffers to store data for local atoms for L-B combination rule
360 * calculations in LJ-PME. lb_buf1 stores either the coefficients
361 * for spreading/gathering (in serial), or the C6 parameters for
362 * local atoms (in parallel). lb_buf2 is only used in parallel,
363 * and stores the sigma values for local atoms. */
364 real *lb_buf1, *lb_buf2;
365 int lb_buf_nalloc; /* Allocation size for the above buffers. */
367 pme_overlap_t overlap[2]; /* Indexed on dimension, 0=x, 1=y */
369 pme_atomcomm_t atc_energy; /* Only for gmx_pme_calc_energy */
371 rvec *bufv; /* Communication buffer */
372 real *bufr; /* Communication buffer */
373 int buf_nalloc; /* The communication buffer size */
375 /* thread local work data for solve_pme */
378 /* Work data for sum_qgrid */
379 real * sum_qgrid_tmp;
380 real * sum_qgrid_dd_tmp;
383 static void calc_interpolation_idx(gmx_pme_t pme, pme_atomcomm_t *atc,
384 int start, int end, int thread)
387 int *idxptr, tix, tiy, tiz;
388 real *xptr, *fptr, tx, ty, tz;
389 real rxx, ryx, ryy, rzx, rzy, rzz;
391 int start_ix, start_iy, start_iz;
392 int *g2tx, *g2ty, *g2tz;
394 int *thread_idx = NULL;
395 thread_plist_t *tpl = NULL;
403 start_ix = pme->pmegrid_start_ix;
404 start_iy = pme->pmegrid_start_iy;
405 start_iz = pme->pmegrid_start_iz;
407 rxx = pme->recipbox[XX][XX];
408 ryx = pme->recipbox[YY][XX];
409 ryy = pme->recipbox[YY][YY];
410 rzx = pme->recipbox[ZZ][XX];
411 rzy = pme->recipbox[ZZ][YY];
412 rzz = pme->recipbox[ZZ][ZZ];
414 g2tx = pme->pmegrid[PME_GRID_QA].g2t[XX];
415 g2ty = pme->pmegrid[PME_GRID_QA].g2t[YY];
416 g2tz = pme->pmegrid[PME_GRID_QA].g2t[ZZ];
418 bThreads = (atc->nthread > 1);
421 thread_idx = atc->thread_idx;
423 tpl = &atc->thread_plist[thread];
425 for (i = 0; i < atc->nthread; i++)
431 for (i = start; i < end; i++)
434 idxptr = atc->idx[i];
435 fptr = atc->fractx[i];
437 /* Fractional coordinates along box vectors, add 2.0 to make 100% sure we are positive for triclinic boxes */
438 tx = nx * ( xptr[XX] * rxx + xptr[YY] * ryx + xptr[ZZ] * rzx + 2.0 );
439 ty = ny * ( xptr[YY] * ryy + xptr[ZZ] * rzy + 2.0 );
440 tz = nz * ( xptr[ZZ] * rzz + 2.0 );
446 /* Because decomposition only occurs in x and y,
447 * we never have a fraction correction in z.
449 fptr[XX] = tx - tix + pme->fshx[tix];
450 fptr[YY] = ty - tiy + pme->fshy[tiy];
453 idxptr[XX] = pme->nnx[tix];
454 idxptr[YY] = pme->nny[tiy];
455 idxptr[ZZ] = pme->nnz[tiz];
458 range_check(idxptr[XX], 0, pme->pmegrid_nx);
459 range_check(idxptr[YY], 0, pme->pmegrid_ny);
460 range_check(idxptr[ZZ], 0, pme->pmegrid_nz);
465 thread_i = g2tx[idxptr[XX]] + g2ty[idxptr[YY]] + g2tz[idxptr[ZZ]];
466 thread_idx[i] = thread_i;
473 /* Make a list of particle indices sorted on thread */
475 /* Get the cumulative count */
476 for (i = 1; i < atc->nthread; i++)
478 tpl_n[i] += tpl_n[i-1];
480 /* The current implementation distributes particles equally
481 * over the threads, so we could actually allocate for that
482 * in pme_realloc_atomcomm_things.
484 if (tpl_n[atc->nthread-1] > tpl->nalloc)
486 tpl->nalloc = over_alloc_large(tpl_n[atc->nthread-1]);
487 srenew(tpl->i, tpl->nalloc);
489 /* Set tpl_n to the cumulative start */
490 for (i = atc->nthread-1; i >= 1; i--)
492 tpl_n[i] = tpl_n[i-1];
496 /* Fill our thread local array with indices sorted on thread */
497 for (i = start; i < end; i++)
499 tpl->i[tpl_n[atc->thread_idx[i]]++] = i;
501 /* Now tpl_n contains the cummulative count again */
505 static void make_thread_local_ind(pme_atomcomm_t *atc,
506 int thread, splinedata_t *spline)
508 int n, t, i, start, end;
511 /* Combine the indices made by each thread into one index */
515 for (t = 0; t < atc->nthread; t++)
517 tpl = &atc->thread_plist[t];
518 /* Copy our part (start - end) from the list of thread t */
521 start = tpl->n[thread-1];
523 end = tpl->n[thread];
524 for (i = start; i < end; i++)
526 spline->ind[n++] = tpl->i[i];
534 static void pme_calc_pidx(int start, int end,
535 matrix recipbox, rvec x[],
536 pme_atomcomm_t *atc, int *count)
541 real rxx, ryx, rzx, ryy, rzy;
544 /* Calculate PME task index (pidx) for each grid index.
545 * Here we always assign equally sized slabs to each node
546 * for load balancing reasons (the PME grid spacing is not used).
552 /* Reset the count */
553 for (i = 0; i < nslab; i++)
558 if (atc->dimind == 0)
560 rxx = recipbox[XX][XX];
561 ryx = recipbox[YY][XX];
562 rzx = recipbox[ZZ][XX];
563 /* Calculate the node index in x-dimension */
564 for (i = start; i < end; i++)
567 /* Fractional coordinates along box vectors */
568 s = nslab*(xptr[XX]*rxx + xptr[YY]*ryx + xptr[ZZ]*rzx);
569 si = (int)(s + 2*nslab) % nslab;
576 ryy = recipbox[YY][YY];
577 rzy = recipbox[ZZ][YY];
578 /* Calculate the node index in y-dimension */
579 for (i = start; i < end; i++)
582 /* Fractional coordinates along box vectors */
583 s = nslab*(xptr[YY]*ryy + xptr[ZZ]*rzy);
584 si = (int)(s + 2*nslab) % nslab;
591 static void pme_calc_pidx_wrapper(int natoms, matrix recipbox, rvec x[],
594 int nthread, thread, slab;
596 nthread = atc->nthread;
598 #pragma omp parallel for num_threads(nthread) schedule(static)
599 for (thread = 0; thread < nthread; thread++)
601 pme_calc_pidx(natoms* thread /nthread,
602 natoms*(thread+1)/nthread,
603 recipbox, x, atc, atc->count_thread[thread]);
605 /* Non-parallel reduction, since nslab is small */
607 for (thread = 1; thread < nthread; thread++)
609 for (slab = 0; slab < atc->nslab; slab++)
611 atc->count_thread[0][slab] += atc->count_thread[thread][slab];
616 static void realloc_splinevec(splinevec th, real **ptr_z, int nalloc)
618 const int padding = 4;
621 srenew(th[XX], nalloc);
622 srenew(th[YY], nalloc);
623 /* In z we add padding, this is only required for the aligned SIMD code */
624 sfree_aligned(*ptr_z);
625 snew_aligned(*ptr_z, nalloc+2*padding, SIMD4_ALIGNMENT);
626 th[ZZ] = *ptr_z + padding;
628 for (i = 0; i < padding; i++)
631 (*ptr_z)[padding+nalloc+i] = 0;
635 static void pme_realloc_splinedata(splinedata_t *spline, pme_atomcomm_t *atc)
639 srenew(spline->ind, atc->nalloc);
640 /* Initialize the index to identity so it works without threads */
641 for (i = 0; i < atc->nalloc; i++)
646 realloc_splinevec(spline->theta, &spline->ptr_theta_z,
647 atc->pme_order*atc->nalloc);
648 realloc_splinevec(spline->dtheta, &spline->ptr_dtheta_z,
649 atc->pme_order*atc->nalloc);
652 static void pme_realloc_atomcomm_things(pme_atomcomm_t *atc)
654 int nalloc_old, i, j, nalloc_tpl;
656 /* We have to avoid a NULL pointer for atc->x to avoid
657 * possible fatal errors in MPI routines.
659 if (atc->n > atc->nalloc || atc->nalloc == 0)
661 nalloc_old = atc->nalloc;
662 atc->nalloc = over_alloc_dd(max(atc->n, 1));
666 srenew(atc->x, atc->nalloc);
667 srenew(atc->q, atc->nalloc);
668 srenew(atc->f, atc->nalloc);
669 for (i = nalloc_old; i < atc->nalloc; i++)
671 clear_rvec(atc->f[i]);
676 srenew(atc->fractx, atc->nalloc);
677 srenew(atc->idx, atc->nalloc);
679 if (atc->nthread > 1)
681 srenew(atc->thread_idx, atc->nalloc);
684 for (i = 0; i < atc->nthread; i++)
686 pme_realloc_splinedata(&atc->spline[i], atc);
692 static void pme_dd_sendrecv(pme_atomcomm_t gmx_unused *atc,
693 gmx_bool gmx_unused bBackward, int gmx_unused shift,
694 void gmx_unused *buf_s, int gmx_unused nbyte_s,
695 void gmx_unused *buf_r, int gmx_unused nbyte_r)
701 if (bBackward == FALSE)
703 dest = atc->node_dest[shift];
704 src = atc->node_src[shift];
708 dest = atc->node_src[shift];
709 src = atc->node_dest[shift];
712 if (nbyte_s > 0 && nbyte_r > 0)
714 MPI_Sendrecv(buf_s, nbyte_s, MPI_BYTE,
716 buf_r, nbyte_r, MPI_BYTE,
718 atc->mpi_comm, &stat);
720 else if (nbyte_s > 0)
722 MPI_Send(buf_s, nbyte_s, MPI_BYTE,
726 else if (nbyte_r > 0)
728 MPI_Recv(buf_r, nbyte_r, MPI_BYTE,
730 atc->mpi_comm, &stat);
735 static void dd_pmeredist_x_q(gmx_pme_t pme,
736 int n, gmx_bool bX, rvec *x, real *charge,
739 int *commnode, *buf_index;
740 int nnodes_comm, i, nsend, local_pos, buf_pos, node, scount, rcount;
742 commnode = atc->node_dest;
743 buf_index = atc->buf_index;
745 nnodes_comm = min(2*atc->maxshift, atc->nslab-1);
748 for (i = 0; i < nnodes_comm; i++)
750 buf_index[commnode[i]] = nsend;
751 nsend += atc->count[commnode[i]];
755 if (atc->count[atc->nodeid] + nsend != n)
757 gmx_fatal(FARGS, "%d particles communicated to PME node %d are more than 2/3 times the cut-off out of the domain decomposition cell of their charge group in dimension %c.\n"
758 "This usually means that your system is not well equilibrated.",
759 n - (atc->count[atc->nodeid] + nsend),
760 pme->nodeid, 'x'+atc->dimind);
763 if (nsend > pme->buf_nalloc)
765 pme->buf_nalloc = over_alloc_dd(nsend);
766 srenew(pme->bufv, pme->buf_nalloc);
767 srenew(pme->bufr, pme->buf_nalloc);
770 atc->n = atc->count[atc->nodeid];
771 for (i = 0; i < nnodes_comm; i++)
773 scount = atc->count[commnode[i]];
774 /* Communicate the count */
777 fprintf(debug, "dimind %d PME node %d send to node %d: %d\n",
778 atc->dimind, atc->nodeid, commnode[i], scount);
780 pme_dd_sendrecv(atc, FALSE, i,
781 &scount, sizeof(int),
782 &atc->rcount[i], sizeof(int));
783 atc->n += atc->rcount[i];
786 pme_realloc_atomcomm_things(atc);
790 for (i = 0; i < n; i++)
793 if (node == atc->nodeid)
795 /* Copy direct to the receive buffer */
798 copy_rvec(x[i], atc->x[local_pos]);
800 atc->q[local_pos] = charge[i];
805 /* Copy to the send buffer */
808 copy_rvec(x[i], pme->bufv[buf_index[node]]);
810 pme->bufr[buf_index[node]] = charge[i];
816 for (i = 0; i < nnodes_comm; i++)
818 scount = atc->count[commnode[i]];
819 rcount = atc->rcount[i];
820 if (scount > 0 || rcount > 0)
824 /* Communicate the coordinates */
825 pme_dd_sendrecv(atc, FALSE, i,
826 pme->bufv[buf_pos], scount*sizeof(rvec),
827 atc->x[local_pos], rcount*sizeof(rvec));
829 /* Communicate the charges */
830 pme_dd_sendrecv(atc, FALSE, i,
831 pme->bufr+buf_pos, scount*sizeof(real),
832 atc->q+local_pos, rcount*sizeof(real));
834 local_pos += atc->rcount[i];
839 static void dd_pmeredist_f(gmx_pme_t pme, pme_atomcomm_t *atc,
843 int *commnode, *buf_index;
844 int nnodes_comm, local_pos, buf_pos, i, scount, rcount, node;
846 commnode = atc->node_dest;
847 buf_index = atc->buf_index;
849 nnodes_comm = min(2*atc->maxshift, atc->nslab-1);
851 local_pos = atc->count[atc->nodeid];
853 for (i = 0; i < nnodes_comm; i++)
855 scount = atc->rcount[i];
856 rcount = atc->count[commnode[i]];
857 if (scount > 0 || rcount > 0)
859 /* Communicate the forces */
860 pme_dd_sendrecv(atc, TRUE, i,
861 atc->f[local_pos], scount*sizeof(rvec),
862 pme->bufv[buf_pos], rcount*sizeof(rvec));
865 buf_index[commnode[i]] = buf_pos;
872 for (i = 0; i < n; i++)
875 if (node == atc->nodeid)
877 /* Add from the local force array */
878 rvec_inc(f[i], atc->f[local_pos]);
883 /* Add from the receive buffer */
884 rvec_inc(f[i], pme->bufv[buf_index[node]]);
891 for (i = 0; i < n; i++)
894 if (node == atc->nodeid)
896 /* Copy from the local force array */
897 copy_rvec(atc->f[local_pos], f[i]);
902 /* Copy from the receive buffer */
903 copy_rvec(pme->bufv[buf_index[node]], f[i]);
911 static void gmx_sum_qgrid_dd(gmx_pme_t pme, real *grid, int direction)
913 pme_overlap_t *overlap;
914 int send_index0, send_nindex;
915 int recv_index0, recv_nindex;
917 int i, j, k, ix, iy, iz, icnt;
918 int ipulse, send_id, recv_id, datasize;
920 real *sendptr, *recvptr;
922 /* Start with minor-rank communication. This is a bit of a pain since it is not contiguous */
923 overlap = &pme->overlap[1];
925 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
927 /* Since we have already (un)wrapped the overlap in the z-dimension,
928 * we only have to communicate 0 to nkz (not pmegrid_nz).
930 if (direction == GMX_SUM_QGRID_FORWARD)
932 send_id = overlap->send_id[ipulse];
933 recv_id = overlap->recv_id[ipulse];
934 send_index0 = overlap->comm_data[ipulse].send_index0;
935 send_nindex = overlap->comm_data[ipulse].send_nindex;
936 recv_index0 = overlap->comm_data[ipulse].recv_index0;
937 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
941 send_id = overlap->recv_id[ipulse];
942 recv_id = overlap->send_id[ipulse];
943 send_index0 = overlap->comm_data[ipulse].recv_index0;
944 send_nindex = overlap->comm_data[ipulse].recv_nindex;
945 recv_index0 = overlap->comm_data[ipulse].send_index0;
946 recv_nindex = overlap->comm_data[ipulse].send_nindex;
949 /* Copy data to contiguous send buffer */
952 fprintf(debug, "PME send node %d %d -> %d grid start %d Communicating %d to %d\n",
953 pme->nodeid, overlap->nodeid, send_id,
954 pme->pmegrid_start_iy,
955 send_index0-pme->pmegrid_start_iy,
956 send_index0-pme->pmegrid_start_iy+send_nindex);
959 for (i = 0; i < pme->pmegrid_nx; i++)
962 for (j = 0; j < send_nindex; j++)
964 iy = j + send_index0 - pme->pmegrid_start_iy;
965 for (k = 0; k < pme->nkz; k++)
968 overlap->sendbuf[icnt++] = grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz];
973 datasize = pme->pmegrid_nx * pme->nkz;
975 MPI_Sendrecv(overlap->sendbuf, send_nindex*datasize, GMX_MPI_REAL,
977 overlap->recvbuf, recv_nindex*datasize, GMX_MPI_REAL,
979 overlap->mpi_comm, &stat);
981 /* Get data from contiguous recv buffer */
984 fprintf(debug, "PME recv node %d %d <- %d grid start %d Communicating %d to %d\n",
985 pme->nodeid, overlap->nodeid, recv_id,
986 pme->pmegrid_start_iy,
987 recv_index0-pme->pmegrid_start_iy,
988 recv_index0-pme->pmegrid_start_iy+recv_nindex);
991 for (i = 0; i < pme->pmegrid_nx; i++)
994 for (j = 0; j < recv_nindex; j++)
996 iy = j + recv_index0 - pme->pmegrid_start_iy;
997 for (k = 0; k < pme->nkz; k++)
1000 if (direction == GMX_SUM_QGRID_FORWARD)
1002 grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz] += overlap->recvbuf[icnt++];
1006 grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz] = overlap->recvbuf[icnt++];
1013 /* Major dimension is easier, no copying required,
1014 * but we might have to sum to separate array.
1015 * Since we don't copy, we have to communicate up to pmegrid_nz,
1016 * not nkz as for the minor direction.
1018 overlap = &pme->overlap[0];
1020 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
1022 if (direction == GMX_SUM_QGRID_FORWARD)
1024 send_id = overlap->send_id[ipulse];
1025 recv_id = overlap->recv_id[ipulse];
1026 send_index0 = overlap->comm_data[ipulse].send_index0;
1027 send_nindex = overlap->comm_data[ipulse].send_nindex;
1028 recv_index0 = overlap->comm_data[ipulse].recv_index0;
1029 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
1030 recvptr = overlap->recvbuf;
1034 send_id = overlap->recv_id[ipulse];
1035 recv_id = overlap->send_id[ipulse];
1036 send_index0 = overlap->comm_data[ipulse].recv_index0;
1037 send_nindex = overlap->comm_data[ipulse].recv_nindex;
1038 recv_index0 = overlap->comm_data[ipulse].send_index0;
1039 recv_nindex = overlap->comm_data[ipulse].send_nindex;
1040 recvptr = grid + (recv_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1043 sendptr = grid + (send_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1044 datasize = pme->pmegrid_ny * pme->pmegrid_nz;
1048 fprintf(debug, "PME send node %d %d -> %d grid start %d Communicating %d to %d\n",
1049 pme->nodeid, overlap->nodeid, send_id,
1050 pme->pmegrid_start_ix,
1051 send_index0-pme->pmegrid_start_ix,
1052 send_index0-pme->pmegrid_start_ix+send_nindex);
1053 fprintf(debug, "PME recv node %d %d <- %d grid start %d Communicating %d to %d\n",
1054 pme->nodeid, overlap->nodeid, recv_id,
1055 pme->pmegrid_start_ix,
1056 recv_index0-pme->pmegrid_start_ix,
1057 recv_index0-pme->pmegrid_start_ix+recv_nindex);
1060 MPI_Sendrecv(sendptr, send_nindex*datasize, GMX_MPI_REAL,
1062 recvptr, recv_nindex*datasize, GMX_MPI_REAL,
1064 overlap->mpi_comm, &stat);
1066 /* ADD data from contiguous recv buffer */
1067 if (direction == GMX_SUM_QGRID_FORWARD)
1069 p = grid + (recv_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1070 for (i = 0; i < recv_nindex*datasize; i++)
1072 p[i] += overlap->recvbuf[i];
1080 static int copy_pmegrid_to_fftgrid(gmx_pme_t pme, real *pmegrid, real *fftgrid, int grid_index)
1082 ivec local_fft_ndata, local_fft_offset, local_fft_size;
1083 ivec local_pme_size;
1087 /* Dimensions should be identical for A/B grid, so we just use A here */
1088 gmx_parallel_3dfft_real_limits(pme->pfft_setup[grid_index],
1093 local_pme_size[0] = pme->pmegrid_nx;
1094 local_pme_size[1] = pme->pmegrid_ny;
1095 local_pme_size[2] = pme->pmegrid_nz;
1097 /* The fftgrid is always 'justified' to the lower-left corner of the PME grid,
1098 the offset is identical, and the PME grid always has more data (due to overlap)
1103 char fn[STRLEN], format[STRLEN];
1105 sprintf(fn, "pmegrid%d.pdb", pme->nodeid);
1106 fp = gmx_ffopen(fn, "w");
1107 sprintf(fn, "pmegrid%d.txt", pme->nodeid);
1108 fp2 = gmx_ffopen(fn, "w");
1109 sprintf(format, "%s%s\n", pdbformat, "%6.2f%6.2f");
1112 for (ix = 0; ix < local_fft_ndata[XX]; ix++)
1114 for (iy = 0; iy < local_fft_ndata[YY]; iy++)
1116 for (iz = 0; iz < local_fft_ndata[ZZ]; iz++)
1118 pmeidx = ix*(local_pme_size[YY]*local_pme_size[ZZ])+iy*(local_pme_size[ZZ])+iz;
1119 fftidx = ix*(local_fft_size[YY]*local_fft_size[ZZ])+iy*(local_fft_size[ZZ])+iz;
1120 fftgrid[fftidx] = pmegrid[pmeidx];
1122 val = 100*pmegrid[pmeidx];
1123 if (pmegrid[pmeidx] != 0)
1125 fprintf(fp, format, "ATOM", pmeidx, "CA", "GLY", ' ', pmeidx, ' ',
1126 5.0*ix, 5.0*iy, 5.0*iz, 1.0, val);
1128 if (pmegrid[pmeidx] != 0)
1130 fprintf(fp2, "%-12s %5d %5d %5d %12.5e\n",
1132 pme->pmegrid_start_ix + ix,
1133 pme->pmegrid_start_iy + iy,
1134 pme->pmegrid_start_iz + iz,
1150 static gmx_cycles_t omp_cyc_start()
1152 return gmx_cycles_read();
1155 static gmx_cycles_t omp_cyc_end(gmx_cycles_t c)
1157 return gmx_cycles_read() - c;
1161 static int copy_fftgrid_to_pmegrid(gmx_pme_t pme, const real *fftgrid, real *pmegrid, int grid_index,
1162 int nthread, int thread)
1164 ivec local_fft_ndata, local_fft_offset, local_fft_size;
1165 ivec local_pme_size;
1166 int ixy0, ixy1, ixy, ix, iy, iz;
1168 #ifdef PME_TIME_THREADS
1170 static double cs1 = 0;
1174 #ifdef PME_TIME_THREADS
1175 c1 = omp_cyc_start();
1177 /* Dimensions should be identical for A/B grid, so we just use A here */
1178 gmx_parallel_3dfft_real_limits(pme->pfft_setup[grid_index],
1183 local_pme_size[0] = pme->pmegrid_nx;
1184 local_pme_size[1] = pme->pmegrid_ny;
1185 local_pme_size[2] = pme->pmegrid_nz;
1187 /* The fftgrid is always 'justified' to the lower-left corner of the PME grid,
1188 the offset is identical, and the PME grid always has more data (due to overlap)
1190 ixy0 = ((thread )*local_fft_ndata[XX]*local_fft_ndata[YY])/nthread;
1191 ixy1 = ((thread+1)*local_fft_ndata[XX]*local_fft_ndata[YY])/nthread;
1193 for (ixy = ixy0; ixy < ixy1; ixy++)
1195 ix = ixy/local_fft_ndata[YY];
1196 iy = ixy - ix*local_fft_ndata[YY];
1198 pmeidx = (ix*local_pme_size[YY] + iy)*local_pme_size[ZZ];
1199 fftidx = (ix*local_fft_size[YY] + iy)*local_fft_size[ZZ];
1200 for (iz = 0; iz < local_fft_ndata[ZZ]; iz++)
1202 pmegrid[pmeidx+iz] = fftgrid[fftidx+iz];
1206 #ifdef PME_TIME_THREADS
1207 c1 = omp_cyc_end(c1);
1212 printf("copy %.2f\n", cs1*1e-9);
1220 static void wrap_periodic_pmegrid(gmx_pme_t pme, real *pmegrid)
1222 int nx, ny, nz, pnx, pny, pnz, ny_x, overlap, ix, iy, iz;
1228 pnx = pme->pmegrid_nx;
1229 pny = pme->pmegrid_ny;
1230 pnz = pme->pmegrid_nz;
1232 overlap = pme->pme_order - 1;
1234 /* Add periodic overlap in z */
1235 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1237 for (iy = 0; iy < pme->pmegrid_ny; iy++)
1239 for (iz = 0; iz < overlap; iz++)
1241 pmegrid[(ix*pny+iy)*pnz+iz] +=
1242 pmegrid[(ix*pny+iy)*pnz+nz+iz];
1247 if (pme->nnodes_minor == 1)
1249 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1251 for (iy = 0; iy < overlap; iy++)
1253 for (iz = 0; iz < nz; iz++)
1255 pmegrid[(ix*pny+iy)*pnz+iz] +=
1256 pmegrid[(ix*pny+ny+iy)*pnz+iz];
1262 if (pme->nnodes_major == 1)
1264 ny_x = (pme->nnodes_minor == 1 ? ny : pme->pmegrid_ny);
1266 for (ix = 0; ix < overlap; ix++)
1268 for (iy = 0; iy < ny_x; iy++)
1270 for (iz = 0; iz < nz; iz++)
1272 pmegrid[(ix*pny+iy)*pnz+iz] +=
1273 pmegrid[((nx+ix)*pny+iy)*pnz+iz];
1281 static void unwrap_periodic_pmegrid(gmx_pme_t pme, real *pmegrid)
1283 int nx, ny, nz, pnx, pny, pnz, ny_x, overlap, ix;
1289 pnx = pme->pmegrid_nx;
1290 pny = pme->pmegrid_ny;
1291 pnz = pme->pmegrid_nz;
1293 overlap = pme->pme_order - 1;
1295 if (pme->nnodes_major == 1)
1297 ny_x = (pme->nnodes_minor == 1 ? ny : pme->pmegrid_ny);
1299 for (ix = 0; ix < overlap; ix++)
1303 for (iy = 0; iy < ny_x; iy++)
1305 for (iz = 0; iz < nz; iz++)
1307 pmegrid[((nx+ix)*pny+iy)*pnz+iz] =
1308 pmegrid[(ix*pny+iy)*pnz+iz];
1314 if (pme->nnodes_minor == 1)
1316 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1317 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1321 for (iy = 0; iy < overlap; iy++)
1323 for (iz = 0; iz < nz; iz++)
1325 pmegrid[(ix*pny+ny+iy)*pnz+iz] =
1326 pmegrid[(ix*pny+iy)*pnz+iz];
1332 /* Copy periodic overlap in z */
1333 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1334 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1338 for (iy = 0; iy < pme->pmegrid_ny; iy++)
1340 for (iz = 0; iz < overlap; iz++)
1342 pmegrid[(ix*pny+iy)*pnz+nz+iz] =
1343 pmegrid[(ix*pny+iy)*pnz+iz];
1350 /* This has to be a macro to enable full compiler optimization with xlC (and probably others too) */
1351 #define DO_BSPLINE(order) \
1352 for (ithx = 0; (ithx < order); ithx++) \
1354 index_x = (i0+ithx)*pny*pnz; \
1355 valx = qn*thx[ithx]; \
1357 for (ithy = 0; (ithy < order); ithy++) \
1359 valxy = valx*thy[ithy]; \
1360 index_xy = index_x+(j0+ithy)*pnz; \
1362 for (ithz = 0; (ithz < order); ithz++) \
1364 index_xyz = index_xy+(k0+ithz); \
1365 grid[index_xyz] += valxy*thz[ithz]; \
1371 static void spread_q_bsplines_thread(pmegrid_t *pmegrid,
1372 pme_atomcomm_t *atc,
1373 splinedata_t *spline,
1374 pme_spline_work_t gmx_unused *work)
1377 /* spread charges from home atoms to local grid */
1380 int b, i, nn, n, ithx, ithy, ithz, i0, j0, k0;
1382 int order, norder, index_x, index_xy, index_xyz;
1383 real valx, valxy, qn;
1384 real *thx, *thy, *thz;
1385 int localsize, bndsize;
1386 int pnx, pny, pnz, ndatatot;
1387 int offx, offy, offz;
1389 #if defined PME_SIMD4_SPREAD_GATHER && !defined PME_SIMD4_UNALIGNED
1390 real thz_buffer[GMX_SIMD4_WIDTH*3], *thz_aligned;
1392 thz_aligned = gmx_simd4_align_r(thz_buffer);
1395 pnx = pmegrid->s[XX];
1396 pny = pmegrid->s[YY];
1397 pnz = pmegrid->s[ZZ];
1399 offx = pmegrid->offset[XX];
1400 offy = pmegrid->offset[YY];
1401 offz = pmegrid->offset[ZZ];
1403 ndatatot = pnx*pny*pnz;
1404 grid = pmegrid->grid;
1405 for (i = 0; i < ndatatot; i++)
1410 order = pmegrid->order;
1412 for (nn = 0; nn < spline->n; nn++)
1414 n = spline->ind[nn];
1419 idxptr = atc->idx[n];
1422 i0 = idxptr[XX] - offx;
1423 j0 = idxptr[YY] - offy;
1424 k0 = idxptr[ZZ] - offz;
1426 thx = spline->theta[XX] + norder;
1427 thy = spline->theta[YY] + norder;
1428 thz = spline->theta[ZZ] + norder;
1433 #ifdef PME_SIMD4_SPREAD_GATHER
1434 #ifdef PME_SIMD4_UNALIGNED
1435 #define PME_SPREAD_SIMD4_ORDER4
1437 #define PME_SPREAD_SIMD4_ALIGNED
1440 #include "pme_simd4.h"
1446 #ifdef PME_SIMD4_SPREAD_GATHER
1447 #define PME_SPREAD_SIMD4_ALIGNED
1449 #include "pme_simd4.h"
1462 static void set_grid_alignment(int gmx_unused *pmegrid_nz, int gmx_unused pme_order)
1464 #ifdef PME_SIMD4_SPREAD_GATHER
1466 #ifndef PME_SIMD4_UNALIGNED
1471 /* Round nz up to a multiple of 4 to ensure alignment */
1472 *pmegrid_nz = ((*pmegrid_nz + 3) & ~3);
1477 static void set_gridsize_alignment(int gmx_unused *gridsize, int gmx_unused pme_order)
1479 #ifdef PME_SIMD4_SPREAD_GATHER
1480 #ifndef PME_SIMD4_UNALIGNED
1483 /* Add extra elements to ensured aligned operations do not go
1484 * beyond the allocated grid size.
1485 * Note that for pme_order=5, the pme grid z-size alignment
1486 * ensures that we will not go beyond the grid size.
1494 static void pmegrid_init(pmegrid_t *grid,
1495 int cx, int cy, int cz,
1496 int x0, int y0, int z0,
1497 int x1, int y1, int z1,
1498 gmx_bool set_alignment,
1507 grid->offset[XX] = x0;
1508 grid->offset[YY] = y0;
1509 grid->offset[ZZ] = z0;
1510 grid->n[XX] = x1 - x0 + pme_order - 1;
1511 grid->n[YY] = y1 - y0 + pme_order - 1;
1512 grid->n[ZZ] = z1 - z0 + pme_order - 1;
1513 copy_ivec(grid->n, grid->s);
1516 set_grid_alignment(&nz, pme_order);
1521 else if (nz != grid->s[ZZ])
1523 gmx_incons("pmegrid_init call with an unaligned z size");
1526 grid->order = pme_order;
1529 gridsize = grid->s[XX]*grid->s[YY]*grid->s[ZZ];
1530 set_gridsize_alignment(&gridsize, pme_order);
1531 snew_aligned(grid->grid, gridsize, SIMD4_ALIGNMENT);
1539 static int div_round_up(int enumerator, int denominator)
1541 return (enumerator + denominator - 1)/denominator;
1544 static void make_subgrid_division(const ivec n, int ovl, int nthread,
1547 int gsize_opt, gsize;
1552 for (nsx = 1; nsx <= nthread; nsx++)
1554 if (nthread % nsx == 0)
1556 for (nsy = 1; nsy <= nthread; nsy++)
1558 if (nsx*nsy <= nthread && nthread % (nsx*nsy) == 0)
1560 nsz = nthread/(nsx*nsy);
1562 /* Determine the number of grid points per thread */
1564 (div_round_up(n[XX], nsx) + ovl)*
1565 (div_round_up(n[YY], nsy) + ovl)*
1566 (div_round_up(n[ZZ], nsz) + ovl);
1568 /* Minimize the number of grids points per thread
1569 * and, secondarily, the number of cuts in minor dimensions.
1571 if (gsize_opt == -1 ||
1572 gsize < gsize_opt ||
1573 (gsize == gsize_opt &&
1574 (nsz < nsub[ZZ] || (nsz == nsub[ZZ] && nsy < nsub[YY]))))
1586 env = getenv("GMX_PME_THREAD_DIVISION");
1589 sscanf(env, "%d %d %d", &nsub[XX], &nsub[YY], &nsub[ZZ]);
1592 if (nsub[XX]*nsub[YY]*nsub[ZZ] != nthread)
1594 gmx_fatal(FARGS, "PME grid thread division (%d x %d x %d) does not match the total number of threads (%d)", nsub[XX], nsub[YY], nsub[ZZ], nthread);
1598 static void pmegrids_init(pmegrids_t *grids,
1599 int nx, int ny, int nz, int nz_base,
1601 gmx_bool bUseThreads,
1606 ivec n, n_base, g0, g1;
1607 int t, x, y, z, d, i, tfac;
1608 int max_comm_lines = -1;
1610 n[XX] = nx - (pme_order - 1);
1611 n[YY] = ny - (pme_order - 1);
1612 n[ZZ] = nz - (pme_order - 1);
1614 copy_ivec(n, n_base);
1615 n_base[ZZ] = nz_base;
1617 pmegrid_init(&grids->grid, 0, 0, 0, 0, 0, 0, n[XX], n[YY], n[ZZ], FALSE, pme_order,
1620 grids->nthread = nthread;
1622 make_subgrid_division(n_base, pme_order-1, grids->nthread, grids->nc);
1629 for (d = 0; d < DIM; d++)
1631 nst[d] = div_round_up(n[d], grids->nc[d]) + pme_order - 1;
1633 set_grid_alignment(&nst[ZZ], pme_order);
1637 fprintf(debug, "pmegrid thread local division: %d x %d x %d\n",
1638 grids->nc[XX], grids->nc[YY], grids->nc[ZZ]);
1639 fprintf(debug, "pmegrid %d %d %d max thread pmegrid %d %d %d\n",
1641 nst[XX], nst[YY], nst[ZZ]);
1644 snew(grids->grid_th, grids->nthread);
1646 gridsize = nst[XX]*nst[YY]*nst[ZZ];
1647 set_gridsize_alignment(&gridsize, pme_order);
1648 snew_aligned(grids->grid_all,
1649 grids->nthread*gridsize+(grids->nthread+1)*GMX_CACHE_SEP,
1652 for (x = 0; x < grids->nc[XX]; x++)
1654 for (y = 0; y < grids->nc[YY]; y++)
1656 for (z = 0; z < grids->nc[ZZ]; z++)
1658 pmegrid_init(&grids->grid_th[t],
1660 (n[XX]*(x ))/grids->nc[XX],
1661 (n[YY]*(y ))/grids->nc[YY],
1662 (n[ZZ]*(z ))/grids->nc[ZZ],
1663 (n[XX]*(x+1))/grids->nc[XX],
1664 (n[YY]*(y+1))/grids->nc[YY],
1665 (n[ZZ]*(z+1))/grids->nc[ZZ],
1668 grids->grid_all+GMX_CACHE_SEP+t*(gridsize+GMX_CACHE_SEP));
1676 grids->grid_th = NULL;
1679 snew(grids->g2t, DIM);
1681 for (d = DIM-1; d >= 0; d--)
1683 snew(grids->g2t[d], n[d]);
1685 for (i = 0; i < n[d]; i++)
1687 /* The second check should match the parameters
1688 * of the pmegrid_init call above.
1690 while (t + 1 < grids->nc[d] && i >= (n[d]*(t+1))/grids->nc[d])
1694 grids->g2t[d][i] = t*tfac;
1697 tfac *= grids->nc[d];
1701 case XX: max_comm_lines = overlap_x; break;
1702 case YY: max_comm_lines = overlap_y; break;
1703 case ZZ: max_comm_lines = pme_order - 1; break;
1705 grids->nthread_comm[d] = 0;
1706 while ((n[d]*grids->nthread_comm[d])/grids->nc[d] < max_comm_lines &&
1707 grids->nthread_comm[d] < grids->nc[d])
1709 grids->nthread_comm[d]++;
1713 fprintf(debug, "pmegrid thread grid communication range in %c: %d\n",
1714 'x'+d, grids->nthread_comm[d]);
1716 /* It should be possible to make grids->nthread_comm[d]==grids->nc[d]
1717 * work, but this is not a problematic restriction.
1719 if (grids->nc[d] > 1 && grids->nthread_comm[d] > grids->nc[d])
1721 gmx_fatal(FARGS, "Too many threads for PME (%d) compared to the number of grid lines, reduce the number of threads doing PME", grids->nthread);
1727 static void pmegrids_destroy(pmegrids_t *grids)
1731 if (grids->grid.grid != NULL)
1733 sfree(grids->grid.grid);
1735 if (grids->nthread > 0)
1737 for (t = 0; t < grids->nthread; t++)
1739 sfree(grids->grid_th[t].grid);
1741 sfree(grids->grid_th);
1747 static void realloc_work(pme_work_t *work, int nkx)
1751 if (nkx > work->nalloc)
1754 srenew(work->mhx, work->nalloc);
1755 srenew(work->mhy, work->nalloc);
1756 srenew(work->mhz, work->nalloc);
1757 srenew(work->m2, work->nalloc);
1758 /* Allocate an aligned pointer for SIMD operations, including extra
1759 * elements at the end for padding.
1761 #ifdef PME_SIMD_SOLVE
1762 simd_width = GMX_SIMD_REAL_WIDTH;
1764 /* We can use any alignment, apart from 0, so we use 4 */
1767 sfree_aligned(work->denom);
1768 sfree_aligned(work->tmp1);
1769 sfree_aligned(work->tmp2);
1770 sfree_aligned(work->eterm);
1771 snew_aligned(work->denom, work->nalloc+simd_width, simd_width*sizeof(real));
1772 snew_aligned(work->tmp1, work->nalloc+simd_width, simd_width*sizeof(real));
1773 snew_aligned(work->tmp2, work->nalloc+simd_width, simd_width*sizeof(real));
1774 snew_aligned(work->eterm, work->nalloc+simd_width, simd_width*sizeof(real));
1775 srenew(work->m2inv, work->nalloc);
1780 static void free_work(pme_work_t *work)
1786 sfree_aligned(work->denom);
1787 sfree_aligned(work->tmp1);
1788 sfree_aligned(work->tmp2);
1789 sfree_aligned(work->eterm);
1794 #if defined PME_SIMD_SOLVE
1795 /* Calculate exponentials through SIMD */
1796 gmx_inline static void calc_exponentials_q(int gmx_unused start, int end, real f, real *d_aligned, real *r_aligned, real *e_aligned)
1799 const gmx_simd_real_t two = gmx_simd_set1_r(2.0);
1800 gmx_simd_real_t f_simd;
1802 gmx_simd_real_t tmp_d1, d_inv, tmp_r, tmp_e;
1804 f_simd = gmx_simd_set1_r(f);
1805 /* We only need to calculate from start. But since start is 0 or 1
1806 * and we want to use aligned loads/stores, we always start from 0.
1808 for (kx = 0; kx < end; kx += GMX_SIMD_REAL_WIDTH)
1810 tmp_d1 = gmx_simd_load_r(d_aligned+kx);
1811 d_inv = gmx_simd_inv_r(tmp_d1);
1812 tmp_r = gmx_simd_load_r(r_aligned+kx);
1813 tmp_r = gmx_simd_exp_r(tmp_r);
1814 tmp_e = gmx_simd_mul_r(f_simd, d_inv);
1815 tmp_e = gmx_simd_mul_r(tmp_e, tmp_r);
1816 gmx_simd_store_r(e_aligned+kx, tmp_e);
1821 gmx_inline static void calc_exponentials_q(int start, int end, real f, real *d, real *r, real *e)
1824 for (kx = start; kx < end; kx++)
1828 for (kx = start; kx < end; kx++)
1832 for (kx = start; kx < end; kx++)
1834 e[kx] = f*r[kx]*d[kx];
1839 #if defined PME_SIMD_SOLVE
1840 /* Calculate exponentials through SIMD */
1841 gmx_inline static void calc_exponentials_lj(int gmx_unused start, int end, real *r_aligned, real *factor_aligned, real *d_aligned)
1843 gmx_simd_real_t tmp_r, tmp_d, tmp_fac, d_inv, tmp_mk;
1844 const gmx_simd_real_t sqr_PI = gmx_simd_sqrt_r(gmx_simd_set1_r(M_PI));
1846 for (kx = 0; kx < end; kx += GMX_SIMD_REAL_WIDTH)
1848 /* We only need to calculate from start. But since start is 0 or 1
1849 * and we want to use aligned loads/stores, we always start from 0.
1851 tmp_d = gmx_simd_load_r(d_aligned+kx);
1852 d_inv = gmx_simd_inv_r(tmp_d);
1853 gmx_simd_store_r(d_aligned+kx, d_inv);
1854 tmp_r = gmx_simd_load_r(r_aligned+kx);
1855 tmp_r = gmx_simd_exp_r(tmp_r);
1856 gmx_simd_store_r(r_aligned+kx, tmp_r);
1857 tmp_mk = gmx_simd_load_r(factor_aligned+kx);
1858 tmp_fac = gmx_simd_mul_r(sqr_PI, gmx_simd_mul_r(tmp_mk, gmx_simd_erfc_r(tmp_mk)));
1859 gmx_simd_store_r(factor_aligned+kx, tmp_fac);
1863 gmx_inline static void calc_exponentials_lj(int start, int end, real *r, real *tmp2, real *d)
1867 for (kx = start; kx < end; kx++)
1872 for (kx = start; kx < end; kx++)
1877 for (kx = start; kx < end; kx++)
1880 tmp2[kx] = sqrt(M_PI)*mk*gmx_erfc(mk);
1885 static int solve_pme_yzx(gmx_pme_t pme, t_complex *grid,
1886 real ewaldcoeff, real vol,
1888 int nthread, int thread)
1890 /* do recip sum over local cells in grid */
1891 /* y major, z middle, x minor or continuous */
1893 int kx, ky, kz, maxkx, maxky, maxkz;
1894 int nx, ny, nz, iyz0, iyz1, iyz, iy, iz, kxstart, kxend;
1896 real factor = M_PI*M_PI/(ewaldcoeff*ewaldcoeff);
1897 real ets2, struct2, vfactor, ets2vf;
1898 real d1, d2, energy = 0;
1900 real virxx = 0, virxy = 0, virxz = 0, viryy = 0, viryz = 0, virzz = 0;
1901 real rxx, ryx, ryy, rzx, rzy, rzz;
1903 real *mhx, *mhy, *mhz, *m2, *denom, *tmp1, *eterm, *m2inv;
1904 real mhxk, mhyk, mhzk, m2k;
1907 ivec local_ndata, local_offset, local_size;
1910 elfac = ONE_4PI_EPS0/pme->epsilon_r;
1916 /* Dimensions should be identical for A/B grid, so we just use A here */
1917 gmx_parallel_3dfft_complex_limits(pme->pfft_setup[PME_GRID_QA],
1923 rxx = pme->recipbox[XX][XX];
1924 ryx = pme->recipbox[YY][XX];
1925 ryy = pme->recipbox[YY][YY];
1926 rzx = pme->recipbox[ZZ][XX];
1927 rzy = pme->recipbox[ZZ][YY];
1928 rzz = pme->recipbox[ZZ][ZZ];
1934 work = &pme->work[thread];
1939 denom = work->denom;
1941 eterm = work->eterm;
1942 m2inv = work->m2inv;
1944 iyz0 = local_ndata[YY]*local_ndata[ZZ]* thread /nthread;
1945 iyz1 = local_ndata[YY]*local_ndata[ZZ]*(thread+1)/nthread;
1947 for (iyz = iyz0; iyz < iyz1; iyz++)
1949 iy = iyz/local_ndata[ZZ];
1950 iz = iyz - iy*local_ndata[ZZ];
1952 ky = iy + local_offset[YY];
1963 by = M_PI*vol*pme->bsp_mod[YY][ky];
1965 kz = iz + local_offset[ZZ];
1969 bz = pme->bsp_mod[ZZ][kz];
1971 /* 0.5 correction for corner points */
1973 if (kz == 0 || kz == (nz+1)/2)
1978 p0 = grid + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
1980 /* We should skip the k-space point (0,0,0) */
1981 /* Note that since here x is the minor index, local_offset[XX]=0 */
1982 if (local_offset[XX] > 0 || ky > 0 || kz > 0)
1984 kxstart = local_offset[XX];
1988 kxstart = local_offset[XX] + 1;
1991 kxend = local_offset[XX] + local_ndata[XX];
1995 /* More expensive inner loop, especially because of the storage
1996 * of the mh elements in array's.
1997 * Because x is the minor grid index, all mh elements
1998 * depend on kx for triclinic unit cells.
2001 /* Two explicit loops to avoid a conditional inside the loop */
2002 for (kx = kxstart; kx < maxkx; kx++)
2007 mhyk = mx * ryx + my * ryy;
2008 mhzk = mx * rzx + my * rzy + mz * rzz;
2009 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2014 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2015 tmp1[kx] = -factor*m2k;
2018 for (kx = maxkx; kx < kxend; kx++)
2023 mhyk = mx * ryx + my * ryy;
2024 mhzk = mx * rzx + my * rzy + mz * rzz;
2025 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2030 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2031 tmp1[kx] = -factor*m2k;
2034 for (kx = kxstart; kx < kxend; kx++)
2036 m2inv[kx] = 1.0/m2[kx];
2039 calc_exponentials_q(kxstart, kxend, elfac, denom, tmp1, eterm);
2041 for (kx = kxstart; kx < kxend; kx++, p0++)
2046 p0->re = d1*eterm[kx];
2047 p0->im = d2*eterm[kx];
2049 struct2 = 2.0*(d1*d1+d2*d2);
2051 tmp1[kx] = eterm[kx]*struct2;
2054 for (kx = kxstart; kx < kxend; kx++)
2056 ets2 = corner_fac*tmp1[kx];
2057 vfactor = (factor*m2[kx] + 1.0)*2.0*m2inv[kx];
2060 ets2vf = ets2*vfactor;
2061 virxx += ets2vf*mhx[kx]*mhx[kx] - ets2;
2062 virxy += ets2vf*mhx[kx]*mhy[kx];
2063 virxz += ets2vf*mhx[kx]*mhz[kx];
2064 viryy += ets2vf*mhy[kx]*mhy[kx] - ets2;
2065 viryz += ets2vf*mhy[kx]*mhz[kx];
2066 virzz += ets2vf*mhz[kx]*mhz[kx] - ets2;
2071 /* We don't need to calculate the energy and the virial.
2072 * In this case the triclinic overhead is small.
2075 /* Two explicit loops to avoid a conditional inside the loop */
2077 for (kx = kxstart; kx < maxkx; kx++)
2082 mhyk = mx * ryx + my * ryy;
2083 mhzk = mx * rzx + my * rzy + mz * rzz;
2084 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2085 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2086 tmp1[kx] = -factor*m2k;
2089 for (kx = maxkx; kx < kxend; kx++)
2094 mhyk = mx * ryx + my * ryy;
2095 mhzk = mx * rzx + my * rzy + mz * rzz;
2096 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2097 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2098 tmp1[kx] = -factor*m2k;
2101 calc_exponentials_q(kxstart, kxend, elfac, denom, tmp1, eterm);
2103 for (kx = kxstart; kx < kxend; kx++, p0++)
2108 p0->re = d1*eterm[kx];
2109 p0->im = d2*eterm[kx];
2116 /* Update virial with local values.
2117 * The virial is symmetric by definition.
2118 * this virial seems ok for isotropic scaling, but I'm
2119 * experiencing problems on semiisotropic membranes.
2120 * IS THAT COMMENT STILL VALID??? (DvdS, 2001/02/07).
2122 work->vir_q[XX][XX] = 0.25*virxx;
2123 work->vir_q[YY][YY] = 0.25*viryy;
2124 work->vir_q[ZZ][ZZ] = 0.25*virzz;
2125 work->vir_q[XX][YY] = work->vir_q[YY][XX] = 0.25*virxy;
2126 work->vir_q[XX][ZZ] = work->vir_q[ZZ][XX] = 0.25*virxz;
2127 work->vir_q[YY][ZZ] = work->vir_q[ZZ][YY] = 0.25*viryz;
2129 /* This energy should be corrected for a charged system */
2130 work->energy_q = 0.5*energy;
2133 /* Return the loop count */
2134 return local_ndata[YY]*local_ndata[XX];
2137 static int solve_pme_lj_yzx(gmx_pme_t pme, t_complex **grid, gmx_bool bLB,
2138 real ewaldcoeff, real vol,
2139 gmx_bool bEnerVir, int nthread, int thread)
2141 /* do recip sum over local cells in grid */
2142 /* y major, z middle, x minor or continuous */
2144 int kx, ky, kz, maxkx, maxky, maxkz;
2145 int nx, ny, nz, iy, iyz0, iyz1, iyz, iz, kxstart, kxend;
2147 real factor = M_PI*M_PI/(ewaldcoeff*ewaldcoeff);
2149 real eterm, vterm, d1, d2, energy = 0;
2151 real virxx = 0, virxy = 0, virxz = 0, viryy = 0, viryz = 0, virzz = 0;
2152 real rxx, ryx, ryy, rzx, rzy, rzz;
2153 real *mhx, *mhy, *mhz, *m2, *denom, *tmp1, *tmp2;
2154 real mhxk, mhyk, mhzk, m2k;
2159 ivec local_ndata, local_offset, local_size;
2164 /* Dimensions should be identical for A/B grid, so we just use A here */
2165 gmx_parallel_3dfft_complex_limits(pme->pfft_setup[PME_GRID_C6A],
2170 rxx = pme->recipbox[XX][XX];
2171 ryx = pme->recipbox[YY][XX];
2172 ryy = pme->recipbox[YY][YY];
2173 rzx = pme->recipbox[ZZ][XX];
2174 rzy = pme->recipbox[ZZ][YY];
2175 rzz = pme->recipbox[ZZ][ZZ];
2181 work = &pme->work[thread];
2186 denom = work->denom;
2190 iyz0 = local_ndata[YY]*local_ndata[ZZ]* thread /nthread;
2191 iyz1 = local_ndata[YY]*local_ndata[ZZ]*(thread+1)/nthread;
2193 for (iyz = iyz0; iyz < iyz1; iyz++)
2195 iy = iyz/local_ndata[ZZ];
2196 iz = iyz - iy*local_ndata[ZZ];
2198 ky = iy + local_offset[YY];
2209 by = 3.0*vol*pme->bsp_mod[YY][ky]
2210 / (M_PI*sqrt(M_PI)*ewaldcoeff*ewaldcoeff*ewaldcoeff);
2212 kz = iz + local_offset[ZZ];
2216 bz = pme->bsp_mod[ZZ][kz];
2218 /* 0.5 correction for corner points */
2220 if (kz == 0 || kz == (nz+1)/2)
2225 kxstart = local_offset[XX];
2226 kxend = local_offset[XX] + local_ndata[XX];
2229 /* More expensive inner loop, especially because of the
2230 * storage of the mh elements in array's. Because x is the
2231 * minor grid index, all mh elements depend on kx for
2232 * triclinic unit cells.
2235 /* Two explicit loops to avoid a conditional inside the loop */
2236 for (kx = kxstart; kx < maxkx; kx++)
2241 mhyk = mx * ryx + my * ryy;
2242 mhzk = mx * rzx + my * rzy + mz * rzz;
2243 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2248 denom[kx] = bz*by*pme->bsp_mod[XX][kx];
2249 tmp1[kx] = -factor*m2k;
2250 tmp2[kx] = sqrt(factor*m2k);
2253 for (kx = maxkx; kx < kxend; kx++)
2258 mhyk = mx * ryx + my * ryy;
2259 mhzk = mx * rzx + my * rzy + mz * rzz;
2260 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2265 denom[kx] = bz*by*pme->bsp_mod[XX][kx];
2266 tmp1[kx] = -factor*m2k;
2267 tmp2[kx] = sqrt(factor*m2k);
2270 calc_exponentials_lj(kxstart, kxend, tmp1, tmp2, denom);
2272 for (kx = kxstart; kx < kxend; kx++)
2274 m2k = factor*m2[kx];
2275 eterm = -((1.0 - 2.0*m2k)*tmp1[kx]
2276 + 2.0*m2k*tmp2[kx]);
2277 vterm = 3.0*(-tmp1[kx] + tmp2[kx]);
2278 tmp1[kx] = eterm*denom[kx];
2279 tmp2[kx] = vterm*denom[kx];
2287 p0 = grid[0] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2288 for (kx = kxstart; kx < kxend; kx++, p0++)
2298 struct2 = 2.0*(d1*d1+d2*d2);
2300 tmp1[kx] = eterm*struct2;
2301 tmp2[kx] = vterm*struct2;
2306 real *struct2 = denom;
2309 for (kx = kxstart; kx < kxend; kx++)
2313 /* Due to symmetry we only need to calculate 4 of the 7 terms */
2314 for (ig = 0; ig <= 3; ++ig)
2319 p0 = grid[ig] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2320 p1 = grid[6-ig] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2321 scale = 2.0*lb_scale_factor_symm[ig];
2322 for (kx = kxstart; kx < kxend; ++kx, ++p0, ++p1)
2324 struct2[kx] += scale*(p0->re*p1->re + p0->im*p1->im);
2328 for (ig = 0; ig <= 6; ++ig)
2332 p0 = grid[ig] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2333 for (kx = kxstart; kx < kxend; kx++, p0++)
2343 for (kx = kxstart; kx < kxend; kx++)
2348 tmp1[kx] = eterm*str2;
2349 tmp2[kx] = vterm*str2;
2353 for (kx = kxstart; kx < kxend; kx++)
2355 ets2 = corner_fac*tmp1[kx];
2356 vterm = 2.0*factor*tmp2[kx];
2358 ets2vf = corner_fac*vterm;
2359 virxx += ets2vf*mhx[kx]*mhx[kx] - ets2;
2360 virxy += ets2vf*mhx[kx]*mhy[kx];
2361 virxz += ets2vf*mhx[kx]*mhz[kx];
2362 viryy += ets2vf*mhy[kx]*mhy[kx] - ets2;
2363 viryz += ets2vf*mhy[kx]*mhz[kx];
2364 virzz += ets2vf*mhz[kx]*mhz[kx] - ets2;
2369 /* We don't need to calculate the energy and the virial.
2370 * In this case the triclinic overhead is small.
2373 /* Two explicit loops to avoid a conditional inside the loop */
2375 for (kx = kxstart; kx < maxkx; kx++)
2380 mhyk = mx * ryx + my * ryy;
2381 mhzk = mx * rzx + my * rzy + mz * rzz;
2382 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2384 denom[kx] = bz*by*pme->bsp_mod[XX][kx];
2385 tmp1[kx] = -factor*m2k;
2386 tmp2[kx] = sqrt(factor*m2k);
2389 for (kx = maxkx; kx < kxend; kx++)
2394 mhyk = mx * ryx + my * ryy;
2395 mhzk = mx * rzx + my * rzy + mz * rzz;
2396 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2398 denom[kx] = bz*by*pme->bsp_mod[XX][kx];
2399 tmp1[kx] = -factor*m2k;
2400 tmp2[kx] = sqrt(factor*m2k);
2403 calc_exponentials_lj(kxstart, kxend, tmp1, tmp2, denom);
2405 for (kx = kxstart; kx < kxend; kx++)
2407 m2k = factor*m2[kx];
2408 eterm = -((1.0 - 2.0*m2k)*tmp1[kx]
2409 + 2.0*m2k*tmp2[kx]);
2410 tmp1[kx] = eterm*denom[kx];
2412 gcount = (bLB ? 7 : 1);
2413 for (ig = 0; ig < gcount; ++ig)
2417 p0 = grid[ig] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2418 for (kx = kxstart; kx < kxend; kx++, p0++)
2433 work->vir_lj[XX][XX] = 0.25*virxx;
2434 work->vir_lj[YY][YY] = 0.25*viryy;
2435 work->vir_lj[ZZ][ZZ] = 0.25*virzz;
2436 work->vir_lj[XX][YY] = work->vir_lj[YY][XX] = 0.25*virxy;
2437 work->vir_lj[XX][ZZ] = work->vir_lj[ZZ][XX] = 0.25*virxz;
2438 work->vir_lj[YY][ZZ] = work->vir_lj[ZZ][YY] = 0.25*viryz;
2439 /* This energy should be corrected for a charged system */
2440 work->energy_lj = 0.5*energy;
2442 /* Return the loop count */
2443 return local_ndata[YY]*local_ndata[XX];
2446 static void get_pme_ener_vir_q(const gmx_pme_t pme, int nthread,
2447 real *mesh_energy, matrix vir)
2449 /* This function sums output over threads and should therefore
2450 * only be called after thread synchronization.
2454 *mesh_energy = pme->work[0].energy_q;
2455 copy_mat(pme->work[0].vir_q, vir);
2457 for (thread = 1; thread < nthread; thread++)
2459 *mesh_energy += pme->work[thread].energy_q;
2460 m_add(vir, pme->work[thread].vir_q, vir);
2464 static void get_pme_ener_vir_lj(const gmx_pme_t pme, int nthread,
2465 real *mesh_energy, matrix vir)
2467 /* This function sums output over threads and should therefore
2468 * only be called after thread synchronization.
2472 *mesh_energy = pme->work[0].energy_lj;
2473 copy_mat(pme->work[0].vir_lj, vir);
2475 for (thread = 1; thread < nthread; thread++)
2477 *mesh_energy += pme->work[thread].energy_lj;
2478 m_add(vir, pme->work[thread].vir_lj, vir);
2483 #define DO_FSPLINE(order) \
2484 for (ithx = 0; (ithx < order); ithx++) \
2486 index_x = (i0+ithx)*pny*pnz; \
2490 for (ithy = 0; (ithy < order); ithy++) \
2492 index_xy = index_x+(j0+ithy)*pnz; \
2497 for (ithz = 0; (ithz < order); ithz++) \
2499 gval = grid[index_xy+(k0+ithz)]; \
2500 fxy1 += thz[ithz]*gval; \
2501 fz1 += dthz[ithz]*gval; \
2510 static void gather_f_bsplines(gmx_pme_t pme, real *grid,
2511 gmx_bool bClearF, pme_atomcomm_t *atc,
2512 splinedata_t *spline,
2515 /* sum forces for local particles */
2516 int nn, n, ithx, ithy, ithz, i0, j0, k0;
2517 int index_x, index_xy;
2518 int nx, ny, nz, pnx, pny, pnz;
2520 real tx, ty, dx, dy, qn;
2521 real fx, fy, fz, gval;
2523 real *thx, *thy, *thz, *dthx, *dthy, *dthz;
2525 real rxx, ryx, ryy, rzx, rzy, rzz;
2528 pme_spline_work_t *work;
2530 #if defined PME_SIMD4_SPREAD_GATHER && !defined PME_SIMD4_UNALIGNED
2531 real thz_buffer[GMX_SIMD4_WIDTH*3], *thz_aligned;
2532 real dthz_buffer[GMX_SIMD4_WIDTH*3], *dthz_aligned;
2534 thz_aligned = gmx_simd4_align_r(thz_buffer);
2535 dthz_aligned = gmx_simd4_align_r(dthz_buffer);
2538 work = pme->spline_work;
2540 order = pme->pme_order;
2541 thx = spline->theta[XX];
2542 thy = spline->theta[YY];
2543 thz = spline->theta[ZZ];
2544 dthx = spline->dtheta[XX];
2545 dthy = spline->dtheta[YY];
2546 dthz = spline->dtheta[ZZ];
2550 pnx = pme->pmegrid_nx;
2551 pny = pme->pmegrid_ny;
2552 pnz = pme->pmegrid_nz;
2554 rxx = pme->recipbox[XX][XX];
2555 ryx = pme->recipbox[YY][XX];
2556 ryy = pme->recipbox[YY][YY];
2557 rzx = pme->recipbox[ZZ][XX];
2558 rzy = pme->recipbox[ZZ][YY];
2559 rzz = pme->recipbox[ZZ][ZZ];
2561 for (nn = 0; nn < spline->n; nn++)
2563 n = spline->ind[nn];
2564 qn = scale*atc->q[n];
2577 idxptr = atc->idx[n];
2584 /* Pointer arithmetic alert, next six statements */
2585 thx = spline->theta[XX] + norder;
2586 thy = spline->theta[YY] + norder;
2587 thz = spline->theta[ZZ] + norder;
2588 dthx = spline->dtheta[XX] + norder;
2589 dthy = spline->dtheta[YY] + norder;
2590 dthz = spline->dtheta[ZZ] + norder;
2595 #ifdef PME_SIMD4_SPREAD_GATHER
2596 #ifdef PME_SIMD4_UNALIGNED
2597 #define PME_GATHER_F_SIMD4_ORDER4
2599 #define PME_GATHER_F_SIMD4_ALIGNED
2602 #include "pme_simd4.h"
2608 #ifdef PME_SIMD4_SPREAD_GATHER
2609 #define PME_GATHER_F_SIMD4_ALIGNED
2611 #include "pme_simd4.h"
2621 atc->f[n][XX] += -qn*( fx*nx*rxx );
2622 atc->f[n][YY] += -qn*( fx*nx*ryx + fy*ny*ryy );
2623 atc->f[n][ZZ] += -qn*( fx*nx*rzx + fy*ny*rzy + fz*nz*rzz );
2626 /* Since the energy and not forces are interpolated
2627 * the net force might not be exactly zero.
2628 * This can be solved by also interpolating F, but
2629 * that comes at a cost.
2630 * A better hack is to remove the net force every
2631 * step, but that must be done at a higher level
2632 * since this routine doesn't see all atoms if running
2633 * in parallel. Don't know how important it is? EL 990726
2638 static real gather_energy_bsplines(gmx_pme_t pme, real *grid,
2639 pme_atomcomm_t *atc)
2641 splinedata_t *spline;
2642 int n, ithx, ithy, ithz, i0, j0, k0;
2643 int index_x, index_xy;
2645 real energy, pot, tx, ty, qn, gval;
2646 real *thx, *thy, *thz;
2650 spline = &atc->spline[0];
2652 order = pme->pme_order;
2655 for (n = 0; (n < atc->n); n++)
2661 idxptr = atc->idx[n];
2668 /* Pointer arithmetic alert, next three statements */
2669 thx = spline->theta[XX] + norder;
2670 thy = spline->theta[YY] + norder;
2671 thz = spline->theta[ZZ] + norder;
2674 for (ithx = 0; (ithx < order); ithx++)
2676 index_x = (i0+ithx)*pme->pmegrid_ny*pme->pmegrid_nz;
2679 for (ithy = 0; (ithy < order); ithy++)
2681 index_xy = index_x+(j0+ithy)*pme->pmegrid_nz;
2684 for (ithz = 0; (ithz < order); ithz++)
2686 gval = grid[index_xy+(k0+ithz)];
2687 pot += tx*ty*thz[ithz]*gval;
2700 /* Macro to force loop unrolling by fixing order.
2701 * This gives a significant performance gain.
2703 #define CALC_SPLINE(order) \
2707 real data[PME_ORDER_MAX]; \
2708 real ddata[PME_ORDER_MAX]; \
2710 for (j = 0; (j < DIM); j++) \
2714 /* dr is relative offset from lower cell limit */ \
2715 data[order-1] = 0; \
2719 for (k = 3; (k < order); k++) \
2721 div = 1.0/(k - 1.0); \
2722 data[k-1] = div*dr*data[k-2]; \
2723 for (l = 1; (l < (k-1)); l++) \
2725 data[k-l-1] = div*((dr+l)*data[k-l-2]+(k-l-dr)* \
2728 data[0] = div*(1-dr)*data[0]; \
2730 /* differentiate */ \
2731 ddata[0] = -data[0]; \
2732 for (k = 1; (k < order); k++) \
2734 ddata[k] = data[k-1] - data[k]; \
2737 div = 1.0/(order - 1); \
2738 data[order-1] = div*dr*data[order-2]; \
2739 for (l = 1; (l < (order-1)); l++) \
2741 data[order-l-1] = div*((dr+l)*data[order-l-2]+ \
2742 (order-l-dr)*data[order-l-1]); \
2744 data[0] = div*(1 - dr)*data[0]; \
2746 for (k = 0; k < order; k++) \
2748 theta[j][i*order+k] = data[k]; \
2749 dtheta[j][i*order+k] = ddata[k]; \
2754 void make_bsplines(splinevec theta, splinevec dtheta, int order,
2755 rvec fractx[], int nr, int ind[], real charge[],
2756 gmx_bool bDoSplines)
2758 /* construct splines for local atoms */
2762 for (i = 0; i < nr; i++)
2764 /* With free energy we do not use the charge check.
2765 * In most cases this will be more efficient than calling make_bsplines
2766 * twice, since usually more than half the particles have charges.
2769 if (bDoSplines || charge[ii] != 0.0)
2774 case 4: CALC_SPLINE(4); break;
2775 case 5: CALC_SPLINE(5); break;
2776 default: CALC_SPLINE(order); break;
2783 void make_dft_mod(real *mod, real *data, int ndata)
2788 for (i = 0; i < ndata; i++)
2791 for (j = 0; j < ndata; j++)
2793 arg = (2.0*M_PI*i*j)/ndata;
2794 sc += data[j]*cos(arg);
2795 ss += data[j]*sin(arg);
2797 mod[i] = sc*sc+ss*ss;
2799 for (i = 0; i < ndata; i++)
2803 mod[i] = (mod[i-1]+mod[i+1])*0.5;
2809 static void make_bspline_moduli(splinevec bsp_mod,
2810 int nx, int ny, int nz, int order)
2812 int nmax = max(nx, max(ny, nz));
2813 real *data, *ddata, *bsp_data;
2819 snew(bsp_data, nmax);
2825 for (k = 3; k < order; k++)
2829 for (l = 1; l < (k-1); l++)
2831 data[k-l-1] = div*(l*data[k-l-2]+(k-l)*data[k-l-1]);
2833 data[0] = div*data[0];
2836 ddata[0] = -data[0];
2837 for (k = 1; k < order; k++)
2839 ddata[k] = data[k-1]-data[k];
2841 div = 1.0/(order-1);
2843 for (l = 1; l < (order-1); l++)
2845 data[order-l-1] = div*(l*data[order-l-2]+(order-l)*data[order-l-1]);
2847 data[0] = div*data[0];
2849 for (i = 0; i < nmax; i++)
2853 for (i = 1; i <= order; i++)
2855 bsp_data[i] = data[i-1];
2858 make_dft_mod(bsp_mod[XX], bsp_data, nx);
2859 make_dft_mod(bsp_mod[YY], bsp_data, ny);
2860 make_dft_mod(bsp_mod[ZZ], bsp_data, nz);
2868 /* Return the P3M optimal influence function */
2869 static double do_p3m_influence(double z, int order)
2876 /* The formula and most constants can be found in:
2877 * Ballenegger et al., JCTC 8, 936 (2012)
2882 return 1.0 - 2.0*z2/3.0;
2885 return 1.0 - z2 + 2.0*z4/15.0;
2888 return 1.0 - 4.0*z2/3.0 + 2.0*z4/5.0 + 4.0*z2*z4/315.0;
2891 return 1.0 - 5.0*z2/3.0 + 7.0*z4/9.0 - 17.0*z2*z4/189.0 + 2.0*z4*z4/2835.0;
2894 return 1.0 - 2.0*z2 + 19.0*z4/15.0 - 256.0*z2*z4/945.0 + 62.0*z4*z4/4725.0 + 4.0*z2*z4*z4/155925.0;
2897 return 1.0 - 7.0*z2/3.0 + 28.0*z4/15.0 - 16.0*z2*z4/27.0 + 26.0*z4*z4/405.0 - 2.0*z2*z4*z4/1485.0 + 4.0*z4*z4*z4/6081075.0;
2899 return 1.0 - 8.0*z2/3.0 + 116.0*z4/45.0 - 344.0*z2*z4/315.0 + 914.0*z4*z4/4725.0 - 248.0*z4*z4*z2/22275.0 + 21844.0*z4*z4*z4/212837625.0 - 8.0*z4*z4*z4*z2/638512875.0;
2906 /* Calculate the P3M B-spline moduli for one dimension */
2907 static void make_p3m_bspline_moduli_dim(real *bsp_mod, int n, int order)
2909 double zarg, zai, sinzai, infl;
2914 gmx_fatal(FARGS, "The current P3M code only supports orders up to 8");
2921 for (i = -maxk; i < 0; i++)
2925 infl = do_p3m_influence(sinzai, order);
2926 bsp_mod[n+i] = infl*infl*pow(sinzai/zai, -2.0*order);
2929 for (i = 1; i < maxk; i++)
2933 infl = do_p3m_influence(sinzai, order);
2934 bsp_mod[i] = infl*infl*pow(sinzai/zai, -2.0*order);
2938 /* Calculate the P3M B-spline moduli */
2939 static void make_p3m_bspline_moduli(splinevec bsp_mod,
2940 int nx, int ny, int nz, int order)
2942 make_p3m_bspline_moduli_dim(bsp_mod[XX], nx, order);
2943 make_p3m_bspline_moduli_dim(bsp_mod[YY], ny, order);
2944 make_p3m_bspline_moduli_dim(bsp_mod[ZZ], nz, order);
2948 static void setup_coordinate_communication(pme_atomcomm_t *atc)
2956 for (i = 1; i <= nslab/2; i++)
2958 fw = (atc->nodeid + i) % nslab;
2959 bw = (atc->nodeid - i + nslab) % nslab;
2962 atc->node_dest[n] = fw;
2963 atc->node_src[n] = bw;
2968 atc->node_dest[n] = bw;
2969 atc->node_src[n] = fw;
2975 int gmx_pme_destroy(FILE *log, gmx_pme_t *pmedata)
2981 fprintf(log, "Destroying PME data structures.\n");
2984 sfree((*pmedata)->nnx);
2985 sfree((*pmedata)->nny);
2986 sfree((*pmedata)->nnz);
2988 for (i = 0; i < (*pmedata)->ngrids; ++i)
2990 pmegrids_destroy(&(*pmedata)->pmegrid[i]);
2991 sfree((*pmedata)->fftgrid[i]);
2992 sfree((*pmedata)->cfftgrid[i]);
2993 gmx_parallel_3dfft_destroy((*pmedata)->pfft_setup[i]);
2996 sfree((*pmedata)->lb_buf1);
2997 sfree((*pmedata)->lb_buf2);
2999 for (thread = 0; thread < (*pmedata)->nthread; thread++)
3001 free_work(&(*pmedata)->work[thread]);
3003 sfree((*pmedata)->work);
3011 static int mult_up(int n, int f)
3013 return ((n + f - 1)/f)*f;
3017 static double pme_load_imbalance(gmx_pme_t pme)
3022 nma = pme->nnodes_major;
3023 nmi = pme->nnodes_minor;
3025 n1 = mult_up(pme->nkx, nma)*mult_up(pme->nky, nmi)*pme->nkz;
3026 n2 = mult_up(pme->nkx, nma)*mult_up(pme->nkz, nmi)*pme->nky;
3027 n3 = mult_up(pme->nky, nma)*mult_up(pme->nkz, nmi)*pme->nkx;
3029 /* pme_solve is roughly double the cost of an fft */
3031 return (n1 + n2 + 3*n3)/(double)(6*pme->nkx*pme->nky*pme->nkz);
3034 static void init_atomcomm(gmx_pme_t pme, pme_atomcomm_t *atc,
3035 int dimind, gmx_bool bSpread)
3037 int nk, k, s, thread;
3039 atc->dimind = dimind;
3044 if (pme->nnodes > 1)
3046 atc->mpi_comm = pme->mpi_comm_d[dimind];
3047 MPI_Comm_size(atc->mpi_comm, &atc->nslab);
3048 MPI_Comm_rank(atc->mpi_comm, &atc->nodeid);
3052 fprintf(debug, "For PME atom communication in dimind %d: nslab %d rank %d\n", atc->dimind, atc->nslab, atc->nodeid);
3056 atc->bSpread = bSpread;
3057 atc->pme_order = pme->pme_order;
3061 snew(atc->node_dest, atc->nslab);
3062 snew(atc->node_src, atc->nslab);
3063 setup_coordinate_communication(atc);
3065 snew(atc->count_thread, pme->nthread);
3066 for (thread = 0; thread < pme->nthread; thread++)
3068 snew(atc->count_thread[thread], atc->nslab);
3070 atc->count = atc->count_thread[0];
3071 snew(atc->rcount, atc->nslab);
3072 snew(atc->buf_index, atc->nslab);
3075 atc->nthread = pme->nthread;
3076 if (atc->nthread > 1)
3078 snew(atc->thread_plist, atc->nthread);
3080 snew(atc->spline, atc->nthread);
3081 for (thread = 0; thread < atc->nthread; thread++)
3083 if (atc->nthread > 1)
3085 snew(atc->thread_plist[thread].n, atc->nthread+2*GMX_CACHE_SEP);
3086 atc->thread_plist[thread].n += GMX_CACHE_SEP;
3088 snew(atc->spline[thread].thread_one, pme->nthread);
3089 atc->spline[thread].thread_one[thread] = 1;
3094 init_overlap_comm(pme_overlap_t * ol,
3104 int lbnd, rbnd, maxlr, b, i;
3107 pme_grid_comm_t *pgc;
3109 int fft_start, fft_end, send_index1, recv_index1;
3113 ol->mpi_comm = comm;
3116 ol->nnodes = nnodes;
3117 ol->nodeid = nodeid;
3119 /* Linear translation of the PME grid won't affect reciprocal space
3120 * calculations, so to optimize we only interpolate "upwards",
3121 * which also means we only have to consider overlap in one direction.
3122 * I.e., particles on this node might also be spread to grid indices
3123 * that belong to higher nodes (modulo nnodes)
3126 snew(ol->s2g0, ol->nnodes+1);
3127 snew(ol->s2g1, ol->nnodes);
3130 fprintf(debug, "PME slab boundaries:");
3132 for (i = 0; i < nnodes; i++)
3134 /* s2g0 the local interpolation grid start.
3135 * s2g1 the local interpolation grid end.
3136 * Because grid overlap communication only goes forward,
3137 * the grid the slabs for fft's should be rounded down.
3139 ol->s2g0[i] = ( i *ndata + 0 )/nnodes;
3140 ol->s2g1[i] = ((i+1)*ndata + nnodes-1)/nnodes + norder - 1;
3144 fprintf(debug, " %3d %3d", ol->s2g0[i], ol->s2g1[i]);
3147 ol->s2g0[nnodes] = ndata;
3150 fprintf(debug, "\n");
3153 /* Determine with how many nodes we need to communicate the grid overlap */
3159 for (i = 0; i < nnodes; i++)
3161 if ((i+b < nnodes && ol->s2g1[i] > ol->s2g0[i+b]) ||
3162 (i+b >= nnodes && ol->s2g1[i] > ol->s2g0[i+b-nnodes] + ndata))
3168 while (bCont && b < nnodes);
3169 ol->noverlap_nodes = b - 1;
3171 snew(ol->send_id, ol->noverlap_nodes);
3172 snew(ol->recv_id, ol->noverlap_nodes);
3173 for (b = 0; b < ol->noverlap_nodes; b++)
3175 ol->send_id[b] = (ol->nodeid + (b + 1)) % ol->nnodes;
3176 ol->recv_id[b] = (ol->nodeid - (b + 1) + ol->nnodes) % ol->nnodes;
3178 snew(ol->comm_data, ol->noverlap_nodes);
3181 for (b = 0; b < ol->noverlap_nodes; b++)
3183 pgc = &ol->comm_data[b];
3185 fft_start = ol->s2g0[ol->send_id[b]];
3186 fft_end = ol->s2g0[ol->send_id[b]+1];
3187 if (ol->send_id[b] < nodeid)
3192 send_index1 = ol->s2g1[nodeid];
3193 send_index1 = min(send_index1, fft_end);
3194 pgc->send_index0 = fft_start;
3195 pgc->send_nindex = max(0, send_index1 - pgc->send_index0);
3196 ol->send_size += pgc->send_nindex;
3198 /* We always start receiving to the first index of our slab */
3199 fft_start = ol->s2g0[ol->nodeid];
3200 fft_end = ol->s2g0[ol->nodeid+1];
3201 recv_index1 = ol->s2g1[ol->recv_id[b]];
3202 if (ol->recv_id[b] > nodeid)
3204 recv_index1 -= ndata;
3206 recv_index1 = min(recv_index1, fft_end);
3207 pgc->recv_index0 = fft_start;
3208 pgc->recv_nindex = max(0, recv_index1 - pgc->recv_index0);
3212 /* Communicate the buffer sizes to receive */
3213 for (b = 0; b < ol->noverlap_nodes; b++)
3215 MPI_Sendrecv(&ol->send_size, 1, MPI_INT, ol->send_id[b], b,
3216 &ol->comm_data[b].recv_size, 1, MPI_INT, ol->recv_id[b], b,
3217 ol->mpi_comm, &stat);
3221 /* For non-divisible grid we need pme_order iso pme_order-1 */
3222 snew(ol->sendbuf, norder*commplainsize);
3223 snew(ol->recvbuf, norder*commplainsize);
3227 make_gridindex5_to_localindex(int n, int local_start, int local_range,
3228 int **global_to_local,
3229 real **fraction_shift)
3237 for (i = 0; (i < 5*n); i++)
3239 /* Determine the global to local grid index */
3240 gtl[i] = (i - local_start + n) % n;
3241 /* For coordinates that fall within the local grid the fraction
3242 * is correct, we don't need to shift it.
3245 if (local_range < n)
3247 /* Due to rounding issues i could be 1 beyond the lower or
3248 * upper boundary of the local grid. Correct the index for this.
3249 * If we shift the index, we need to shift the fraction by
3250 * the same amount in the other direction to not affect
3252 * Note that due to this shifting the weights at the end of
3253 * the spline might change, but that will only involve values
3254 * between zero and values close to the precision of a real,
3255 * which is anyhow the accuracy of the whole mesh calculation.
3257 /* With local_range=0 we should not change i=local_start */
3258 if (i % n != local_start)
3265 else if (gtl[i] == local_range)
3267 gtl[i] = local_range - 1;
3274 *global_to_local = gtl;
3275 *fraction_shift = fsh;
3278 static pme_spline_work_t *make_pme_spline_work(int gmx_unused order)
3280 pme_spline_work_t *work;
3282 #ifdef PME_SIMD4_SPREAD_GATHER
3283 real tmp[GMX_SIMD4_WIDTH*3], *tmp_aligned;
3284 gmx_simd4_real_t zero_S;
3285 gmx_simd4_real_t real_mask_S0, real_mask_S1;
3288 snew_aligned(work, 1, SIMD4_ALIGNMENT);
3290 tmp_aligned = gmx_simd4_align_r(tmp);
3292 zero_S = gmx_simd4_setzero_r();
3294 /* Generate bit masks to mask out the unused grid entries,
3295 * as we only operate on order of the 8 grid entries that are
3296 * load into 2 SIMD registers.
3298 for (of = 0; of < 2*GMX_SIMD4_WIDTH-(order-1); of++)
3300 for (i = 0; i < 2*GMX_SIMD4_WIDTH; i++)
3302 tmp_aligned[i] = (i >= of && i < of+order ? -1.0 : 1.0);
3304 real_mask_S0 = gmx_simd4_load_r(tmp_aligned);
3305 real_mask_S1 = gmx_simd4_load_r(tmp_aligned+GMX_SIMD4_WIDTH);
3306 work->mask_S0[of] = gmx_simd4_cmplt_r(real_mask_S0, zero_S);
3307 work->mask_S1[of] = gmx_simd4_cmplt_r(real_mask_S1, zero_S);
3316 void gmx_pme_check_restrictions(int pme_order,
3317 int nkx, int nky, int nkz,
3320 gmx_bool bUseThreads,
3322 gmx_bool *bValidSettings)
3324 if (pme_order > PME_ORDER_MAX)
3328 *bValidSettings = FALSE;
3331 gmx_fatal(FARGS, "pme_order (%d) is larger than the maximum allowed value (%d). Modify and recompile the code if you really need such a high order.",
3332 pme_order, PME_ORDER_MAX);
3335 if (nkx <= pme_order*(nnodes_major > 1 ? 2 : 1) ||
3336 nky <= pme_order*(nnodes_minor > 1 ? 2 : 1) ||
3341 *bValidSettings = FALSE;
3344 gmx_fatal(FARGS, "The PME grid sizes need to be larger than pme_order (%d) and for dimensions with domain decomposition larger than 2*pme_order",
3348 /* Check for a limitation of the (current) sum_fftgrid_dd code.
3349 * We only allow multiple communication pulses in dim 1, not in dim 0.
3351 if (bUseThreads && (nkx < nnodes_major*pme_order &&
3352 nkx != nnodes_major*(pme_order - 1)))
3356 *bValidSettings = FALSE;
3359 gmx_fatal(FARGS, "The number of PME grid lines per node along x is %g. But when using OpenMP threads, the number of grid lines per node along x should be >= pme_order (%d) or = pmeorder-1. To resolve this issue, use less nodes along x (and possibly more along y and/or z) by specifying -dd manually.",
3360 nkx/(double)nnodes_major, pme_order);
3363 if (bValidSettings != NULL)
3365 *bValidSettings = TRUE;
3371 int gmx_pme_init(gmx_pme_t * pmedata,
3377 gmx_bool bFreeEnergy_q,
3378 gmx_bool bFreeEnergy_lj,
3379 gmx_bool bReproducible,
3382 gmx_pme_t pme = NULL;
3384 int use_threads, sum_use_threads, i;
3389 fprintf(debug, "Creating PME data structures.\n");
3393 pme->sum_qgrid_tmp = NULL;
3394 pme->sum_qgrid_dd_tmp = NULL;
3395 pme->buf_nalloc = 0;
3398 pme->bPPnode = TRUE;
3400 pme->nnodes_major = nnodes_major;
3401 pme->nnodes_minor = nnodes_minor;
3404 if (nnodes_major*nnodes_minor > 1)
3406 pme->mpi_comm = cr->mpi_comm_mygroup;
3408 MPI_Comm_rank(pme->mpi_comm, &pme->nodeid);
3409 MPI_Comm_size(pme->mpi_comm, &pme->nnodes);
3410 if (pme->nnodes != nnodes_major*nnodes_minor)
3412 gmx_incons("PME node count mismatch");
3417 pme->mpi_comm = MPI_COMM_NULL;
3421 if (pme->nnodes == 1)
3424 pme->mpi_comm_d[0] = MPI_COMM_NULL;
3425 pme->mpi_comm_d[1] = MPI_COMM_NULL;
3427 pme->ndecompdim = 0;
3428 pme->nodeid_major = 0;
3429 pme->nodeid_minor = 0;
3431 pme->mpi_comm_d[0] = pme->mpi_comm_d[1] = MPI_COMM_NULL;
3436 if (nnodes_minor == 1)
3439 pme->mpi_comm_d[0] = pme->mpi_comm;
3440 pme->mpi_comm_d[1] = MPI_COMM_NULL;
3442 pme->ndecompdim = 1;
3443 pme->nodeid_major = pme->nodeid;
3444 pme->nodeid_minor = 0;
3447 else if (nnodes_major == 1)
3450 pme->mpi_comm_d[0] = MPI_COMM_NULL;
3451 pme->mpi_comm_d[1] = pme->mpi_comm;
3453 pme->ndecompdim = 1;
3454 pme->nodeid_major = 0;
3455 pme->nodeid_minor = pme->nodeid;
3459 if (pme->nnodes % nnodes_major != 0)
3461 gmx_incons("For 2D PME decomposition, #PME nodes must be divisible by the number of nodes in the major dimension");
3463 pme->ndecompdim = 2;
3466 MPI_Comm_split(pme->mpi_comm, pme->nodeid % nnodes_minor,
3467 pme->nodeid, &pme->mpi_comm_d[0]); /* My communicator along major dimension */
3468 MPI_Comm_split(pme->mpi_comm, pme->nodeid/nnodes_minor,
3469 pme->nodeid, &pme->mpi_comm_d[1]); /* My communicator along minor dimension */
3471 MPI_Comm_rank(pme->mpi_comm_d[0], &pme->nodeid_major);
3472 MPI_Comm_size(pme->mpi_comm_d[0], &pme->nnodes_major);
3473 MPI_Comm_rank(pme->mpi_comm_d[1], &pme->nodeid_minor);
3474 MPI_Comm_size(pme->mpi_comm_d[1], &pme->nnodes_minor);
3477 pme->bPPnode = (cr->duty & DUTY_PP);
3480 pme->nthread = nthread;
3482 /* Check if any of the PME MPI ranks uses threads */
3483 use_threads = (pme->nthread > 1 ? 1 : 0);
3485 if (pme->nnodes > 1)
3487 MPI_Allreduce(&use_threads, &sum_use_threads, 1, MPI_INT,
3488 MPI_SUM, pme->mpi_comm);
3493 sum_use_threads = use_threads;
3495 pme->bUseThreads = (sum_use_threads > 0);
3497 if (ir->ePBC == epbcSCREW)
3499 gmx_fatal(FARGS, "pme does not (yet) work with pbc = screw");
3502 pme->bFEP_q = ((ir->efep != efepNO) && bFreeEnergy_q);
3503 pme->bFEP_lj = ((ir->efep != efepNO) && bFreeEnergy_lj);
3504 pme->bFEP = (pme->bFEP_q || pme->bFEP_lj);
3508 pme->bP3M = (ir->coulombtype == eelP3M_AD || getenv("GMX_PME_P3M") != NULL);
3509 pme->pme_order = ir->pme_order;
3511 /* Always constant electrostatics parameters */
3512 pme->epsilon_r = ir->epsilon_r;
3514 /* Always constant LJ parameters */
3515 pme->ljpme_combination_rule = ir->ljpme_combination_rule;
3517 /* If we violate restrictions, generate a fatal error here */
3518 gmx_pme_check_restrictions(pme->pme_order,
3519 pme->nkx, pme->nky, pme->nkz,
3526 if (pme->nnodes > 1)
3531 MPI_Type_contiguous(DIM, mpi_type, &(pme->rvec_mpi));
3532 MPI_Type_commit(&(pme->rvec_mpi));
3535 /* Note that the charge spreading and force gathering, which usually
3536 * takes about the same amount of time as FFT+solve_pme,
3537 * is always fully load balanced
3538 * (unless the charge distribution is inhomogeneous).
3541 imbal = pme_load_imbalance(pme);
3542 if (imbal >= 1.2 && pme->nodeid_major == 0 && pme->nodeid_minor == 0)
3546 "NOTE: The load imbalance in PME FFT and solve is %d%%.\n"
3547 " For optimal PME load balancing\n"
3548 " PME grid_x (%d) and grid_y (%d) should be divisible by #PME_nodes_x (%d)\n"
3549 " and PME grid_y (%d) and grid_z (%d) should be divisible by #PME_nodes_y (%d)\n"
3551 (int)((imbal-1)*100 + 0.5),
3552 pme->nkx, pme->nky, pme->nnodes_major,
3553 pme->nky, pme->nkz, pme->nnodes_minor);
3557 /* For non-divisible grid we need pme_order iso pme_order-1 */
3558 /* In sum_qgrid_dd x overlap is copied in place: take padding into account.
3559 * y is always copied through a buffer: we don't need padding in z,
3560 * but we do need the overlap in x because of the communication order.
3562 init_overlap_comm(&pme->overlap[0], pme->pme_order,
3566 pme->nnodes_major, pme->nodeid_major,
3568 (div_round_up(pme->nky, pme->nnodes_minor)+pme->pme_order)*(pme->nkz+pme->pme_order-1));
3570 /* Along overlap dim 1 we can send in multiple pulses in sum_fftgrid_dd.
3571 * We do this with an offset buffer of equal size, so we need to allocate
3572 * extra for the offset. That's what the (+1)*pme->nkz is for.
3574 init_overlap_comm(&pme->overlap[1], pme->pme_order,
3578 pme->nnodes_minor, pme->nodeid_minor,
3580 (div_round_up(pme->nkx, pme->nnodes_major)+pme->pme_order+1)*pme->nkz);
3582 /* Double-check for a limitation of the (current) sum_fftgrid_dd code.
3583 * Note that gmx_pme_check_restrictions checked for this already.
3585 if (pme->bUseThreads && pme->overlap[0].noverlap_nodes > 1)
3587 gmx_incons("More than one communication pulse required for grid overlap communication along the major dimension while using threads");
3590 snew(pme->bsp_mod[XX], pme->nkx);
3591 snew(pme->bsp_mod[YY], pme->nky);
3592 snew(pme->bsp_mod[ZZ], pme->nkz);
3594 /* The required size of the interpolation grid, including overlap.
3595 * The allocated size (pmegrid_n?) might be slightly larger.
3597 pme->pmegrid_nx = pme->overlap[0].s2g1[pme->nodeid_major] -
3598 pme->overlap[0].s2g0[pme->nodeid_major];
3599 pme->pmegrid_ny = pme->overlap[1].s2g1[pme->nodeid_minor] -
3600 pme->overlap[1].s2g0[pme->nodeid_minor];
3601 pme->pmegrid_nz_base = pme->nkz;
3602 pme->pmegrid_nz = pme->pmegrid_nz_base + pme->pme_order - 1;
3603 set_grid_alignment(&pme->pmegrid_nz, pme->pme_order);
3605 pme->pmegrid_start_ix = pme->overlap[0].s2g0[pme->nodeid_major];
3606 pme->pmegrid_start_iy = pme->overlap[1].s2g0[pme->nodeid_minor];
3607 pme->pmegrid_start_iz = 0;
3609 make_gridindex5_to_localindex(pme->nkx,
3610 pme->pmegrid_start_ix,
3611 pme->pmegrid_nx - (pme->pme_order-1),
3612 &pme->nnx, &pme->fshx);
3613 make_gridindex5_to_localindex(pme->nky,
3614 pme->pmegrid_start_iy,
3615 pme->pmegrid_ny - (pme->pme_order-1),
3616 &pme->nny, &pme->fshy);
3617 make_gridindex5_to_localindex(pme->nkz,
3618 pme->pmegrid_start_iz,
3619 pme->pmegrid_nz_base,
3620 &pme->nnz, &pme->fshz);
3622 pme->spline_work = make_pme_spline_work(pme->pme_order);
3624 ndata[0] = pme->nkx;
3625 ndata[1] = pme->nky;
3626 ndata[2] = pme->nkz;
3627 /* It doesn't matter if we allocate too many grids here,
3628 * we only allocate and use the ones we need.
3630 if (EVDW_PME(ir->vdwtype))
3632 pme->ngrids = ((ir->ljpme_combination_rule == eljpmeLB) ? DO_Q_AND_LJ_LB : DO_Q_AND_LJ);
3638 snew(pme->fftgrid, pme->ngrids);
3639 snew(pme->cfftgrid, pme->ngrids);
3640 snew(pme->pfft_setup, pme->ngrids);
3642 for (i = 0; i < pme->ngrids; ++i)
3644 if ((i < DO_Q && EEL_PME(ir->coulombtype) && (i == 0 ||
3646 (i >= DO_Q && EVDW_PME(ir->vdwtype) && (i == 2 ||
3648 ir->ljpme_combination_rule == eljpmeLB)))
3650 pmegrids_init(&pme->pmegrid[i],
3651 pme->pmegrid_nx, pme->pmegrid_ny, pme->pmegrid_nz,
3652 pme->pmegrid_nz_base,
3656 pme->overlap[0].s2g1[pme->nodeid_major]-pme->overlap[0].s2g0[pme->nodeid_major+1],
3657 pme->overlap[1].s2g1[pme->nodeid_minor]-pme->overlap[1].s2g0[pme->nodeid_minor+1]);
3658 /* This routine will allocate the grid data to fit the FFTs */
3659 gmx_parallel_3dfft_init(&pme->pfft_setup[i], ndata,
3660 &pme->fftgrid[i], &pme->cfftgrid[i],
3662 bReproducible, pme->nthread);
3669 /* Use plain SPME B-spline interpolation */
3670 make_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
3674 /* Use the P3M grid-optimized influence function */
3675 make_p3m_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
3678 /* Use atc[0] for spreading */
3679 init_atomcomm(pme, &pme->atc[0], nnodes_major > 1 ? 0 : 1, TRUE);
3680 if (pme->ndecompdim >= 2)
3682 init_atomcomm(pme, &pme->atc[1], 1, FALSE);
3685 if (pme->nnodes == 1)
3687 pme->atc[0].n = homenr;
3688 pme_realloc_atomcomm_things(&pme->atc[0]);
3691 pme->lb_buf1 = NULL;
3692 pme->lb_buf2 = NULL;
3693 pme->lb_buf_nalloc = 0;
3698 /* Use fft5d, order after FFT is y major, z, x minor */
3700 snew(pme->work, pme->nthread);
3701 for (thread = 0; thread < pme->nthread; thread++)
3703 realloc_work(&pme->work[thread], pme->nkx);
3712 static void reuse_pmegrids(const pmegrids_t *old, pmegrids_t *new)
3716 for (d = 0; d < DIM; d++)
3718 if (new->grid.n[d] > old->grid.n[d])
3724 sfree_aligned(new->grid.grid);
3725 new->grid.grid = old->grid.grid;
3727 if (new->grid_th != NULL && new->nthread == old->nthread)
3729 sfree_aligned(new->grid_all);
3730 for (t = 0; t < new->nthread; t++)
3732 new->grid_th[t].grid = old->grid_th[t].grid;
3737 int gmx_pme_reinit(gmx_pme_t * pmedata,
3740 const t_inputrec * ir,
3748 irc.nkx = grid_size[XX];
3749 irc.nky = grid_size[YY];
3750 irc.nkz = grid_size[ZZ];
3752 if (pme_src->nnodes == 1)
3754 homenr = pme_src->atc[0].n;
3761 ret = gmx_pme_init(pmedata, cr, pme_src->nnodes_major, pme_src->nnodes_minor,
3762 &irc, homenr, pme_src->bFEP_q, pme_src->bFEP_lj, FALSE, pme_src->nthread);
3766 /* We can easily reuse the allocated pme grids in pme_src */
3767 reuse_pmegrids(&pme_src->pmegrid[PME_GRID_QA], &(*pmedata)->pmegrid[PME_GRID_QA]);
3768 /* We would like to reuse the fft grids, but that's harder */
3775 static void copy_local_grid(gmx_pme_t pme,
3776 pmegrids_t *pmegrids, int thread, real *fftgrid)
3778 ivec local_fft_ndata, local_fft_offset, local_fft_size;
3782 int offx, offy, offz, x, y, z, i0, i0t;
3787 gmx_parallel_3dfft_real_limits(pme->pfft_setup[PME_GRID_QA],
3791 fft_my = local_fft_size[YY];
3792 fft_mz = local_fft_size[ZZ];
3794 pmegrid = &pmegrids->grid_th[thread];
3796 nsx = pmegrid->s[XX];
3797 nsy = pmegrid->s[YY];
3798 nsz = pmegrid->s[ZZ];
3800 for (d = 0; d < DIM; d++)
3802 nf[d] = min(pmegrid->n[d] - (pmegrid->order - 1),
3803 local_fft_ndata[d] - pmegrid->offset[d]);
3806 offx = pmegrid->offset[XX];
3807 offy = pmegrid->offset[YY];
3808 offz = pmegrid->offset[ZZ];
3810 /* Directly copy the non-overlapping parts of the local grids.
3811 * This also initializes the full grid.
3813 grid_th = pmegrid->grid;
3814 for (x = 0; x < nf[XX]; x++)
3816 for (y = 0; y < nf[YY]; y++)
3818 i0 = ((offx + x)*fft_my + (offy + y))*fft_mz + offz;
3819 i0t = (x*nsy + y)*nsz;
3820 for (z = 0; z < nf[ZZ]; z++)
3822 fftgrid[i0+z] = grid_th[i0t+z];
3829 reduce_threadgrid_overlap(gmx_pme_t pme,
3830 const pmegrids_t *pmegrids, int thread,
3831 real *fftgrid, real *commbuf_x, real *commbuf_y)
3833 ivec local_fft_ndata, local_fft_offset, local_fft_size;
3834 int fft_nx, fft_ny, fft_nz;
3839 int offx, offy, offz, x, y, z, i0, i0t;
3840 int sx, sy, sz, fx, fy, fz, tx1, ty1, tz1, ox, oy, oz;
3841 gmx_bool bClearBufX, bClearBufY, bClearBufXY, bClearBuf;
3842 gmx_bool bCommX, bCommY;
3845 const pmegrid_t *pmegrid, *pmegrid_g, *pmegrid_f;
3846 const real *grid_th;
3847 real *commbuf = NULL;
3849 gmx_parallel_3dfft_real_limits(pme->pfft_setup[PME_GRID_QA],
3853 fft_nx = local_fft_ndata[XX];
3854 fft_ny = local_fft_ndata[YY];
3855 fft_nz = local_fft_ndata[ZZ];
3857 fft_my = local_fft_size[YY];
3858 fft_mz = local_fft_size[ZZ];
3860 /* This routine is called when all thread have finished spreading.
3861 * Here each thread sums grid contributions calculated by other threads
3862 * to the thread local grid volume.
3863 * To minimize the number of grid copying operations,
3864 * this routines sums immediately from the pmegrid to the fftgrid.
3867 /* Determine which part of the full node grid we should operate on,
3868 * this is our thread local part of the full grid.
3870 pmegrid = &pmegrids->grid_th[thread];
3872 for (d = 0; d < DIM; d++)
3874 ne[d] = min(pmegrid->offset[d]+pmegrid->n[d]-(pmegrid->order-1),
3875 local_fft_ndata[d]);
3878 offx = pmegrid->offset[XX];
3879 offy = pmegrid->offset[YY];
3880 offz = pmegrid->offset[ZZ];
3887 /* Now loop over all the thread data blocks that contribute
3888 * to the grid region we (our thread) are operating on.
3890 /* Note that ffy_nx/y is equal to the number of grid points
3891 * between the first point of our node grid and the one of the next node.
3893 for (sx = 0; sx >= -pmegrids->nthread_comm[XX]; sx--)
3895 fx = pmegrid->ci[XX] + sx;
3900 fx += pmegrids->nc[XX];
3902 bCommX = (pme->nnodes_major > 1);
3904 pmegrid_g = &pmegrids->grid_th[fx*pmegrids->nc[YY]*pmegrids->nc[ZZ]];
3905 ox += pmegrid_g->offset[XX];
3908 tx1 = min(ox + pmegrid_g->n[XX], ne[XX]);
3912 tx1 = min(ox + pmegrid_g->n[XX], pme->pme_order);
3915 for (sy = 0; sy >= -pmegrids->nthread_comm[YY]; sy--)
3917 fy = pmegrid->ci[YY] + sy;
3922 fy += pmegrids->nc[YY];
3924 bCommY = (pme->nnodes_minor > 1);
3926 pmegrid_g = &pmegrids->grid_th[fy*pmegrids->nc[ZZ]];
3927 oy += pmegrid_g->offset[YY];
3930 ty1 = min(oy + pmegrid_g->n[YY], ne[YY]);
3934 ty1 = min(oy + pmegrid_g->n[YY], pme->pme_order);
3937 for (sz = 0; sz >= -pmegrids->nthread_comm[ZZ]; sz--)
3939 fz = pmegrid->ci[ZZ] + sz;
3943 fz += pmegrids->nc[ZZ];
3946 pmegrid_g = &pmegrids->grid_th[fz];
3947 oz += pmegrid_g->offset[ZZ];
3948 tz1 = min(oz + pmegrid_g->n[ZZ], ne[ZZ]);
3950 if (sx == 0 && sy == 0 && sz == 0)
3952 /* We have already added our local contribution
3953 * before calling this routine, so skip it here.
3958 thread_f = (fx*pmegrids->nc[YY] + fy)*pmegrids->nc[ZZ] + fz;
3960 pmegrid_f = &pmegrids->grid_th[thread_f];
3962 grid_th = pmegrid_f->grid;
3964 nsx = pmegrid_f->s[XX];
3965 nsy = pmegrid_f->s[YY];
3966 nsz = pmegrid_f->s[ZZ];
3968 #ifdef DEBUG_PME_REDUCE
3969 printf("n%d t%d add %d %2d %2d %2d %2d %2d %2d %2d-%2d %2d-%2d, %2d-%2d %2d-%2d, %2d-%2d %2d-%2d\n",
3970 pme->nodeid, thread, thread_f,
3971 pme->pmegrid_start_ix,
3972 pme->pmegrid_start_iy,
3973 pme->pmegrid_start_iz,
3975 offx-ox, tx1-ox, offx, tx1,
3976 offy-oy, ty1-oy, offy, ty1,
3977 offz-oz, tz1-oz, offz, tz1);
3980 if (!(bCommX || bCommY))
3982 /* Copy from the thread local grid to the node grid */
3983 for (x = offx; x < tx1; x++)
3985 for (y = offy; y < ty1; y++)
3987 i0 = (x*fft_my + y)*fft_mz;
3988 i0t = ((x - ox)*nsy + (y - oy))*nsz - oz;
3989 for (z = offz; z < tz1; z++)
3991 fftgrid[i0+z] += grid_th[i0t+z];
3998 /* The order of this conditional decides
3999 * where the corner volume gets stored with x+y decomp.
4003 commbuf = commbuf_y;
4004 buf_my = ty1 - offy;
4007 /* We index commbuf modulo the local grid size */
4008 commbuf += buf_my*fft_nx*fft_nz;
4010 bClearBuf = bClearBufXY;
4011 bClearBufXY = FALSE;
4015 bClearBuf = bClearBufY;
4021 commbuf = commbuf_x;
4023 bClearBuf = bClearBufX;
4027 /* Copy to the communication buffer */
4028 for (x = offx; x < tx1; x++)
4030 for (y = offy; y < ty1; y++)
4032 i0 = (x*buf_my + y)*fft_nz;
4033 i0t = ((x - ox)*nsy + (y - oy))*nsz - oz;
4037 /* First access of commbuf, initialize it */
4038 for (z = offz; z < tz1; z++)
4040 commbuf[i0+z] = grid_th[i0t+z];
4045 for (z = offz; z < tz1; z++)
4047 commbuf[i0+z] += grid_th[i0t+z];
4059 static void sum_fftgrid_dd(gmx_pme_t pme, real *fftgrid)
4061 ivec local_fft_ndata, local_fft_offset, local_fft_size;
4062 pme_overlap_t *overlap;
4063 int send_index0, send_nindex;
4068 int send_size_y, recv_size_y;
4069 int ipulse, send_id, recv_id, datasize, gridsize, size_yx;
4070 real *sendptr, *recvptr;
4071 int x, y, z, indg, indb;
4073 /* Note that this routine is only used for forward communication.
4074 * Since the force gathering, unlike the charge spreading,
4075 * can be trivially parallelized over the particles,
4076 * the backwards process is much simpler and can use the "old"
4077 * communication setup.
4080 gmx_parallel_3dfft_real_limits(pme->pfft_setup[PME_GRID_QA],
4085 if (pme->nnodes_minor > 1)
4087 /* Major dimension */
4088 overlap = &pme->overlap[1];
4090 if (pme->nnodes_major > 1)
4092 size_yx = pme->overlap[0].comm_data[0].send_nindex;
4098 datasize = (local_fft_ndata[XX] + size_yx)*local_fft_ndata[ZZ];
4100 send_size_y = overlap->send_size;
4102 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
4104 send_id = overlap->send_id[ipulse];
4105 recv_id = overlap->recv_id[ipulse];
4107 overlap->comm_data[ipulse].send_index0 -
4108 overlap->comm_data[0].send_index0;
4109 send_nindex = overlap->comm_data[ipulse].send_nindex;
4110 /* We don't use recv_index0, as we always receive starting at 0 */
4111 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
4112 recv_size_y = overlap->comm_data[ipulse].recv_size;
4114 sendptr = overlap->sendbuf + send_index0*local_fft_ndata[ZZ];
4115 recvptr = overlap->recvbuf;
4118 MPI_Sendrecv(sendptr, send_size_y*datasize, GMX_MPI_REAL,
4120 recvptr, recv_size_y*datasize, GMX_MPI_REAL,
4122 overlap->mpi_comm, &stat);
4125 for (x = 0; x < local_fft_ndata[XX]; x++)
4127 for (y = 0; y < recv_nindex; y++)
4129 indg = (x*local_fft_size[YY] + y)*local_fft_size[ZZ];
4130 indb = (x*recv_size_y + y)*local_fft_ndata[ZZ];
4131 for (z = 0; z < local_fft_ndata[ZZ]; z++)
4133 fftgrid[indg+z] += recvptr[indb+z];
4138 if (pme->nnodes_major > 1)
4140 /* Copy from the received buffer to the send buffer for dim 0 */
4141 sendptr = pme->overlap[0].sendbuf;
4142 for (x = 0; x < size_yx; x++)
4144 for (y = 0; y < recv_nindex; y++)
4146 indg = (x*local_fft_ndata[YY] + y)*local_fft_ndata[ZZ];
4147 indb = ((local_fft_ndata[XX] + x)*recv_size_y + y)*local_fft_ndata[ZZ];
4148 for (z = 0; z < local_fft_ndata[ZZ]; z++)
4150 sendptr[indg+z] += recvptr[indb+z];
4158 /* We only support a single pulse here.
4159 * This is not a severe limitation, as this code is only used
4160 * with OpenMP and with OpenMP the (PME) domains can be larger.
4162 if (pme->nnodes_major > 1)
4164 /* Major dimension */
4165 overlap = &pme->overlap[0];
4167 datasize = local_fft_ndata[YY]*local_fft_ndata[ZZ];
4168 gridsize = local_fft_size[YY] *local_fft_size[ZZ];
4172 send_id = overlap->send_id[ipulse];
4173 recv_id = overlap->recv_id[ipulse];
4174 send_nindex = overlap->comm_data[ipulse].send_nindex;
4175 /* We don't use recv_index0, as we always receive starting at 0 */
4176 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
4178 sendptr = overlap->sendbuf;
4179 recvptr = overlap->recvbuf;
4183 fprintf(debug, "PME fftgrid comm %2d x %2d x %2d\n",
4184 send_nindex, local_fft_ndata[YY], local_fft_ndata[ZZ]);
4188 MPI_Sendrecv(sendptr, send_nindex*datasize, GMX_MPI_REAL,
4190 recvptr, recv_nindex*datasize, GMX_MPI_REAL,
4192 overlap->mpi_comm, &stat);
4195 for (x = 0; x < recv_nindex; x++)
4197 for (y = 0; y < local_fft_ndata[YY]; y++)
4199 indg = (x*local_fft_size[YY] + y)*local_fft_size[ZZ];
4200 indb = (x*local_fft_ndata[YY] + y)*local_fft_ndata[ZZ];
4201 for (z = 0; z < local_fft_ndata[ZZ]; z++)
4203 fftgrid[indg+z] += recvptr[indb+z];
4211 static void spread_on_grid(gmx_pme_t pme,
4212 pme_atomcomm_t *atc, pmegrids_t *grids,
4213 gmx_bool bCalcSplines, gmx_bool bSpread,
4214 real *fftgrid, gmx_bool bDoSplines )
4216 int nthread, thread;
4217 #ifdef PME_TIME_THREADS
4218 gmx_cycles_t c1, c2, c3, ct1a, ct1b, ct1c;
4219 static double cs1 = 0, cs2 = 0, cs3 = 0;
4220 static double cs1a[6] = {0, 0, 0, 0, 0, 0};
4224 nthread = pme->nthread;
4225 assert(nthread > 0);
4227 #ifdef PME_TIME_THREADS
4228 c1 = omp_cyc_start();
4232 #pragma omp parallel for num_threads(nthread) schedule(static)
4233 for (thread = 0; thread < nthread; thread++)
4237 start = atc->n* thread /nthread;
4238 end = atc->n*(thread+1)/nthread;
4240 /* Compute fftgrid index for all atoms,
4241 * with help of some extra variables.
4243 calc_interpolation_idx(pme, atc, start, end, thread);
4246 #ifdef PME_TIME_THREADS
4247 c1 = omp_cyc_end(c1);
4251 #ifdef PME_TIME_THREADS
4252 c2 = omp_cyc_start();
4254 #pragma omp parallel for num_threads(nthread) schedule(static)
4255 for (thread = 0; thread < nthread; thread++)
4257 splinedata_t *spline;
4258 pmegrid_t *grid = NULL;
4260 /* make local bsplines */
4261 if (grids == NULL || !pme->bUseThreads)
4263 spline = &atc->spline[0];
4269 grid = &grids->grid;
4274 spline = &atc->spline[thread];
4276 if (grids->nthread == 1)
4278 /* One thread, we operate on all charges */
4283 /* Get the indices our thread should operate on */
4284 make_thread_local_ind(atc, thread, spline);
4287 grid = &grids->grid_th[thread];
4292 make_bsplines(spline->theta, spline->dtheta, pme->pme_order,
4293 atc->fractx, spline->n, spline->ind, atc->q, bDoSplines);
4298 /* put local atoms on grid. */
4299 #ifdef PME_TIME_SPREAD
4300 ct1a = omp_cyc_start();
4302 spread_q_bsplines_thread(grid, atc, spline, pme->spline_work);
4304 if (pme->bUseThreads)
4306 copy_local_grid(pme, grids, thread, fftgrid);
4308 #ifdef PME_TIME_SPREAD
4309 ct1a = omp_cyc_end(ct1a);
4310 cs1a[thread] += (double)ct1a;
4314 #ifdef PME_TIME_THREADS
4315 c2 = omp_cyc_end(c2);
4319 if (bSpread && pme->bUseThreads)
4321 #ifdef PME_TIME_THREADS
4322 c3 = omp_cyc_start();
4324 #pragma omp parallel for num_threads(grids->nthread) schedule(static)
4325 for (thread = 0; thread < grids->nthread; thread++)
4327 reduce_threadgrid_overlap(pme, grids, thread,
4329 pme->overlap[0].sendbuf,
4330 pme->overlap[1].sendbuf);
4332 #ifdef PME_TIME_THREADS
4333 c3 = omp_cyc_end(c3);
4337 if (pme->nnodes > 1)
4339 /* Communicate the overlapping part of the fftgrid.
4340 * For this communication call we need to check pme->bUseThreads
4341 * to have all ranks communicate here, regardless of pme->nthread.
4343 sum_fftgrid_dd(pme, fftgrid);
4347 #ifdef PME_TIME_THREADS
4351 printf("idx %.2f spread %.2f red %.2f",
4352 cs1*1e-9, cs2*1e-9, cs3*1e-9);
4353 #ifdef PME_TIME_SPREAD
4354 for (thread = 0; thread < nthread; thread++)
4356 printf(" %.2f", cs1a[thread]*1e-9);
4365 static void dump_grid(FILE *fp,
4366 int sx, int sy, int sz, int nx, int ny, int nz,
4367 int my, int mz, const real *g)
4371 for (x = 0; x < nx; x++)
4373 for (y = 0; y < ny; y++)
4375 for (z = 0; z < nz; z++)
4377 fprintf(fp, "%2d %2d %2d %6.3f\n",
4378 sx+x, sy+y, sz+z, g[(x*my + y)*mz + z]);
4384 static void dump_local_fftgrid(gmx_pme_t pme, const real *fftgrid)
4386 ivec local_fft_ndata, local_fft_offset, local_fft_size;
4388 gmx_parallel_3dfft_real_limits(pme->pfft_setup[PME_GRID_QA],
4394 pme->pmegrid_start_ix,
4395 pme->pmegrid_start_iy,
4396 pme->pmegrid_start_iz,
4397 pme->pmegrid_nx-pme->pme_order+1,
4398 pme->pmegrid_ny-pme->pme_order+1,
4399 pme->pmegrid_nz-pme->pme_order+1,
4406 void gmx_pme_calc_energy(gmx_pme_t pme, int n, rvec *x, real *q, real *V)
4408 pme_atomcomm_t *atc;
4411 if (pme->nnodes > 1)
4413 gmx_incons("gmx_pme_calc_energy called in parallel");
4415 if (pme->bFEP_q > 1)
4417 gmx_incons("gmx_pme_calc_energy with free energy");
4420 atc = &pme->atc_energy;
4422 if (atc->spline == NULL)
4424 snew(atc->spline, atc->nthread);
4427 atc->bSpread = TRUE;
4428 atc->pme_order = pme->pme_order;
4430 pme_realloc_atomcomm_things(atc);
4434 /* We only use the A-charges grid */
4435 grid = &pme->pmegrid[PME_GRID_QA];
4437 /* Only calculate the spline coefficients, don't actually spread */
4438 spread_on_grid(pme, atc, NULL, TRUE, FALSE, pme->fftgrid[PME_GRID_QA], FALSE);
4440 *V = gather_energy_bsplines(pme, grid->grid.grid, atc);
4444 static void reset_pmeonly_counters(gmx_wallcycle_t wcycle,
4445 gmx_walltime_accounting_t walltime_accounting,
4446 t_nrnb *nrnb, t_inputrec *ir,
4449 /* Reset all the counters related to performance over the run */
4450 wallcycle_stop(wcycle, ewcRUN);
4451 wallcycle_reset_all(wcycle);
4453 if (ir->nsteps >= 0)
4455 /* ir->nsteps is not used here, but we update it for consistency */
4456 ir->nsteps -= step - ir->init_step;
4458 ir->init_step = step;
4459 wallcycle_start(wcycle, ewcRUN);
4460 walltime_accounting_start(walltime_accounting);
4464 static void gmx_pmeonly_switch(int *npmedata, gmx_pme_t **pmedata,
4466 t_commrec *cr, t_inputrec *ir,
4470 gmx_pme_t pme = NULL;
4473 while (ind < *npmedata)
4475 pme = (*pmedata)[ind];
4476 if (pme->nkx == grid_size[XX] &&
4477 pme->nky == grid_size[YY] &&
4478 pme->nkz == grid_size[ZZ])
4489 srenew(*pmedata, *npmedata);
4491 /* Generate a new PME data structure, copying part of the old pointers */
4492 gmx_pme_reinit(&((*pmedata)[ind]), cr, pme, ir, grid_size);
4494 *pme_ret = (*pmedata)[ind];
4497 int gmx_pmeonly(gmx_pme_t pme,
4498 t_commrec *cr, t_nrnb *mynrnb,
4499 gmx_wallcycle_t wcycle,
4500 gmx_walltime_accounting_t walltime_accounting,
4501 real ewaldcoeff_q, real ewaldcoeff_lj,
4506 gmx_pme_pp_t pme_pp;
4510 rvec *x_pp = NULL, *f_pp = NULL;
4511 real *chargeA = NULL, *chargeB = NULL;
4512 real *c6A = NULL, *c6B = NULL;
4513 real *sigmaA = NULL, *sigmaB = NULL;
4516 int maxshift_x = 0, maxshift_y = 0;
4517 real energy_q, energy_lj, dvdlambda_q, dvdlambda_lj;
4518 matrix vir_q, vir_lj;
4523 gmx_int64_t step, step_rel;
4526 /* This data will only use with PME tuning, i.e. switching PME grids */
4528 snew(pmedata, npmedata);
4531 pme_pp = gmx_pme_pp_init(cr);
4536 do /****** this is a quasi-loop over time steps! */
4538 /* The reason for having a loop here is PME grid tuning/switching */
4541 /* Domain decomposition */
4542 ret = gmx_pme_recv_params_coords(pme_pp,
4548 &maxshift_x, &maxshift_y,
4549 &pme->bFEP_q, &pme->bFEP_lj,
4550 &lambda_q, &lambda_lj,
4554 grid_switch, &ewaldcoeff_q, &ewaldcoeff_lj);
4556 if (ret == pmerecvqxSWITCHGRID)
4558 /* Switch the PME grid to grid_switch */
4559 gmx_pmeonly_switch(&npmedata, &pmedata, grid_switch, cr, ir, &pme);
4562 if (ret == pmerecvqxRESETCOUNTERS)
4564 /* Reset the cycle and flop counters */
4565 reset_pmeonly_counters(wcycle, walltime_accounting, mynrnb, ir, step);
4568 while (ret == pmerecvqxSWITCHGRID || ret == pmerecvqxRESETCOUNTERS);
4570 if (ret == pmerecvqxFINISH)
4572 /* We should stop: break out of the loop */
4576 step_rel = step - ir->init_step;
4580 wallcycle_start(wcycle, ewcRUN);
4581 walltime_accounting_start(walltime_accounting);
4584 wallcycle_start(wcycle, ewcPMEMESH);
4591 gmx_pme_do(pme, 0, natoms, x_pp, f_pp,
4592 chargeA, chargeB, c6A, c6B, sigmaA, sigmaB, box,
4593 cr, maxshift_x, maxshift_y, mynrnb, wcycle,
4594 vir_q, ewaldcoeff_q, vir_lj, ewaldcoeff_lj,
4595 &energy_q, &energy_lj, lambda_q, lambda_lj, &dvdlambda_q, &dvdlambda_lj,
4596 pme_flags | GMX_PME_DO_ALL_F | (bEnerVir ? GMX_PME_CALC_ENER_VIR : 0));
4598 cycles = wallcycle_stop(wcycle, ewcPMEMESH);
4600 gmx_pme_send_force_vir_ener(pme_pp,
4601 f_pp, vir_q, energy_q, vir_lj, energy_lj,
4602 dvdlambda_q, dvdlambda_lj, cycles);
4605 } /***** end of quasi-loop, we stop with the break above */
4608 walltime_accounting_end(walltime_accounting);
4614 calc_initial_lb_coeffs(gmx_pme_t pme, real *local_c6, real *local_sigma)
4618 for (i = 0; i < pme->atc[0].n; ++i)
4622 sigma4 = local_sigma[i];
4623 sigma4 = sigma4*sigma4;
4624 sigma4 = sigma4*sigma4;
4625 pme->atc[0].q[i] = local_c6[i] / sigma4;
4630 calc_next_lb_coeffs(gmx_pme_t pme, real *local_sigma)
4634 for (i = 0; i < pme->atc[0].n; ++i)
4636 pme->atc[0].q[i] *= local_sigma[i];
4641 do_redist_x_q(gmx_pme_t pme, t_commrec *cr, int start, int homenr,
4642 gmx_bool bFirst, rvec x[], real *data)
4645 pme_atomcomm_t *atc;
4648 for (d = pme->ndecompdim - 1; d >= 0; d--)
4654 if (d == pme->ndecompdim - 1)
4662 n_d = pme->atc[d + 1].n;
4668 if (atc->npd > atc->pd_nalloc)
4670 atc->pd_nalloc = over_alloc_dd(atc->npd);
4671 srenew(atc->pd, atc->pd_nalloc);
4673 pme_calc_pidx_wrapper(n_d, pme->recipbox, x_d, atc);
4675 /* Redistribute x (only once) and qA/c6A or qB/c6B */
4676 if (DOMAINDECOMP(cr))
4678 dd_pmeredist_x_q(pme, n_d, bFirst, x_d, q_d, atc);
4683 /* TODO: when adding free-energy support, sigmaB will no longer be
4685 int gmx_pme_do(gmx_pme_t pme,
4686 int start, int homenr,
4688 real *chargeA, real *chargeB,
4689 real *c6A, real *c6B,
4690 real *sigmaA, real gmx_unused *sigmaB,
4691 matrix box, t_commrec *cr,
4692 int maxshift_x, int maxshift_y,
4693 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
4694 matrix vir_q, real ewaldcoeff_q,
4695 matrix vir_lj, real ewaldcoeff_lj,
4696 real *energy_q, real *energy_lj,
4697 real lambda_q, real lambda_lj,
4698 real *dvdlambda_q, real *dvdlambda_lj,
4701 int q, d, i, j, ntot, npme, qmax;
4704 pme_atomcomm_t *atc = NULL;
4705 pmegrids_t *pmegrid = NULL;
4709 real *charge = NULL, *q_d;
4714 gmx_parallel_3dfft_t pfft_setup;
4716 t_complex * cfftgrid;
4718 gmx_bool bFirst, bDoSplines;
4719 const gmx_bool bCalcEnerVir = flags & GMX_PME_CALC_ENER_VIR;
4720 const gmx_bool bCalcF = flags & GMX_PME_CALC_F;
4722 assert(pme->nnodes > 0);
4723 assert(pme->nnodes == 1 || pme->ndecompdim > 0);
4725 if (pme->nnodes > 1)
4729 if (atc->npd > atc->pd_nalloc)
4731 atc->pd_nalloc = over_alloc_dd(atc->npd);
4732 srenew(atc->pd, atc->pd_nalloc);
4734 for (d = pme->ndecompdim-1; d >= 0; d--)
4737 atc->maxshift = (atc->dimind == 0 ? maxshift_x : maxshift_y);
4743 /* This could be necessary for TPI */
4744 pme->atc[0].n = homenr;
4745 if (DOMAINDECOMP(cr))
4747 pme_realloc_atomcomm_things(atc);
4753 m_inv_ur0(box, pme->recipbox);
4756 /* For simplicity, we construct the splines for all particles if
4757 * more than one PME calculations is needed. Some optimization
4758 * could be done by keeping track of which atoms have splines
4759 * constructed, and construct new splines on each pass for atoms
4760 * that don't yet have them.
4763 bDoSplines = pme->bFEP || ((flags & GMX_PME_DO_COULOMB) && (flags & GMX_PME_DO_LJ));
4765 /* We need a maximum of four separate PME calculations:
4766 * q=0: Coulomb PME with charges from state A
4767 * q=1: Coulomb PME with charges from state B
4768 * q=2: LJ PME with C6 from state A
4769 * q=3: LJ PME with C6 from state B
4770 * For Lorentz-Berthelot combination rules, a separate loop is used to
4771 * calculate all the terms
4774 /* If we are doing LJ-PME with LB, we only do Q here */
4775 qmax = (pme->ljpme_combination_rule == eljpmeLB) ? DO_Q : DO_Q_AND_LJ;
4777 for (q = 0; q < qmax; ++q)
4779 /* Check if we should do calculations at this q
4780 * If q is odd we should be doing FEP
4781 * If q < 2 we should be doing electrostatic PME
4782 * If q >= 2 we should be doing LJ-PME
4784 if ((q < DO_Q && (!(flags & GMX_PME_DO_COULOMB) ||
4785 (q == 1 && !pme->bFEP_q))) ||
4786 (q >= DO_Q && (!(flags & GMX_PME_DO_LJ) ||
4787 (q == 3 && !pme->bFEP_lj))))
4791 /* Unpack structure */
4792 pmegrid = &pme->pmegrid[q];
4793 fftgrid = pme->fftgrid[q];
4794 cfftgrid = pme->cfftgrid[q];
4795 pfft_setup = pme->pfft_setup[q];
4798 case 0: charge = chargeA + start; break;
4799 case 1: charge = chargeB + start; break;
4800 case 2: charge = c6A + start; break;
4801 case 3: charge = c6B + start; break;
4804 grid = pmegrid->grid.grid;
4808 fprintf(debug, "PME: nnodes = %d, nodeid = %d\n",
4809 cr->nnodes, cr->nodeid);
4810 fprintf(debug, "Grid = %p\n", (void*)grid);
4813 gmx_fatal(FARGS, "No grid!");
4818 if (pme->nnodes == 1)
4824 wallcycle_start(wcycle, ewcPME_REDISTXF);
4825 do_redist_x_q(pme, cr, start, homenr, bFirst, x, charge);
4828 wallcycle_stop(wcycle, ewcPME_REDISTXF);
4833 fprintf(debug, "Node= %6d, pme local particles=%6d\n",
4834 cr->nodeid, atc->n);
4837 if (flags & GMX_PME_SPREAD_Q)
4839 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
4841 /* Spread the charges on a grid */
4842 spread_on_grid(pme, &pme->atc[0], pmegrid, bFirst, TRUE, fftgrid, bDoSplines);
4846 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
4848 inc_nrnb(nrnb, eNR_SPREADQBSP,
4849 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
4851 if (!pme->bUseThreads)
4853 wrap_periodic_pmegrid(pme, grid);
4855 /* sum contributions to local grid from other nodes */
4857 if (pme->nnodes > 1)
4859 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_QGRID_FORWARD);
4864 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, q);
4867 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
4870 dump_local_fftgrid(pme,fftgrid);
4875 /* Here we start a large thread parallel region */
4876 #pragma omp parallel num_threads(pme->nthread) private(thread)
4878 thread = gmx_omp_get_thread_num();
4879 if (flags & GMX_PME_SOLVE)
4886 wallcycle_start(wcycle, ewcPME_FFT);
4888 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
4892 wallcycle_stop(wcycle, ewcPME_FFT);
4896 /* solve in k-space for our local cells */
4899 wallcycle_start(wcycle, (q < DO_Q ? ewcPME_SOLVE : ewcLJPME));
4904 solve_pme_yzx(pme, cfftgrid, ewaldcoeff_q,
4905 box[XX][XX]*box[YY][YY]*box[ZZ][ZZ],
4907 pme->nthread, thread);
4912 solve_pme_lj_yzx(pme, &cfftgrid, FALSE, ewaldcoeff_lj,
4913 box[XX][XX]*box[YY][YY]*box[ZZ][ZZ],
4915 pme->nthread, thread);
4920 wallcycle_stop(wcycle, (q < DO_Q ? ewcPME_SOLVE : ewcLJPME));
4922 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
4932 wallcycle_start(wcycle, ewcPME_FFT);
4934 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
4938 wallcycle_stop(wcycle, ewcPME_FFT);
4942 if (pme->nodeid == 0)
4944 ntot = pme->nkx*pme->nky*pme->nkz;
4945 npme = ntot*log((real)ntot)/log(2.0);
4946 inc_nrnb(nrnb, eNR_FFT, 2*npme);
4949 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
4952 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, q, pme->nthread, thread);
4955 /* End of thread parallel section.
4956 * With MPI we have to synchronize here before gmx_sum_qgrid_dd.
4961 /* distribute local grid to all nodes */
4963 if (pme->nnodes > 1)
4965 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_QGRID_BACKWARD);
4970 unwrap_periodic_pmegrid(pme, grid);
4972 /* interpolate forces for our local atoms */
4976 /* If we are running without parallelization,
4977 * atc->f is the actual force array, not a buffer,
4978 * therefore we should not clear it.
4980 lambda = q < DO_Q ? lambda_q : lambda_lj;
4981 bClearF = (bFirst && PAR(cr));
4982 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
4983 for (thread = 0; thread < pme->nthread; thread++)
4985 gather_f_bsplines(pme, grid, bClearF, atc,
4986 &atc->spline[thread],
4987 pme->bFEP ? (q % 2 == 0 ? 1.0-lambda : lambda) : 1.0);
4992 inc_nrnb(nrnb, eNR_GATHERFBSP,
4993 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
4994 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
4999 /* This should only be called on the master thread
5000 * and after the threads have synchronized.
5004 get_pme_ener_vir_q(pme, pme->nthread, &energy_AB[q], vir_AB[q]);
5008 get_pme_ener_vir_lj(pme, pme->nthread, &energy_AB[q], vir_AB[q]);
5014 /* For Lorentz-Berthelot combination rules in LJ-PME, we need to calculate
5017 if ((flags & GMX_PME_DO_LJ) && pme->ljpme_combination_rule == eljpmeLB)
5019 real *local_c6 = NULL, *local_sigma = NULL;
5020 if (pme->nnodes == 1)
5022 if (pme->lb_buf1 == NULL)
5024 pme->lb_buf_nalloc = pme->atc[0].n;
5025 snew(pme->lb_buf1, pme->lb_buf_nalloc);
5027 pme->atc[0].q = pme->lb_buf1;
5029 local_sigma = sigmaA;
5035 wallcycle_start(wcycle, ewcPME_REDISTXF);
5037 do_redist_x_q(pme, cr, start, homenr, bFirst, x, c6A);
5038 if (pme->lb_buf_nalloc < atc->n)
5040 pme->lb_buf_nalloc = atc->nalloc;
5041 srenew(pme->lb_buf1, pme->lb_buf_nalloc);
5042 srenew(pme->lb_buf2, pme->lb_buf_nalloc);
5044 local_c6 = pme->lb_buf1;
5045 for (i = 0; i < atc->n; ++i)
5047 local_c6[i] = atc->q[i];
5051 do_redist_x_q(pme, cr, start, homenr, FALSE, x, sigmaA);
5052 local_sigma = pme->lb_buf2;
5053 for (i = 0; i < atc->n; ++i)
5055 local_sigma[i] = atc->q[i];
5059 wallcycle_stop(wcycle, ewcPME_REDISTXF);
5061 calc_initial_lb_coeffs(pme, local_c6, local_sigma);
5063 /*Seven terms in LJ-PME with LB, q < 2 reserved for electrostatics*/
5064 for (q = 2; q < 9; ++q)
5066 /* Unpack structure */
5067 pmegrid = &pme->pmegrid[q];
5068 fftgrid = pme->fftgrid[q];
5069 cfftgrid = pme->cfftgrid[q];
5070 pfft_setup = pme->pfft_setup[q];
5071 calc_next_lb_coeffs(pme, local_sigma);
5072 grid = pmegrid->grid.grid;
5075 if (flags & GMX_PME_SPREAD_Q)
5077 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
5078 /* Spread the charges on a grid */
5079 spread_on_grid(pme, &pme->atc[0], pmegrid, bFirst, TRUE, fftgrid, bDoSplines);
5083 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
5085 inc_nrnb(nrnb, eNR_SPREADQBSP,
5086 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
5087 if (pme->nthread == 1)
5089 wrap_periodic_pmegrid(pme, grid);
5091 /* sum contributions to local grid from other nodes */
5093 if (pme->nnodes > 1)
5095 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_QGRID_FORWARD);
5099 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, q);
5102 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
5105 /*Here we start a large thread parallel region*/
5106 #pragma omp parallel num_threads(pme->nthread) private(thread)
5108 thread = gmx_omp_get_thread_num();
5109 if (flags & GMX_PME_SOLVE)
5114 wallcycle_start(wcycle, ewcPME_FFT);
5117 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
5121 wallcycle_stop(wcycle, ewcPME_FFT);
5128 if (flags & GMX_PME_SOLVE)
5131 /* solve in k-space for our local cells */
5132 #pragma omp parallel num_threads(pme->nthread) private(thread)
5134 thread = gmx_omp_get_thread_num();
5137 wallcycle_start(wcycle, ewcLJPME);
5141 solve_pme_lj_yzx(pme, &pme->cfftgrid[2], TRUE, ewaldcoeff_lj,
5142 box[XX][XX]*box[YY][YY]*box[ZZ][ZZ],
5144 pme->nthread, thread);
5147 wallcycle_stop(wcycle, ewcLJPME);
5149 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
5156 /* This should only be called on the master thread and
5157 * after the threads have synchronized.
5159 get_pme_ener_vir_lj(pme, pme->nthread, &energy_AB[2], vir_AB[2]);
5164 bFirst = !(flags & GMX_PME_DO_COULOMB);
5165 calc_initial_lb_coeffs(pme, local_c6, local_sigma);
5166 for (q = 8; q >= 2; --q)
5168 /* Unpack structure */
5169 pmegrid = &pme->pmegrid[q];
5170 fftgrid = pme->fftgrid[q];
5171 cfftgrid = pme->cfftgrid[q];
5172 pfft_setup = pme->pfft_setup[q];
5173 grid = pmegrid->grid.grid;
5174 calc_next_lb_coeffs(pme, local_sigma);
5176 #pragma omp parallel num_threads(pme->nthread) private(thread)
5178 thread = gmx_omp_get_thread_num();
5183 wallcycle_start(wcycle, ewcPME_FFT);
5186 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
5190 wallcycle_stop(wcycle, ewcPME_FFT);
5194 if (pme->nodeid == 0)
5196 ntot = pme->nkx*pme->nky*pme->nkz;
5197 npme = ntot*log((real)ntot)/log(2.0);
5198 inc_nrnb(nrnb, eNR_FFT, 2*npme);
5200 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
5203 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, q, pme->nthread, thread);
5205 } /*#pragma omp parallel*/
5207 /* distribute local grid to all nodes */
5209 if (pme->nnodes > 1)
5211 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_QGRID_BACKWARD);
5216 unwrap_periodic_pmegrid(pme, grid);
5218 /* interpolate forces for our local atoms */
5220 bClearF = (bFirst && PAR(cr));
5221 scale = pme->bFEP ? (q < 9 ? 1.0-lambda_lj : lambda_lj) : 1.0;
5222 scale *= lb_scale_factor[q-2];
5223 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
5224 for (thread = 0; thread < pme->nthread; thread++)
5226 gather_f_bsplines(pme, grid, bClearF, &pme->atc[0],
5227 &pme->atc[0].spline[thread],
5232 inc_nrnb(nrnb, eNR_GATHERFBSP,
5233 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
5234 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
5237 } /*for (q = 8; q >= 2; --q)*/
5239 } /*if ((flags & GMX_PME_DO_LJ) && pme->ljpme_combination_rule == eljpmeLB) */
5241 if (bCalcF && pme->nnodes > 1)
5243 wallcycle_start(wcycle, ewcPME_REDISTXF);
5244 for (d = 0; d < pme->ndecompdim; d++)
5247 if (d == pme->ndecompdim - 1)
5254 n_d = pme->atc[d+1].n;
5255 f_d = pme->atc[d+1].f;
5257 if (DOMAINDECOMP(cr))
5259 dd_pmeredist_f(pme, atc, n_d, f_d,
5260 d == pme->ndecompdim-1 && pme->bPPnode);
5264 wallcycle_stop(wcycle, ewcPME_REDISTXF);
5270 if (flags & GMX_PME_DO_COULOMB)
5274 *energy_q = energy_AB[0];
5275 m_add(vir_q, vir_AB[0], vir_q);
5279 *energy_q = (1.0-lambda_q)*energy_AB[0] + lambda_q*energy_AB[1];
5280 *dvdlambda_q += energy_AB[1] - energy_AB[0];
5281 for (i = 0; i < DIM; i++)
5283 for (j = 0; j < DIM; j++)
5285 vir_q[i][j] += (1.0-lambda_q)*vir_AB[0][i][j] +
5286 lambda_q*vir_AB[1][i][j];
5292 fprintf(debug, "Electrostatic PME mesh energy: %g\n", *energy_q);
5300 if (flags & GMX_PME_DO_LJ)
5304 *energy_lj = energy_AB[2];
5305 m_add(vir_lj, vir_AB[2], vir_lj);
5309 *energy_lj = (1.0-lambda_lj)*energy_AB[2] + lambda_lj*energy_AB[3];
5310 *dvdlambda_lj += energy_AB[3] - energy_AB[2];
5311 for (i = 0; i < DIM; i++)
5313 for (j = 0; j < DIM; j++)
5315 vir_lj[i][j] += (1.0-lambda_lj)*vir_AB[2][i][j] + lambda_lj*vir_AB[3][i][j];
5321 fprintf(debug, "Lennard-Jones PME mesh energy: %g\n", *energy_lj);