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.
71 #include "gromacs/utility/smalloc.h"
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 coefficient */
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 */
307 gmx_bool bFEP; /* Compute Free energy contribution */
310 int nkx, nky, nkz; /* Grid dimensions */
311 gmx_bool bP3M; /* Do P3M: optimize the influence function */
315 int ljpme_combination_rule; /* Type of combination rule in LJ-PME */
317 int ngrids; /* number of grids we maintain for pmegrid, (c)fftgrid and pfft_setups*/
319 pmegrids_t pmegrid[DO_Q_AND_LJ_LB]; /* Grids on which we do spreading/interpolation,
320 * includes overlap Grid indices are ordered as
322 * 0: Coloumb PME, state A
323 * 1: Coloumb PME, state B
325 * This can probably be done in a better way
326 * but this simple hack works for now
328 /* The PME coefficient spreading grid sizes/strides, includes pme_order-1 */
329 int pmegrid_nx, pmegrid_ny, pmegrid_nz;
330 /* pmegrid_nz might be larger than strictly necessary to ensure
331 * memory alignment, pmegrid_nz_base gives the real base size.
334 /* The local PME grid starting indices */
335 int pmegrid_start_ix, pmegrid_start_iy, pmegrid_start_iz;
337 /* Work data for spreading and gathering */
338 pme_spline_work_t *spline_work;
340 real **fftgrid; /* Grids for FFT. With 1D FFT decomposition this can be a pointer */
341 /* inside the interpolation grid, but separate for 2D PME decomp. */
342 int fftgrid_nx, fftgrid_ny, fftgrid_nz;
344 t_complex **cfftgrid; /* Grids for complex FFT data */
346 int cfftgrid_nx, cfftgrid_ny, cfftgrid_nz;
348 gmx_parallel_3dfft_t *pfft_setup;
350 int *nnx, *nny, *nnz;
351 real *fshx, *fshy, *fshz;
353 pme_atomcomm_t atc[2]; /* Indexed on decomposition index */
356 /* Buffers to store data for local atoms for L-B combination rule
357 * calculations in LJ-PME. lb_buf1 stores either the coefficients
358 * for spreading/gathering (in serial), or the C6 coefficient for
359 * local atoms (in parallel). lb_buf2 is only used in parallel,
360 * and stores the sigma values for local atoms. */
361 real *lb_buf1, *lb_buf2;
362 int lb_buf_nalloc; /* Allocation size for the above buffers. */
364 pme_overlap_t overlap[2]; /* Indexed on dimension, 0=x, 1=y */
366 pme_atomcomm_t atc_energy; /* Only for gmx_pme_calc_energy */
368 rvec *bufv; /* Communication buffer */
369 real *bufr; /* Communication buffer */
370 int buf_nalloc; /* The communication buffer size */
372 /* thread local work data for solve_pme */
375 /* Work data for sum_qgrid */
376 real * sum_qgrid_tmp;
377 real * sum_qgrid_dd_tmp;
380 static void calc_interpolation_idx(gmx_pme_t pme, pme_atomcomm_t *atc,
381 int start, int grid_index, int end, int thread)
384 int *idxptr, tix, tiy, tiz;
385 real *xptr, *fptr, tx, ty, tz;
386 real rxx, ryx, ryy, rzx, rzy, rzz;
388 int start_ix, start_iy, start_iz;
389 int *g2tx, *g2ty, *g2tz;
391 int *thread_idx = NULL;
392 thread_plist_t *tpl = NULL;
400 start_ix = pme->pmegrid_start_ix;
401 start_iy = pme->pmegrid_start_iy;
402 start_iz = pme->pmegrid_start_iz;
404 rxx = pme->recipbox[XX][XX];
405 ryx = pme->recipbox[YY][XX];
406 ryy = pme->recipbox[YY][YY];
407 rzx = pme->recipbox[ZZ][XX];
408 rzy = pme->recipbox[ZZ][YY];
409 rzz = pme->recipbox[ZZ][ZZ];
411 g2tx = pme->pmegrid[grid_index].g2t[XX];
412 g2ty = pme->pmegrid[grid_index].g2t[YY];
413 g2tz = pme->pmegrid[grid_index].g2t[ZZ];
415 bThreads = (atc->nthread > 1);
418 thread_idx = atc->thread_idx;
420 tpl = &atc->thread_plist[thread];
422 for (i = 0; i < atc->nthread; i++)
428 for (i = start; i < end; i++)
431 idxptr = atc->idx[i];
432 fptr = atc->fractx[i];
434 /* Fractional coordinates along box vectors, add 2.0 to make 100% sure we are positive for triclinic boxes */
435 tx = nx * ( xptr[XX] * rxx + xptr[YY] * ryx + xptr[ZZ] * rzx + 2.0 );
436 ty = ny * ( xptr[YY] * ryy + xptr[ZZ] * rzy + 2.0 );
437 tz = nz * ( xptr[ZZ] * rzz + 2.0 );
443 /* Because decomposition only occurs in x and y,
444 * we never have a fraction correction in z.
446 fptr[XX] = tx - tix + pme->fshx[tix];
447 fptr[YY] = ty - tiy + pme->fshy[tiy];
450 idxptr[XX] = pme->nnx[tix];
451 idxptr[YY] = pme->nny[tiy];
452 idxptr[ZZ] = pme->nnz[tiz];
455 range_check(idxptr[XX], 0, pme->pmegrid_nx);
456 range_check(idxptr[YY], 0, pme->pmegrid_ny);
457 range_check(idxptr[ZZ], 0, pme->pmegrid_nz);
462 thread_i = g2tx[idxptr[XX]] + g2ty[idxptr[YY]] + g2tz[idxptr[ZZ]];
463 thread_idx[i] = thread_i;
470 /* Make a list of particle indices sorted on thread */
472 /* Get the cumulative count */
473 for (i = 1; i < atc->nthread; i++)
475 tpl_n[i] += tpl_n[i-1];
477 /* The current implementation distributes particles equally
478 * over the threads, so we could actually allocate for that
479 * in pme_realloc_atomcomm_things.
481 if (tpl_n[atc->nthread-1] > tpl->nalloc)
483 tpl->nalloc = over_alloc_large(tpl_n[atc->nthread-1]);
484 srenew(tpl->i, tpl->nalloc);
486 /* Set tpl_n to the cumulative start */
487 for (i = atc->nthread-1; i >= 1; i--)
489 tpl_n[i] = tpl_n[i-1];
493 /* Fill our thread local array with indices sorted on thread */
494 for (i = start; i < end; i++)
496 tpl->i[tpl_n[atc->thread_idx[i]]++] = i;
498 /* Now tpl_n contains the cummulative count again */
502 static void make_thread_local_ind(pme_atomcomm_t *atc,
503 int thread, splinedata_t *spline)
505 int n, t, i, start, end;
508 /* Combine the indices made by each thread into one index */
512 for (t = 0; t < atc->nthread; t++)
514 tpl = &atc->thread_plist[t];
515 /* Copy our part (start - end) from the list of thread t */
518 start = tpl->n[thread-1];
520 end = tpl->n[thread];
521 for (i = start; i < end; i++)
523 spline->ind[n++] = tpl->i[i];
531 static void pme_calc_pidx(int start, int end,
532 matrix recipbox, rvec x[],
533 pme_atomcomm_t *atc, int *count)
538 real rxx, ryx, rzx, ryy, rzy;
541 /* Calculate PME task index (pidx) for each grid index.
542 * Here we always assign equally sized slabs to each node
543 * for load balancing reasons (the PME grid spacing is not used).
549 /* Reset the count */
550 for (i = 0; i < nslab; i++)
555 if (atc->dimind == 0)
557 rxx = recipbox[XX][XX];
558 ryx = recipbox[YY][XX];
559 rzx = recipbox[ZZ][XX];
560 /* Calculate the node index in x-dimension */
561 for (i = start; i < end; i++)
564 /* Fractional coordinates along box vectors */
565 s = nslab*(xptr[XX]*rxx + xptr[YY]*ryx + xptr[ZZ]*rzx);
566 si = (int)(s + 2*nslab) % nslab;
573 ryy = recipbox[YY][YY];
574 rzy = recipbox[ZZ][YY];
575 /* Calculate the node index in y-dimension */
576 for (i = start; i < end; i++)
579 /* Fractional coordinates along box vectors */
580 s = nslab*(xptr[YY]*ryy + xptr[ZZ]*rzy);
581 si = (int)(s + 2*nslab) % nslab;
588 static void pme_calc_pidx_wrapper(int natoms, matrix recipbox, rvec x[],
591 int nthread, thread, slab;
593 nthread = atc->nthread;
595 #pragma omp parallel for num_threads(nthread) schedule(static)
596 for (thread = 0; thread < nthread; thread++)
598 pme_calc_pidx(natoms* thread /nthread,
599 natoms*(thread+1)/nthread,
600 recipbox, x, atc, atc->count_thread[thread]);
602 /* Non-parallel reduction, since nslab is small */
604 for (thread = 1; thread < nthread; thread++)
606 for (slab = 0; slab < atc->nslab; slab++)
608 atc->count_thread[0][slab] += atc->count_thread[thread][slab];
613 static void realloc_splinevec(splinevec th, real **ptr_z, int nalloc)
615 const int padding = 4;
618 srenew(th[XX], nalloc);
619 srenew(th[YY], nalloc);
620 /* In z we add padding, this is only required for the aligned SIMD code */
621 sfree_aligned(*ptr_z);
622 snew_aligned(*ptr_z, nalloc+2*padding, SIMD4_ALIGNMENT);
623 th[ZZ] = *ptr_z + padding;
625 for (i = 0; i < padding; i++)
628 (*ptr_z)[padding+nalloc+i] = 0;
632 static void pme_realloc_splinedata(splinedata_t *spline, pme_atomcomm_t *atc)
636 srenew(spline->ind, atc->nalloc);
637 /* Initialize the index to identity so it works without threads */
638 for (i = 0; i < atc->nalloc; i++)
643 realloc_splinevec(spline->theta, &spline->ptr_theta_z,
644 atc->pme_order*atc->nalloc);
645 realloc_splinevec(spline->dtheta, &spline->ptr_dtheta_z,
646 atc->pme_order*atc->nalloc);
649 static void pme_realloc_atomcomm_things(pme_atomcomm_t *atc)
651 int nalloc_old, i, j, nalloc_tpl;
653 /* We have to avoid a NULL pointer for atc->x to avoid
654 * possible fatal errors in MPI routines.
656 if (atc->n > atc->nalloc || atc->nalloc == 0)
658 nalloc_old = atc->nalloc;
659 atc->nalloc = over_alloc_dd(max(atc->n, 1));
663 srenew(atc->x, atc->nalloc);
664 srenew(atc->coefficient, atc->nalloc);
665 srenew(atc->f, atc->nalloc);
666 for (i = nalloc_old; i < atc->nalloc; i++)
668 clear_rvec(atc->f[i]);
673 srenew(atc->fractx, atc->nalloc);
674 srenew(atc->idx, atc->nalloc);
676 if (atc->nthread > 1)
678 srenew(atc->thread_idx, atc->nalloc);
681 for (i = 0; i < atc->nthread; i++)
683 pme_realloc_splinedata(&atc->spline[i], atc);
689 static void pme_dd_sendrecv(pme_atomcomm_t gmx_unused *atc,
690 gmx_bool gmx_unused bBackward, int gmx_unused shift,
691 void gmx_unused *buf_s, int gmx_unused nbyte_s,
692 void gmx_unused *buf_r, int gmx_unused nbyte_r)
698 if (bBackward == FALSE)
700 dest = atc->node_dest[shift];
701 src = atc->node_src[shift];
705 dest = atc->node_src[shift];
706 src = atc->node_dest[shift];
709 if (nbyte_s > 0 && nbyte_r > 0)
711 MPI_Sendrecv(buf_s, nbyte_s, MPI_BYTE,
713 buf_r, nbyte_r, MPI_BYTE,
715 atc->mpi_comm, &stat);
717 else if (nbyte_s > 0)
719 MPI_Send(buf_s, nbyte_s, MPI_BYTE,
723 else if (nbyte_r > 0)
725 MPI_Recv(buf_r, nbyte_r, MPI_BYTE,
727 atc->mpi_comm, &stat);
732 static void dd_pmeredist_pos_coeffs(gmx_pme_t pme,
733 int n, gmx_bool bX, rvec *x, real *data,
736 int *commnode, *buf_index;
737 int nnodes_comm, i, nsend, local_pos, buf_pos, node, scount, rcount;
739 commnode = atc->node_dest;
740 buf_index = atc->buf_index;
742 nnodes_comm = min(2*atc->maxshift, atc->nslab-1);
745 for (i = 0; i < nnodes_comm; i++)
747 buf_index[commnode[i]] = nsend;
748 nsend += atc->count[commnode[i]];
752 if (atc->count[atc->nodeid] + nsend != n)
754 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"
755 "This usually means that your system is not well equilibrated.",
756 n - (atc->count[atc->nodeid] + nsend),
757 pme->nodeid, 'x'+atc->dimind);
760 if (nsend > pme->buf_nalloc)
762 pme->buf_nalloc = over_alloc_dd(nsend);
763 srenew(pme->bufv, pme->buf_nalloc);
764 srenew(pme->bufr, pme->buf_nalloc);
767 atc->n = atc->count[atc->nodeid];
768 for (i = 0; i < nnodes_comm; i++)
770 scount = atc->count[commnode[i]];
771 /* Communicate the count */
774 fprintf(debug, "dimind %d PME node %d send to node %d: %d\n",
775 atc->dimind, atc->nodeid, commnode[i], scount);
777 pme_dd_sendrecv(atc, FALSE, i,
778 &scount, sizeof(int),
779 &atc->rcount[i], sizeof(int));
780 atc->n += atc->rcount[i];
783 pme_realloc_atomcomm_things(atc);
787 for (i = 0; i < n; i++)
790 if (node == atc->nodeid)
792 /* Copy direct to the receive buffer */
795 copy_rvec(x[i], atc->x[local_pos]);
797 atc->coefficient[local_pos] = data[i];
802 /* Copy to the send buffer */
805 copy_rvec(x[i], pme->bufv[buf_index[node]]);
807 pme->bufr[buf_index[node]] = data[i];
813 for (i = 0; i < nnodes_comm; i++)
815 scount = atc->count[commnode[i]];
816 rcount = atc->rcount[i];
817 if (scount > 0 || rcount > 0)
821 /* Communicate the coordinates */
822 pme_dd_sendrecv(atc, FALSE, i,
823 pme->bufv[buf_pos], scount*sizeof(rvec),
824 atc->x[local_pos], rcount*sizeof(rvec));
826 /* Communicate the coefficients */
827 pme_dd_sendrecv(atc, FALSE, i,
828 pme->bufr+buf_pos, scount*sizeof(real),
829 atc->coefficient+local_pos, rcount*sizeof(real));
831 local_pos += atc->rcount[i];
836 static void dd_pmeredist_f(gmx_pme_t pme, pme_atomcomm_t *atc,
840 int *commnode, *buf_index;
841 int nnodes_comm, local_pos, buf_pos, i, scount, rcount, node;
843 commnode = atc->node_dest;
844 buf_index = atc->buf_index;
846 nnodes_comm = min(2*atc->maxshift, atc->nslab-1);
848 local_pos = atc->count[atc->nodeid];
850 for (i = 0; i < nnodes_comm; i++)
852 scount = atc->rcount[i];
853 rcount = atc->count[commnode[i]];
854 if (scount > 0 || rcount > 0)
856 /* Communicate the forces */
857 pme_dd_sendrecv(atc, TRUE, i,
858 atc->f[local_pos], scount*sizeof(rvec),
859 pme->bufv[buf_pos], rcount*sizeof(rvec));
862 buf_index[commnode[i]] = buf_pos;
869 for (i = 0; i < n; i++)
872 if (node == atc->nodeid)
874 /* Add from the local force array */
875 rvec_inc(f[i], atc->f[local_pos]);
880 /* Add from the receive buffer */
881 rvec_inc(f[i], pme->bufv[buf_index[node]]);
888 for (i = 0; i < n; i++)
891 if (node == atc->nodeid)
893 /* Copy from the local force array */
894 copy_rvec(atc->f[local_pos], f[i]);
899 /* Copy from the receive buffer */
900 copy_rvec(pme->bufv[buf_index[node]], f[i]);
908 static void gmx_sum_qgrid_dd(gmx_pme_t pme, real *grid, int direction)
910 pme_overlap_t *overlap;
911 int send_index0, send_nindex;
912 int recv_index0, recv_nindex;
914 int i, j, k, ix, iy, iz, icnt;
915 int ipulse, send_id, recv_id, datasize;
917 real *sendptr, *recvptr;
919 /* Start with minor-rank communication. This is a bit of a pain since it is not contiguous */
920 overlap = &pme->overlap[1];
922 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
924 /* Since we have already (un)wrapped the overlap in the z-dimension,
925 * we only have to communicate 0 to nkz (not pmegrid_nz).
927 if (direction == GMX_SUM_GRID_FORWARD)
929 send_id = overlap->send_id[ipulse];
930 recv_id = overlap->recv_id[ipulse];
931 send_index0 = overlap->comm_data[ipulse].send_index0;
932 send_nindex = overlap->comm_data[ipulse].send_nindex;
933 recv_index0 = overlap->comm_data[ipulse].recv_index0;
934 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
938 send_id = overlap->recv_id[ipulse];
939 recv_id = overlap->send_id[ipulse];
940 send_index0 = overlap->comm_data[ipulse].recv_index0;
941 send_nindex = overlap->comm_data[ipulse].recv_nindex;
942 recv_index0 = overlap->comm_data[ipulse].send_index0;
943 recv_nindex = overlap->comm_data[ipulse].send_nindex;
946 /* Copy data to contiguous send buffer */
949 fprintf(debug, "PME send node %d %d -> %d grid start %d Communicating %d to %d\n",
950 pme->nodeid, overlap->nodeid, send_id,
951 pme->pmegrid_start_iy,
952 send_index0-pme->pmegrid_start_iy,
953 send_index0-pme->pmegrid_start_iy+send_nindex);
956 for (i = 0; i < pme->pmegrid_nx; i++)
959 for (j = 0; j < send_nindex; j++)
961 iy = j + send_index0 - pme->pmegrid_start_iy;
962 for (k = 0; k < pme->nkz; k++)
965 overlap->sendbuf[icnt++] = grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz];
970 datasize = pme->pmegrid_nx * pme->nkz;
972 MPI_Sendrecv(overlap->sendbuf, send_nindex*datasize, GMX_MPI_REAL,
974 overlap->recvbuf, recv_nindex*datasize, GMX_MPI_REAL,
976 overlap->mpi_comm, &stat);
978 /* Get data from contiguous recv buffer */
981 fprintf(debug, "PME recv node %d %d <- %d grid start %d Communicating %d to %d\n",
982 pme->nodeid, overlap->nodeid, recv_id,
983 pme->pmegrid_start_iy,
984 recv_index0-pme->pmegrid_start_iy,
985 recv_index0-pme->pmegrid_start_iy+recv_nindex);
988 for (i = 0; i < pme->pmegrid_nx; i++)
991 for (j = 0; j < recv_nindex; j++)
993 iy = j + recv_index0 - pme->pmegrid_start_iy;
994 for (k = 0; k < pme->nkz; k++)
997 if (direction == GMX_SUM_GRID_FORWARD)
999 grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz] += overlap->recvbuf[icnt++];
1003 grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz] = overlap->recvbuf[icnt++];
1010 /* Major dimension is easier, no copying required,
1011 * but we might have to sum to separate array.
1012 * Since we don't copy, we have to communicate up to pmegrid_nz,
1013 * not nkz as for the minor direction.
1015 overlap = &pme->overlap[0];
1017 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
1019 if (direction == GMX_SUM_GRID_FORWARD)
1021 send_id = overlap->send_id[ipulse];
1022 recv_id = overlap->recv_id[ipulse];
1023 send_index0 = overlap->comm_data[ipulse].send_index0;
1024 send_nindex = overlap->comm_data[ipulse].send_nindex;
1025 recv_index0 = overlap->comm_data[ipulse].recv_index0;
1026 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
1027 recvptr = overlap->recvbuf;
1031 send_id = overlap->recv_id[ipulse];
1032 recv_id = overlap->send_id[ipulse];
1033 send_index0 = overlap->comm_data[ipulse].recv_index0;
1034 send_nindex = overlap->comm_data[ipulse].recv_nindex;
1035 recv_index0 = overlap->comm_data[ipulse].send_index0;
1036 recv_nindex = overlap->comm_data[ipulse].send_nindex;
1037 recvptr = grid + (recv_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1040 sendptr = grid + (send_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1041 datasize = pme->pmegrid_ny * pme->pmegrid_nz;
1045 fprintf(debug, "PME send node %d %d -> %d grid start %d Communicating %d to %d\n",
1046 pme->nodeid, overlap->nodeid, send_id,
1047 pme->pmegrid_start_ix,
1048 send_index0-pme->pmegrid_start_ix,
1049 send_index0-pme->pmegrid_start_ix+send_nindex);
1050 fprintf(debug, "PME recv node %d %d <- %d grid start %d Communicating %d to %d\n",
1051 pme->nodeid, overlap->nodeid, recv_id,
1052 pme->pmegrid_start_ix,
1053 recv_index0-pme->pmegrid_start_ix,
1054 recv_index0-pme->pmegrid_start_ix+recv_nindex);
1057 MPI_Sendrecv(sendptr, send_nindex*datasize, GMX_MPI_REAL,
1059 recvptr, recv_nindex*datasize, GMX_MPI_REAL,
1061 overlap->mpi_comm, &stat);
1063 /* ADD data from contiguous recv buffer */
1064 if (direction == GMX_SUM_GRID_FORWARD)
1066 p = grid + (recv_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1067 for (i = 0; i < recv_nindex*datasize; i++)
1069 p[i] += overlap->recvbuf[i];
1077 static int copy_pmegrid_to_fftgrid(gmx_pme_t pme, real *pmegrid, real *fftgrid, int grid_index)
1079 ivec local_fft_ndata, local_fft_offset, local_fft_size;
1080 ivec local_pme_size;
1084 /* Dimensions should be identical for A/B grid, so we just use A here */
1085 gmx_parallel_3dfft_real_limits(pme->pfft_setup[grid_index],
1090 local_pme_size[0] = pme->pmegrid_nx;
1091 local_pme_size[1] = pme->pmegrid_ny;
1092 local_pme_size[2] = pme->pmegrid_nz;
1094 /* The fftgrid is always 'justified' to the lower-left corner of the PME grid,
1095 the offset is identical, and the PME grid always has more data (due to overlap)
1102 sprintf(fn, "pmegrid%d.pdb", pme->nodeid);
1103 fp = gmx_ffopen(fn, "w");
1104 sprintf(fn, "pmegrid%d.txt", pme->nodeid);
1105 fp2 = gmx_ffopen(fn, "w");
1108 for (ix = 0; ix < local_fft_ndata[XX]; ix++)
1110 for (iy = 0; iy < local_fft_ndata[YY]; iy++)
1112 for (iz = 0; iz < local_fft_ndata[ZZ]; iz++)
1114 pmeidx = ix*(local_pme_size[YY]*local_pme_size[ZZ])+iy*(local_pme_size[ZZ])+iz;
1115 fftidx = ix*(local_fft_size[YY]*local_fft_size[ZZ])+iy*(local_fft_size[ZZ])+iz;
1116 fftgrid[fftidx] = pmegrid[pmeidx];
1118 val = 100*pmegrid[pmeidx];
1119 if (pmegrid[pmeidx] != 0)
1121 gmx_fprintf_pdb_atomline(fp, epdbATOM, pmeidx, "CA", ' ', "GLY", ' ', pmeidx, ' ',
1122 5.0*ix, 5.0*iy, 5.0*iz, 1.0, val, "");
1124 if (pmegrid[pmeidx] != 0)
1126 fprintf(fp2, "%-12s %5d %5d %5d %12.5e\n",
1128 pme->pmegrid_start_ix + ix,
1129 pme->pmegrid_start_iy + iy,
1130 pme->pmegrid_start_iz + iz,
1146 static gmx_cycles_t omp_cyc_start()
1148 return gmx_cycles_read();
1151 static gmx_cycles_t omp_cyc_end(gmx_cycles_t c)
1153 return gmx_cycles_read() - c;
1157 static int copy_fftgrid_to_pmegrid(gmx_pme_t pme, const real *fftgrid, real *pmegrid, int grid_index,
1158 int nthread, int thread)
1160 ivec local_fft_ndata, local_fft_offset, local_fft_size;
1161 ivec local_pme_size;
1162 int ixy0, ixy1, ixy, ix, iy, iz;
1164 #ifdef PME_TIME_THREADS
1166 static double cs1 = 0;
1170 #ifdef PME_TIME_THREADS
1171 c1 = omp_cyc_start();
1173 /* Dimensions should be identical for A/B grid, so we just use A here */
1174 gmx_parallel_3dfft_real_limits(pme->pfft_setup[grid_index],
1179 local_pme_size[0] = pme->pmegrid_nx;
1180 local_pme_size[1] = pme->pmegrid_ny;
1181 local_pme_size[2] = pme->pmegrid_nz;
1183 /* The fftgrid is always 'justified' to the lower-left corner of the PME grid,
1184 the offset is identical, and the PME grid always has more data (due to overlap)
1186 ixy0 = ((thread )*local_fft_ndata[XX]*local_fft_ndata[YY])/nthread;
1187 ixy1 = ((thread+1)*local_fft_ndata[XX]*local_fft_ndata[YY])/nthread;
1189 for (ixy = ixy0; ixy < ixy1; ixy++)
1191 ix = ixy/local_fft_ndata[YY];
1192 iy = ixy - ix*local_fft_ndata[YY];
1194 pmeidx = (ix*local_pme_size[YY] + iy)*local_pme_size[ZZ];
1195 fftidx = (ix*local_fft_size[YY] + iy)*local_fft_size[ZZ];
1196 for (iz = 0; iz < local_fft_ndata[ZZ]; iz++)
1198 pmegrid[pmeidx+iz] = fftgrid[fftidx+iz];
1202 #ifdef PME_TIME_THREADS
1203 c1 = omp_cyc_end(c1);
1208 printf("copy %.2f\n", cs1*1e-9);
1216 static void wrap_periodic_pmegrid(gmx_pme_t pme, real *pmegrid)
1218 int nx, ny, nz, pnx, pny, pnz, ny_x, overlap, ix, iy, iz;
1224 pnx = pme->pmegrid_nx;
1225 pny = pme->pmegrid_ny;
1226 pnz = pme->pmegrid_nz;
1228 overlap = pme->pme_order - 1;
1230 /* Add periodic overlap in z */
1231 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1233 for (iy = 0; iy < pme->pmegrid_ny; iy++)
1235 for (iz = 0; iz < overlap; iz++)
1237 pmegrid[(ix*pny+iy)*pnz+iz] +=
1238 pmegrid[(ix*pny+iy)*pnz+nz+iz];
1243 if (pme->nnodes_minor == 1)
1245 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1247 for (iy = 0; iy < overlap; iy++)
1249 for (iz = 0; iz < nz; iz++)
1251 pmegrid[(ix*pny+iy)*pnz+iz] +=
1252 pmegrid[(ix*pny+ny+iy)*pnz+iz];
1258 if (pme->nnodes_major == 1)
1260 ny_x = (pme->nnodes_minor == 1 ? ny : pme->pmegrid_ny);
1262 for (ix = 0; ix < overlap; ix++)
1264 for (iy = 0; iy < ny_x; iy++)
1266 for (iz = 0; iz < nz; iz++)
1268 pmegrid[(ix*pny+iy)*pnz+iz] +=
1269 pmegrid[((nx+ix)*pny+iy)*pnz+iz];
1277 static void unwrap_periodic_pmegrid(gmx_pme_t pme, real *pmegrid)
1279 int nx, ny, nz, pnx, pny, pnz, ny_x, overlap, ix;
1285 pnx = pme->pmegrid_nx;
1286 pny = pme->pmegrid_ny;
1287 pnz = pme->pmegrid_nz;
1289 overlap = pme->pme_order - 1;
1291 if (pme->nnodes_major == 1)
1293 ny_x = (pme->nnodes_minor == 1 ? ny : pme->pmegrid_ny);
1295 for (ix = 0; ix < overlap; ix++)
1299 for (iy = 0; iy < ny_x; iy++)
1301 for (iz = 0; iz < nz; iz++)
1303 pmegrid[((nx+ix)*pny+iy)*pnz+iz] =
1304 pmegrid[(ix*pny+iy)*pnz+iz];
1310 if (pme->nnodes_minor == 1)
1312 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1313 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1317 for (iy = 0; iy < overlap; iy++)
1319 for (iz = 0; iz < nz; iz++)
1321 pmegrid[(ix*pny+ny+iy)*pnz+iz] =
1322 pmegrid[(ix*pny+iy)*pnz+iz];
1328 /* Copy periodic overlap in z */
1329 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1330 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1334 for (iy = 0; iy < pme->pmegrid_ny; iy++)
1336 for (iz = 0; iz < overlap; iz++)
1338 pmegrid[(ix*pny+iy)*pnz+nz+iz] =
1339 pmegrid[(ix*pny+iy)*pnz+iz];
1346 /* This has to be a macro to enable full compiler optimization with xlC (and probably others too) */
1347 #define DO_BSPLINE(order) \
1348 for (ithx = 0; (ithx < order); ithx++) \
1350 index_x = (i0+ithx)*pny*pnz; \
1351 valx = coefficient*thx[ithx]; \
1353 for (ithy = 0; (ithy < order); ithy++) \
1355 valxy = valx*thy[ithy]; \
1356 index_xy = index_x+(j0+ithy)*pnz; \
1358 for (ithz = 0; (ithz < order); ithz++) \
1360 index_xyz = index_xy+(k0+ithz); \
1361 grid[index_xyz] += valxy*thz[ithz]; \
1367 static void spread_coefficients_bsplines_thread(pmegrid_t *pmegrid,
1368 pme_atomcomm_t *atc,
1369 splinedata_t *spline,
1370 pme_spline_work_t gmx_unused *work)
1373 /* spread coefficients from home atoms to local grid */
1376 int b, i, nn, n, ithx, ithy, ithz, i0, j0, k0;
1378 int order, norder, index_x, index_xy, index_xyz;
1379 real valx, valxy, coefficient;
1380 real *thx, *thy, *thz;
1381 int localsize, bndsize;
1382 int pnx, pny, pnz, ndatatot;
1383 int offx, offy, offz;
1385 #if defined PME_SIMD4_SPREAD_GATHER && !defined PME_SIMD4_UNALIGNED
1386 real thz_buffer[GMX_SIMD4_WIDTH*3], *thz_aligned;
1388 thz_aligned = gmx_simd4_align_r(thz_buffer);
1391 pnx = pmegrid->s[XX];
1392 pny = pmegrid->s[YY];
1393 pnz = pmegrid->s[ZZ];
1395 offx = pmegrid->offset[XX];
1396 offy = pmegrid->offset[YY];
1397 offz = pmegrid->offset[ZZ];
1399 ndatatot = pnx*pny*pnz;
1400 grid = pmegrid->grid;
1401 for (i = 0; i < ndatatot; i++)
1406 order = pmegrid->order;
1408 for (nn = 0; nn < spline->n; nn++)
1410 n = spline->ind[nn];
1411 coefficient = atc->coefficient[n];
1413 if (coefficient != 0)
1415 idxptr = atc->idx[n];
1418 i0 = idxptr[XX] - offx;
1419 j0 = idxptr[YY] - offy;
1420 k0 = idxptr[ZZ] - offz;
1422 thx = spline->theta[XX] + norder;
1423 thy = spline->theta[YY] + norder;
1424 thz = spline->theta[ZZ] + norder;
1429 #ifdef PME_SIMD4_SPREAD_GATHER
1430 #ifdef PME_SIMD4_UNALIGNED
1431 #define PME_SPREAD_SIMD4_ORDER4
1433 #define PME_SPREAD_SIMD4_ALIGNED
1436 #include "pme_simd4.h"
1442 #ifdef PME_SIMD4_SPREAD_GATHER
1443 #define PME_SPREAD_SIMD4_ALIGNED
1445 #include "pme_simd4.h"
1458 static void set_grid_alignment(int gmx_unused *pmegrid_nz, int gmx_unused pme_order)
1460 #ifdef PME_SIMD4_SPREAD_GATHER
1462 #ifndef PME_SIMD4_UNALIGNED
1467 /* Round nz up to a multiple of 4 to ensure alignment */
1468 *pmegrid_nz = ((*pmegrid_nz + 3) & ~3);
1473 static void set_gridsize_alignment(int gmx_unused *gridsize, int gmx_unused pme_order)
1475 #ifdef PME_SIMD4_SPREAD_GATHER
1476 #ifndef PME_SIMD4_UNALIGNED
1479 /* Add extra elements to ensured aligned operations do not go
1480 * beyond the allocated grid size.
1481 * Note that for pme_order=5, the pme grid z-size alignment
1482 * ensures that we will not go beyond the grid size.
1490 static void pmegrid_init(pmegrid_t *grid,
1491 int cx, int cy, int cz,
1492 int x0, int y0, int z0,
1493 int x1, int y1, int z1,
1494 gmx_bool set_alignment,
1503 grid->offset[XX] = x0;
1504 grid->offset[YY] = y0;
1505 grid->offset[ZZ] = z0;
1506 grid->n[XX] = x1 - x0 + pme_order - 1;
1507 grid->n[YY] = y1 - y0 + pme_order - 1;
1508 grid->n[ZZ] = z1 - z0 + pme_order - 1;
1509 copy_ivec(grid->n, grid->s);
1512 set_grid_alignment(&nz, pme_order);
1517 else if (nz != grid->s[ZZ])
1519 gmx_incons("pmegrid_init call with an unaligned z size");
1522 grid->order = pme_order;
1525 gridsize = grid->s[XX]*grid->s[YY]*grid->s[ZZ];
1526 set_gridsize_alignment(&gridsize, pme_order);
1527 snew_aligned(grid->grid, gridsize, SIMD4_ALIGNMENT);
1535 static int div_round_up(int enumerator, int denominator)
1537 return (enumerator + denominator - 1)/denominator;
1540 static void make_subgrid_division(const ivec n, int ovl, int nthread,
1543 int gsize_opt, gsize;
1548 for (nsx = 1; nsx <= nthread; nsx++)
1550 if (nthread % nsx == 0)
1552 for (nsy = 1; nsy <= nthread; nsy++)
1554 if (nsx*nsy <= nthread && nthread % (nsx*nsy) == 0)
1556 nsz = nthread/(nsx*nsy);
1558 /* Determine the number of grid points per thread */
1560 (div_round_up(n[XX], nsx) + ovl)*
1561 (div_round_up(n[YY], nsy) + ovl)*
1562 (div_round_up(n[ZZ], nsz) + ovl);
1564 /* Minimize the number of grids points per thread
1565 * and, secondarily, the number of cuts in minor dimensions.
1567 if (gsize_opt == -1 ||
1568 gsize < gsize_opt ||
1569 (gsize == gsize_opt &&
1570 (nsz < nsub[ZZ] || (nsz == nsub[ZZ] && nsy < nsub[YY]))))
1582 env = getenv("GMX_PME_THREAD_DIVISION");
1585 sscanf(env, "%d %d %d", &nsub[XX], &nsub[YY], &nsub[ZZ]);
1588 if (nsub[XX]*nsub[YY]*nsub[ZZ] != nthread)
1590 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);
1594 static void pmegrids_init(pmegrids_t *grids,
1595 int nx, int ny, int nz, int nz_base,
1597 gmx_bool bUseThreads,
1602 ivec n, n_base, g0, g1;
1603 int t, x, y, z, d, i, tfac;
1604 int max_comm_lines = -1;
1606 n[XX] = nx - (pme_order - 1);
1607 n[YY] = ny - (pme_order - 1);
1608 n[ZZ] = nz - (pme_order - 1);
1610 copy_ivec(n, n_base);
1611 n_base[ZZ] = nz_base;
1613 pmegrid_init(&grids->grid, 0, 0, 0, 0, 0, 0, n[XX], n[YY], n[ZZ], FALSE, pme_order,
1616 grids->nthread = nthread;
1618 make_subgrid_division(n_base, pme_order-1, grids->nthread, grids->nc);
1625 for (d = 0; d < DIM; d++)
1627 nst[d] = div_round_up(n[d], grids->nc[d]) + pme_order - 1;
1629 set_grid_alignment(&nst[ZZ], pme_order);
1633 fprintf(debug, "pmegrid thread local division: %d x %d x %d\n",
1634 grids->nc[XX], grids->nc[YY], grids->nc[ZZ]);
1635 fprintf(debug, "pmegrid %d %d %d max thread pmegrid %d %d %d\n",
1637 nst[XX], nst[YY], nst[ZZ]);
1640 snew(grids->grid_th, grids->nthread);
1642 gridsize = nst[XX]*nst[YY]*nst[ZZ];
1643 set_gridsize_alignment(&gridsize, pme_order);
1644 snew_aligned(grids->grid_all,
1645 grids->nthread*gridsize+(grids->nthread+1)*GMX_CACHE_SEP,
1648 for (x = 0; x < grids->nc[XX]; x++)
1650 for (y = 0; y < grids->nc[YY]; y++)
1652 for (z = 0; z < grids->nc[ZZ]; z++)
1654 pmegrid_init(&grids->grid_th[t],
1656 (n[XX]*(x ))/grids->nc[XX],
1657 (n[YY]*(y ))/grids->nc[YY],
1658 (n[ZZ]*(z ))/grids->nc[ZZ],
1659 (n[XX]*(x+1))/grids->nc[XX],
1660 (n[YY]*(y+1))/grids->nc[YY],
1661 (n[ZZ]*(z+1))/grids->nc[ZZ],
1664 grids->grid_all+GMX_CACHE_SEP+t*(gridsize+GMX_CACHE_SEP));
1672 grids->grid_th = NULL;
1675 snew(grids->g2t, DIM);
1677 for (d = DIM-1; d >= 0; d--)
1679 snew(grids->g2t[d], n[d]);
1681 for (i = 0; i < n[d]; i++)
1683 /* The second check should match the parameters
1684 * of the pmegrid_init call above.
1686 while (t + 1 < grids->nc[d] && i >= (n[d]*(t+1))/grids->nc[d])
1690 grids->g2t[d][i] = t*tfac;
1693 tfac *= grids->nc[d];
1697 case XX: max_comm_lines = overlap_x; break;
1698 case YY: max_comm_lines = overlap_y; break;
1699 case ZZ: max_comm_lines = pme_order - 1; break;
1701 grids->nthread_comm[d] = 0;
1702 while ((n[d]*grids->nthread_comm[d])/grids->nc[d] < max_comm_lines &&
1703 grids->nthread_comm[d] < grids->nc[d])
1705 grids->nthread_comm[d]++;
1709 fprintf(debug, "pmegrid thread grid communication range in %c: %d\n",
1710 'x'+d, grids->nthread_comm[d]);
1712 /* It should be possible to make grids->nthread_comm[d]==grids->nc[d]
1713 * work, but this is not a problematic restriction.
1715 if (grids->nc[d] > 1 && grids->nthread_comm[d] > grids->nc[d])
1717 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);
1723 static void pmegrids_destroy(pmegrids_t *grids)
1727 if (grids->grid.grid != NULL)
1729 sfree(grids->grid.grid);
1731 if (grids->nthread > 0)
1733 for (t = 0; t < grids->nthread; t++)
1735 sfree(grids->grid_th[t].grid);
1737 sfree(grids->grid_th);
1743 static void realloc_work(pme_work_t *work, int nkx)
1747 if (nkx > work->nalloc)
1750 srenew(work->mhx, work->nalloc);
1751 srenew(work->mhy, work->nalloc);
1752 srenew(work->mhz, work->nalloc);
1753 srenew(work->m2, work->nalloc);
1754 /* Allocate an aligned pointer for SIMD operations, including extra
1755 * elements at the end for padding.
1757 #ifdef PME_SIMD_SOLVE
1758 simd_width = GMX_SIMD_REAL_WIDTH;
1760 /* We can use any alignment, apart from 0, so we use 4 */
1763 sfree_aligned(work->denom);
1764 sfree_aligned(work->tmp1);
1765 sfree_aligned(work->tmp2);
1766 sfree_aligned(work->eterm);
1767 snew_aligned(work->denom, work->nalloc+simd_width, simd_width*sizeof(real));
1768 snew_aligned(work->tmp1, work->nalloc+simd_width, simd_width*sizeof(real));
1769 snew_aligned(work->tmp2, work->nalloc+simd_width, simd_width*sizeof(real));
1770 snew_aligned(work->eterm, work->nalloc+simd_width, simd_width*sizeof(real));
1771 srenew(work->m2inv, work->nalloc);
1776 static void free_work(pme_work_t *work)
1782 sfree_aligned(work->denom);
1783 sfree_aligned(work->tmp1);
1784 sfree_aligned(work->tmp2);
1785 sfree_aligned(work->eterm);
1790 #if defined PME_SIMD_SOLVE
1791 /* Calculate exponentials through SIMD */
1792 gmx_inline static void calc_exponentials_q(int gmx_unused start, int end, real f, real *d_aligned, real *r_aligned, real *e_aligned)
1795 const gmx_simd_real_t two = gmx_simd_set1_r(2.0);
1796 gmx_simd_real_t f_simd;
1798 gmx_simd_real_t tmp_d1, d_inv, tmp_r, tmp_e;
1800 f_simd = gmx_simd_set1_r(f);
1801 /* We only need to calculate from start. But since start is 0 or 1
1802 * and we want to use aligned loads/stores, we always start from 0.
1804 for (kx = 0; kx < end; kx += GMX_SIMD_REAL_WIDTH)
1806 tmp_d1 = gmx_simd_load_r(d_aligned+kx);
1807 d_inv = gmx_simd_inv_r(tmp_d1);
1808 tmp_r = gmx_simd_load_r(r_aligned+kx);
1809 tmp_r = gmx_simd_exp_r(tmp_r);
1810 tmp_e = gmx_simd_mul_r(f_simd, d_inv);
1811 tmp_e = gmx_simd_mul_r(tmp_e, tmp_r);
1812 gmx_simd_store_r(e_aligned+kx, tmp_e);
1817 gmx_inline static void calc_exponentials_q(int start, int end, real f, real *d, real *r, real *e)
1820 for (kx = start; kx < end; kx++)
1824 for (kx = start; kx < end; kx++)
1828 for (kx = start; kx < end; kx++)
1830 e[kx] = f*r[kx]*d[kx];
1835 #if defined PME_SIMD_SOLVE
1836 /* Calculate exponentials through SIMD */
1837 gmx_inline static void calc_exponentials_lj(int gmx_unused start, int end, real *r_aligned, real *factor_aligned, real *d_aligned)
1839 gmx_simd_real_t tmp_r, tmp_d, tmp_fac, d_inv, tmp_mk;
1840 const gmx_simd_real_t sqr_PI = gmx_simd_sqrt_r(gmx_simd_set1_r(M_PI));
1842 for (kx = 0; kx < end; kx += GMX_SIMD_REAL_WIDTH)
1844 /* We only need to calculate from start. But since start is 0 or 1
1845 * and we want to use aligned loads/stores, we always start from 0.
1847 tmp_d = gmx_simd_load_r(d_aligned+kx);
1848 d_inv = gmx_simd_inv_r(tmp_d);
1849 gmx_simd_store_r(d_aligned+kx, d_inv);
1850 tmp_r = gmx_simd_load_r(r_aligned+kx);
1851 tmp_r = gmx_simd_exp_r(tmp_r);
1852 gmx_simd_store_r(r_aligned+kx, tmp_r);
1853 tmp_mk = gmx_simd_load_r(factor_aligned+kx);
1854 tmp_fac = gmx_simd_mul_r(sqr_PI, gmx_simd_mul_r(tmp_mk, gmx_simd_erfc_r(tmp_mk)));
1855 gmx_simd_store_r(factor_aligned+kx, tmp_fac);
1859 gmx_inline static void calc_exponentials_lj(int start, int end, real *r, real *tmp2, real *d)
1863 for (kx = start; kx < end; kx++)
1868 for (kx = start; kx < end; kx++)
1873 for (kx = start; kx < end; kx++)
1876 tmp2[kx] = sqrt(M_PI)*mk*gmx_erfc(mk);
1881 static int solve_pme_yzx(gmx_pme_t pme, t_complex *grid,
1882 real ewaldcoeff, real vol,
1884 int nthread, int thread)
1886 /* do recip sum over local cells in grid */
1887 /* y major, z middle, x minor or continuous */
1889 int kx, ky, kz, maxkx, maxky, maxkz;
1890 int nx, ny, nz, iyz0, iyz1, iyz, iy, iz, kxstart, kxend;
1892 real factor = M_PI*M_PI/(ewaldcoeff*ewaldcoeff);
1893 real ets2, struct2, vfactor, ets2vf;
1894 real d1, d2, energy = 0;
1896 real virxx = 0, virxy = 0, virxz = 0, viryy = 0, viryz = 0, virzz = 0;
1897 real rxx, ryx, ryy, rzx, rzy, rzz;
1899 real *mhx, *mhy, *mhz, *m2, *denom, *tmp1, *eterm, *m2inv;
1900 real mhxk, mhyk, mhzk, m2k;
1903 ivec local_ndata, local_offset, local_size;
1906 elfac = ONE_4PI_EPS0/pme->epsilon_r;
1912 /* Dimensions should be identical for A/B grid, so we just use A here */
1913 gmx_parallel_3dfft_complex_limits(pme->pfft_setup[PME_GRID_QA],
1919 rxx = pme->recipbox[XX][XX];
1920 ryx = pme->recipbox[YY][XX];
1921 ryy = pme->recipbox[YY][YY];
1922 rzx = pme->recipbox[ZZ][XX];
1923 rzy = pme->recipbox[ZZ][YY];
1924 rzz = pme->recipbox[ZZ][ZZ];
1930 work = &pme->work[thread];
1935 denom = work->denom;
1937 eterm = work->eterm;
1938 m2inv = work->m2inv;
1940 iyz0 = local_ndata[YY]*local_ndata[ZZ]* thread /nthread;
1941 iyz1 = local_ndata[YY]*local_ndata[ZZ]*(thread+1)/nthread;
1943 for (iyz = iyz0; iyz < iyz1; iyz++)
1945 iy = iyz/local_ndata[ZZ];
1946 iz = iyz - iy*local_ndata[ZZ];
1948 ky = iy + local_offset[YY];
1959 by = M_PI*vol*pme->bsp_mod[YY][ky];
1961 kz = iz + local_offset[ZZ];
1965 bz = pme->bsp_mod[ZZ][kz];
1967 /* 0.5 correction for corner points */
1969 if (kz == 0 || kz == (nz+1)/2)
1974 p0 = grid + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
1976 /* We should skip the k-space point (0,0,0) */
1977 /* Note that since here x is the minor index, local_offset[XX]=0 */
1978 if (local_offset[XX] > 0 || ky > 0 || kz > 0)
1980 kxstart = local_offset[XX];
1984 kxstart = local_offset[XX] + 1;
1987 kxend = local_offset[XX] + local_ndata[XX];
1991 /* More expensive inner loop, especially because of the storage
1992 * of the mh elements in array's.
1993 * Because x is the minor grid index, all mh elements
1994 * depend on kx for triclinic unit cells.
1997 /* Two explicit loops to avoid a conditional inside the loop */
1998 for (kx = kxstart; kx < maxkx; kx++)
2003 mhyk = mx * ryx + my * ryy;
2004 mhzk = mx * rzx + my * rzy + mz * rzz;
2005 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2010 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2011 tmp1[kx] = -factor*m2k;
2014 for (kx = maxkx; kx < kxend; kx++)
2019 mhyk = mx * ryx + my * ryy;
2020 mhzk = mx * rzx + my * rzy + mz * rzz;
2021 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2026 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2027 tmp1[kx] = -factor*m2k;
2030 for (kx = kxstart; kx < kxend; kx++)
2032 m2inv[kx] = 1.0/m2[kx];
2035 calc_exponentials_q(kxstart, kxend, elfac, denom, tmp1, eterm);
2037 for (kx = kxstart; kx < kxend; kx++, p0++)
2042 p0->re = d1*eterm[kx];
2043 p0->im = d2*eterm[kx];
2045 struct2 = 2.0*(d1*d1+d2*d2);
2047 tmp1[kx] = eterm[kx]*struct2;
2050 for (kx = kxstart; kx < kxend; kx++)
2052 ets2 = corner_fac*tmp1[kx];
2053 vfactor = (factor*m2[kx] + 1.0)*2.0*m2inv[kx];
2056 ets2vf = ets2*vfactor;
2057 virxx += ets2vf*mhx[kx]*mhx[kx] - ets2;
2058 virxy += ets2vf*mhx[kx]*mhy[kx];
2059 virxz += ets2vf*mhx[kx]*mhz[kx];
2060 viryy += ets2vf*mhy[kx]*mhy[kx] - ets2;
2061 viryz += ets2vf*mhy[kx]*mhz[kx];
2062 virzz += ets2vf*mhz[kx]*mhz[kx] - ets2;
2067 /* We don't need to calculate the energy and the virial.
2068 * In this case the triclinic overhead is small.
2071 /* Two explicit loops to avoid a conditional inside the loop */
2073 for (kx = kxstart; kx < maxkx; kx++)
2078 mhyk = mx * ryx + my * ryy;
2079 mhzk = mx * rzx + my * rzy + mz * rzz;
2080 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2081 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2082 tmp1[kx] = -factor*m2k;
2085 for (kx = maxkx; kx < kxend; kx++)
2090 mhyk = mx * ryx + my * ryy;
2091 mhzk = mx * rzx + my * rzy + mz * rzz;
2092 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2093 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2094 tmp1[kx] = -factor*m2k;
2097 calc_exponentials_q(kxstart, kxend, elfac, denom, tmp1, eterm);
2099 for (kx = kxstart; kx < kxend; kx++, p0++)
2104 p0->re = d1*eterm[kx];
2105 p0->im = d2*eterm[kx];
2112 /* Update virial with local values.
2113 * The virial is symmetric by definition.
2114 * this virial seems ok for isotropic scaling, but I'm
2115 * experiencing problems on semiisotropic membranes.
2116 * IS THAT COMMENT STILL VALID??? (DvdS, 2001/02/07).
2118 work->vir_q[XX][XX] = 0.25*virxx;
2119 work->vir_q[YY][YY] = 0.25*viryy;
2120 work->vir_q[ZZ][ZZ] = 0.25*virzz;
2121 work->vir_q[XX][YY] = work->vir_q[YY][XX] = 0.25*virxy;
2122 work->vir_q[XX][ZZ] = work->vir_q[ZZ][XX] = 0.25*virxz;
2123 work->vir_q[YY][ZZ] = work->vir_q[ZZ][YY] = 0.25*viryz;
2125 /* This energy should be corrected for a charged system */
2126 work->energy_q = 0.5*energy;
2129 /* Return the loop count */
2130 return local_ndata[YY]*local_ndata[XX];
2133 static int solve_pme_lj_yzx(gmx_pme_t pme, t_complex **grid, gmx_bool bLB,
2134 real ewaldcoeff, real vol,
2135 gmx_bool bEnerVir, int nthread, int thread)
2137 /* do recip sum over local cells in grid */
2138 /* y major, z middle, x minor or continuous */
2140 int kx, ky, kz, maxkx, maxky, maxkz;
2141 int nx, ny, nz, iy, iyz0, iyz1, iyz, iz, kxstart, kxend;
2143 real factor = M_PI*M_PI/(ewaldcoeff*ewaldcoeff);
2145 real eterm, vterm, d1, d2, energy = 0;
2147 real virxx = 0, virxy = 0, virxz = 0, viryy = 0, viryz = 0, virzz = 0;
2148 real rxx, ryx, ryy, rzx, rzy, rzz;
2149 real *mhx, *mhy, *mhz, *m2, *denom, *tmp1, *tmp2;
2150 real mhxk, mhyk, mhzk, m2k;
2155 ivec local_ndata, local_offset, local_size;
2160 /* Dimensions should be identical for A/B grid, so we just use A here */
2161 gmx_parallel_3dfft_complex_limits(pme->pfft_setup[PME_GRID_C6A],
2166 rxx = pme->recipbox[XX][XX];
2167 ryx = pme->recipbox[YY][XX];
2168 ryy = pme->recipbox[YY][YY];
2169 rzx = pme->recipbox[ZZ][XX];
2170 rzy = pme->recipbox[ZZ][YY];
2171 rzz = pme->recipbox[ZZ][ZZ];
2177 work = &pme->work[thread];
2182 denom = work->denom;
2186 iyz0 = local_ndata[YY]*local_ndata[ZZ]* thread /nthread;
2187 iyz1 = local_ndata[YY]*local_ndata[ZZ]*(thread+1)/nthread;
2189 for (iyz = iyz0; iyz < iyz1; iyz++)
2191 iy = iyz/local_ndata[ZZ];
2192 iz = iyz - iy*local_ndata[ZZ];
2194 ky = iy + local_offset[YY];
2205 by = 3.0*vol*pme->bsp_mod[YY][ky]
2206 / (M_PI*sqrt(M_PI)*ewaldcoeff*ewaldcoeff*ewaldcoeff);
2208 kz = iz + local_offset[ZZ];
2212 bz = pme->bsp_mod[ZZ][kz];
2214 /* 0.5 correction for corner points */
2216 if (kz == 0 || kz == (nz+1)/2)
2221 kxstart = local_offset[XX];
2222 kxend = local_offset[XX] + local_ndata[XX];
2225 /* More expensive inner loop, especially because of the
2226 * storage of the mh elements in array's. Because x is the
2227 * minor grid index, all mh elements depend on kx for
2228 * triclinic unit cells.
2231 /* Two explicit loops to avoid a conditional inside the loop */
2232 for (kx = kxstart; kx < maxkx; kx++)
2237 mhyk = mx * ryx + my * ryy;
2238 mhzk = mx * rzx + my * rzy + mz * rzz;
2239 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2244 denom[kx] = bz*by*pme->bsp_mod[XX][kx];
2245 tmp1[kx] = -factor*m2k;
2246 tmp2[kx] = sqrt(factor*m2k);
2249 for (kx = maxkx; kx < kxend; kx++)
2254 mhyk = mx * ryx + my * ryy;
2255 mhzk = mx * rzx + my * rzy + mz * rzz;
2256 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2261 denom[kx] = bz*by*pme->bsp_mod[XX][kx];
2262 tmp1[kx] = -factor*m2k;
2263 tmp2[kx] = sqrt(factor*m2k);
2266 calc_exponentials_lj(kxstart, kxend, tmp1, tmp2, denom);
2268 for (kx = kxstart; kx < kxend; kx++)
2270 m2k = factor*m2[kx];
2271 eterm = -((1.0 - 2.0*m2k)*tmp1[kx]
2272 + 2.0*m2k*tmp2[kx]);
2273 vterm = 3.0*(-tmp1[kx] + tmp2[kx]);
2274 tmp1[kx] = eterm*denom[kx];
2275 tmp2[kx] = vterm*denom[kx];
2283 p0 = grid[0] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2284 for (kx = kxstart; kx < kxend; kx++, p0++)
2294 struct2 = 2.0*(d1*d1+d2*d2);
2296 tmp1[kx] = eterm*struct2;
2297 tmp2[kx] = vterm*struct2;
2302 real *struct2 = denom;
2305 for (kx = kxstart; kx < kxend; kx++)
2309 /* Due to symmetry we only need to calculate 4 of the 7 terms */
2310 for (ig = 0; ig <= 3; ++ig)
2315 p0 = grid[ig] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2316 p1 = grid[6-ig] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2317 scale = 2.0*lb_scale_factor_symm[ig];
2318 for (kx = kxstart; kx < kxend; ++kx, ++p0, ++p1)
2320 struct2[kx] += scale*(p0->re*p1->re + p0->im*p1->im);
2324 for (ig = 0; ig <= 6; ++ig)
2328 p0 = grid[ig] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2329 for (kx = kxstart; kx < kxend; kx++, p0++)
2339 for (kx = kxstart; kx < kxend; kx++)
2344 tmp1[kx] = eterm*str2;
2345 tmp2[kx] = vterm*str2;
2349 for (kx = kxstart; kx < kxend; kx++)
2351 ets2 = corner_fac*tmp1[kx];
2352 vterm = 2.0*factor*tmp2[kx];
2354 ets2vf = corner_fac*vterm;
2355 virxx += ets2vf*mhx[kx]*mhx[kx] - ets2;
2356 virxy += ets2vf*mhx[kx]*mhy[kx];
2357 virxz += ets2vf*mhx[kx]*mhz[kx];
2358 viryy += ets2vf*mhy[kx]*mhy[kx] - ets2;
2359 viryz += ets2vf*mhy[kx]*mhz[kx];
2360 virzz += ets2vf*mhz[kx]*mhz[kx] - ets2;
2365 /* We don't need to calculate the energy and the virial.
2366 * In this case the triclinic overhead is small.
2369 /* Two explicit loops to avoid a conditional inside the loop */
2371 for (kx = kxstart; kx < maxkx; kx++)
2376 mhyk = mx * ryx + my * ryy;
2377 mhzk = mx * rzx + my * rzy + mz * rzz;
2378 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2380 denom[kx] = bz*by*pme->bsp_mod[XX][kx];
2381 tmp1[kx] = -factor*m2k;
2382 tmp2[kx] = sqrt(factor*m2k);
2385 for (kx = maxkx; kx < kxend; kx++)
2390 mhyk = mx * ryx + my * ryy;
2391 mhzk = mx * rzx + my * rzy + mz * rzz;
2392 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2394 denom[kx] = bz*by*pme->bsp_mod[XX][kx];
2395 tmp1[kx] = -factor*m2k;
2396 tmp2[kx] = sqrt(factor*m2k);
2399 calc_exponentials_lj(kxstart, kxend, tmp1, tmp2, denom);
2401 for (kx = kxstart; kx < kxend; kx++)
2403 m2k = factor*m2[kx];
2404 eterm = -((1.0 - 2.0*m2k)*tmp1[kx]
2405 + 2.0*m2k*tmp2[kx]);
2406 tmp1[kx] = eterm*denom[kx];
2408 gcount = (bLB ? 7 : 1);
2409 for (ig = 0; ig < gcount; ++ig)
2413 p0 = grid[ig] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2414 for (kx = kxstart; kx < kxend; kx++, p0++)
2429 work->vir_lj[XX][XX] = 0.25*virxx;
2430 work->vir_lj[YY][YY] = 0.25*viryy;
2431 work->vir_lj[ZZ][ZZ] = 0.25*virzz;
2432 work->vir_lj[XX][YY] = work->vir_lj[YY][XX] = 0.25*virxy;
2433 work->vir_lj[XX][ZZ] = work->vir_lj[ZZ][XX] = 0.25*virxz;
2434 work->vir_lj[YY][ZZ] = work->vir_lj[ZZ][YY] = 0.25*viryz;
2436 /* This energy should be corrected for a charged system */
2437 work->energy_lj = 0.5*energy;
2439 /* Return the loop count */
2440 return local_ndata[YY]*local_ndata[XX];
2443 static void get_pme_ener_vir_q(const gmx_pme_t pme, int nthread,
2444 real *mesh_energy, matrix vir)
2446 /* This function sums output over threads and should therefore
2447 * only be called after thread synchronization.
2451 *mesh_energy = pme->work[0].energy_q;
2452 copy_mat(pme->work[0].vir_q, vir);
2454 for (thread = 1; thread < nthread; thread++)
2456 *mesh_energy += pme->work[thread].energy_q;
2457 m_add(vir, pme->work[thread].vir_q, vir);
2461 static void get_pme_ener_vir_lj(const gmx_pme_t pme, int nthread,
2462 real *mesh_energy, matrix vir)
2464 /* This function sums output over threads and should therefore
2465 * only be called after thread synchronization.
2469 *mesh_energy = pme->work[0].energy_lj;
2470 copy_mat(pme->work[0].vir_lj, vir);
2472 for (thread = 1; thread < nthread; thread++)
2474 *mesh_energy += pme->work[thread].energy_lj;
2475 m_add(vir, pme->work[thread].vir_lj, vir);
2480 #define DO_FSPLINE(order) \
2481 for (ithx = 0; (ithx < order); ithx++) \
2483 index_x = (i0+ithx)*pny*pnz; \
2487 for (ithy = 0; (ithy < order); ithy++) \
2489 index_xy = index_x+(j0+ithy)*pnz; \
2494 for (ithz = 0; (ithz < order); ithz++) \
2496 gval = grid[index_xy+(k0+ithz)]; \
2497 fxy1 += thz[ithz]*gval; \
2498 fz1 += dthz[ithz]*gval; \
2507 static void gather_f_bsplines(gmx_pme_t pme, real *grid,
2508 gmx_bool bClearF, pme_atomcomm_t *atc,
2509 splinedata_t *spline,
2512 /* sum forces for local particles */
2513 int nn, n, ithx, ithy, ithz, i0, j0, k0;
2514 int index_x, index_xy;
2515 int nx, ny, nz, pnx, pny, pnz;
2517 real tx, ty, dx, dy, coefficient;
2518 real fx, fy, fz, gval;
2520 real *thx, *thy, *thz, *dthx, *dthy, *dthz;
2522 real rxx, ryx, ryy, rzx, rzy, rzz;
2525 pme_spline_work_t *work;
2527 #if defined PME_SIMD4_SPREAD_GATHER && !defined PME_SIMD4_UNALIGNED
2528 real thz_buffer[GMX_SIMD4_WIDTH*3], *thz_aligned;
2529 real dthz_buffer[GMX_SIMD4_WIDTH*3], *dthz_aligned;
2531 thz_aligned = gmx_simd4_align_r(thz_buffer);
2532 dthz_aligned = gmx_simd4_align_r(dthz_buffer);
2535 work = pme->spline_work;
2537 order = pme->pme_order;
2538 thx = spline->theta[XX];
2539 thy = spline->theta[YY];
2540 thz = spline->theta[ZZ];
2541 dthx = spline->dtheta[XX];
2542 dthy = spline->dtheta[YY];
2543 dthz = spline->dtheta[ZZ];
2547 pnx = pme->pmegrid_nx;
2548 pny = pme->pmegrid_ny;
2549 pnz = pme->pmegrid_nz;
2551 rxx = pme->recipbox[XX][XX];
2552 ryx = pme->recipbox[YY][XX];
2553 ryy = pme->recipbox[YY][YY];
2554 rzx = pme->recipbox[ZZ][XX];
2555 rzy = pme->recipbox[ZZ][YY];
2556 rzz = pme->recipbox[ZZ][ZZ];
2558 for (nn = 0; nn < spline->n; nn++)
2560 n = spline->ind[nn];
2561 coefficient = scale*atc->coefficient[n];
2569 if (coefficient != 0)
2574 idxptr = atc->idx[n];
2581 /* Pointer arithmetic alert, next six statements */
2582 thx = spline->theta[XX] + norder;
2583 thy = spline->theta[YY] + norder;
2584 thz = spline->theta[ZZ] + norder;
2585 dthx = spline->dtheta[XX] + norder;
2586 dthy = spline->dtheta[YY] + norder;
2587 dthz = spline->dtheta[ZZ] + norder;
2592 #ifdef PME_SIMD4_SPREAD_GATHER
2593 #ifdef PME_SIMD4_UNALIGNED
2594 #define PME_GATHER_F_SIMD4_ORDER4
2596 #define PME_GATHER_F_SIMD4_ALIGNED
2599 #include "pme_simd4.h"
2605 #ifdef PME_SIMD4_SPREAD_GATHER
2606 #define PME_GATHER_F_SIMD4_ALIGNED
2608 #include "pme_simd4.h"
2618 atc->f[n][XX] += -coefficient*( fx*nx*rxx );
2619 atc->f[n][YY] += -coefficient*( fx*nx*ryx + fy*ny*ryy );
2620 atc->f[n][ZZ] += -coefficient*( fx*nx*rzx + fy*ny*rzy + fz*nz*rzz );
2623 /* Since the energy and not forces are interpolated
2624 * the net force might not be exactly zero.
2625 * This can be solved by also interpolating F, but
2626 * that comes at a cost.
2627 * A better hack is to remove the net force every
2628 * step, but that must be done at a higher level
2629 * since this routine doesn't see all atoms if running
2630 * in parallel. Don't know how important it is? EL 990726
2635 static real gather_energy_bsplines(gmx_pme_t pme, real *grid,
2636 pme_atomcomm_t *atc)
2638 splinedata_t *spline;
2639 int n, ithx, ithy, ithz, i0, j0, k0;
2640 int index_x, index_xy;
2642 real energy, pot, tx, ty, coefficient, gval;
2643 real *thx, *thy, *thz;
2647 spline = &atc->spline[0];
2649 order = pme->pme_order;
2652 for (n = 0; (n < atc->n); n++)
2654 coefficient = atc->coefficient[n];
2656 if (coefficient != 0)
2658 idxptr = atc->idx[n];
2665 /* Pointer arithmetic alert, next three statements */
2666 thx = spline->theta[XX] + norder;
2667 thy = spline->theta[YY] + norder;
2668 thz = spline->theta[ZZ] + norder;
2671 for (ithx = 0; (ithx < order); ithx++)
2673 index_x = (i0+ithx)*pme->pmegrid_ny*pme->pmegrid_nz;
2676 for (ithy = 0; (ithy < order); ithy++)
2678 index_xy = index_x+(j0+ithy)*pme->pmegrid_nz;
2681 for (ithz = 0; (ithz < order); ithz++)
2683 gval = grid[index_xy+(k0+ithz)];
2684 pot += tx*ty*thz[ithz]*gval;
2690 energy += pot*coefficient;
2697 /* Macro to force loop unrolling by fixing order.
2698 * This gives a significant performance gain.
2700 #define CALC_SPLINE(order) \
2704 real data[PME_ORDER_MAX]; \
2705 real ddata[PME_ORDER_MAX]; \
2707 for (j = 0; (j < DIM); j++) \
2711 /* dr is relative offset from lower cell limit */ \
2712 data[order-1] = 0; \
2716 for (k = 3; (k < order); k++) \
2718 div = 1.0/(k - 1.0); \
2719 data[k-1] = div*dr*data[k-2]; \
2720 for (l = 1; (l < (k-1)); l++) \
2722 data[k-l-1] = div*((dr+l)*data[k-l-2]+(k-l-dr)* \
2725 data[0] = div*(1-dr)*data[0]; \
2727 /* differentiate */ \
2728 ddata[0] = -data[0]; \
2729 for (k = 1; (k < order); k++) \
2731 ddata[k] = data[k-1] - data[k]; \
2734 div = 1.0/(order - 1); \
2735 data[order-1] = div*dr*data[order-2]; \
2736 for (l = 1; (l < (order-1)); l++) \
2738 data[order-l-1] = div*((dr+l)*data[order-l-2]+ \
2739 (order-l-dr)*data[order-l-1]); \
2741 data[0] = div*(1 - dr)*data[0]; \
2743 for (k = 0; k < order; k++) \
2745 theta[j][i*order+k] = data[k]; \
2746 dtheta[j][i*order+k] = ddata[k]; \
2751 void make_bsplines(splinevec theta, splinevec dtheta, int order,
2752 rvec fractx[], int nr, int ind[], real coefficient[],
2753 gmx_bool bDoSplines)
2755 /* construct splines for local atoms */
2759 for (i = 0; i < nr; i++)
2761 /* With free energy we do not use the coefficient check.
2762 * In most cases this will be more efficient than calling make_bsplines
2763 * twice, since usually more than half the particles have non-zero coefficients.
2766 if (bDoSplines || coefficient[ii] != 0.0)
2771 case 4: CALC_SPLINE(4); break;
2772 case 5: CALC_SPLINE(5); break;
2773 default: CALC_SPLINE(order); break;
2780 void make_dft_mod(real *mod, real *data, int ndata)
2785 for (i = 0; i < ndata; i++)
2788 for (j = 0; j < ndata; j++)
2790 arg = (2.0*M_PI*i*j)/ndata;
2791 sc += data[j]*cos(arg);
2792 ss += data[j]*sin(arg);
2794 mod[i] = sc*sc+ss*ss;
2796 for (i = 0; i < ndata; i++)
2800 mod[i] = (mod[i-1]+mod[i+1])*0.5;
2806 static void make_bspline_moduli(splinevec bsp_mod,
2807 int nx, int ny, int nz, int order)
2809 int nmax = max(nx, max(ny, nz));
2810 real *data, *ddata, *bsp_data;
2816 snew(bsp_data, nmax);
2822 for (k = 3; k < order; k++)
2826 for (l = 1; l < (k-1); l++)
2828 data[k-l-1] = div*(l*data[k-l-2]+(k-l)*data[k-l-1]);
2830 data[0] = div*data[0];
2833 ddata[0] = -data[0];
2834 for (k = 1; k < order; k++)
2836 ddata[k] = data[k-1]-data[k];
2838 div = 1.0/(order-1);
2840 for (l = 1; l < (order-1); l++)
2842 data[order-l-1] = div*(l*data[order-l-2]+(order-l)*data[order-l-1]);
2844 data[0] = div*data[0];
2846 for (i = 0; i < nmax; i++)
2850 for (i = 1; i <= order; i++)
2852 bsp_data[i] = data[i-1];
2855 make_dft_mod(bsp_mod[XX], bsp_data, nx);
2856 make_dft_mod(bsp_mod[YY], bsp_data, ny);
2857 make_dft_mod(bsp_mod[ZZ], bsp_data, nz);
2865 /* Return the P3M optimal influence function */
2866 static double do_p3m_influence(double z, int order)
2873 /* The formula and most constants can be found in:
2874 * Ballenegger et al., JCTC 8, 936 (2012)
2879 return 1.0 - 2.0*z2/3.0;
2882 return 1.0 - z2 + 2.0*z4/15.0;
2885 return 1.0 - 4.0*z2/3.0 + 2.0*z4/5.0 + 4.0*z2*z4/315.0;
2888 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;
2891 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;
2894 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;
2896 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;
2903 /* Calculate the P3M B-spline moduli for one dimension */
2904 static void make_p3m_bspline_moduli_dim(real *bsp_mod, int n, int order)
2906 double zarg, zai, sinzai, infl;
2911 gmx_fatal(FARGS, "The current P3M code only supports orders up to 8");
2918 for (i = -maxk; i < 0; i++)
2922 infl = do_p3m_influence(sinzai, order);
2923 bsp_mod[n+i] = infl*infl*pow(sinzai/zai, -2.0*order);
2926 for (i = 1; i < maxk; i++)
2930 infl = do_p3m_influence(sinzai, order);
2931 bsp_mod[i] = infl*infl*pow(sinzai/zai, -2.0*order);
2935 /* Calculate the P3M B-spline moduli */
2936 static void make_p3m_bspline_moduli(splinevec bsp_mod,
2937 int nx, int ny, int nz, int order)
2939 make_p3m_bspline_moduli_dim(bsp_mod[XX], nx, order);
2940 make_p3m_bspline_moduli_dim(bsp_mod[YY], ny, order);
2941 make_p3m_bspline_moduli_dim(bsp_mod[ZZ], nz, order);
2945 static void setup_coordinate_communication(pme_atomcomm_t *atc)
2953 for (i = 1; i <= nslab/2; i++)
2955 fw = (atc->nodeid + i) % nslab;
2956 bw = (atc->nodeid - i + nslab) % nslab;
2959 atc->node_dest[n] = fw;
2960 atc->node_src[n] = bw;
2965 atc->node_dest[n] = bw;
2966 atc->node_src[n] = fw;
2972 int gmx_pme_destroy(FILE *log, gmx_pme_t *pmedata)
2978 fprintf(log, "Destroying PME data structures.\n");
2981 sfree((*pmedata)->nnx);
2982 sfree((*pmedata)->nny);
2983 sfree((*pmedata)->nnz);
2985 for (i = 0; i < (*pmedata)->ngrids; ++i)
2987 pmegrids_destroy(&(*pmedata)->pmegrid[i]);
2988 sfree((*pmedata)->fftgrid[i]);
2989 sfree((*pmedata)->cfftgrid[i]);
2990 gmx_parallel_3dfft_destroy((*pmedata)->pfft_setup[i]);
2993 sfree((*pmedata)->lb_buf1);
2994 sfree((*pmedata)->lb_buf2);
2996 for (thread = 0; thread < (*pmedata)->nthread; thread++)
2998 free_work(&(*pmedata)->work[thread]);
3000 sfree((*pmedata)->work);
3008 static int mult_up(int n, int f)
3010 return ((n + f - 1)/f)*f;
3014 static double pme_load_imbalance(gmx_pme_t pme)
3019 nma = pme->nnodes_major;
3020 nmi = pme->nnodes_minor;
3022 n1 = mult_up(pme->nkx, nma)*mult_up(pme->nky, nmi)*pme->nkz;
3023 n2 = mult_up(pme->nkx, nma)*mult_up(pme->nkz, nmi)*pme->nky;
3024 n3 = mult_up(pme->nky, nma)*mult_up(pme->nkz, nmi)*pme->nkx;
3026 /* pme_solve is roughly double the cost of an fft */
3028 return (n1 + n2 + 3*n3)/(double)(6*pme->nkx*pme->nky*pme->nkz);
3031 static void init_atomcomm(gmx_pme_t pme, pme_atomcomm_t *atc,
3032 int dimind, gmx_bool bSpread)
3034 int nk, k, s, thread;
3036 atc->dimind = dimind;
3041 if (pme->nnodes > 1)
3043 atc->mpi_comm = pme->mpi_comm_d[dimind];
3044 MPI_Comm_size(atc->mpi_comm, &atc->nslab);
3045 MPI_Comm_rank(atc->mpi_comm, &atc->nodeid);
3049 fprintf(debug, "For PME atom communication in dimind %d: nslab %d rank %d\n", atc->dimind, atc->nslab, atc->nodeid);
3053 atc->bSpread = bSpread;
3054 atc->pme_order = pme->pme_order;
3058 snew(atc->node_dest, atc->nslab);
3059 snew(atc->node_src, atc->nslab);
3060 setup_coordinate_communication(atc);
3062 snew(atc->count_thread, pme->nthread);
3063 for (thread = 0; thread < pme->nthread; thread++)
3065 snew(atc->count_thread[thread], atc->nslab);
3067 atc->count = atc->count_thread[0];
3068 snew(atc->rcount, atc->nslab);
3069 snew(atc->buf_index, atc->nslab);
3072 atc->nthread = pme->nthread;
3073 if (atc->nthread > 1)
3075 snew(atc->thread_plist, atc->nthread);
3077 snew(atc->spline, atc->nthread);
3078 for (thread = 0; thread < atc->nthread; thread++)
3080 if (atc->nthread > 1)
3082 snew(atc->thread_plist[thread].n, atc->nthread+2*GMX_CACHE_SEP);
3083 atc->thread_plist[thread].n += GMX_CACHE_SEP;
3085 snew(atc->spline[thread].thread_one, pme->nthread);
3086 atc->spline[thread].thread_one[thread] = 1;
3091 init_overlap_comm(pme_overlap_t * ol,
3101 int lbnd, rbnd, maxlr, b, i;
3104 pme_grid_comm_t *pgc;
3106 int fft_start, fft_end, send_index1, recv_index1;
3110 ol->mpi_comm = comm;
3113 ol->nnodes = nnodes;
3114 ol->nodeid = nodeid;
3116 /* Linear translation of the PME grid won't affect reciprocal space
3117 * calculations, so to optimize we only interpolate "upwards",
3118 * which also means we only have to consider overlap in one direction.
3119 * I.e., particles on this node might also be spread to grid indices
3120 * that belong to higher nodes (modulo nnodes)
3123 snew(ol->s2g0, ol->nnodes+1);
3124 snew(ol->s2g1, ol->nnodes);
3127 fprintf(debug, "PME slab boundaries:");
3129 for (i = 0; i < nnodes; i++)
3131 /* s2g0 the local interpolation grid start.
3132 * s2g1 the local interpolation grid end.
3133 * Because grid overlap communication only goes forward,
3134 * the grid the slabs for fft's should be rounded down.
3136 ol->s2g0[i] = ( i *ndata + 0 )/nnodes;
3137 ol->s2g1[i] = ((i+1)*ndata + nnodes-1)/nnodes + norder - 1;
3141 fprintf(debug, " %3d %3d", ol->s2g0[i], ol->s2g1[i]);
3144 ol->s2g0[nnodes] = ndata;
3147 fprintf(debug, "\n");
3150 /* Determine with how many nodes we need to communicate the grid overlap */
3156 for (i = 0; i < nnodes; i++)
3158 if ((i+b < nnodes && ol->s2g1[i] > ol->s2g0[i+b]) ||
3159 (i+b >= nnodes && ol->s2g1[i] > ol->s2g0[i+b-nnodes] + ndata))
3165 while (bCont && b < nnodes);
3166 ol->noverlap_nodes = b - 1;
3168 snew(ol->send_id, ol->noverlap_nodes);
3169 snew(ol->recv_id, ol->noverlap_nodes);
3170 for (b = 0; b < ol->noverlap_nodes; b++)
3172 ol->send_id[b] = (ol->nodeid + (b + 1)) % ol->nnodes;
3173 ol->recv_id[b] = (ol->nodeid - (b + 1) + ol->nnodes) % ol->nnodes;
3175 snew(ol->comm_data, ol->noverlap_nodes);
3178 for (b = 0; b < ol->noverlap_nodes; b++)
3180 pgc = &ol->comm_data[b];
3182 fft_start = ol->s2g0[ol->send_id[b]];
3183 fft_end = ol->s2g0[ol->send_id[b]+1];
3184 if (ol->send_id[b] < nodeid)
3189 send_index1 = ol->s2g1[nodeid];
3190 send_index1 = min(send_index1, fft_end);
3191 pgc->send_index0 = fft_start;
3192 pgc->send_nindex = max(0, send_index1 - pgc->send_index0);
3193 ol->send_size += pgc->send_nindex;
3195 /* We always start receiving to the first index of our slab */
3196 fft_start = ol->s2g0[ol->nodeid];
3197 fft_end = ol->s2g0[ol->nodeid+1];
3198 recv_index1 = ol->s2g1[ol->recv_id[b]];
3199 if (ol->recv_id[b] > nodeid)
3201 recv_index1 -= ndata;
3203 recv_index1 = min(recv_index1, fft_end);
3204 pgc->recv_index0 = fft_start;
3205 pgc->recv_nindex = max(0, recv_index1 - pgc->recv_index0);
3209 /* Communicate the buffer sizes to receive */
3210 for (b = 0; b < ol->noverlap_nodes; b++)
3212 MPI_Sendrecv(&ol->send_size, 1, MPI_INT, ol->send_id[b], b,
3213 &ol->comm_data[b].recv_size, 1, MPI_INT, ol->recv_id[b], b,
3214 ol->mpi_comm, &stat);
3218 /* For non-divisible grid we need pme_order iso pme_order-1 */
3219 snew(ol->sendbuf, norder*commplainsize);
3220 snew(ol->recvbuf, norder*commplainsize);
3224 make_gridindex5_to_localindex(int n, int local_start, int local_range,
3225 int **global_to_local,
3226 real **fraction_shift)
3234 for (i = 0; (i < 5*n); i++)
3236 /* Determine the global to local grid index */
3237 gtl[i] = (i - local_start + n) % n;
3238 /* For coordinates that fall within the local grid the fraction
3239 * is correct, we don't need to shift it.
3242 if (local_range < n)
3244 /* Due to rounding issues i could be 1 beyond the lower or
3245 * upper boundary of the local grid. Correct the index for this.
3246 * If we shift the index, we need to shift the fraction by
3247 * the same amount in the other direction to not affect
3249 * Note that due to this shifting the weights at the end of
3250 * the spline might change, but that will only involve values
3251 * between zero and values close to the precision of a real,
3252 * which is anyhow the accuracy of the whole mesh calculation.
3254 /* With local_range=0 we should not change i=local_start */
3255 if (i % n != local_start)
3262 else if (gtl[i] == local_range)
3264 gtl[i] = local_range - 1;
3271 *global_to_local = gtl;
3272 *fraction_shift = fsh;
3275 static pme_spline_work_t *make_pme_spline_work(int gmx_unused order)
3277 pme_spline_work_t *work;
3279 #ifdef PME_SIMD4_SPREAD_GATHER
3280 real tmp[GMX_SIMD4_WIDTH*3], *tmp_aligned;
3281 gmx_simd4_real_t zero_S;
3282 gmx_simd4_real_t real_mask_S0, real_mask_S1;
3285 snew_aligned(work, 1, SIMD4_ALIGNMENT);
3287 tmp_aligned = gmx_simd4_align_r(tmp);
3289 zero_S = gmx_simd4_setzero_r();
3291 /* Generate bit masks to mask out the unused grid entries,
3292 * as we only operate on order of the 8 grid entries that are
3293 * load into 2 SIMD registers.
3295 for (of = 0; of < 2*GMX_SIMD4_WIDTH-(order-1); of++)
3297 for (i = 0; i < 2*GMX_SIMD4_WIDTH; i++)
3299 tmp_aligned[i] = (i >= of && i < of+order ? -1.0 : 1.0);
3301 real_mask_S0 = gmx_simd4_load_r(tmp_aligned);
3302 real_mask_S1 = gmx_simd4_load_r(tmp_aligned+GMX_SIMD4_WIDTH);
3303 work->mask_S0[of] = gmx_simd4_cmplt_r(real_mask_S0, zero_S);
3304 work->mask_S1[of] = gmx_simd4_cmplt_r(real_mask_S1, zero_S);
3313 void gmx_pme_check_restrictions(int pme_order,
3314 int nkx, int nky, int nkz,
3317 gmx_bool bUseThreads,
3319 gmx_bool *bValidSettings)
3321 if (pme_order > PME_ORDER_MAX)
3325 *bValidSettings = FALSE;
3328 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.",
3329 pme_order, PME_ORDER_MAX);
3332 if (nkx <= pme_order*(nnodes_major > 1 ? 2 : 1) ||
3333 nky <= pme_order*(nnodes_minor > 1 ? 2 : 1) ||
3338 *bValidSettings = FALSE;
3341 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",
3345 /* Check for a limitation of the (current) sum_fftgrid_dd code.
3346 * We only allow multiple communication pulses in dim 1, not in dim 0.
3348 if (bUseThreads && (nkx < nnodes_major*pme_order &&
3349 nkx != nnodes_major*(pme_order - 1)))
3353 *bValidSettings = FALSE;
3356 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.",
3357 nkx/(double)nnodes_major, pme_order);
3360 if (bValidSettings != NULL)
3362 *bValidSettings = TRUE;
3368 int gmx_pme_init(gmx_pme_t * pmedata,
3374 gmx_bool bFreeEnergy_q,
3375 gmx_bool bFreeEnergy_lj,
3376 gmx_bool bReproducible,
3379 gmx_pme_t pme = NULL;
3381 int use_threads, sum_use_threads, i;
3386 fprintf(debug, "Creating PME data structures.\n");
3390 pme->sum_qgrid_tmp = NULL;
3391 pme->sum_qgrid_dd_tmp = NULL;
3392 pme->buf_nalloc = 0;
3395 pme->bPPnode = TRUE;
3397 pme->nnodes_major = nnodes_major;
3398 pme->nnodes_minor = nnodes_minor;
3401 if (nnodes_major*nnodes_minor > 1)
3403 pme->mpi_comm = cr->mpi_comm_mygroup;
3405 MPI_Comm_rank(pme->mpi_comm, &pme->nodeid);
3406 MPI_Comm_size(pme->mpi_comm, &pme->nnodes);
3407 if (pme->nnodes != nnodes_major*nnodes_minor)
3409 gmx_incons("PME node count mismatch");
3414 pme->mpi_comm = MPI_COMM_NULL;
3418 if (pme->nnodes == 1)
3421 pme->mpi_comm_d[0] = MPI_COMM_NULL;
3422 pme->mpi_comm_d[1] = MPI_COMM_NULL;
3424 pme->ndecompdim = 0;
3425 pme->nodeid_major = 0;
3426 pme->nodeid_minor = 0;
3428 pme->mpi_comm_d[0] = pme->mpi_comm_d[1] = MPI_COMM_NULL;
3433 if (nnodes_minor == 1)
3436 pme->mpi_comm_d[0] = pme->mpi_comm;
3437 pme->mpi_comm_d[1] = MPI_COMM_NULL;
3439 pme->ndecompdim = 1;
3440 pme->nodeid_major = pme->nodeid;
3441 pme->nodeid_minor = 0;
3444 else if (nnodes_major == 1)
3447 pme->mpi_comm_d[0] = MPI_COMM_NULL;
3448 pme->mpi_comm_d[1] = pme->mpi_comm;
3450 pme->ndecompdim = 1;
3451 pme->nodeid_major = 0;
3452 pme->nodeid_minor = pme->nodeid;
3456 if (pme->nnodes % nnodes_major != 0)
3458 gmx_incons("For 2D PME decomposition, #PME nodes must be divisible by the number of nodes in the major dimension");
3460 pme->ndecompdim = 2;
3463 MPI_Comm_split(pme->mpi_comm, pme->nodeid % nnodes_minor,
3464 pme->nodeid, &pme->mpi_comm_d[0]); /* My communicator along major dimension */
3465 MPI_Comm_split(pme->mpi_comm, pme->nodeid/nnodes_minor,
3466 pme->nodeid, &pme->mpi_comm_d[1]); /* My communicator along minor dimension */
3468 MPI_Comm_rank(pme->mpi_comm_d[0], &pme->nodeid_major);
3469 MPI_Comm_size(pme->mpi_comm_d[0], &pme->nnodes_major);
3470 MPI_Comm_rank(pme->mpi_comm_d[1], &pme->nodeid_minor);
3471 MPI_Comm_size(pme->mpi_comm_d[1], &pme->nnodes_minor);
3474 pme->bPPnode = (cr->duty & DUTY_PP);
3477 pme->nthread = nthread;
3479 /* Check if any of the PME MPI ranks uses threads */
3480 use_threads = (pme->nthread > 1 ? 1 : 0);
3482 if (pme->nnodes > 1)
3484 MPI_Allreduce(&use_threads, &sum_use_threads, 1, MPI_INT,
3485 MPI_SUM, pme->mpi_comm);
3490 sum_use_threads = use_threads;
3492 pme->bUseThreads = (sum_use_threads > 0);
3494 if (ir->ePBC == epbcSCREW)
3496 gmx_fatal(FARGS, "pme does not (yet) work with pbc = screw");
3499 pme->bFEP_q = ((ir->efep != efepNO) && bFreeEnergy_q);
3500 pme->bFEP_lj = ((ir->efep != efepNO) && bFreeEnergy_lj);
3501 pme->bFEP = (pme->bFEP_q || pme->bFEP_lj);
3505 pme->bP3M = (ir->coulombtype == eelP3M_AD || getenv("GMX_PME_P3M") != NULL);
3506 pme->pme_order = ir->pme_order;
3508 /* Always constant electrostatics coefficients */
3509 pme->epsilon_r = ir->epsilon_r;
3511 /* Always constant LJ coefficients */
3512 pme->ljpme_combination_rule = ir->ljpme_combination_rule;
3514 /* If we violate restrictions, generate a fatal error here */
3515 gmx_pme_check_restrictions(pme->pme_order,
3516 pme->nkx, pme->nky, pme->nkz,
3523 if (pme->nnodes > 1)
3528 MPI_Type_contiguous(DIM, mpi_type, &(pme->rvec_mpi));
3529 MPI_Type_commit(&(pme->rvec_mpi));
3532 /* Note that the coefficient spreading and force gathering, which usually
3533 * takes about the same amount of time as FFT+solve_pme,
3534 * is always fully load balanced
3535 * (unless the coefficient distribution is inhomogeneous).
3538 imbal = pme_load_imbalance(pme);
3539 if (imbal >= 1.2 && pme->nodeid_major == 0 && pme->nodeid_minor == 0)
3543 "NOTE: The load imbalance in PME FFT and solve is %d%%.\n"
3544 " For optimal PME load balancing\n"
3545 " PME grid_x (%d) and grid_y (%d) should be divisible by #PME_nodes_x (%d)\n"
3546 " and PME grid_y (%d) and grid_z (%d) should be divisible by #PME_nodes_y (%d)\n"
3548 (int)((imbal-1)*100 + 0.5),
3549 pme->nkx, pme->nky, pme->nnodes_major,
3550 pme->nky, pme->nkz, pme->nnodes_minor);
3554 /* For non-divisible grid we need pme_order iso pme_order-1 */
3555 /* In sum_qgrid_dd x overlap is copied in place: take padding into account.
3556 * y is always copied through a buffer: we don't need padding in z,
3557 * but we do need the overlap in x because of the communication order.
3559 init_overlap_comm(&pme->overlap[0], pme->pme_order,
3563 pme->nnodes_major, pme->nodeid_major,
3565 (div_round_up(pme->nky, pme->nnodes_minor)+pme->pme_order)*(pme->nkz+pme->pme_order-1));
3567 /* Along overlap dim 1 we can send in multiple pulses in sum_fftgrid_dd.
3568 * We do this with an offset buffer of equal size, so we need to allocate
3569 * extra for the offset. That's what the (+1)*pme->nkz is for.
3571 init_overlap_comm(&pme->overlap[1], pme->pme_order,
3575 pme->nnodes_minor, pme->nodeid_minor,
3577 (div_round_up(pme->nkx, pme->nnodes_major)+pme->pme_order+1)*pme->nkz);
3579 /* Double-check for a limitation of the (current) sum_fftgrid_dd code.
3580 * Note that gmx_pme_check_restrictions checked for this already.
3582 if (pme->bUseThreads && pme->overlap[0].noverlap_nodes > 1)
3584 gmx_incons("More than one communication pulse required for grid overlap communication along the major dimension while using threads");
3587 snew(pme->bsp_mod[XX], pme->nkx);
3588 snew(pme->bsp_mod[YY], pme->nky);
3589 snew(pme->bsp_mod[ZZ], pme->nkz);
3591 /* The required size of the interpolation grid, including overlap.
3592 * The allocated size (pmegrid_n?) might be slightly larger.
3594 pme->pmegrid_nx = pme->overlap[0].s2g1[pme->nodeid_major] -
3595 pme->overlap[0].s2g0[pme->nodeid_major];
3596 pme->pmegrid_ny = pme->overlap[1].s2g1[pme->nodeid_minor] -
3597 pme->overlap[1].s2g0[pme->nodeid_minor];
3598 pme->pmegrid_nz_base = pme->nkz;
3599 pme->pmegrid_nz = pme->pmegrid_nz_base + pme->pme_order - 1;
3600 set_grid_alignment(&pme->pmegrid_nz, pme->pme_order);
3602 pme->pmegrid_start_ix = pme->overlap[0].s2g0[pme->nodeid_major];
3603 pme->pmegrid_start_iy = pme->overlap[1].s2g0[pme->nodeid_minor];
3604 pme->pmegrid_start_iz = 0;
3606 make_gridindex5_to_localindex(pme->nkx,
3607 pme->pmegrid_start_ix,
3608 pme->pmegrid_nx - (pme->pme_order-1),
3609 &pme->nnx, &pme->fshx);
3610 make_gridindex5_to_localindex(pme->nky,
3611 pme->pmegrid_start_iy,
3612 pme->pmegrid_ny - (pme->pme_order-1),
3613 &pme->nny, &pme->fshy);
3614 make_gridindex5_to_localindex(pme->nkz,
3615 pme->pmegrid_start_iz,
3616 pme->pmegrid_nz_base,
3617 &pme->nnz, &pme->fshz);
3619 pme->spline_work = make_pme_spline_work(pme->pme_order);
3621 ndata[0] = pme->nkx;
3622 ndata[1] = pme->nky;
3623 ndata[2] = pme->nkz;
3624 /* It doesn't matter if we allocate too many grids here,
3625 * we only allocate and use the ones we need.
3627 if (EVDW_PME(ir->vdwtype))
3629 pme->ngrids = ((ir->ljpme_combination_rule == eljpmeLB) ? DO_Q_AND_LJ_LB : DO_Q_AND_LJ);
3635 snew(pme->fftgrid, pme->ngrids);
3636 snew(pme->cfftgrid, pme->ngrids);
3637 snew(pme->pfft_setup, pme->ngrids);
3639 for (i = 0; i < pme->ngrids; ++i)
3641 if ((i < DO_Q && EEL_PME(ir->coulombtype) && (i == 0 ||
3643 (i >= DO_Q && EVDW_PME(ir->vdwtype) && (i == 2 ||
3645 ir->ljpme_combination_rule == eljpmeLB)))
3647 pmegrids_init(&pme->pmegrid[i],
3648 pme->pmegrid_nx, pme->pmegrid_ny, pme->pmegrid_nz,
3649 pme->pmegrid_nz_base,
3653 pme->overlap[0].s2g1[pme->nodeid_major]-pme->overlap[0].s2g0[pme->nodeid_major+1],
3654 pme->overlap[1].s2g1[pme->nodeid_minor]-pme->overlap[1].s2g0[pme->nodeid_minor+1]);
3655 /* This routine will allocate the grid data to fit the FFTs */
3656 gmx_parallel_3dfft_init(&pme->pfft_setup[i], ndata,
3657 &pme->fftgrid[i], &pme->cfftgrid[i],
3659 bReproducible, pme->nthread);
3666 /* Use plain SPME B-spline interpolation */
3667 make_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
3671 /* Use the P3M grid-optimized influence function */
3672 make_p3m_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
3675 /* Use atc[0] for spreading */
3676 init_atomcomm(pme, &pme->atc[0], nnodes_major > 1 ? 0 : 1, TRUE);
3677 if (pme->ndecompdim >= 2)
3679 init_atomcomm(pme, &pme->atc[1], 1, FALSE);
3682 if (pme->nnodes == 1)
3684 pme->atc[0].n = homenr;
3685 pme_realloc_atomcomm_things(&pme->atc[0]);
3688 pme->lb_buf1 = NULL;
3689 pme->lb_buf2 = NULL;
3690 pme->lb_buf_nalloc = 0;
3695 /* Use fft5d, order after FFT is y major, z, x minor */
3697 snew(pme->work, pme->nthread);
3698 for (thread = 0; thread < pme->nthread; thread++)
3700 realloc_work(&pme->work[thread], pme->nkx);
3709 static void reuse_pmegrids(const pmegrids_t *old, pmegrids_t *new)
3713 for (d = 0; d < DIM; d++)
3715 if (new->grid.n[d] > old->grid.n[d])
3721 sfree_aligned(new->grid.grid);
3722 new->grid.grid = old->grid.grid;
3724 if (new->grid_th != NULL && new->nthread == old->nthread)
3726 sfree_aligned(new->grid_all);
3727 for (t = 0; t < new->nthread; t++)
3729 new->grid_th[t].grid = old->grid_th[t].grid;
3734 int gmx_pme_reinit(gmx_pme_t * pmedata,
3737 const t_inputrec * ir,
3745 irc.nkx = grid_size[XX];
3746 irc.nky = grid_size[YY];
3747 irc.nkz = grid_size[ZZ];
3749 if (pme_src->nnodes == 1)
3751 homenr = pme_src->atc[0].n;
3758 ret = gmx_pme_init(pmedata, cr, pme_src->nnodes_major, pme_src->nnodes_minor,
3759 &irc, homenr, pme_src->bFEP_q, pme_src->bFEP_lj, FALSE, pme_src->nthread);
3763 /* We can easily reuse the allocated pme grids in pme_src */
3764 reuse_pmegrids(&pme_src->pmegrid[PME_GRID_QA], &(*pmedata)->pmegrid[PME_GRID_QA]);
3765 /* We would like to reuse the fft grids, but that's harder */
3772 static void copy_local_grid(gmx_pme_t pme, pmegrids_t *pmegrids,
3773 int grid_index, int thread, real *fftgrid)
3775 ivec local_fft_ndata, local_fft_offset, local_fft_size;
3779 int offx, offy, offz, x, y, z, i0, i0t;
3784 gmx_parallel_3dfft_real_limits(pme->pfft_setup[grid_index],
3788 fft_my = local_fft_size[YY];
3789 fft_mz = local_fft_size[ZZ];
3791 pmegrid = &pmegrids->grid_th[thread];
3793 nsx = pmegrid->s[XX];
3794 nsy = pmegrid->s[YY];
3795 nsz = pmegrid->s[ZZ];
3797 for (d = 0; d < DIM; d++)
3799 nf[d] = min(pmegrid->n[d] - (pmegrid->order - 1),
3800 local_fft_ndata[d] - pmegrid->offset[d]);
3803 offx = pmegrid->offset[XX];
3804 offy = pmegrid->offset[YY];
3805 offz = pmegrid->offset[ZZ];
3807 /* Directly copy the non-overlapping parts of the local grids.
3808 * This also initializes the full grid.
3810 grid_th = pmegrid->grid;
3811 for (x = 0; x < nf[XX]; x++)
3813 for (y = 0; y < nf[YY]; y++)
3815 i0 = ((offx + x)*fft_my + (offy + y))*fft_mz + offz;
3816 i0t = (x*nsy + y)*nsz;
3817 for (z = 0; z < nf[ZZ]; z++)
3819 fftgrid[i0+z] = grid_th[i0t+z];
3826 reduce_threadgrid_overlap(gmx_pme_t pme,
3827 const pmegrids_t *pmegrids, int thread,
3828 real *fftgrid, real *commbuf_x, real *commbuf_y,
3831 ivec local_fft_ndata, local_fft_offset, local_fft_size;
3832 int fft_nx, fft_ny, fft_nz;
3837 int offx, offy, offz, x, y, z, i0, i0t;
3838 int sx, sy, sz, fx, fy, fz, tx1, ty1, tz1, ox, oy, oz;
3839 gmx_bool bClearBufX, bClearBufY, bClearBufXY, bClearBuf;
3840 gmx_bool bCommX, bCommY;
3843 const pmegrid_t *pmegrid, *pmegrid_g, *pmegrid_f;
3844 const real *grid_th;
3845 real *commbuf = NULL;
3847 gmx_parallel_3dfft_real_limits(pme->pfft_setup[grid_index],
3851 fft_nx = local_fft_ndata[XX];
3852 fft_ny = local_fft_ndata[YY];
3853 fft_nz = local_fft_ndata[ZZ];
3855 fft_my = local_fft_size[YY];
3856 fft_mz = local_fft_size[ZZ];
3858 /* This routine is called when all thread have finished spreading.
3859 * Here each thread sums grid contributions calculated by other threads
3860 * to the thread local grid volume.
3861 * To minimize the number of grid copying operations,
3862 * this routines sums immediately from the pmegrid to the fftgrid.
3865 /* Determine which part of the full node grid we should operate on,
3866 * this is our thread local part of the full grid.
3868 pmegrid = &pmegrids->grid_th[thread];
3870 for (d = 0; d < DIM; d++)
3872 ne[d] = min(pmegrid->offset[d]+pmegrid->n[d]-(pmegrid->order-1),
3873 local_fft_ndata[d]);
3876 offx = pmegrid->offset[XX];
3877 offy = pmegrid->offset[YY];
3878 offz = pmegrid->offset[ZZ];
3885 /* Now loop over all the thread data blocks that contribute
3886 * to the grid region we (our thread) are operating on.
3888 /* Note that ffy_nx/y is equal to the number of grid points
3889 * between the first point of our node grid and the one of the next node.
3891 for (sx = 0; sx >= -pmegrids->nthread_comm[XX]; sx--)
3893 fx = pmegrid->ci[XX] + sx;
3898 fx += pmegrids->nc[XX];
3900 bCommX = (pme->nnodes_major > 1);
3902 pmegrid_g = &pmegrids->grid_th[fx*pmegrids->nc[YY]*pmegrids->nc[ZZ]];
3903 ox += pmegrid_g->offset[XX];
3906 tx1 = min(ox + pmegrid_g->n[XX], ne[XX]);
3910 tx1 = min(ox + pmegrid_g->n[XX], pme->pme_order);
3913 for (sy = 0; sy >= -pmegrids->nthread_comm[YY]; sy--)
3915 fy = pmegrid->ci[YY] + sy;
3920 fy += pmegrids->nc[YY];
3922 bCommY = (pme->nnodes_minor > 1);
3924 pmegrid_g = &pmegrids->grid_th[fy*pmegrids->nc[ZZ]];
3925 oy += pmegrid_g->offset[YY];
3928 ty1 = min(oy + pmegrid_g->n[YY], ne[YY]);
3932 ty1 = min(oy + pmegrid_g->n[YY], pme->pme_order);
3935 for (sz = 0; sz >= -pmegrids->nthread_comm[ZZ]; sz--)
3937 fz = pmegrid->ci[ZZ] + sz;
3941 fz += pmegrids->nc[ZZ];
3944 pmegrid_g = &pmegrids->grid_th[fz];
3945 oz += pmegrid_g->offset[ZZ];
3946 tz1 = min(oz + pmegrid_g->n[ZZ], ne[ZZ]);
3948 if (sx == 0 && sy == 0 && sz == 0)
3950 /* We have already added our local contribution
3951 * before calling this routine, so skip it here.
3956 thread_f = (fx*pmegrids->nc[YY] + fy)*pmegrids->nc[ZZ] + fz;
3958 pmegrid_f = &pmegrids->grid_th[thread_f];
3960 grid_th = pmegrid_f->grid;
3962 nsx = pmegrid_f->s[XX];
3963 nsy = pmegrid_f->s[YY];
3964 nsz = pmegrid_f->s[ZZ];
3966 #ifdef DEBUG_PME_REDUCE
3967 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",
3968 pme->nodeid, thread, thread_f,
3969 pme->pmegrid_start_ix,
3970 pme->pmegrid_start_iy,
3971 pme->pmegrid_start_iz,
3973 offx-ox, tx1-ox, offx, tx1,
3974 offy-oy, ty1-oy, offy, ty1,
3975 offz-oz, tz1-oz, offz, tz1);
3978 if (!(bCommX || bCommY))
3980 /* Copy from the thread local grid to the node grid */
3981 for (x = offx; x < tx1; x++)
3983 for (y = offy; y < ty1; y++)
3985 i0 = (x*fft_my + y)*fft_mz;
3986 i0t = ((x - ox)*nsy + (y - oy))*nsz - oz;
3987 for (z = offz; z < tz1; z++)
3989 fftgrid[i0+z] += grid_th[i0t+z];
3996 /* The order of this conditional decides
3997 * where the corner volume gets stored with x+y decomp.
4001 commbuf = commbuf_y;
4002 buf_my = ty1 - offy;
4005 /* We index commbuf modulo the local grid size */
4006 commbuf += buf_my*fft_nx*fft_nz;
4008 bClearBuf = bClearBufXY;
4009 bClearBufXY = FALSE;
4013 bClearBuf = bClearBufY;
4019 commbuf = commbuf_x;
4021 bClearBuf = bClearBufX;
4025 /* Copy to the communication buffer */
4026 for (x = offx; x < tx1; x++)
4028 for (y = offy; y < ty1; y++)
4030 i0 = (x*buf_my + y)*fft_nz;
4031 i0t = ((x - ox)*nsy + (y - oy))*nsz - oz;
4035 /* First access of commbuf, initialize it */
4036 for (z = offz; z < tz1; z++)
4038 commbuf[i0+z] = grid_th[i0t+z];
4043 for (z = offz; z < tz1; z++)
4045 commbuf[i0+z] += grid_th[i0t+z];
4057 static void sum_fftgrid_dd(gmx_pme_t pme, real *fftgrid, int grid_index)
4059 ivec local_fft_ndata, local_fft_offset, local_fft_size;
4060 pme_overlap_t *overlap;
4061 int send_index0, send_nindex;
4066 int send_size_y, recv_size_y;
4067 int ipulse, send_id, recv_id, datasize, gridsize, size_yx;
4068 real *sendptr, *recvptr;
4069 int x, y, z, indg, indb;
4071 /* Note that this routine is only used for forward communication.
4072 * Since the force gathering, unlike the coefficient spreading,
4073 * can be trivially parallelized over the particles,
4074 * the backwards process is much simpler and can use the "old"
4075 * communication setup.
4078 gmx_parallel_3dfft_real_limits(pme->pfft_setup[grid_index],
4083 if (pme->nnodes_minor > 1)
4085 /* Major dimension */
4086 overlap = &pme->overlap[1];
4088 if (pme->nnodes_major > 1)
4090 size_yx = pme->overlap[0].comm_data[0].send_nindex;
4096 datasize = (local_fft_ndata[XX] + size_yx)*local_fft_ndata[ZZ];
4098 send_size_y = overlap->send_size;
4100 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
4102 send_id = overlap->send_id[ipulse];
4103 recv_id = overlap->recv_id[ipulse];
4105 overlap->comm_data[ipulse].send_index0 -
4106 overlap->comm_data[0].send_index0;
4107 send_nindex = overlap->comm_data[ipulse].send_nindex;
4108 /* We don't use recv_index0, as we always receive starting at 0 */
4109 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
4110 recv_size_y = overlap->comm_data[ipulse].recv_size;
4112 sendptr = overlap->sendbuf + send_index0*local_fft_ndata[ZZ];
4113 recvptr = overlap->recvbuf;
4116 MPI_Sendrecv(sendptr, send_size_y*datasize, GMX_MPI_REAL,
4118 recvptr, recv_size_y*datasize, GMX_MPI_REAL,
4120 overlap->mpi_comm, &stat);
4123 for (x = 0; x < local_fft_ndata[XX]; x++)
4125 for (y = 0; y < recv_nindex; y++)
4127 indg = (x*local_fft_size[YY] + y)*local_fft_size[ZZ];
4128 indb = (x*recv_size_y + y)*local_fft_ndata[ZZ];
4129 for (z = 0; z < local_fft_ndata[ZZ]; z++)
4131 fftgrid[indg+z] += recvptr[indb+z];
4136 if (pme->nnodes_major > 1)
4138 /* Copy from the received buffer to the send buffer for dim 0 */
4139 sendptr = pme->overlap[0].sendbuf;
4140 for (x = 0; x < size_yx; x++)
4142 for (y = 0; y < recv_nindex; y++)
4144 indg = (x*local_fft_ndata[YY] + y)*local_fft_ndata[ZZ];
4145 indb = ((local_fft_ndata[XX] + x)*recv_size_y + y)*local_fft_ndata[ZZ];
4146 for (z = 0; z < local_fft_ndata[ZZ]; z++)
4148 sendptr[indg+z] += recvptr[indb+z];
4156 /* We only support a single pulse here.
4157 * This is not a severe limitation, as this code is only used
4158 * with OpenMP and with OpenMP the (PME) domains can be larger.
4160 if (pme->nnodes_major > 1)
4162 /* Major dimension */
4163 overlap = &pme->overlap[0];
4165 datasize = local_fft_ndata[YY]*local_fft_ndata[ZZ];
4166 gridsize = local_fft_size[YY] *local_fft_size[ZZ];
4170 send_id = overlap->send_id[ipulse];
4171 recv_id = overlap->recv_id[ipulse];
4172 send_nindex = overlap->comm_data[ipulse].send_nindex;
4173 /* We don't use recv_index0, as we always receive starting at 0 */
4174 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
4176 sendptr = overlap->sendbuf;
4177 recvptr = overlap->recvbuf;
4181 fprintf(debug, "PME fftgrid comm %2d x %2d x %2d\n",
4182 send_nindex, local_fft_ndata[YY], local_fft_ndata[ZZ]);
4186 MPI_Sendrecv(sendptr, send_nindex*datasize, GMX_MPI_REAL,
4188 recvptr, recv_nindex*datasize, GMX_MPI_REAL,
4190 overlap->mpi_comm, &stat);
4193 for (x = 0; x < recv_nindex; x++)
4195 for (y = 0; y < local_fft_ndata[YY]; y++)
4197 indg = (x*local_fft_size[YY] + y)*local_fft_size[ZZ];
4198 indb = (x*local_fft_ndata[YY] + y)*local_fft_ndata[ZZ];
4199 for (z = 0; z < local_fft_ndata[ZZ]; z++)
4201 fftgrid[indg+z] += recvptr[indb+z];
4209 static void spread_on_grid(gmx_pme_t pme,
4210 pme_atomcomm_t *atc, pmegrids_t *grids,
4211 gmx_bool bCalcSplines, gmx_bool bSpread,
4212 real *fftgrid, gmx_bool bDoSplines, int grid_index)
4214 int nthread, thread;
4215 #ifdef PME_TIME_THREADS
4216 gmx_cycles_t c1, c2, c3, ct1a, ct1b, ct1c;
4217 static double cs1 = 0, cs2 = 0, cs3 = 0;
4218 static double cs1a[6] = {0, 0, 0, 0, 0, 0};
4222 nthread = pme->nthread;
4223 assert(nthread > 0);
4225 #ifdef PME_TIME_THREADS
4226 c1 = omp_cyc_start();
4230 #pragma omp parallel for num_threads(nthread) schedule(static)
4231 for (thread = 0; thread < nthread; thread++)
4235 start = atc->n* thread /nthread;
4236 end = atc->n*(thread+1)/nthread;
4238 /* Compute fftgrid index for all atoms,
4239 * with help of some extra variables.
4241 calc_interpolation_idx(pme, atc, start, grid_index, end, thread);
4244 #ifdef PME_TIME_THREADS
4245 c1 = omp_cyc_end(c1);
4249 #ifdef PME_TIME_THREADS
4250 c2 = omp_cyc_start();
4252 #pragma omp parallel for num_threads(nthread) schedule(static)
4253 for (thread = 0; thread < nthread; thread++)
4255 splinedata_t *spline;
4256 pmegrid_t *grid = NULL;
4258 /* make local bsplines */
4259 if (grids == NULL || !pme->bUseThreads)
4261 spline = &atc->spline[0];
4267 grid = &grids->grid;
4272 spline = &atc->spline[thread];
4274 if (grids->nthread == 1)
4276 /* One thread, we operate on all coefficients */
4281 /* Get the indices our thread should operate on */
4282 make_thread_local_ind(atc, thread, spline);
4285 grid = &grids->grid_th[thread];
4290 make_bsplines(spline->theta, spline->dtheta, pme->pme_order,
4291 atc->fractx, spline->n, spline->ind, atc->coefficient, bDoSplines);
4296 /* put local atoms on grid. */
4297 #ifdef PME_TIME_SPREAD
4298 ct1a = omp_cyc_start();
4300 spread_coefficients_bsplines_thread(grid, atc, spline, pme->spline_work);
4302 if (pme->bUseThreads)
4304 copy_local_grid(pme, grids, grid_index, thread, fftgrid);
4306 #ifdef PME_TIME_SPREAD
4307 ct1a = omp_cyc_end(ct1a);
4308 cs1a[thread] += (double)ct1a;
4312 #ifdef PME_TIME_THREADS
4313 c2 = omp_cyc_end(c2);
4317 if (bSpread && pme->bUseThreads)
4319 #ifdef PME_TIME_THREADS
4320 c3 = omp_cyc_start();
4322 #pragma omp parallel for num_threads(grids->nthread) schedule(static)
4323 for (thread = 0; thread < grids->nthread; thread++)
4325 reduce_threadgrid_overlap(pme, grids, thread,
4327 pme->overlap[0].sendbuf,
4328 pme->overlap[1].sendbuf,
4331 #ifdef PME_TIME_THREADS
4332 c3 = omp_cyc_end(c3);
4336 if (pme->nnodes > 1)
4338 /* Communicate the overlapping part of the fftgrid.
4339 * For this communication call we need to check pme->bUseThreads
4340 * to have all ranks communicate here, regardless of pme->nthread.
4342 sum_fftgrid_dd(pme, fftgrid, grid_index);
4346 #ifdef PME_TIME_THREADS
4350 printf("idx %.2f spread %.2f red %.2f",
4351 cs1*1e-9, cs2*1e-9, cs3*1e-9);
4352 #ifdef PME_TIME_SPREAD
4353 for (thread = 0; thread < nthread; thread++)
4355 printf(" %.2f", cs1a[thread]*1e-9);
4364 static void dump_grid(FILE *fp,
4365 int sx, int sy, int sz, int nx, int ny, int nz,
4366 int my, int mz, const real *g)
4370 for (x = 0; x < nx; x++)
4372 for (y = 0; y < ny; y++)
4374 for (z = 0; z < nz; z++)
4376 fprintf(fp, "%2d %2d %2d %6.3f\n",
4377 sx+x, sy+y, sz+z, g[(x*my + y)*mz + z]);
4383 static void dump_local_fftgrid(gmx_pme_t pme, const real *fftgrid)
4385 ivec local_fft_ndata, local_fft_offset, local_fft_size;
4387 gmx_parallel_3dfft_real_limits(pme->pfft_setup[PME_GRID_QA],
4393 pme->pmegrid_start_ix,
4394 pme->pmegrid_start_iy,
4395 pme->pmegrid_start_iz,
4396 pme->pmegrid_nx-pme->pme_order+1,
4397 pme->pmegrid_ny-pme->pme_order+1,
4398 pme->pmegrid_nz-pme->pme_order+1,
4405 void gmx_pme_calc_energy(gmx_pme_t pme, int n, rvec *x, real *q, real *V)
4407 pme_atomcomm_t *atc;
4410 if (pme->nnodes > 1)
4412 gmx_incons("gmx_pme_calc_energy called in parallel");
4414 if (pme->bFEP_q > 1)
4416 gmx_incons("gmx_pme_calc_energy with free energy");
4419 atc = &pme->atc_energy;
4421 if (atc->spline == NULL)
4423 snew(atc->spline, atc->nthread);
4426 atc->bSpread = TRUE;
4427 atc->pme_order = pme->pme_order;
4429 pme_realloc_atomcomm_things(atc);
4431 atc->coefficient = q;
4433 /* We only use the A-charges grid */
4434 grid = &pme->pmegrid[PME_GRID_QA];
4436 /* Only calculate the spline coefficients, don't actually spread */
4437 spread_on_grid(pme, atc, NULL, TRUE, FALSE, pme->fftgrid[PME_GRID_QA], FALSE, PME_GRID_QA);
4439 *V = gather_energy_bsplines(pme, grid->grid.grid, atc);
4443 static void reset_pmeonly_counters(gmx_wallcycle_t wcycle,
4444 gmx_walltime_accounting_t walltime_accounting,
4445 t_nrnb *nrnb, t_inputrec *ir,
4448 /* Reset all the counters related to performance over the run */
4449 wallcycle_stop(wcycle, ewcRUN);
4450 wallcycle_reset_all(wcycle);
4452 if (ir->nsteps >= 0)
4454 /* ir->nsteps is not used here, but we update it for consistency */
4455 ir->nsteps -= step - ir->init_step;
4457 ir->init_step = step;
4458 wallcycle_start(wcycle, ewcRUN);
4459 walltime_accounting_start(walltime_accounting);
4463 static void gmx_pmeonly_switch(int *npmedata, gmx_pme_t **pmedata,
4465 t_commrec *cr, t_inputrec *ir,
4469 gmx_pme_t pme = NULL;
4472 while (ind < *npmedata)
4474 pme = (*pmedata)[ind];
4475 if (pme->nkx == grid_size[XX] &&
4476 pme->nky == grid_size[YY] &&
4477 pme->nkz == grid_size[ZZ])
4488 srenew(*pmedata, *npmedata);
4490 /* Generate a new PME data structure, copying part of the old pointers */
4491 gmx_pme_reinit(&((*pmedata)[ind]), cr, pme, ir, grid_size);
4493 *pme_ret = (*pmedata)[ind];
4496 int gmx_pmeonly(gmx_pme_t pme,
4497 t_commrec *cr, t_nrnb *mynrnb,
4498 gmx_wallcycle_t wcycle,
4499 gmx_walltime_accounting_t walltime_accounting,
4500 real ewaldcoeff_q, real ewaldcoeff_lj,
4505 gmx_pme_pp_t pme_pp;
4509 rvec *x_pp = NULL, *f_pp = NULL;
4510 real *chargeA = NULL, *chargeB = NULL;
4511 real *c6A = NULL, *c6B = NULL;
4512 real *sigmaA = NULL, *sigmaB = NULL;
4515 int maxshift_x = 0, maxshift_y = 0;
4516 real energy_q, energy_lj, dvdlambda_q, dvdlambda_lj;
4517 matrix vir_q, vir_lj;
4522 gmx_int64_t step, step_rel;
4525 /* This data will only use with PME tuning, i.e. switching PME grids */
4527 snew(pmedata, npmedata);
4530 pme_pp = gmx_pme_pp_init(cr);
4535 do /****** this is a quasi-loop over time steps! */
4537 /* The reason for having a loop here is PME grid tuning/switching */
4540 /* Domain decomposition */
4541 ret = gmx_pme_recv_coeffs_coords(pme_pp,
4547 &maxshift_x, &maxshift_y,
4548 &pme->bFEP_q, &pme->bFEP_lj,
4549 &lambda_q, &lambda_lj,
4553 grid_switch, &ewaldcoeff_q, &ewaldcoeff_lj);
4555 if (ret == pmerecvqxSWITCHGRID)
4557 /* Switch the PME grid to grid_switch */
4558 gmx_pmeonly_switch(&npmedata, &pmedata, grid_switch, cr, ir, &pme);
4561 if (ret == pmerecvqxRESETCOUNTERS)
4563 /* Reset the cycle and flop counters */
4564 reset_pmeonly_counters(wcycle, walltime_accounting, mynrnb, ir, step);
4567 while (ret == pmerecvqxSWITCHGRID || ret == pmerecvqxRESETCOUNTERS);
4569 if (ret == pmerecvqxFINISH)
4571 /* We should stop: break out of the loop */
4575 step_rel = step - ir->init_step;
4579 wallcycle_start(wcycle, ewcRUN);
4580 walltime_accounting_start(walltime_accounting);
4583 wallcycle_start(wcycle, ewcPMEMESH);
4590 gmx_pme_do(pme, 0, natoms, x_pp, f_pp,
4591 chargeA, chargeB, c6A, c6B, sigmaA, sigmaB, box,
4592 cr, maxshift_x, maxshift_y, mynrnb, wcycle,
4593 vir_q, ewaldcoeff_q, vir_lj, ewaldcoeff_lj,
4594 &energy_q, &energy_lj, lambda_q, lambda_lj, &dvdlambda_q, &dvdlambda_lj,
4595 pme_flags | GMX_PME_DO_ALL_F | (bEnerVir ? GMX_PME_CALC_ENER_VIR : 0));
4597 cycles = wallcycle_stop(wcycle, ewcPMEMESH);
4599 gmx_pme_send_force_vir_ener(pme_pp,
4600 f_pp, vir_q, energy_q, vir_lj, energy_lj,
4601 dvdlambda_q, dvdlambda_lj, cycles);
4604 } /***** end of quasi-loop, we stop with the break above */
4607 walltime_accounting_end(walltime_accounting);
4613 calc_initial_lb_coeffs(gmx_pme_t pme, real *local_c6, real *local_sigma)
4617 for (i = 0; i < pme->atc[0].n; ++i)
4621 sigma4 = local_sigma[i];
4622 sigma4 = sigma4*sigma4;
4623 sigma4 = sigma4*sigma4;
4624 pme->atc[0].coefficient[i] = local_c6[i] / sigma4;
4629 calc_next_lb_coeffs(gmx_pme_t pme, real *local_sigma)
4633 for (i = 0; i < pme->atc[0].n; ++i)
4635 pme->atc[0].coefficient[i] *= local_sigma[i];
4640 do_redist_pos_coeffs(gmx_pme_t pme, t_commrec *cr, int start, int homenr,
4641 gmx_bool bFirst, rvec x[], real *data)
4644 pme_atomcomm_t *atc;
4647 for (d = pme->ndecompdim - 1; d >= 0; d--)
4653 if (d == pme->ndecompdim - 1)
4661 n_d = pme->atc[d + 1].n;
4663 param_d = atc->coefficient;
4667 if (atc->npd > atc->pd_nalloc)
4669 atc->pd_nalloc = over_alloc_dd(atc->npd);
4670 srenew(atc->pd, atc->pd_nalloc);
4672 pme_calc_pidx_wrapper(n_d, pme->recipbox, x_d, atc);
4674 /* Redistribute x (only once) and qA/c6A or qB/c6B */
4675 if (DOMAINDECOMP(cr))
4677 dd_pmeredist_pos_coeffs(pme, n_d, bFirst, x_d, param_d, atc);
4682 int gmx_pme_do(gmx_pme_t pme,
4683 int start, int homenr,
4685 real *chargeA, real *chargeB,
4686 real *c6A, real *c6B,
4687 real *sigmaA, real *sigmaB,
4688 matrix box, t_commrec *cr,
4689 int maxshift_x, int maxshift_y,
4690 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
4691 matrix vir_q, real ewaldcoeff_q,
4692 matrix vir_lj, real ewaldcoeff_lj,
4693 real *energy_q, real *energy_lj,
4694 real lambda_q, real lambda_lj,
4695 real *dvdlambda_q, real *dvdlambda_lj,
4698 int d, i, j, k, ntot, npme, grid_index, max_grid_index;
4701 pme_atomcomm_t *atc = NULL;
4702 pmegrids_t *pmegrid = NULL;
4706 real *coefficient = NULL;
4711 gmx_parallel_3dfft_t pfft_setup;
4713 t_complex * cfftgrid;
4715 gmx_bool bFirst, bDoSplines;
4717 int fep_states_lj = pme->bFEP_lj ? 2 : 1;
4718 const gmx_bool bCalcEnerVir = flags & GMX_PME_CALC_ENER_VIR;
4719 const gmx_bool bCalcF = flags & GMX_PME_CALC_F;
4721 assert(pme->nnodes > 0);
4722 assert(pme->nnodes == 1 || pme->ndecompdim > 0);
4724 if (pme->nnodes > 1)
4728 if (atc->npd > atc->pd_nalloc)
4730 atc->pd_nalloc = over_alloc_dd(atc->npd);
4731 srenew(atc->pd, atc->pd_nalloc);
4733 for (d = pme->ndecompdim-1; d >= 0; d--)
4736 atc->maxshift = (atc->dimind == 0 ? maxshift_x : maxshift_y);
4742 /* This could be necessary for TPI */
4743 pme->atc[0].n = homenr;
4744 if (DOMAINDECOMP(cr))
4746 pme_realloc_atomcomm_things(atc);
4752 m_inv_ur0(box, pme->recipbox);
4755 /* For simplicity, we construct the splines for all particles if
4756 * more than one PME calculations is needed. Some optimization
4757 * could be done by keeping track of which atoms have splines
4758 * constructed, and construct new splines on each pass for atoms
4759 * that don't yet have them.
4762 bDoSplines = pme->bFEP || ((flags & GMX_PME_DO_COULOMB) && (flags & GMX_PME_DO_LJ));
4764 /* We need a maximum of four separate PME calculations:
4765 * grid_index=0: Coulomb PME with charges from state A
4766 * grid_index=1: Coulomb PME with charges from state B
4767 * grid_index=2: LJ PME with C6 from state A
4768 * grid_index=3: LJ PME with C6 from state B
4769 * For Lorentz-Berthelot combination rules, a separate loop is used to
4770 * calculate all the terms
4773 /* If we are doing LJ-PME with LB, we only do Q here */
4774 max_grid_index = (pme->ljpme_combination_rule == eljpmeLB) ? DO_Q : DO_Q_AND_LJ;
4776 for (grid_index = 0; grid_index < max_grid_index; ++grid_index)
4778 /* Check if we should do calculations at this grid_index
4779 * If grid_index is odd we should be doing FEP
4780 * If grid_index < 2 we should be doing electrostatic PME
4781 * If grid_index >= 2 we should be doing LJ-PME
4783 if ((grid_index < DO_Q && (!(flags & GMX_PME_DO_COULOMB) ||
4784 (grid_index == 1 && !pme->bFEP_q))) ||
4785 (grid_index >= DO_Q && (!(flags & GMX_PME_DO_LJ) ||
4786 (grid_index == 3 && !pme->bFEP_lj))))
4790 /* Unpack structure */
4791 pmegrid = &pme->pmegrid[grid_index];
4792 fftgrid = pme->fftgrid[grid_index];
4793 cfftgrid = pme->cfftgrid[grid_index];
4794 pfft_setup = pme->pfft_setup[grid_index];
4797 case 0: coefficient = chargeA + start; break;
4798 case 1: coefficient = chargeB + start; break;
4799 case 2: coefficient = c6A + start; break;
4800 case 3: coefficient = c6B + start; break;
4803 grid = pmegrid->grid.grid;
4807 fprintf(debug, "PME: nnodes = %d, nodeid = %d\n",
4808 cr->nnodes, cr->nodeid);
4809 fprintf(debug, "Grid = %p\n", (void*)grid);
4812 gmx_fatal(FARGS, "No grid!");
4817 if (pme->nnodes == 1)
4819 atc->coefficient = coefficient;
4823 wallcycle_start(wcycle, ewcPME_REDISTXF);
4824 do_redist_pos_coeffs(pme, cr, start, homenr, bFirst, x, coefficient);
4827 wallcycle_stop(wcycle, ewcPME_REDISTXF);
4832 fprintf(debug, "Node= %6d, pme local particles=%6d\n",
4833 cr->nodeid, atc->n);
4836 if (flags & GMX_PME_SPREAD)
4838 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
4840 /* Spread the coefficients on a grid */
4841 spread_on_grid(pme, &pme->atc[0], pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
4845 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
4847 inc_nrnb(nrnb, eNR_SPREADBSP,
4848 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
4850 if (!pme->bUseThreads)
4852 wrap_periodic_pmegrid(pme, grid);
4854 /* sum contributions to local grid from other nodes */
4856 if (pme->nnodes > 1)
4858 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
4863 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
4866 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
4869 dump_local_fftgrid(pme,fftgrid);
4874 /* Here we start a large thread parallel region */
4875 #pragma omp parallel num_threads(pme->nthread) private(thread)
4877 thread = gmx_omp_get_thread_num();
4878 if (flags & GMX_PME_SOLVE)
4885 wallcycle_start(wcycle, ewcPME_FFT);
4887 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
4891 wallcycle_stop(wcycle, ewcPME_FFT);
4895 /* solve in k-space for our local cells */
4898 wallcycle_start(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
4900 if (grid_index < DO_Q)
4903 solve_pme_yzx(pme, cfftgrid, ewaldcoeff_q,
4904 box[XX][XX]*box[YY][YY]*box[ZZ][ZZ],
4906 pme->nthread, thread);
4911 solve_pme_lj_yzx(pme, &cfftgrid, FALSE, ewaldcoeff_lj,
4912 box[XX][XX]*box[YY][YY]*box[ZZ][ZZ],
4914 pme->nthread, thread);
4919 wallcycle_stop(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
4921 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
4931 wallcycle_start(wcycle, ewcPME_FFT);
4933 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
4937 wallcycle_stop(wcycle, ewcPME_FFT);
4941 if (pme->nodeid == 0)
4943 ntot = pme->nkx*pme->nky*pme->nkz;
4944 npme = ntot*log((real)ntot)/log(2.0);
4945 inc_nrnb(nrnb, eNR_FFT, 2*npme);
4948 /* Note: this wallcycle region is closed below
4949 outside an OpenMP region, so take care if
4950 refactoring code here. */
4951 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
4954 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
4957 /* End of thread parallel section.
4958 * With MPI we have to synchronize here before gmx_sum_qgrid_dd.
4963 /* distribute local grid to all nodes */
4965 if (pme->nnodes > 1)
4967 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
4972 unwrap_periodic_pmegrid(pme, grid);
4974 /* interpolate forces for our local atoms */
4978 /* If we are running without parallelization,
4979 * atc->f is the actual force array, not a buffer,
4980 * therefore we should not clear it.
4982 lambda = grid_index < DO_Q ? lambda_q : lambda_lj;
4983 bClearF = (bFirst && PAR(cr));
4984 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
4985 for (thread = 0; thread < pme->nthread; thread++)
4987 gather_f_bsplines(pme, grid, bClearF, atc,
4988 &atc->spline[thread],
4989 pme->bFEP ? (grid_index % 2 == 0 ? 1.0-lambda : lambda) : 1.0);
4994 inc_nrnb(nrnb, eNR_GATHERFBSP,
4995 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
4996 /* Note: this wallcycle region is opened above inside an OpenMP
4997 region, so take care if refactoring code here. */
4998 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
5003 /* This should only be called on the master thread
5004 * and after the threads have synchronized.
5008 get_pme_ener_vir_q(pme, pme->nthread, &energy_AB[grid_index], vir_AB[grid_index]);
5012 get_pme_ener_vir_lj(pme, pme->nthread, &energy_AB[grid_index], vir_AB[grid_index]);
5016 } /* of grid_index-loop */
5018 /* For Lorentz-Berthelot combination rules in LJ-PME, we need to calculate
5021 if ((flags & GMX_PME_DO_LJ) && pme->ljpme_combination_rule == eljpmeLB)
5023 /* Loop over A- and B-state if we are doing FEP */
5024 for (fep_state = 0; fep_state < fep_states_lj; ++fep_state)
5026 real *local_c6 = NULL, *local_sigma = NULL, *RedistC6 = NULL, *RedistSigma = NULL;
5027 if (pme->nnodes == 1)
5029 if (pme->lb_buf1 == NULL)
5031 pme->lb_buf_nalloc = pme->atc[0].n;
5032 snew(pme->lb_buf1, pme->lb_buf_nalloc);
5034 pme->atc[0].coefficient = pme->lb_buf1;
5039 local_sigma = sigmaA;
5043 local_sigma = sigmaB;
5046 gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
5056 RedistSigma = sigmaA;
5060 RedistSigma = sigmaB;
5063 gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
5065 wallcycle_start(wcycle, ewcPME_REDISTXF);
5067 do_redist_pos_coeffs(pme, cr, start, homenr, bFirst, x, RedistC6);
5068 if (pme->lb_buf_nalloc < atc->n)
5070 pme->lb_buf_nalloc = atc->nalloc;
5071 srenew(pme->lb_buf1, pme->lb_buf_nalloc);
5072 srenew(pme->lb_buf2, pme->lb_buf_nalloc);
5074 local_c6 = pme->lb_buf1;
5075 for (i = 0; i < atc->n; ++i)
5077 local_c6[i] = atc->coefficient[i];
5081 do_redist_pos_coeffs(pme, cr, start, homenr, FALSE, x, RedistSigma);
5082 local_sigma = pme->lb_buf2;
5083 for (i = 0; i < atc->n; ++i)
5085 local_sigma[i] = atc->coefficient[i];
5089 wallcycle_stop(wcycle, ewcPME_REDISTXF);
5091 calc_initial_lb_coeffs(pme, local_c6, local_sigma);
5093 /*Seven terms in LJ-PME with LB, grid_index < 2 reserved for electrostatics*/
5094 for (grid_index = 2; grid_index < 9; ++grid_index)
5096 /* Unpack structure */
5097 pmegrid = &pme->pmegrid[grid_index];
5098 fftgrid = pme->fftgrid[grid_index];
5099 cfftgrid = pme->cfftgrid[grid_index];
5100 pfft_setup = pme->pfft_setup[grid_index];
5101 calc_next_lb_coeffs(pme, local_sigma);
5102 grid = pmegrid->grid.grid;
5105 if (flags & GMX_PME_SPREAD)
5107 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
5108 /* Spread the c6 on a grid */
5109 spread_on_grid(pme, &pme->atc[0], pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
5113 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
5116 inc_nrnb(nrnb, eNR_SPREADBSP,
5117 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
5118 if (pme->nthread == 1)
5120 wrap_periodic_pmegrid(pme, grid);
5121 /* sum contributions to local grid from other nodes */
5123 if (pme->nnodes > 1)
5125 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
5129 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
5131 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
5133 /*Here we start a large thread parallel region*/
5134 #pragma omp parallel num_threads(pme->nthread) private(thread)
5136 thread = gmx_omp_get_thread_num();
5137 if (flags & GMX_PME_SOLVE)
5142 wallcycle_start(wcycle, ewcPME_FFT);
5145 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
5149 wallcycle_stop(wcycle, ewcPME_FFT);
5156 if (flags & GMX_PME_SOLVE)
5158 /* solve in k-space for our local cells */
5159 #pragma omp parallel num_threads(pme->nthread) private(thread)
5162 thread = gmx_omp_get_thread_num();
5165 wallcycle_start(wcycle, ewcLJPME);
5169 solve_pme_lj_yzx(pme, &pme->cfftgrid[2], TRUE, ewaldcoeff_lj,
5170 box[XX][XX]*box[YY][YY]*box[ZZ][ZZ],
5172 pme->nthread, thread);
5175 wallcycle_stop(wcycle, ewcLJPME);
5177 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
5184 /* This should only be called on the master thread and
5185 * after the threads have synchronized.
5187 get_pme_ener_vir_lj(pme, pme->nthread, &energy_AB[2+fep_state], vir_AB[2+fep_state]);
5192 bFirst = !(flags & GMX_PME_DO_COULOMB);
5193 calc_initial_lb_coeffs(pme, local_c6, local_sigma);
5194 for (grid_index = 8; grid_index >= 2; --grid_index)
5196 /* Unpack structure */
5197 pmegrid = &pme->pmegrid[grid_index];
5198 fftgrid = pme->fftgrid[grid_index];
5199 cfftgrid = pme->cfftgrid[grid_index];
5200 pfft_setup = pme->pfft_setup[grid_index];
5201 grid = pmegrid->grid.grid;
5202 calc_next_lb_coeffs(pme, local_sigma);
5204 #pragma omp parallel num_threads(pme->nthread) private(thread)
5206 thread = gmx_omp_get_thread_num();
5211 wallcycle_start(wcycle, ewcPME_FFT);
5214 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
5218 wallcycle_stop(wcycle, ewcPME_FFT);
5222 if (pme->nodeid == 0)
5224 ntot = pme->nkx*pme->nky*pme->nkz;
5225 npme = ntot*log((real)ntot)/log(2.0);
5226 inc_nrnb(nrnb, eNR_FFT, 2*npme);
5228 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
5231 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
5233 } /*#pragma omp parallel*/
5235 /* distribute local grid to all nodes */
5237 if (pme->nnodes > 1)
5239 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
5244 unwrap_periodic_pmegrid(pme, grid);
5246 /* interpolate forces for our local atoms */
5248 bClearF = (bFirst && PAR(cr));
5249 scale = pme->bFEP ? (fep_state < 1 ? 1.0-lambda_lj : lambda_lj) : 1.0;
5250 scale *= lb_scale_factor[grid_index-2];
5251 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
5252 for (thread = 0; thread < pme->nthread; thread++)
5254 gather_f_bsplines(pme, grid, bClearF, &pme->atc[0],
5255 &pme->atc[0].spline[thread],
5260 inc_nrnb(nrnb, eNR_GATHERFBSP,
5261 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
5262 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
5265 } /* for (grid_index = 8; grid_index >= 2; --grid_index) */
5267 } /* for (fep_state = 0; fep_state < fep_states_lj; ++fep_state) */
5268 } /* if ((flags & GMX_PME_DO_LJ) && pme->ljpme_combination_rule == eljpmeLB) */
5270 if (bCalcF && pme->nnodes > 1)
5272 wallcycle_start(wcycle, ewcPME_REDISTXF);
5273 for (d = 0; d < pme->ndecompdim; d++)
5276 if (d == pme->ndecompdim - 1)
5283 n_d = pme->atc[d+1].n;
5284 f_d = pme->atc[d+1].f;
5286 if (DOMAINDECOMP(cr))
5288 dd_pmeredist_f(pme, atc, n_d, f_d,
5289 d == pme->ndecompdim-1 && pme->bPPnode);
5293 wallcycle_stop(wcycle, ewcPME_REDISTXF);
5299 if (flags & GMX_PME_DO_COULOMB)
5303 *energy_q = energy_AB[0];
5304 m_add(vir_q, vir_AB[0], vir_q);
5308 *energy_q = (1.0-lambda_q)*energy_AB[0] + lambda_q*energy_AB[1];
5309 *dvdlambda_q += energy_AB[1] - energy_AB[0];
5310 for (i = 0; i < DIM; i++)
5312 for (j = 0; j < DIM; j++)
5314 vir_q[i][j] += (1.0-lambda_q)*vir_AB[0][i][j] +
5315 lambda_q*vir_AB[1][i][j];
5321 fprintf(debug, "Electrostatic PME mesh energy: %g\n", *energy_q);
5329 if (flags & GMX_PME_DO_LJ)
5333 *energy_lj = energy_AB[2];
5334 m_add(vir_lj, vir_AB[2], vir_lj);
5338 *energy_lj = (1.0-lambda_lj)*energy_AB[2] + lambda_lj*energy_AB[3];
5339 *dvdlambda_lj += energy_AB[3] - energy_AB[2];
5340 for (i = 0; i < DIM; i++)
5342 for (j = 0; j < DIM; j++)
5344 vir_lj[i][j] += (1.0-lambda_lj)*vir_AB[2][i][j] + lambda_lj*vir_AB[3][i][j];
5350 fprintf(debug, "Lennard-Jones PME mesh energy: %g\n", *energy_lj);