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/ewald/pme-internal.h"
72 #include "gromacs/fft/fft.h"
73 #include "gromacs/fft/parallel_3dfft.h"
74 #include "gromacs/legacyheaders/macros.h"
75 #include "gromacs/legacyheaders/nrnb.h"
76 #include "gromacs/legacyheaders/types/commrec.h"
77 #include "gromacs/legacyheaders/types/enums.h"
78 #include "gromacs/legacyheaders/types/forcerec.h"
79 #include "gromacs/legacyheaders/types/inputrec.h"
80 #include "gromacs/legacyheaders/types/nrnb.h"
81 #include "gromacs/math/gmxcomplex.h"
82 #include "gromacs/math/units.h"
83 #include "gromacs/math/vec.h"
84 #include "gromacs/math/vectypes.h"
85 /* Include the SIMD macro file and then check for support */
86 #include "gromacs/simd/simd.h"
87 #include "gromacs/simd/simd_math.h"
88 #include "gromacs/timing/cyclecounter.h"
89 #include "gromacs/timing/wallcycle.h"
90 #include "gromacs/timing/walltime_accounting.h"
91 #include "gromacs/utility/basedefinitions.h"
92 #include "gromacs/utility/fatalerror.h"
93 #include "gromacs/utility/gmxmpi.h"
94 #include "gromacs/utility/gmxomp.h"
95 #include "gromacs/utility/real.h"
96 #include "gromacs/utility/smalloc.h"
99 #include "gromacs/fileio/pdbio.h"
100 #include "gromacs/utility/cstringutil.h"
101 #include "gromacs/utility/futil.h"
104 #ifdef GMX_SIMD_HAVE_REAL
105 /* Turn on arbitrary width SIMD intrinsics for PME solve */
106 # define PME_SIMD_SOLVE
109 #define PME_GRID_QA 0 /* Gridindex for A-state for Q */
110 #define PME_GRID_C6A 2 /* Gridindex for A-state for LJ */
111 #define DO_Q 2 /* Electrostatic grids have index q<2 */
112 #define DO_Q_AND_LJ 4 /* non-LB LJ grids have index 2 <= q < 4 */
113 #define DO_Q_AND_LJ_LB 9 /* With LB rules we need a total of 2+7 grids */
115 /* Pascal triangle coefficients scaled with (1/2)^6 for LJ-PME with LB-rules */
116 const real lb_scale_factor[] = {
117 1.0/64, 6.0/64, 15.0/64, 20.0/64,
118 15.0/64, 6.0/64, 1.0/64
121 /* Pascal triangle coefficients used in solve_pme_lj_yzx, only need to do 4 calculations due to symmetry */
122 const real lb_scale_factor_symm[] = { 2.0/64, 12.0/64, 30.0/64, 20.0/64 };
124 /* Check if we have 4-wide SIMD macro support */
125 #if (defined GMX_SIMD4_HAVE_REAL)
126 /* Do PME spread and gather with 4-wide SIMD.
127 * NOTE: SIMD is only used with PME order 4 and 5 (which are the most common).
129 # define PME_SIMD4_SPREAD_GATHER
131 # if (defined GMX_SIMD_HAVE_LOADU) && (defined GMX_SIMD_HAVE_STOREU)
132 /* With PME-order=4 on x86, unaligned load+store is slightly faster
133 * than doubling all SIMD operations when using aligned load+store.
135 # define PME_SIMD4_UNALIGNED
140 /* #define PRT_FORCE */
141 /* conditions for on the fly time-measurement */
142 /* #define TAKETIME (step > 1 && timesteps < 10) */
143 #define TAKETIME FALSE
145 /* #define PME_TIME_THREADS */
148 #define mpi_type MPI_DOUBLE
150 #define mpi_type MPI_FLOAT
153 #ifdef PME_SIMD4_SPREAD_GATHER
154 # define SIMD4_ALIGNMENT (GMX_SIMD4_WIDTH*sizeof(real))
156 /* We can use any alignment, apart from 0, so we use 4 reals */
157 # define SIMD4_ALIGNMENT (4*sizeof(real))
160 /* GMX_CACHE_SEP should be a multiple of the SIMD and SIMD4 register size
161 * to preserve alignment.
163 #define GMX_CACHE_SEP 64
165 /* We only define a maximum to be able to use local arrays without allocation.
166 * An order larger than 12 should never be needed, even for test cases.
167 * If needed it can be changed here.
169 #define PME_ORDER_MAX 12
171 /* Internal datastructures */
177 int recv_size; /* Receive buffer width, used with OpenMP */
180 typedef real *splinevec[DIM];
190 int *send_id, *recv_id;
191 int send_size; /* Send buffer width, used with OpenMP */
192 pme_grid_comm_t *comm_data;
198 int *n; /* Cumulative counts of the number of particles per thread */
199 int nalloc; /* Allocation size of i */
200 int *i; /* Particle indices ordered on thread index (n) */
214 int dimind; /* The index of the dimension, 0=x, 1=y */
221 int *node_dest; /* The nodes to send x and q to with DD */
222 int *node_src; /* The nodes to receive x and q from with DD */
223 int *buf_index; /* Index for commnode into the buffers */
230 int *count; /* The number of atoms to send to each node */
232 int *rcount; /* The number of atoms to receive */
239 gmx_bool bSpread; /* These coordinates are used for spreading */
242 rvec *fractx; /* Fractional coordinate relative to
243 * the lower cell boundary
246 int *thread_idx; /* Which thread should spread which coefficient */
247 thread_plist_t *thread_plist;
248 splinedata_t *spline;
255 ivec ci; /* The spatial location of this grid */
256 ivec n; /* The used size of *grid, including order-1 */
257 ivec offset; /* The grid offset from the full node grid */
258 int order; /* PME spreading order */
259 ivec s; /* The allocated size of *grid, s >= n */
260 real *grid; /* The grid local thread, size n */
264 pmegrid_t grid; /* The full node grid (non thread-local) */
265 int nthread; /* The number of threads operating on this grid */
266 ivec nc; /* The local spatial decomposition over the threads */
267 pmegrid_t *grid_th; /* Array of grids for each thread */
268 real *grid_all; /* Allocated array for the grids in *grid_th */
269 int **g2t; /* The grid to thread index */
270 ivec nthread_comm; /* The number of threads to communicate with */
274 #ifdef PME_SIMD4_SPREAD_GATHER
275 /* Masks for 4-wide SIMD aligned spreading and gathering */
276 gmx_simd4_bool_t mask_S0[6], mask_S1[6];
278 int dummy; /* C89 requires that struct has at least one member */
283 /* work data for solve_pme */
302 typedef struct gmx_pme {
303 int ndecompdim; /* The number of decomposition dimensions */
304 int nodeid; /* Our nodeid in mpi->mpi_comm */
307 int nnodes; /* The number of nodes doing PME */
312 MPI_Comm mpi_comm_d[2]; /* Indexed on dimension, 0=x, 1=y */
314 MPI_Datatype rvec_mpi; /* the pme vector's MPI type */
317 gmx_bool bUseThreads; /* Does any of the PME ranks have nthread>1 ? */
318 int nthread; /* The number of threads doing PME on our rank */
320 gmx_bool bPPnode; /* Node also does particle-particle forces */
321 gmx_bool bFEP; /* Compute Free energy contribution */
324 int nkx, nky, nkz; /* Grid dimensions */
325 gmx_bool bP3M; /* Do P3M: optimize the influence function */
329 int ljpme_combination_rule; /* Type of combination rule in LJ-PME */
331 int ngrids; /* number of grids we maintain for pmegrid, (c)fftgrid and pfft_setups*/
333 pmegrids_t pmegrid[DO_Q_AND_LJ_LB]; /* Grids on which we do spreading/interpolation,
334 * includes overlap Grid indices are ordered as
336 * 0: Coloumb PME, state A
337 * 1: Coloumb PME, state B
339 * This can probably be done in a better way
340 * but this simple hack works for now
342 /* The PME coefficient spreading grid sizes/strides, includes pme_order-1 */
343 int pmegrid_nx, pmegrid_ny, pmegrid_nz;
344 /* pmegrid_nz might be larger than strictly necessary to ensure
345 * memory alignment, pmegrid_nz_base gives the real base size.
348 /* The local PME grid starting indices */
349 int pmegrid_start_ix, pmegrid_start_iy, pmegrid_start_iz;
351 /* Work data for spreading and gathering */
352 pme_spline_work_t *spline_work;
354 real **fftgrid; /* Grids for FFT. With 1D FFT decomposition this can be a pointer */
355 /* inside the interpolation grid, but separate for 2D PME decomp. */
356 int fftgrid_nx, fftgrid_ny, fftgrid_nz;
358 t_complex **cfftgrid; /* Grids for complex FFT data */
360 int cfftgrid_nx, cfftgrid_ny, cfftgrid_nz;
362 gmx_parallel_3dfft_t *pfft_setup;
364 int *nnx, *nny, *nnz;
365 real *fshx, *fshy, *fshz;
367 pme_atomcomm_t atc[2]; /* Indexed on decomposition index */
370 /* Buffers to store data for local atoms for L-B combination rule
371 * calculations in LJ-PME. lb_buf1 stores either the coefficients
372 * for spreading/gathering (in serial), or the C6 coefficient for
373 * local atoms (in parallel). lb_buf2 is only used in parallel,
374 * and stores the sigma values for local atoms. */
375 real *lb_buf1, *lb_buf2;
376 int lb_buf_nalloc; /* Allocation size for the above buffers. */
378 pme_overlap_t overlap[2]; /* Indexed on dimension, 0=x, 1=y */
380 pme_atomcomm_t atc_energy; /* Only for gmx_pme_calc_energy */
382 rvec *bufv; /* Communication buffer */
383 real *bufr; /* Communication buffer */
384 int buf_nalloc; /* The communication buffer size */
386 /* thread local work data for solve_pme */
389 /* Work data for sum_qgrid */
390 real * sum_qgrid_tmp;
391 real * sum_qgrid_dd_tmp;
394 static void calc_interpolation_idx(gmx_pme_t pme, pme_atomcomm_t *atc,
395 int start, int grid_index, int end, int thread)
398 int *idxptr, tix, tiy, tiz;
399 real *xptr, *fptr, tx, ty, tz;
400 real rxx, ryx, ryy, rzx, rzy, rzz;
402 int start_ix, start_iy, start_iz;
403 int *g2tx, *g2ty, *g2tz;
405 int *thread_idx = NULL;
406 thread_plist_t *tpl = NULL;
414 start_ix = pme->pmegrid_start_ix;
415 start_iy = pme->pmegrid_start_iy;
416 start_iz = pme->pmegrid_start_iz;
418 rxx = pme->recipbox[XX][XX];
419 ryx = pme->recipbox[YY][XX];
420 ryy = pme->recipbox[YY][YY];
421 rzx = pme->recipbox[ZZ][XX];
422 rzy = pme->recipbox[ZZ][YY];
423 rzz = pme->recipbox[ZZ][ZZ];
425 g2tx = pme->pmegrid[grid_index].g2t[XX];
426 g2ty = pme->pmegrid[grid_index].g2t[YY];
427 g2tz = pme->pmegrid[grid_index].g2t[ZZ];
429 bThreads = (atc->nthread > 1);
432 thread_idx = atc->thread_idx;
434 tpl = &atc->thread_plist[thread];
436 for (i = 0; i < atc->nthread; i++)
442 for (i = start; i < end; i++)
445 idxptr = atc->idx[i];
446 fptr = atc->fractx[i];
448 /* Fractional coordinates along box vectors, add 2.0 to make 100% sure we are positive for triclinic boxes */
449 tx = nx * ( xptr[XX] * rxx + xptr[YY] * ryx + xptr[ZZ] * rzx + 2.0 );
450 ty = ny * ( xptr[YY] * ryy + xptr[ZZ] * rzy + 2.0 );
451 tz = nz * ( xptr[ZZ] * rzz + 2.0 );
457 /* Because decomposition only occurs in x and y,
458 * we never have a fraction correction in z.
460 fptr[XX] = tx - tix + pme->fshx[tix];
461 fptr[YY] = ty - tiy + pme->fshy[tiy];
464 idxptr[XX] = pme->nnx[tix];
465 idxptr[YY] = pme->nny[tiy];
466 idxptr[ZZ] = pme->nnz[tiz];
469 range_check(idxptr[XX], 0, pme->pmegrid_nx);
470 range_check(idxptr[YY], 0, pme->pmegrid_ny);
471 range_check(idxptr[ZZ], 0, pme->pmegrid_nz);
476 thread_i = g2tx[idxptr[XX]] + g2ty[idxptr[YY]] + g2tz[idxptr[ZZ]];
477 thread_idx[i] = thread_i;
484 /* Make a list of particle indices sorted on thread */
486 /* Get the cumulative count */
487 for (i = 1; i < atc->nthread; i++)
489 tpl_n[i] += tpl_n[i-1];
491 /* The current implementation distributes particles equally
492 * over the threads, so we could actually allocate for that
493 * in pme_realloc_atomcomm_things.
495 if (tpl_n[atc->nthread-1] > tpl->nalloc)
497 tpl->nalloc = over_alloc_large(tpl_n[atc->nthread-1]);
498 srenew(tpl->i, tpl->nalloc);
500 /* Set tpl_n to the cumulative start */
501 for (i = atc->nthread-1; i >= 1; i--)
503 tpl_n[i] = tpl_n[i-1];
507 /* Fill our thread local array with indices sorted on thread */
508 for (i = start; i < end; i++)
510 tpl->i[tpl_n[atc->thread_idx[i]]++] = i;
512 /* Now tpl_n contains the cummulative count again */
516 static void make_thread_local_ind(pme_atomcomm_t *atc,
517 int thread, splinedata_t *spline)
519 int n, t, i, start, end;
522 /* Combine the indices made by each thread into one index */
526 for (t = 0; t < atc->nthread; t++)
528 tpl = &atc->thread_plist[t];
529 /* Copy our part (start - end) from the list of thread t */
532 start = tpl->n[thread-1];
534 end = tpl->n[thread];
535 for (i = start; i < end; i++)
537 spline->ind[n++] = tpl->i[i];
545 static void pme_calc_pidx(int start, int end,
546 matrix recipbox, rvec x[],
547 pme_atomcomm_t *atc, int *count)
552 real rxx, ryx, rzx, ryy, rzy;
555 /* Calculate PME task index (pidx) for each grid index.
556 * Here we always assign equally sized slabs to each node
557 * for load balancing reasons (the PME grid spacing is not used).
563 /* Reset the count */
564 for (i = 0; i < nslab; i++)
569 if (atc->dimind == 0)
571 rxx = recipbox[XX][XX];
572 ryx = recipbox[YY][XX];
573 rzx = recipbox[ZZ][XX];
574 /* Calculate the node index in x-dimension */
575 for (i = start; i < end; i++)
578 /* Fractional coordinates along box vectors */
579 s = nslab*(xptr[XX]*rxx + xptr[YY]*ryx + xptr[ZZ]*rzx);
580 si = (int)(s + 2*nslab) % nslab;
587 ryy = recipbox[YY][YY];
588 rzy = recipbox[ZZ][YY];
589 /* Calculate the node index in y-dimension */
590 for (i = start; i < end; i++)
593 /* Fractional coordinates along box vectors */
594 s = nslab*(xptr[YY]*ryy + xptr[ZZ]*rzy);
595 si = (int)(s + 2*nslab) % nslab;
602 static void pme_calc_pidx_wrapper(int natoms, matrix recipbox, rvec x[],
605 int nthread, thread, slab;
607 nthread = atc->nthread;
609 #pragma omp parallel for num_threads(nthread) schedule(static)
610 for (thread = 0; thread < nthread; thread++)
612 pme_calc_pidx(natoms* thread /nthread,
613 natoms*(thread+1)/nthread,
614 recipbox, x, atc, atc->count_thread[thread]);
616 /* Non-parallel reduction, since nslab is small */
618 for (thread = 1; thread < nthread; thread++)
620 for (slab = 0; slab < atc->nslab; slab++)
622 atc->count_thread[0][slab] += atc->count_thread[thread][slab];
627 static void realloc_splinevec(splinevec th, real **ptr_z, int nalloc)
629 const int padding = 4;
632 srenew(th[XX], nalloc);
633 srenew(th[YY], nalloc);
634 /* In z we add padding, this is only required for the aligned SIMD code */
635 sfree_aligned(*ptr_z);
636 snew_aligned(*ptr_z, nalloc+2*padding, SIMD4_ALIGNMENT);
637 th[ZZ] = *ptr_z + padding;
639 for (i = 0; i < padding; i++)
642 (*ptr_z)[padding+nalloc+i] = 0;
646 static void pme_realloc_splinedata(splinedata_t *spline, pme_atomcomm_t *atc)
650 srenew(spline->ind, atc->nalloc);
651 /* Initialize the index to identity so it works without threads */
652 for (i = 0; i < atc->nalloc; i++)
657 realloc_splinevec(spline->theta, &spline->ptr_theta_z,
658 atc->pme_order*atc->nalloc);
659 realloc_splinevec(spline->dtheta, &spline->ptr_dtheta_z,
660 atc->pme_order*atc->nalloc);
663 static void pme_realloc_atomcomm_things(pme_atomcomm_t *atc)
665 int nalloc_old, i, j, nalloc_tpl;
667 /* We have to avoid a NULL pointer for atc->x to avoid
668 * possible fatal errors in MPI routines.
670 if (atc->n > atc->nalloc || atc->nalloc == 0)
672 nalloc_old = atc->nalloc;
673 atc->nalloc = over_alloc_dd(max(atc->n, 1));
677 srenew(atc->x, atc->nalloc);
678 srenew(atc->coefficient, atc->nalloc);
679 srenew(atc->f, atc->nalloc);
680 for (i = nalloc_old; i < atc->nalloc; i++)
682 clear_rvec(atc->f[i]);
687 srenew(atc->fractx, atc->nalloc);
688 srenew(atc->idx, atc->nalloc);
690 if (atc->nthread > 1)
692 srenew(atc->thread_idx, atc->nalloc);
695 for (i = 0; i < atc->nthread; i++)
697 pme_realloc_splinedata(&atc->spline[i], atc);
703 static void pme_dd_sendrecv(pme_atomcomm_t gmx_unused *atc,
704 gmx_bool gmx_unused bBackward, int gmx_unused shift,
705 void gmx_unused *buf_s, int gmx_unused nbyte_s,
706 void gmx_unused *buf_r, int gmx_unused nbyte_r)
712 if (bBackward == FALSE)
714 dest = atc->node_dest[shift];
715 src = atc->node_src[shift];
719 dest = atc->node_src[shift];
720 src = atc->node_dest[shift];
723 if (nbyte_s > 0 && nbyte_r > 0)
725 MPI_Sendrecv(buf_s, nbyte_s, MPI_BYTE,
727 buf_r, nbyte_r, MPI_BYTE,
729 atc->mpi_comm, &stat);
731 else if (nbyte_s > 0)
733 MPI_Send(buf_s, nbyte_s, MPI_BYTE,
737 else if (nbyte_r > 0)
739 MPI_Recv(buf_r, nbyte_r, MPI_BYTE,
741 atc->mpi_comm, &stat);
746 static void dd_pmeredist_pos_coeffs(gmx_pme_t pme,
747 int n, gmx_bool bX, rvec *x, real *data,
750 int *commnode, *buf_index;
751 int nnodes_comm, i, nsend, local_pos, buf_pos, node, scount, rcount;
753 commnode = atc->node_dest;
754 buf_index = atc->buf_index;
756 nnodes_comm = min(2*atc->maxshift, atc->nslab-1);
759 for (i = 0; i < nnodes_comm; i++)
761 buf_index[commnode[i]] = nsend;
762 nsend += atc->count[commnode[i]];
766 if (atc->count[atc->nodeid] + nsend != n)
768 gmx_fatal(FARGS, "%d particles communicated to PME rank %d are more than 2/3 times the cut-off out of the domain decomposition cell of their charge group in dimension %c.\n"
769 "This usually means that your system is not well equilibrated.",
770 n - (atc->count[atc->nodeid] + nsend),
771 pme->nodeid, 'x'+atc->dimind);
774 if (nsend > pme->buf_nalloc)
776 pme->buf_nalloc = over_alloc_dd(nsend);
777 srenew(pme->bufv, pme->buf_nalloc);
778 srenew(pme->bufr, pme->buf_nalloc);
781 atc->n = atc->count[atc->nodeid];
782 for (i = 0; i < nnodes_comm; i++)
784 scount = atc->count[commnode[i]];
785 /* Communicate the count */
788 fprintf(debug, "dimind %d PME rank %d send to rank %d: %d\n",
789 atc->dimind, atc->nodeid, commnode[i], scount);
791 pme_dd_sendrecv(atc, FALSE, i,
792 &scount, sizeof(int),
793 &atc->rcount[i], sizeof(int));
794 atc->n += atc->rcount[i];
797 pme_realloc_atomcomm_things(atc);
801 for (i = 0; i < n; i++)
804 if (node == atc->nodeid)
806 /* Copy direct to the receive buffer */
809 copy_rvec(x[i], atc->x[local_pos]);
811 atc->coefficient[local_pos] = data[i];
816 /* Copy to the send buffer */
819 copy_rvec(x[i], pme->bufv[buf_index[node]]);
821 pme->bufr[buf_index[node]] = data[i];
827 for (i = 0; i < nnodes_comm; i++)
829 scount = atc->count[commnode[i]];
830 rcount = atc->rcount[i];
831 if (scount > 0 || rcount > 0)
835 /* Communicate the coordinates */
836 pme_dd_sendrecv(atc, FALSE, i,
837 pme->bufv[buf_pos], scount*sizeof(rvec),
838 atc->x[local_pos], rcount*sizeof(rvec));
840 /* Communicate the coefficients */
841 pme_dd_sendrecv(atc, FALSE, i,
842 pme->bufr+buf_pos, scount*sizeof(real),
843 atc->coefficient+local_pos, rcount*sizeof(real));
845 local_pos += atc->rcount[i];
850 static void dd_pmeredist_f(gmx_pme_t pme, pme_atomcomm_t *atc,
854 int *commnode, *buf_index;
855 int nnodes_comm, local_pos, buf_pos, i, scount, rcount, node;
857 commnode = atc->node_dest;
858 buf_index = atc->buf_index;
860 nnodes_comm = min(2*atc->maxshift, atc->nslab-1);
862 local_pos = atc->count[atc->nodeid];
864 for (i = 0; i < nnodes_comm; i++)
866 scount = atc->rcount[i];
867 rcount = atc->count[commnode[i]];
868 if (scount > 0 || rcount > 0)
870 /* Communicate the forces */
871 pme_dd_sendrecv(atc, TRUE, i,
872 atc->f[local_pos], scount*sizeof(rvec),
873 pme->bufv[buf_pos], rcount*sizeof(rvec));
876 buf_index[commnode[i]] = buf_pos;
883 for (i = 0; i < n; i++)
886 if (node == atc->nodeid)
888 /* Add from the local force array */
889 rvec_inc(f[i], atc->f[local_pos]);
894 /* Add from the receive buffer */
895 rvec_inc(f[i], pme->bufv[buf_index[node]]);
902 for (i = 0; i < n; i++)
905 if (node == atc->nodeid)
907 /* Copy from the local force array */
908 copy_rvec(atc->f[local_pos], f[i]);
913 /* Copy from the receive buffer */
914 copy_rvec(pme->bufv[buf_index[node]], f[i]);
922 static void gmx_sum_qgrid_dd(gmx_pme_t pme, real *grid, int direction)
924 pme_overlap_t *overlap;
925 int send_index0, send_nindex;
926 int recv_index0, recv_nindex;
928 int i, j, k, ix, iy, iz, icnt;
929 int ipulse, send_id, recv_id, datasize;
931 real *sendptr, *recvptr;
933 /* Start with minor-rank communication. This is a bit of a pain since it is not contiguous */
934 overlap = &pme->overlap[1];
936 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
938 /* Since we have already (un)wrapped the overlap in the z-dimension,
939 * we only have to communicate 0 to nkz (not pmegrid_nz).
941 if (direction == GMX_SUM_GRID_FORWARD)
943 send_id = overlap->send_id[ipulse];
944 recv_id = overlap->recv_id[ipulse];
945 send_index0 = overlap->comm_data[ipulse].send_index0;
946 send_nindex = overlap->comm_data[ipulse].send_nindex;
947 recv_index0 = overlap->comm_data[ipulse].recv_index0;
948 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
952 send_id = overlap->recv_id[ipulse];
953 recv_id = overlap->send_id[ipulse];
954 send_index0 = overlap->comm_data[ipulse].recv_index0;
955 send_nindex = overlap->comm_data[ipulse].recv_nindex;
956 recv_index0 = overlap->comm_data[ipulse].send_index0;
957 recv_nindex = overlap->comm_data[ipulse].send_nindex;
960 /* Copy data to contiguous send buffer */
963 fprintf(debug, "PME send rank %d %d -> %d grid start %d Communicating %d to %d\n",
964 pme->nodeid, overlap->nodeid, send_id,
965 pme->pmegrid_start_iy,
966 send_index0-pme->pmegrid_start_iy,
967 send_index0-pme->pmegrid_start_iy+send_nindex);
970 for (i = 0; i < pme->pmegrid_nx; i++)
973 for (j = 0; j < send_nindex; j++)
975 iy = j + send_index0 - pme->pmegrid_start_iy;
976 for (k = 0; k < pme->nkz; k++)
979 overlap->sendbuf[icnt++] = grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz];
984 datasize = pme->pmegrid_nx * pme->nkz;
986 MPI_Sendrecv(overlap->sendbuf, send_nindex*datasize, GMX_MPI_REAL,
988 overlap->recvbuf, recv_nindex*datasize, GMX_MPI_REAL,
990 overlap->mpi_comm, &stat);
992 /* Get data from contiguous recv buffer */
995 fprintf(debug, "PME recv rank %d %d <- %d grid start %d Communicating %d to %d\n",
996 pme->nodeid, overlap->nodeid, recv_id,
997 pme->pmegrid_start_iy,
998 recv_index0-pme->pmegrid_start_iy,
999 recv_index0-pme->pmegrid_start_iy+recv_nindex);
1002 for (i = 0; i < pme->pmegrid_nx; i++)
1005 for (j = 0; j < recv_nindex; j++)
1007 iy = j + recv_index0 - pme->pmegrid_start_iy;
1008 for (k = 0; k < pme->nkz; k++)
1011 if (direction == GMX_SUM_GRID_FORWARD)
1013 grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz] += overlap->recvbuf[icnt++];
1017 grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz] = overlap->recvbuf[icnt++];
1024 /* Major dimension is easier, no copying required,
1025 * but we might have to sum to separate array.
1026 * Since we don't copy, we have to communicate up to pmegrid_nz,
1027 * not nkz as for the minor direction.
1029 overlap = &pme->overlap[0];
1031 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
1033 if (direction == GMX_SUM_GRID_FORWARD)
1035 send_id = overlap->send_id[ipulse];
1036 recv_id = overlap->recv_id[ipulse];
1037 send_index0 = overlap->comm_data[ipulse].send_index0;
1038 send_nindex = overlap->comm_data[ipulse].send_nindex;
1039 recv_index0 = overlap->comm_data[ipulse].recv_index0;
1040 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
1041 recvptr = overlap->recvbuf;
1045 send_id = overlap->recv_id[ipulse];
1046 recv_id = overlap->send_id[ipulse];
1047 send_index0 = overlap->comm_data[ipulse].recv_index0;
1048 send_nindex = overlap->comm_data[ipulse].recv_nindex;
1049 recv_index0 = overlap->comm_data[ipulse].send_index0;
1050 recv_nindex = overlap->comm_data[ipulse].send_nindex;
1051 recvptr = grid + (recv_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1054 sendptr = grid + (send_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1055 datasize = pme->pmegrid_ny * pme->pmegrid_nz;
1059 fprintf(debug, "PME send rank %d %d -> %d grid start %d Communicating %d to %d\n",
1060 pme->nodeid, overlap->nodeid, send_id,
1061 pme->pmegrid_start_ix,
1062 send_index0-pme->pmegrid_start_ix,
1063 send_index0-pme->pmegrid_start_ix+send_nindex);
1064 fprintf(debug, "PME recv rank %d %d <- %d grid start %d Communicating %d to %d\n",
1065 pme->nodeid, overlap->nodeid, recv_id,
1066 pme->pmegrid_start_ix,
1067 recv_index0-pme->pmegrid_start_ix,
1068 recv_index0-pme->pmegrid_start_ix+recv_nindex);
1071 MPI_Sendrecv(sendptr, send_nindex*datasize, GMX_MPI_REAL,
1073 recvptr, recv_nindex*datasize, GMX_MPI_REAL,
1075 overlap->mpi_comm, &stat);
1077 /* ADD data from contiguous recv buffer */
1078 if (direction == GMX_SUM_GRID_FORWARD)
1080 p = grid + (recv_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1081 for (i = 0; i < recv_nindex*datasize; i++)
1083 p[i] += overlap->recvbuf[i];
1091 static int copy_pmegrid_to_fftgrid(gmx_pme_t pme, real *pmegrid, real *fftgrid, int grid_index)
1093 ivec local_fft_ndata, local_fft_offset, local_fft_size;
1094 ivec local_pme_size;
1098 /* Dimensions should be identical for A/B grid, so we just use A here */
1099 gmx_parallel_3dfft_real_limits(pme->pfft_setup[grid_index],
1104 local_pme_size[0] = pme->pmegrid_nx;
1105 local_pme_size[1] = pme->pmegrid_ny;
1106 local_pme_size[2] = pme->pmegrid_nz;
1108 /* The fftgrid is always 'justified' to the lower-left corner of the PME grid,
1109 the offset is identical, and the PME grid always has more data (due to overlap)
1116 sprintf(fn, "pmegrid%d.pdb", pme->nodeid);
1117 fp = gmx_ffopen(fn, "w");
1118 sprintf(fn, "pmegrid%d.txt", pme->nodeid);
1119 fp2 = gmx_ffopen(fn, "w");
1122 for (ix = 0; ix < local_fft_ndata[XX]; ix++)
1124 for (iy = 0; iy < local_fft_ndata[YY]; iy++)
1126 for (iz = 0; iz < local_fft_ndata[ZZ]; iz++)
1128 pmeidx = ix*(local_pme_size[YY]*local_pme_size[ZZ])+iy*(local_pme_size[ZZ])+iz;
1129 fftidx = ix*(local_fft_size[YY]*local_fft_size[ZZ])+iy*(local_fft_size[ZZ])+iz;
1130 fftgrid[fftidx] = pmegrid[pmeidx];
1132 val = 100*pmegrid[pmeidx];
1133 if (pmegrid[pmeidx] != 0)
1135 gmx_fprintf_pdb_atomline(fp, epdbATOM, pmeidx, "CA", ' ', "GLY", ' ', pmeidx, ' ',
1136 5.0*ix, 5.0*iy, 5.0*iz, 1.0, val, "");
1138 if (pmegrid[pmeidx] != 0)
1140 fprintf(fp2, "%-12s %5d %5d %5d %12.5e\n",
1142 pme->pmegrid_start_ix + ix,
1143 pme->pmegrid_start_iy + iy,
1144 pme->pmegrid_start_iz + iz,
1160 static gmx_cycles_t omp_cyc_start()
1162 return gmx_cycles_read();
1165 static gmx_cycles_t omp_cyc_end(gmx_cycles_t c)
1167 return gmx_cycles_read() - c;
1171 static int copy_fftgrid_to_pmegrid(gmx_pme_t pme, const real *fftgrid, real *pmegrid, int grid_index,
1172 int nthread, int thread)
1174 ivec local_fft_ndata, local_fft_offset, local_fft_size;
1175 ivec local_pme_size;
1176 int ixy0, ixy1, ixy, ix, iy, iz;
1178 #ifdef PME_TIME_THREADS
1180 static double cs1 = 0;
1184 #ifdef PME_TIME_THREADS
1185 c1 = omp_cyc_start();
1187 /* Dimensions should be identical for A/B grid, so we just use A here */
1188 gmx_parallel_3dfft_real_limits(pme->pfft_setup[grid_index],
1193 local_pme_size[0] = pme->pmegrid_nx;
1194 local_pme_size[1] = pme->pmegrid_ny;
1195 local_pme_size[2] = pme->pmegrid_nz;
1197 /* The fftgrid is always 'justified' to the lower-left corner of the PME grid,
1198 the offset is identical, and the PME grid always has more data (due to overlap)
1200 ixy0 = ((thread )*local_fft_ndata[XX]*local_fft_ndata[YY])/nthread;
1201 ixy1 = ((thread+1)*local_fft_ndata[XX]*local_fft_ndata[YY])/nthread;
1203 for (ixy = ixy0; ixy < ixy1; ixy++)
1205 ix = ixy/local_fft_ndata[YY];
1206 iy = ixy - ix*local_fft_ndata[YY];
1208 pmeidx = (ix*local_pme_size[YY] + iy)*local_pme_size[ZZ];
1209 fftidx = (ix*local_fft_size[YY] + iy)*local_fft_size[ZZ];
1210 for (iz = 0; iz < local_fft_ndata[ZZ]; iz++)
1212 pmegrid[pmeidx+iz] = fftgrid[fftidx+iz];
1216 #ifdef PME_TIME_THREADS
1217 c1 = omp_cyc_end(c1);
1222 printf("copy %.2f\n", cs1*1e-9);
1230 static void wrap_periodic_pmegrid(gmx_pme_t pme, real *pmegrid)
1232 int nx, ny, nz, pnx, pny, pnz, ny_x, overlap, ix, iy, iz;
1238 pnx = pme->pmegrid_nx;
1239 pny = pme->pmegrid_ny;
1240 pnz = pme->pmegrid_nz;
1242 overlap = pme->pme_order - 1;
1244 /* Add periodic overlap in z */
1245 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1247 for (iy = 0; iy < pme->pmegrid_ny; iy++)
1249 for (iz = 0; iz < overlap; iz++)
1251 pmegrid[(ix*pny+iy)*pnz+iz] +=
1252 pmegrid[(ix*pny+iy)*pnz+nz+iz];
1257 if (pme->nnodes_minor == 1)
1259 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1261 for (iy = 0; iy < overlap; iy++)
1263 for (iz = 0; iz < nz; iz++)
1265 pmegrid[(ix*pny+iy)*pnz+iz] +=
1266 pmegrid[(ix*pny+ny+iy)*pnz+iz];
1272 if (pme->nnodes_major == 1)
1274 ny_x = (pme->nnodes_minor == 1 ? ny : pme->pmegrid_ny);
1276 for (ix = 0; ix < overlap; ix++)
1278 for (iy = 0; iy < ny_x; iy++)
1280 for (iz = 0; iz < nz; iz++)
1282 pmegrid[(ix*pny+iy)*pnz+iz] +=
1283 pmegrid[((nx+ix)*pny+iy)*pnz+iz];
1291 static void unwrap_periodic_pmegrid(gmx_pme_t pme, real *pmegrid)
1293 int nx, ny, nz, pnx, pny, pnz, ny_x, overlap, ix;
1299 pnx = pme->pmegrid_nx;
1300 pny = pme->pmegrid_ny;
1301 pnz = pme->pmegrid_nz;
1303 overlap = pme->pme_order - 1;
1305 if (pme->nnodes_major == 1)
1307 ny_x = (pme->nnodes_minor == 1 ? ny : pme->pmegrid_ny);
1309 for (ix = 0; ix < overlap; ix++)
1313 for (iy = 0; iy < ny_x; iy++)
1315 for (iz = 0; iz < nz; iz++)
1317 pmegrid[((nx+ix)*pny+iy)*pnz+iz] =
1318 pmegrid[(ix*pny+iy)*pnz+iz];
1324 if (pme->nnodes_minor == 1)
1326 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1327 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1331 for (iy = 0; iy < overlap; iy++)
1333 for (iz = 0; iz < nz; iz++)
1335 pmegrid[(ix*pny+ny+iy)*pnz+iz] =
1336 pmegrid[(ix*pny+iy)*pnz+iz];
1342 /* Copy periodic overlap in z */
1343 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1344 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1348 for (iy = 0; iy < pme->pmegrid_ny; iy++)
1350 for (iz = 0; iz < overlap; iz++)
1352 pmegrid[(ix*pny+iy)*pnz+nz+iz] =
1353 pmegrid[(ix*pny+iy)*pnz+iz];
1360 /* This has to be a macro to enable full compiler optimization with xlC (and probably others too) */
1361 #define DO_BSPLINE(order) \
1362 for (ithx = 0; (ithx < order); ithx++) \
1364 index_x = (i0+ithx)*pny*pnz; \
1365 valx = coefficient*thx[ithx]; \
1367 for (ithy = 0; (ithy < order); ithy++) \
1369 valxy = valx*thy[ithy]; \
1370 index_xy = index_x+(j0+ithy)*pnz; \
1372 for (ithz = 0; (ithz < order); ithz++) \
1374 index_xyz = index_xy+(k0+ithz); \
1375 grid[index_xyz] += valxy*thz[ithz]; \
1381 static void spread_coefficients_bsplines_thread(pmegrid_t *pmegrid,
1382 pme_atomcomm_t *atc,
1383 splinedata_t *spline,
1384 pme_spline_work_t gmx_unused *work)
1387 /* spread coefficients from home atoms to local grid */
1390 int b, i, nn, n, ithx, ithy, ithz, i0, j0, k0;
1392 int order, norder, index_x, index_xy, index_xyz;
1393 real valx, valxy, coefficient;
1394 real *thx, *thy, *thz;
1395 int localsize, bndsize;
1396 int pnx, pny, pnz, ndatatot;
1397 int offx, offy, offz;
1399 #if defined PME_SIMD4_SPREAD_GATHER && !defined PME_SIMD4_UNALIGNED
1400 real thz_buffer[GMX_SIMD4_WIDTH*3], *thz_aligned;
1402 thz_aligned = gmx_simd4_align_r(thz_buffer);
1405 pnx = pmegrid->s[XX];
1406 pny = pmegrid->s[YY];
1407 pnz = pmegrid->s[ZZ];
1409 offx = pmegrid->offset[XX];
1410 offy = pmegrid->offset[YY];
1411 offz = pmegrid->offset[ZZ];
1413 ndatatot = pnx*pny*pnz;
1414 grid = pmegrid->grid;
1415 for (i = 0; i < ndatatot; i++)
1420 order = pmegrid->order;
1422 for (nn = 0; nn < spline->n; nn++)
1424 n = spline->ind[nn];
1425 coefficient = atc->coefficient[n];
1427 if (coefficient != 0)
1429 idxptr = atc->idx[n];
1432 i0 = idxptr[XX] - offx;
1433 j0 = idxptr[YY] - offy;
1434 k0 = idxptr[ZZ] - offz;
1436 thx = spline->theta[XX] + norder;
1437 thy = spline->theta[YY] + norder;
1438 thz = spline->theta[ZZ] + norder;
1443 #ifdef PME_SIMD4_SPREAD_GATHER
1444 #ifdef PME_SIMD4_UNALIGNED
1445 #define PME_SPREAD_SIMD4_ORDER4
1447 #define PME_SPREAD_SIMD4_ALIGNED
1450 #include "gromacs/ewald/pme-simd4.h" /* IWYU pragma: keep */
1456 #ifdef PME_SIMD4_SPREAD_GATHER
1457 #define PME_SPREAD_SIMD4_ALIGNED
1459 #include "gromacs/ewald/pme-simd4.h" /* IWYU pragma: keep */
1472 static void set_grid_alignment(int gmx_unused *pmegrid_nz, int gmx_unused pme_order)
1474 #ifdef PME_SIMD4_SPREAD_GATHER
1476 #ifndef PME_SIMD4_UNALIGNED
1481 /* Round nz up to a multiple of 4 to ensure alignment */
1482 *pmegrid_nz = ((*pmegrid_nz + 3) & ~3);
1487 static void set_gridsize_alignment(int gmx_unused *gridsize, int gmx_unused pme_order)
1489 #ifdef PME_SIMD4_SPREAD_GATHER
1490 #ifndef PME_SIMD4_UNALIGNED
1493 /* Add extra elements to ensured aligned operations do not go
1494 * beyond the allocated grid size.
1495 * Note that for pme_order=5, the pme grid z-size alignment
1496 * ensures that we will not go beyond the grid size.
1504 static void pmegrid_init(pmegrid_t *grid,
1505 int cx, int cy, int cz,
1506 int x0, int y0, int z0,
1507 int x1, int y1, int z1,
1508 gmx_bool set_alignment,
1517 grid->offset[XX] = x0;
1518 grid->offset[YY] = y0;
1519 grid->offset[ZZ] = z0;
1520 grid->n[XX] = x1 - x0 + pme_order - 1;
1521 grid->n[YY] = y1 - y0 + pme_order - 1;
1522 grid->n[ZZ] = z1 - z0 + pme_order - 1;
1523 copy_ivec(grid->n, grid->s);
1526 set_grid_alignment(&nz, pme_order);
1531 else if (nz != grid->s[ZZ])
1533 gmx_incons("pmegrid_init call with an unaligned z size");
1536 grid->order = pme_order;
1539 gridsize = grid->s[XX]*grid->s[YY]*grid->s[ZZ];
1540 set_gridsize_alignment(&gridsize, pme_order);
1541 snew_aligned(grid->grid, gridsize, SIMD4_ALIGNMENT);
1549 static int div_round_up(int enumerator, int denominator)
1551 return (enumerator + denominator - 1)/denominator;
1554 static void make_subgrid_division(const ivec n, int ovl, int nthread,
1557 int gsize_opt, gsize;
1562 for (nsx = 1; nsx <= nthread; nsx++)
1564 if (nthread % nsx == 0)
1566 for (nsy = 1; nsy <= nthread; nsy++)
1568 if (nsx*nsy <= nthread && nthread % (nsx*nsy) == 0)
1570 nsz = nthread/(nsx*nsy);
1572 /* Determine the number of grid points per thread */
1574 (div_round_up(n[XX], nsx) + ovl)*
1575 (div_round_up(n[YY], nsy) + ovl)*
1576 (div_round_up(n[ZZ], nsz) + ovl);
1578 /* Minimize the number of grids points per thread
1579 * and, secondarily, the number of cuts in minor dimensions.
1581 if (gsize_opt == -1 ||
1582 gsize < gsize_opt ||
1583 (gsize == gsize_opt &&
1584 (nsz < nsub[ZZ] || (nsz == nsub[ZZ] && nsy < nsub[YY]))))
1596 env = getenv("GMX_PME_THREAD_DIVISION");
1599 sscanf(env, "%d %d %d", &nsub[XX], &nsub[YY], &nsub[ZZ]);
1602 if (nsub[XX]*nsub[YY]*nsub[ZZ] != nthread)
1604 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);
1608 static void pmegrids_init(pmegrids_t *grids,
1609 int nx, int ny, int nz, int nz_base,
1611 gmx_bool bUseThreads,
1616 ivec n, n_base, g0, g1;
1617 int t, x, y, z, d, i, tfac;
1618 int max_comm_lines = -1;
1620 n[XX] = nx - (pme_order - 1);
1621 n[YY] = ny - (pme_order - 1);
1622 n[ZZ] = nz - (pme_order - 1);
1624 copy_ivec(n, n_base);
1625 n_base[ZZ] = nz_base;
1627 pmegrid_init(&grids->grid, 0, 0, 0, 0, 0, 0, n[XX], n[YY], n[ZZ], FALSE, pme_order,
1630 grids->nthread = nthread;
1632 make_subgrid_division(n_base, pme_order-1, grids->nthread, grids->nc);
1639 for (d = 0; d < DIM; d++)
1641 nst[d] = div_round_up(n[d], grids->nc[d]) + pme_order - 1;
1643 set_grid_alignment(&nst[ZZ], pme_order);
1647 fprintf(debug, "pmegrid thread local division: %d x %d x %d\n",
1648 grids->nc[XX], grids->nc[YY], grids->nc[ZZ]);
1649 fprintf(debug, "pmegrid %d %d %d max thread pmegrid %d %d %d\n",
1651 nst[XX], nst[YY], nst[ZZ]);
1654 snew(grids->grid_th, grids->nthread);
1656 gridsize = nst[XX]*nst[YY]*nst[ZZ];
1657 set_gridsize_alignment(&gridsize, pme_order);
1658 snew_aligned(grids->grid_all,
1659 grids->nthread*gridsize+(grids->nthread+1)*GMX_CACHE_SEP,
1662 for (x = 0; x < grids->nc[XX]; x++)
1664 for (y = 0; y < grids->nc[YY]; y++)
1666 for (z = 0; z < grids->nc[ZZ]; z++)
1668 pmegrid_init(&grids->grid_th[t],
1670 (n[XX]*(x ))/grids->nc[XX],
1671 (n[YY]*(y ))/grids->nc[YY],
1672 (n[ZZ]*(z ))/grids->nc[ZZ],
1673 (n[XX]*(x+1))/grids->nc[XX],
1674 (n[YY]*(y+1))/grids->nc[YY],
1675 (n[ZZ]*(z+1))/grids->nc[ZZ],
1678 grids->grid_all+GMX_CACHE_SEP+t*(gridsize+GMX_CACHE_SEP));
1686 grids->grid_th = NULL;
1689 snew(grids->g2t, DIM);
1691 for (d = DIM-1; d >= 0; d--)
1693 snew(grids->g2t[d], n[d]);
1695 for (i = 0; i < n[d]; i++)
1697 /* The second check should match the parameters
1698 * of the pmegrid_init call above.
1700 while (t + 1 < grids->nc[d] && i >= (n[d]*(t+1))/grids->nc[d])
1704 grids->g2t[d][i] = t*tfac;
1707 tfac *= grids->nc[d];
1711 case XX: max_comm_lines = overlap_x; break;
1712 case YY: max_comm_lines = overlap_y; break;
1713 case ZZ: max_comm_lines = pme_order - 1; break;
1715 grids->nthread_comm[d] = 0;
1716 while ((n[d]*grids->nthread_comm[d])/grids->nc[d] < max_comm_lines &&
1717 grids->nthread_comm[d] < grids->nc[d])
1719 grids->nthread_comm[d]++;
1723 fprintf(debug, "pmegrid thread grid communication range in %c: %d\n",
1724 'x'+d, grids->nthread_comm[d]);
1726 /* It should be possible to make grids->nthread_comm[d]==grids->nc[d]
1727 * work, but this is not a problematic restriction.
1729 if (grids->nc[d] > 1 && grids->nthread_comm[d] > grids->nc[d])
1731 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);
1737 static void pmegrids_destroy(pmegrids_t *grids)
1741 if (grids->grid.grid != NULL)
1743 sfree(grids->grid.grid);
1745 if (grids->nthread > 0)
1747 for (t = 0; t < grids->nthread; t++)
1749 sfree(grids->grid_th[t].grid);
1751 sfree(grids->grid_th);
1757 static void realloc_work(pme_work_t *work, int nkx)
1761 if (nkx > work->nalloc)
1764 srenew(work->mhx, work->nalloc);
1765 srenew(work->mhy, work->nalloc);
1766 srenew(work->mhz, work->nalloc);
1767 srenew(work->m2, work->nalloc);
1768 /* Allocate an aligned pointer for SIMD operations, including extra
1769 * elements at the end for padding.
1771 #ifdef PME_SIMD_SOLVE
1772 simd_width = GMX_SIMD_REAL_WIDTH;
1774 /* We can use any alignment, apart from 0, so we use 4 */
1777 sfree_aligned(work->denom);
1778 sfree_aligned(work->tmp1);
1779 sfree_aligned(work->tmp2);
1780 sfree_aligned(work->eterm);
1781 snew_aligned(work->denom, work->nalloc+simd_width, simd_width*sizeof(real));
1782 snew_aligned(work->tmp1, work->nalloc+simd_width, simd_width*sizeof(real));
1783 snew_aligned(work->tmp2, work->nalloc+simd_width, simd_width*sizeof(real));
1784 snew_aligned(work->eterm, work->nalloc+simd_width, simd_width*sizeof(real));
1785 srenew(work->m2inv, work->nalloc);
1790 static void free_work(pme_work_t *work)
1796 sfree_aligned(work->denom);
1797 sfree_aligned(work->tmp1);
1798 sfree_aligned(work->tmp2);
1799 sfree_aligned(work->eterm);
1804 #if defined PME_SIMD_SOLVE
1805 /* Calculate exponentials through SIMD */
1806 gmx_inline static void calc_exponentials_q(int gmx_unused start, int end, real f, real *d_aligned, real *r_aligned, real *e_aligned)
1809 const gmx_simd_real_t two = gmx_simd_set1_r(2.0);
1810 gmx_simd_real_t f_simd;
1812 gmx_simd_real_t tmp_d1, d_inv, tmp_r, tmp_e;
1814 f_simd = gmx_simd_set1_r(f);
1815 /* We only need to calculate from start. But since start is 0 or 1
1816 * and we want to use aligned loads/stores, we always start from 0.
1818 for (kx = 0; kx < end; kx += GMX_SIMD_REAL_WIDTH)
1820 tmp_d1 = gmx_simd_load_r(d_aligned+kx);
1821 d_inv = gmx_simd_inv_r(tmp_d1);
1822 tmp_r = gmx_simd_load_r(r_aligned+kx);
1823 tmp_r = gmx_simd_exp_r(tmp_r);
1824 tmp_e = gmx_simd_mul_r(f_simd, d_inv);
1825 tmp_e = gmx_simd_mul_r(tmp_e, tmp_r);
1826 gmx_simd_store_r(e_aligned+kx, tmp_e);
1831 gmx_inline static void calc_exponentials_q(int start, int end, real f, real *d, real *r, real *e)
1834 for (kx = start; kx < end; kx++)
1838 for (kx = start; kx < end; kx++)
1842 for (kx = start; kx < end; kx++)
1844 e[kx] = f*r[kx]*d[kx];
1849 #if defined PME_SIMD_SOLVE
1850 /* Calculate exponentials through SIMD */
1851 gmx_inline static void calc_exponentials_lj(int gmx_unused start, int end, real *r_aligned, real *factor_aligned, real *d_aligned)
1853 gmx_simd_real_t tmp_r, tmp_d, tmp_fac, d_inv, tmp_mk;
1854 const gmx_simd_real_t sqr_PI = gmx_simd_sqrt_r(gmx_simd_set1_r(M_PI));
1856 for (kx = 0; kx < end; kx += GMX_SIMD_REAL_WIDTH)
1858 /* We only need to calculate from start. But since start is 0 or 1
1859 * and we want to use aligned loads/stores, we always start from 0.
1861 tmp_d = gmx_simd_load_r(d_aligned+kx);
1862 d_inv = gmx_simd_inv_r(tmp_d);
1863 gmx_simd_store_r(d_aligned+kx, d_inv);
1864 tmp_r = gmx_simd_load_r(r_aligned+kx);
1865 tmp_r = gmx_simd_exp_r(tmp_r);
1866 gmx_simd_store_r(r_aligned+kx, tmp_r);
1867 tmp_mk = gmx_simd_load_r(factor_aligned+kx);
1868 tmp_fac = gmx_simd_mul_r(sqr_PI, gmx_simd_mul_r(tmp_mk, gmx_simd_erfc_r(tmp_mk)));
1869 gmx_simd_store_r(factor_aligned+kx, tmp_fac);
1873 gmx_inline static void calc_exponentials_lj(int start, int end, real *r, real *tmp2, real *d)
1877 for (kx = start; kx < end; kx++)
1882 for (kx = start; kx < end; kx++)
1887 for (kx = start; kx < end; kx++)
1890 tmp2[kx] = sqrt(M_PI)*mk*gmx_erfc(mk);
1895 static int solve_pme_yzx(gmx_pme_t pme, t_complex *grid,
1896 real ewaldcoeff, real vol,
1898 int nthread, int thread)
1900 /* do recip sum over local cells in grid */
1901 /* y major, z middle, x minor or continuous */
1903 int kx, ky, kz, maxkx, maxky, maxkz;
1904 int nx, ny, nz, iyz0, iyz1, iyz, iy, iz, kxstart, kxend;
1906 real factor = M_PI*M_PI/(ewaldcoeff*ewaldcoeff);
1907 real ets2, struct2, vfactor, ets2vf;
1908 real d1, d2, energy = 0;
1910 real virxx = 0, virxy = 0, virxz = 0, viryy = 0, viryz = 0, virzz = 0;
1911 real rxx, ryx, ryy, rzx, rzy, rzz;
1913 real *mhx, *mhy, *mhz, *m2, *denom, *tmp1, *eterm, *m2inv;
1914 real mhxk, mhyk, mhzk, m2k;
1917 ivec local_ndata, local_offset, local_size;
1920 elfac = ONE_4PI_EPS0/pme->epsilon_r;
1926 /* Dimensions should be identical for A/B grid, so we just use A here */
1927 gmx_parallel_3dfft_complex_limits(pme->pfft_setup[PME_GRID_QA],
1933 rxx = pme->recipbox[XX][XX];
1934 ryx = pme->recipbox[YY][XX];
1935 ryy = pme->recipbox[YY][YY];
1936 rzx = pme->recipbox[ZZ][XX];
1937 rzy = pme->recipbox[ZZ][YY];
1938 rzz = pme->recipbox[ZZ][ZZ];
1944 work = &pme->work[thread];
1949 denom = work->denom;
1951 eterm = work->eterm;
1952 m2inv = work->m2inv;
1954 iyz0 = local_ndata[YY]*local_ndata[ZZ]* thread /nthread;
1955 iyz1 = local_ndata[YY]*local_ndata[ZZ]*(thread+1)/nthread;
1957 for (iyz = iyz0; iyz < iyz1; iyz++)
1959 iy = iyz/local_ndata[ZZ];
1960 iz = iyz - iy*local_ndata[ZZ];
1962 ky = iy + local_offset[YY];
1973 by = M_PI*vol*pme->bsp_mod[YY][ky];
1975 kz = iz + local_offset[ZZ];
1979 bz = pme->bsp_mod[ZZ][kz];
1981 /* 0.5 correction for corner points */
1983 if (kz == 0 || kz == (nz+1)/2)
1988 p0 = grid + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
1990 /* We should skip the k-space point (0,0,0) */
1991 /* Note that since here x is the minor index, local_offset[XX]=0 */
1992 if (local_offset[XX] > 0 || ky > 0 || kz > 0)
1994 kxstart = local_offset[XX];
1998 kxstart = local_offset[XX] + 1;
2001 kxend = local_offset[XX] + local_ndata[XX];
2005 /* More expensive inner loop, especially because of the storage
2006 * of the mh elements in array's.
2007 * Because x is the minor grid index, all mh elements
2008 * depend on kx for triclinic unit cells.
2011 /* Two explicit loops to avoid a conditional inside the loop */
2012 for (kx = kxstart; kx < maxkx; kx++)
2017 mhyk = mx * ryx + my * ryy;
2018 mhzk = mx * rzx + my * rzy + mz * rzz;
2019 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2024 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2025 tmp1[kx] = -factor*m2k;
2028 for (kx = maxkx; kx < kxend; kx++)
2033 mhyk = mx * ryx + my * ryy;
2034 mhzk = mx * rzx + my * rzy + mz * rzz;
2035 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2040 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2041 tmp1[kx] = -factor*m2k;
2044 for (kx = kxstart; kx < kxend; kx++)
2046 m2inv[kx] = 1.0/m2[kx];
2049 calc_exponentials_q(kxstart, kxend, elfac, denom, tmp1, eterm);
2051 for (kx = kxstart; kx < kxend; kx++, p0++)
2056 p0->re = d1*eterm[kx];
2057 p0->im = d2*eterm[kx];
2059 struct2 = 2.0*(d1*d1+d2*d2);
2061 tmp1[kx] = eterm[kx]*struct2;
2064 for (kx = kxstart; kx < kxend; kx++)
2066 ets2 = corner_fac*tmp1[kx];
2067 vfactor = (factor*m2[kx] + 1.0)*2.0*m2inv[kx];
2070 ets2vf = ets2*vfactor;
2071 virxx += ets2vf*mhx[kx]*mhx[kx] - ets2;
2072 virxy += ets2vf*mhx[kx]*mhy[kx];
2073 virxz += ets2vf*mhx[kx]*mhz[kx];
2074 viryy += ets2vf*mhy[kx]*mhy[kx] - ets2;
2075 viryz += ets2vf*mhy[kx]*mhz[kx];
2076 virzz += ets2vf*mhz[kx]*mhz[kx] - ets2;
2081 /* We don't need to calculate the energy and the virial.
2082 * In this case the triclinic overhead is small.
2085 /* Two explicit loops to avoid a conditional inside the loop */
2087 for (kx = kxstart; kx < maxkx; kx++)
2092 mhyk = mx * ryx + my * ryy;
2093 mhzk = mx * rzx + my * rzy + mz * rzz;
2094 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2095 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2096 tmp1[kx] = -factor*m2k;
2099 for (kx = maxkx; kx < kxend; kx++)
2104 mhyk = mx * ryx + my * ryy;
2105 mhzk = mx * rzx + my * rzy + mz * rzz;
2106 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2107 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2108 tmp1[kx] = -factor*m2k;
2111 calc_exponentials_q(kxstart, kxend, elfac, denom, tmp1, eterm);
2113 for (kx = kxstart; kx < kxend; kx++, p0++)
2118 p0->re = d1*eterm[kx];
2119 p0->im = d2*eterm[kx];
2126 /* Update virial with local values.
2127 * The virial is symmetric by definition.
2128 * this virial seems ok for isotropic scaling, but I'm
2129 * experiencing problems on semiisotropic membranes.
2130 * IS THAT COMMENT STILL VALID??? (DvdS, 2001/02/07).
2132 work->vir_q[XX][XX] = 0.25*virxx;
2133 work->vir_q[YY][YY] = 0.25*viryy;
2134 work->vir_q[ZZ][ZZ] = 0.25*virzz;
2135 work->vir_q[XX][YY] = work->vir_q[YY][XX] = 0.25*virxy;
2136 work->vir_q[XX][ZZ] = work->vir_q[ZZ][XX] = 0.25*virxz;
2137 work->vir_q[YY][ZZ] = work->vir_q[ZZ][YY] = 0.25*viryz;
2139 /* This energy should be corrected for a charged system */
2140 work->energy_q = 0.5*energy;
2143 /* Return the loop count */
2144 return local_ndata[YY]*local_ndata[XX];
2147 static int solve_pme_lj_yzx(gmx_pme_t pme, t_complex **grid, gmx_bool bLB,
2148 real ewaldcoeff, real vol,
2149 gmx_bool bEnerVir, int nthread, int thread)
2151 /* do recip sum over local cells in grid */
2152 /* y major, z middle, x minor or continuous */
2154 int kx, ky, kz, maxkx, maxky, maxkz;
2155 int nx, ny, nz, iy, iyz0, iyz1, iyz, iz, kxstart, kxend;
2157 real factor = M_PI*M_PI/(ewaldcoeff*ewaldcoeff);
2159 real eterm, vterm, d1, d2, energy = 0;
2161 real virxx = 0, virxy = 0, virxz = 0, viryy = 0, viryz = 0, virzz = 0;
2162 real rxx, ryx, ryy, rzx, rzy, rzz;
2163 real *mhx, *mhy, *mhz, *m2, *denom, *tmp1, *tmp2;
2164 real mhxk, mhyk, mhzk, m2k;
2169 ivec local_ndata, local_offset, local_size;
2174 /* Dimensions should be identical for A/B grid, so we just use A here */
2175 gmx_parallel_3dfft_complex_limits(pme->pfft_setup[PME_GRID_C6A],
2180 rxx = pme->recipbox[XX][XX];
2181 ryx = pme->recipbox[YY][XX];
2182 ryy = pme->recipbox[YY][YY];
2183 rzx = pme->recipbox[ZZ][XX];
2184 rzy = pme->recipbox[ZZ][YY];
2185 rzz = pme->recipbox[ZZ][ZZ];
2191 work = &pme->work[thread];
2196 denom = work->denom;
2200 iyz0 = local_ndata[YY]*local_ndata[ZZ]* thread /nthread;
2201 iyz1 = local_ndata[YY]*local_ndata[ZZ]*(thread+1)/nthread;
2203 for (iyz = iyz0; iyz < iyz1; iyz++)
2205 iy = iyz/local_ndata[ZZ];
2206 iz = iyz - iy*local_ndata[ZZ];
2208 ky = iy + local_offset[YY];
2219 by = 3.0*vol*pme->bsp_mod[YY][ky]
2220 / (M_PI*sqrt(M_PI)*ewaldcoeff*ewaldcoeff*ewaldcoeff);
2222 kz = iz + local_offset[ZZ];
2226 bz = pme->bsp_mod[ZZ][kz];
2228 /* 0.5 correction for corner points */
2230 if (kz == 0 || kz == (nz+1)/2)
2235 kxstart = local_offset[XX];
2236 kxend = local_offset[XX] + local_ndata[XX];
2239 /* More expensive inner loop, especially because of the
2240 * storage of the mh elements in array's. Because x is the
2241 * minor grid index, all mh elements depend on kx for
2242 * triclinic unit cells.
2245 /* Two explicit loops to avoid a conditional inside the loop */
2246 for (kx = kxstart; kx < maxkx; kx++)
2251 mhyk = mx * ryx + my * ryy;
2252 mhzk = mx * rzx + my * rzy + mz * rzz;
2253 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2258 denom[kx] = bz*by*pme->bsp_mod[XX][kx];
2259 tmp1[kx] = -factor*m2k;
2260 tmp2[kx] = sqrt(factor*m2k);
2263 for (kx = maxkx; kx < kxend; kx++)
2268 mhyk = mx * ryx + my * ryy;
2269 mhzk = mx * rzx + my * rzy + mz * rzz;
2270 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2275 denom[kx] = bz*by*pme->bsp_mod[XX][kx];
2276 tmp1[kx] = -factor*m2k;
2277 tmp2[kx] = sqrt(factor*m2k);
2280 calc_exponentials_lj(kxstart, kxend, tmp1, tmp2, denom);
2282 for (kx = kxstart; kx < kxend; kx++)
2284 m2k = factor*m2[kx];
2285 eterm = -((1.0 - 2.0*m2k)*tmp1[kx]
2286 + 2.0*m2k*tmp2[kx]);
2287 vterm = 3.0*(-tmp1[kx] + tmp2[kx]);
2288 tmp1[kx] = eterm*denom[kx];
2289 tmp2[kx] = vterm*denom[kx];
2297 p0 = grid[0] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2298 for (kx = kxstart; kx < kxend; kx++, p0++)
2308 struct2 = 2.0*(d1*d1+d2*d2);
2310 tmp1[kx] = eterm*struct2;
2311 tmp2[kx] = vterm*struct2;
2316 real *struct2 = denom;
2319 for (kx = kxstart; kx < kxend; kx++)
2323 /* Due to symmetry we only need to calculate 4 of the 7 terms */
2324 for (ig = 0; ig <= 3; ++ig)
2329 p0 = grid[ig] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2330 p1 = grid[6-ig] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2331 scale = 2.0*lb_scale_factor_symm[ig];
2332 for (kx = kxstart; kx < kxend; ++kx, ++p0, ++p1)
2334 struct2[kx] += scale*(p0->re*p1->re + p0->im*p1->im);
2338 for (ig = 0; ig <= 6; ++ig)
2342 p0 = grid[ig] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2343 for (kx = kxstart; kx < kxend; kx++, p0++)
2353 for (kx = kxstart; kx < kxend; kx++)
2358 tmp1[kx] = eterm*str2;
2359 tmp2[kx] = vterm*str2;
2363 for (kx = kxstart; kx < kxend; kx++)
2365 ets2 = corner_fac*tmp1[kx];
2366 vterm = 2.0*factor*tmp2[kx];
2368 ets2vf = corner_fac*vterm;
2369 virxx += ets2vf*mhx[kx]*mhx[kx] - ets2;
2370 virxy += ets2vf*mhx[kx]*mhy[kx];
2371 virxz += ets2vf*mhx[kx]*mhz[kx];
2372 viryy += ets2vf*mhy[kx]*mhy[kx] - ets2;
2373 viryz += ets2vf*mhy[kx]*mhz[kx];
2374 virzz += ets2vf*mhz[kx]*mhz[kx] - ets2;
2379 /* We don't need to calculate the energy and the virial.
2380 * In this case the triclinic overhead is small.
2383 /* Two explicit loops to avoid a conditional inside the loop */
2385 for (kx = kxstart; kx < maxkx; 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 for (kx = maxkx; kx < kxend; kx++)
2404 mhyk = mx * ryx + my * ryy;
2405 mhzk = mx * rzx + my * rzy + mz * rzz;
2406 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2408 denom[kx] = bz*by*pme->bsp_mod[XX][kx];
2409 tmp1[kx] = -factor*m2k;
2410 tmp2[kx] = sqrt(factor*m2k);
2413 calc_exponentials_lj(kxstart, kxend, tmp1, tmp2, denom);
2415 for (kx = kxstart; kx < kxend; kx++)
2417 m2k = factor*m2[kx];
2418 eterm = -((1.0 - 2.0*m2k)*tmp1[kx]
2419 + 2.0*m2k*tmp2[kx]);
2420 tmp1[kx] = eterm*denom[kx];
2422 gcount = (bLB ? 7 : 1);
2423 for (ig = 0; ig < gcount; ++ig)
2427 p0 = grid[ig] + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2428 for (kx = kxstart; kx < kxend; kx++, p0++)
2443 work->vir_lj[XX][XX] = 0.25*virxx;
2444 work->vir_lj[YY][YY] = 0.25*viryy;
2445 work->vir_lj[ZZ][ZZ] = 0.25*virzz;
2446 work->vir_lj[XX][YY] = work->vir_lj[YY][XX] = 0.25*virxy;
2447 work->vir_lj[XX][ZZ] = work->vir_lj[ZZ][XX] = 0.25*virxz;
2448 work->vir_lj[YY][ZZ] = work->vir_lj[ZZ][YY] = 0.25*viryz;
2450 /* This energy should be corrected for a charged system */
2451 work->energy_lj = 0.5*energy;
2453 /* Return the loop count */
2454 return local_ndata[YY]*local_ndata[XX];
2457 static void get_pme_ener_vir_q(const gmx_pme_t pme, int nthread,
2458 real *mesh_energy, matrix vir)
2460 /* This function sums output over threads and should therefore
2461 * only be called after thread synchronization.
2465 *mesh_energy = pme->work[0].energy_q;
2466 copy_mat(pme->work[0].vir_q, vir);
2468 for (thread = 1; thread < nthread; thread++)
2470 *mesh_energy += pme->work[thread].energy_q;
2471 m_add(vir, pme->work[thread].vir_q, vir);
2475 static void get_pme_ener_vir_lj(const gmx_pme_t pme, int nthread,
2476 real *mesh_energy, matrix vir)
2478 /* This function sums output over threads and should therefore
2479 * only be called after thread synchronization.
2483 *mesh_energy = pme->work[0].energy_lj;
2484 copy_mat(pme->work[0].vir_lj, vir);
2486 for (thread = 1; thread < nthread; thread++)
2488 *mesh_energy += pme->work[thread].energy_lj;
2489 m_add(vir, pme->work[thread].vir_lj, vir);
2494 #define DO_FSPLINE(order) \
2495 for (ithx = 0; (ithx < order); ithx++) \
2497 index_x = (i0+ithx)*pny*pnz; \
2501 for (ithy = 0; (ithy < order); ithy++) \
2503 index_xy = index_x+(j0+ithy)*pnz; \
2508 for (ithz = 0; (ithz < order); ithz++) \
2510 gval = grid[index_xy+(k0+ithz)]; \
2511 fxy1 += thz[ithz]*gval; \
2512 fz1 += dthz[ithz]*gval; \
2521 static void gather_f_bsplines(gmx_pme_t pme, real *grid,
2522 gmx_bool bClearF, pme_atomcomm_t *atc,
2523 splinedata_t *spline,
2526 /* sum forces for local particles */
2527 int nn, n, ithx, ithy, ithz, i0, j0, k0;
2528 int index_x, index_xy;
2529 int nx, ny, nz, pnx, pny, pnz;
2531 real tx, ty, dx, dy, coefficient;
2532 real fx, fy, fz, gval;
2534 real *thx, *thy, *thz, *dthx, *dthy, *dthz;
2536 real rxx, ryx, ryy, rzx, rzy, rzz;
2539 pme_spline_work_t *work;
2541 #if defined PME_SIMD4_SPREAD_GATHER && !defined PME_SIMD4_UNALIGNED
2542 real thz_buffer[GMX_SIMD4_WIDTH*3], *thz_aligned;
2543 real dthz_buffer[GMX_SIMD4_WIDTH*3], *dthz_aligned;
2545 thz_aligned = gmx_simd4_align_r(thz_buffer);
2546 dthz_aligned = gmx_simd4_align_r(dthz_buffer);
2549 work = pme->spline_work;
2551 order = pme->pme_order;
2552 thx = spline->theta[XX];
2553 thy = spline->theta[YY];
2554 thz = spline->theta[ZZ];
2555 dthx = spline->dtheta[XX];
2556 dthy = spline->dtheta[YY];
2557 dthz = spline->dtheta[ZZ];
2561 pnx = pme->pmegrid_nx;
2562 pny = pme->pmegrid_ny;
2563 pnz = pme->pmegrid_nz;
2565 rxx = pme->recipbox[XX][XX];
2566 ryx = pme->recipbox[YY][XX];
2567 ryy = pme->recipbox[YY][YY];
2568 rzx = pme->recipbox[ZZ][XX];
2569 rzy = pme->recipbox[ZZ][YY];
2570 rzz = pme->recipbox[ZZ][ZZ];
2572 for (nn = 0; nn < spline->n; nn++)
2574 n = spline->ind[nn];
2575 coefficient = scale*atc->coefficient[n];
2583 if (coefficient != 0)
2588 idxptr = atc->idx[n];
2595 /* Pointer arithmetic alert, next six statements */
2596 thx = spline->theta[XX] + norder;
2597 thy = spline->theta[YY] + norder;
2598 thz = spline->theta[ZZ] + norder;
2599 dthx = spline->dtheta[XX] + norder;
2600 dthy = spline->dtheta[YY] + norder;
2601 dthz = spline->dtheta[ZZ] + norder;
2606 #ifdef PME_SIMD4_SPREAD_GATHER
2607 #ifdef PME_SIMD4_UNALIGNED
2608 #define PME_GATHER_F_SIMD4_ORDER4
2610 #define PME_GATHER_F_SIMD4_ALIGNED
2613 #include "gromacs/ewald/pme-simd4.h" /* IWYU pragma: keep */
2619 #ifdef PME_SIMD4_SPREAD_GATHER
2620 #define PME_GATHER_F_SIMD4_ALIGNED
2622 #include "gromacs/ewald/pme-simd4.h" /* IWYU pragma: keep */
2632 atc->f[n][XX] += -coefficient*( fx*nx*rxx );
2633 atc->f[n][YY] += -coefficient*( fx*nx*ryx + fy*ny*ryy );
2634 atc->f[n][ZZ] += -coefficient*( fx*nx*rzx + fy*ny*rzy + fz*nz*rzz );
2637 /* Since the energy and not forces are interpolated
2638 * the net force might not be exactly zero.
2639 * This can be solved by also interpolating F, but
2640 * that comes at a cost.
2641 * A better hack is to remove the net force every
2642 * step, but that must be done at a higher level
2643 * since this routine doesn't see all atoms if running
2644 * in parallel. Don't know how important it is? EL 990726
2649 static real gather_energy_bsplines(gmx_pme_t pme, real *grid,
2650 pme_atomcomm_t *atc)
2652 splinedata_t *spline;
2653 int n, ithx, ithy, ithz, i0, j0, k0;
2654 int index_x, index_xy;
2656 real energy, pot, tx, ty, coefficient, gval;
2657 real *thx, *thy, *thz;
2661 spline = &atc->spline[0];
2663 order = pme->pme_order;
2666 for (n = 0; (n < atc->n); n++)
2668 coefficient = atc->coefficient[n];
2670 if (coefficient != 0)
2672 idxptr = atc->idx[n];
2679 /* Pointer arithmetic alert, next three statements */
2680 thx = spline->theta[XX] + norder;
2681 thy = spline->theta[YY] + norder;
2682 thz = spline->theta[ZZ] + norder;
2685 for (ithx = 0; (ithx < order); ithx++)
2687 index_x = (i0+ithx)*pme->pmegrid_ny*pme->pmegrid_nz;
2690 for (ithy = 0; (ithy < order); ithy++)
2692 index_xy = index_x+(j0+ithy)*pme->pmegrid_nz;
2695 for (ithz = 0; (ithz < order); ithz++)
2697 gval = grid[index_xy+(k0+ithz)];
2698 pot += tx*ty*thz[ithz]*gval;
2704 energy += pot*coefficient;
2711 /* Macro to force loop unrolling by fixing order.
2712 * This gives a significant performance gain.
2714 #define CALC_SPLINE(order) \
2718 real data[PME_ORDER_MAX]; \
2719 real ddata[PME_ORDER_MAX]; \
2721 for (j = 0; (j < DIM); j++) \
2725 /* dr is relative offset from lower cell limit */ \
2726 data[order-1] = 0; \
2730 for (k = 3; (k < order); k++) \
2732 div = 1.0/(k - 1.0); \
2733 data[k-1] = div*dr*data[k-2]; \
2734 for (l = 1; (l < (k-1)); l++) \
2736 data[k-l-1] = div*((dr+l)*data[k-l-2]+(k-l-dr)* \
2739 data[0] = div*(1-dr)*data[0]; \
2741 /* differentiate */ \
2742 ddata[0] = -data[0]; \
2743 for (k = 1; (k < order); k++) \
2745 ddata[k] = data[k-1] - data[k]; \
2748 div = 1.0/(order - 1); \
2749 data[order-1] = div*dr*data[order-2]; \
2750 for (l = 1; (l < (order-1)); l++) \
2752 data[order-l-1] = div*((dr+l)*data[order-l-2]+ \
2753 (order-l-dr)*data[order-l-1]); \
2755 data[0] = div*(1 - dr)*data[0]; \
2757 for (k = 0; k < order; k++) \
2759 theta[j][i*order+k] = data[k]; \
2760 dtheta[j][i*order+k] = ddata[k]; \
2765 void make_bsplines(splinevec theta, splinevec dtheta, int order,
2766 rvec fractx[], int nr, int ind[], real coefficient[],
2767 gmx_bool bDoSplines)
2769 /* construct splines for local atoms */
2773 for (i = 0; i < nr; i++)
2775 /* With free energy we do not use the coefficient check.
2776 * In most cases this will be more efficient than calling make_bsplines
2777 * twice, since usually more than half the particles have non-zero coefficients.
2780 if (bDoSplines || coefficient[ii] != 0.0)
2785 case 4: CALC_SPLINE(4); break;
2786 case 5: CALC_SPLINE(5); break;
2787 default: CALC_SPLINE(order); break;
2794 void make_dft_mod(real *mod, real *data, int ndata)
2799 for (i = 0; i < ndata; i++)
2802 for (j = 0; j < ndata; j++)
2804 arg = (2.0*M_PI*i*j)/ndata;
2805 sc += data[j]*cos(arg);
2806 ss += data[j]*sin(arg);
2808 mod[i] = sc*sc+ss*ss;
2810 for (i = 0; i < ndata; i++)
2814 mod[i] = (mod[i-1]+mod[i+1])*0.5;
2820 static void make_bspline_moduli(splinevec bsp_mod,
2821 int nx, int ny, int nz, int order)
2823 int nmax = max(nx, max(ny, nz));
2824 real *data, *ddata, *bsp_data;
2830 snew(bsp_data, nmax);
2836 for (k = 3; k < order; k++)
2840 for (l = 1; l < (k-1); l++)
2842 data[k-l-1] = div*(l*data[k-l-2]+(k-l)*data[k-l-1]);
2844 data[0] = div*data[0];
2847 ddata[0] = -data[0];
2848 for (k = 1; k < order; k++)
2850 ddata[k] = data[k-1]-data[k];
2852 div = 1.0/(order-1);
2854 for (l = 1; l < (order-1); l++)
2856 data[order-l-1] = div*(l*data[order-l-2]+(order-l)*data[order-l-1]);
2858 data[0] = div*data[0];
2860 for (i = 0; i < nmax; i++)
2864 for (i = 1; i <= order; i++)
2866 bsp_data[i] = data[i-1];
2869 make_dft_mod(bsp_mod[XX], bsp_data, nx);
2870 make_dft_mod(bsp_mod[YY], bsp_data, ny);
2871 make_dft_mod(bsp_mod[ZZ], bsp_data, nz);
2879 /* Return the P3M optimal influence function */
2880 static double do_p3m_influence(double z, int order)
2887 /* The formula and most constants can be found in:
2888 * Ballenegger et al., JCTC 8, 936 (2012)
2893 return 1.0 - 2.0*z2/3.0;
2896 return 1.0 - z2 + 2.0*z4/15.0;
2899 return 1.0 - 4.0*z2/3.0 + 2.0*z4/5.0 + 4.0*z2*z4/315.0;
2902 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;
2905 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;
2908 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;
2910 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;
2917 /* Calculate the P3M B-spline moduli for one dimension */
2918 static void make_p3m_bspline_moduli_dim(real *bsp_mod, int n, int order)
2920 double zarg, zai, sinzai, infl;
2925 gmx_fatal(FARGS, "The current P3M code only supports orders up to 8");
2932 for (i = -maxk; i < 0; i++)
2936 infl = do_p3m_influence(sinzai, order);
2937 bsp_mod[n+i] = infl*infl*pow(sinzai/zai, -2.0*order);
2940 for (i = 1; i < maxk; i++)
2944 infl = do_p3m_influence(sinzai, order);
2945 bsp_mod[i] = infl*infl*pow(sinzai/zai, -2.0*order);
2949 /* Calculate the P3M B-spline moduli */
2950 static void make_p3m_bspline_moduli(splinevec bsp_mod,
2951 int nx, int ny, int nz, int order)
2953 make_p3m_bspline_moduli_dim(bsp_mod[XX], nx, order);
2954 make_p3m_bspline_moduli_dim(bsp_mod[YY], ny, order);
2955 make_p3m_bspline_moduli_dim(bsp_mod[ZZ], nz, order);
2959 static void setup_coordinate_communication(pme_atomcomm_t *atc)
2967 for (i = 1; i <= nslab/2; i++)
2969 fw = (atc->nodeid + i) % nslab;
2970 bw = (atc->nodeid - i + nslab) % nslab;
2973 atc->node_dest[n] = fw;
2974 atc->node_src[n] = bw;
2979 atc->node_dest[n] = bw;
2980 atc->node_src[n] = fw;
2986 int gmx_pme_destroy(FILE *log, gmx_pme_t *pmedata)
2992 fprintf(log, "Destroying PME data structures.\n");
2995 sfree((*pmedata)->nnx);
2996 sfree((*pmedata)->nny);
2997 sfree((*pmedata)->nnz);
2999 for (i = 0; i < (*pmedata)->ngrids; ++i)
3001 pmegrids_destroy(&(*pmedata)->pmegrid[i]);
3002 sfree((*pmedata)->fftgrid[i]);
3003 sfree((*pmedata)->cfftgrid[i]);
3004 gmx_parallel_3dfft_destroy((*pmedata)->pfft_setup[i]);
3007 sfree((*pmedata)->lb_buf1);
3008 sfree((*pmedata)->lb_buf2);
3010 for (thread = 0; thread < (*pmedata)->nthread; thread++)
3012 free_work(&(*pmedata)->work[thread]);
3014 sfree((*pmedata)->work);
3022 static int mult_up(int n, int f)
3024 return ((n + f - 1)/f)*f;
3028 static double pme_load_imbalance(gmx_pme_t pme)
3033 nma = pme->nnodes_major;
3034 nmi = pme->nnodes_minor;
3036 n1 = mult_up(pme->nkx, nma)*mult_up(pme->nky, nmi)*pme->nkz;
3037 n2 = mult_up(pme->nkx, nma)*mult_up(pme->nkz, nmi)*pme->nky;
3038 n3 = mult_up(pme->nky, nma)*mult_up(pme->nkz, nmi)*pme->nkx;
3040 /* pme_solve is roughly double the cost of an fft */
3042 return (n1 + n2 + 3*n3)/(double)(6*pme->nkx*pme->nky*pme->nkz);
3045 static void init_atomcomm(gmx_pme_t pme, pme_atomcomm_t *atc,
3046 int dimind, gmx_bool bSpread)
3048 int nk, k, s, thread;
3050 atc->dimind = dimind;
3055 if (pme->nnodes > 1)
3057 atc->mpi_comm = pme->mpi_comm_d[dimind];
3058 MPI_Comm_size(atc->mpi_comm, &atc->nslab);
3059 MPI_Comm_rank(atc->mpi_comm, &atc->nodeid);
3063 fprintf(debug, "For PME atom communication in dimind %d: nslab %d rank %d\n", atc->dimind, atc->nslab, atc->nodeid);
3067 atc->bSpread = bSpread;
3068 atc->pme_order = pme->pme_order;
3072 snew(atc->node_dest, atc->nslab);
3073 snew(atc->node_src, atc->nslab);
3074 setup_coordinate_communication(atc);
3076 snew(atc->count_thread, pme->nthread);
3077 for (thread = 0; thread < pme->nthread; thread++)
3079 snew(atc->count_thread[thread], atc->nslab);
3081 atc->count = atc->count_thread[0];
3082 snew(atc->rcount, atc->nslab);
3083 snew(atc->buf_index, atc->nslab);
3086 atc->nthread = pme->nthread;
3087 if (atc->nthread > 1)
3089 snew(atc->thread_plist, atc->nthread);
3091 snew(atc->spline, atc->nthread);
3092 for (thread = 0; thread < atc->nthread; thread++)
3094 if (atc->nthread > 1)
3096 snew(atc->thread_plist[thread].n, atc->nthread+2*GMX_CACHE_SEP);
3097 atc->thread_plist[thread].n += GMX_CACHE_SEP;
3099 snew(atc->spline[thread].thread_one, pme->nthread);
3100 atc->spline[thread].thread_one[thread] = 1;
3105 init_overlap_comm(pme_overlap_t * ol,
3115 int lbnd, rbnd, maxlr, b, i;
3118 pme_grid_comm_t *pgc;
3120 int fft_start, fft_end, send_index1, recv_index1;
3124 ol->mpi_comm = comm;
3127 ol->nnodes = nnodes;
3128 ol->nodeid = nodeid;
3130 /* Linear translation of the PME grid won't affect reciprocal space
3131 * calculations, so to optimize we only interpolate "upwards",
3132 * which also means we only have to consider overlap in one direction.
3133 * I.e., particles on this node might also be spread to grid indices
3134 * that belong to higher nodes (modulo nnodes)
3137 snew(ol->s2g0, ol->nnodes+1);
3138 snew(ol->s2g1, ol->nnodes);
3141 fprintf(debug, "PME slab boundaries:");
3143 for (i = 0; i < nnodes; i++)
3145 /* s2g0 the local interpolation grid start.
3146 * s2g1 the local interpolation grid end.
3147 * Since in calc_pidx we divide particles, and not grid lines,
3148 * spatially uniform along dimension x or y, we need to round
3149 * s2g0 down and s2g1 up.
3151 ol->s2g0[i] = ( i *ndata + 0 )/nnodes;
3152 ol->s2g1[i] = ((i+1)*ndata + nnodes-1)/nnodes + norder - 1;
3156 fprintf(debug, " %3d %3d", ol->s2g0[i], ol->s2g1[i]);
3159 ol->s2g0[nnodes] = ndata;
3162 fprintf(debug, "\n");
3165 /* Determine with how many nodes we need to communicate the grid overlap */
3171 for (i = 0; i < nnodes; i++)
3173 if ((i+b < nnodes && ol->s2g1[i] > ol->s2g0[i+b]) ||
3174 (i+b >= nnodes && ol->s2g1[i] > ol->s2g0[i+b-nnodes] + ndata))
3180 while (bCont && b < nnodes);
3181 ol->noverlap_nodes = b - 1;
3183 snew(ol->send_id, ol->noverlap_nodes);
3184 snew(ol->recv_id, ol->noverlap_nodes);
3185 for (b = 0; b < ol->noverlap_nodes; b++)
3187 ol->send_id[b] = (ol->nodeid + (b + 1)) % ol->nnodes;
3188 ol->recv_id[b] = (ol->nodeid - (b + 1) + ol->nnodes) % ol->nnodes;
3190 snew(ol->comm_data, ol->noverlap_nodes);
3193 for (b = 0; b < ol->noverlap_nodes; b++)
3195 pgc = &ol->comm_data[b];
3197 fft_start = ol->s2g0[ol->send_id[b]];
3198 fft_end = ol->s2g0[ol->send_id[b]+1];
3199 if (ol->send_id[b] < nodeid)
3204 send_index1 = ol->s2g1[nodeid];
3205 send_index1 = min(send_index1, fft_end);
3206 pgc->send_index0 = fft_start;
3207 pgc->send_nindex = max(0, send_index1 - pgc->send_index0);
3208 ol->send_size += pgc->send_nindex;
3210 /* We always start receiving to the first index of our slab */
3211 fft_start = ol->s2g0[ol->nodeid];
3212 fft_end = ol->s2g0[ol->nodeid+1];
3213 recv_index1 = ol->s2g1[ol->recv_id[b]];
3214 if (ol->recv_id[b] > nodeid)
3216 recv_index1 -= ndata;
3218 recv_index1 = min(recv_index1, fft_end);
3219 pgc->recv_index0 = fft_start;
3220 pgc->recv_nindex = max(0, recv_index1 - pgc->recv_index0);
3224 /* Communicate the buffer sizes to receive */
3225 for (b = 0; b < ol->noverlap_nodes; b++)
3227 MPI_Sendrecv(&ol->send_size, 1, MPI_INT, ol->send_id[b], b,
3228 &ol->comm_data[b].recv_size, 1, MPI_INT, ol->recv_id[b], b,
3229 ol->mpi_comm, &stat);
3233 /* For non-divisible grid we need pme_order iso pme_order-1 */
3234 snew(ol->sendbuf, norder*commplainsize);
3235 snew(ol->recvbuf, norder*commplainsize);
3239 make_gridindex5_to_localindex(int n, int local_start, int local_range,
3240 int **global_to_local,
3241 real **fraction_shift)
3249 for (i = 0; (i < 5*n); i++)
3251 /* Determine the global to local grid index */
3252 gtl[i] = (i - local_start + n) % n;
3253 /* For coordinates that fall within the local grid the fraction
3254 * is correct, we don't need to shift it.
3257 if (local_range < n)
3259 /* Due to rounding issues i could be 1 beyond the lower or
3260 * upper boundary of the local grid. Correct the index for this.
3261 * If we shift the index, we need to shift the fraction by
3262 * the same amount in the other direction to not affect
3264 * Note that due to this shifting the weights at the end of
3265 * the spline might change, but that will only involve values
3266 * between zero and values close to the precision of a real,
3267 * which is anyhow the accuracy of the whole mesh calculation.
3269 /* With local_range=0 we should not change i=local_start */
3270 if (i % n != local_start)
3277 else if (gtl[i] == local_range)
3279 gtl[i] = local_range - 1;
3286 *global_to_local = gtl;
3287 *fraction_shift = fsh;
3290 static pme_spline_work_t *make_pme_spline_work(int gmx_unused order)
3292 pme_spline_work_t *work;
3294 #ifdef PME_SIMD4_SPREAD_GATHER
3295 real tmp[GMX_SIMD4_WIDTH*3], *tmp_aligned;
3296 gmx_simd4_real_t zero_S;
3297 gmx_simd4_real_t real_mask_S0, real_mask_S1;
3300 snew_aligned(work, 1, SIMD4_ALIGNMENT);
3302 tmp_aligned = gmx_simd4_align_r(tmp);
3304 zero_S = gmx_simd4_setzero_r();
3306 /* Generate bit masks to mask out the unused grid entries,
3307 * as we only operate on order of the 8 grid entries that are
3308 * load into 2 SIMD registers.
3310 for (of = 0; of < 2*GMX_SIMD4_WIDTH-(order-1); of++)
3312 for (i = 0; i < 2*GMX_SIMD4_WIDTH; i++)
3314 tmp_aligned[i] = (i >= of && i < of+order ? -1.0 : 1.0);
3316 real_mask_S0 = gmx_simd4_load_r(tmp_aligned);
3317 real_mask_S1 = gmx_simd4_load_r(tmp_aligned+GMX_SIMD4_WIDTH);
3318 work->mask_S0[of] = gmx_simd4_cmplt_r(real_mask_S0, zero_S);
3319 work->mask_S1[of] = gmx_simd4_cmplt_r(real_mask_S1, zero_S);
3328 void gmx_pme_check_restrictions(int pme_order,
3329 int nkx, int nky, int nkz,
3332 gmx_bool bUseThreads,
3334 gmx_bool *bValidSettings)
3336 if (pme_order > PME_ORDER_MAX)
3340 *bValidSettings = FALSE;
3343 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.",
3344 pme_order, PME_ORDER_MAX);
3347 if (nkx <= pme_order*(nnodes_major > 1 ? 2 : 1) ||
3348 nky <= pme_order*(nnodes_minor > 1 ? 2 : 1) ||
3353 *bValidSettings = FALSE;
3356 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",
3360 /* Check for a limitation of the (current) sum_fftgrid_dd code.
3361 * We only allow multiple communication pulses in dim 1, not in dim 0.
3363 if (bUseThreads && (nkx < nnodes_major*pme_order &&
3364 nkx != nnodes_major*(pme_order - 1)))
3368 *bValidSettings = FALSE;
3371 gmx_fatal(FARGS, "The number of PME grid lines per rank along x is %g. But when using OpenMP threads, the number of grid lines per rank along x should be >= pme_order (%d) or = pmeorder-1. To resolve this issue, use fewer ranks along x (and possibly more along y and/or z) by specifying -dd manually.",
3372 nkx/(double)nnodes_major, pme_order);
3375 if (bValidSettings != NULL)
3377 *bValidSettings = TRUE;
3383 int gmx_pme_init(gmx_pme_t * pmedata,
3389 gmx_bool bFreeEnergy_q,
3390 gmx_bool bFreeEnergy_lj,
3391 gmx_bool bReproducible,
3394 gmx_pme_t pme = NULL;
3396 int use_threads, sum_use_threads, i;
3401 fprintf(debug, "Creating PME data structures.\n");
3405 pme->sum_qgrid_tmp = NULL;
3406 pme->sum_qgrid_dd_tmp = NULL;
3407 pme->buf_nalloc = 0;
3410 pme->bPPnode = TRUE;
3412 pme->nnodes_major = nnodes_major;
3413 pme->nnodes_minor = nnodes_minor;
3416 if (nnodes_major*nnodes_minor > 1)
3418 pme->mpi_comm = cr->mpi_comm_mygroup;
3420 MPI_Comm_rank(pme->mpi_comm, &pme->nodeid);
3421 MPI_Comm_size(pme->mpi_comm, &pme->nnodes);
3422 if (pme->nnodes != nnodes_major*nnodes_minor)
3424 gmx_incons("PME rank count mismatch");
3429 pme->mpi_comm = MPI_COMM_NULL;
3433 if (pme->nnodes == 1)
3436 pme->mpi_comm_d[0] = MPI_COMM_NULL;
3437 pme->mpi_comm_d[1] = MPI_COMM_NULL;
3439 pme->ndecompdim = 0;
3440 pme->nodeid_major = 0;
3441 pme->nodeid_minor = 0;
3443 pme->mpi_comm_d[0] = pme->mpi_comm_d[1] = MPI_COMM_NULL;
3448 if (nnodes_minor == 1)
3451 pme->mpi_comm_d[0] = pme->mpi_comm;
3452 pme->mpi_comm_d[1] = MPI_COMM_NULL;
3454 pme->ndecompdim = 1;
3455 pme->nodeid_major = pme->nodeid;
3456 pme->nodeid_minor = 0;
3459 else if (nnodes_major == 1)
3462 pme->mpi_comm_d[0] = MPI_COMM_NULL;
3463 pme->mpi_comm_d[1] = pme->mpi_comm;
3465 pme->ndecompdim = 1;
3466 pme->nodeid_major = 0;
3467 pme->nodeid_minor = pme->nodeid;
3471 if (pme->nnodes % nnodes_major != 0)
3473 gmx_incons("For 2D PME decomposition, #PME ranks must be divisible by the number of ranks in the major dimension");
3475 pme->ndecompdim = 2;
3478 MPI_Comm_split(pme->mpi_comm, pme->nodeid % nnodes_minor,
3479 pme->nodeid, &pme->mpi_comm_d[0]); /* My communicator along major dimension */
3480 MPI_Comm_split(pme->mpi_comm, pme->nodeid/nnodes_minor,
3481 pme->nodeid, &pme->mpi_comm_d[1]); /* My communicator along minor dimension */
3483 MPI_Comm_rank(pme->mpi_comm_d[0], &pme->nodeid_major);
3484 MPI_Comm_size(pme->mpi_comm_d[0], &pme->nnodes_major);
3485 MPI_Comm_rank(pme->mpi_comm_d[1], &pme->nodeid_minor);
3486 MPI_Comm_size(pme->mpi_comm_d[1], &pme->nnodes_minor);
3489 pme->bPPnode = (cr->duty & DUTY_PP);
3492 pme->nthread = nthread;
3494 /* Check if any of the PME MPI ranks uses threads */
3495 use_threads = (pme->nthread > 1 ? 1 : 0);
3497 if (pme->nnodes > 1)
3499 MPI_Allreduce(&use_threads, &sum_use_threads, 1, MPI_INT,
3500 MPI_SUM, pme->mpi_comm);
3505 sum_use_threads = use_threads;
3507 pme->bUseThreads = (sum_use_threads > 0);
3509 if (ir->ePBC == epbcSCREW)
3511 gmx_fatal(FARGS, "pme does not (yet) work with pbc = screw");
3514 pme->bFEP_q = ((ir->efep != efepNO) && bFreeEnergy_q);
3515 pme->bFEP_lj = ((ir->efep != efepNO) && bFreeEnergy_lj);
3516 pme->bFEP = (pme->bFEP_q || pme->bFEP_lj);
3520 pme->bP3M = (ir->coulombtype == eelP3M_AD || getenv("GMX_PME_P3M") != NULL);
3521 pme->pme_order = ir->pme_order;
3523 /* Always constant electrostatics coefficients */
3524 pme->epsilon_r = ir->epsilon_r;
3526 /* Always constant LJ coefficients */
3527 pme->ljpme_combination_rule = ir->ljpme_combination_rule;
3529 /* If we violate restrictions, generate a fatal error here */
3530 gmx_pme_check_restrictions(pme->pme_order,
3531 pme->nkx, pme->nky, pme->nkz,
3538 if (pme->nnodes > 1)
3543 MPI_Type_contiguous(DIM, mpi_type, &(pme->rvec_mpi));
3544 MPI_Type_commit(&(pme->rvec_mpi));
3547 /* Note that the coefficient spreading and force gathering, which usually
3548 * takes about the same amount of time as FFT+solve_pme,
3549 * is always fully load balanced
3550 * (unless the coefficient distribution is inhomogeneous).
3553 imbal = pme_load_imbalance(pme);
3554 if (imbal >= 1.2 && pme->nodeid_major == 0 && pme->nodeid_minor == 0)
3558 "NOTE: The load imbalance in PME FFT and solve is %d%%.\n"
3559 " For optimal PME load balancing\n"
3560 " PME grid_x (%d) and grid_y (%d) should be divisible by #PME_ranks_x (%d)\n"
3561 " and PME grid_y (%d) and grid_z (%d) should be divisible by #PME_ranks_y (%d)\n"
3563 (int)((imbal-1)*100 + 0.5),
3564 pme->nkx, pme->nky, pme->nnodes_major,
3565 pme->nky, pme->nkz, pme->nnodes_minor);
3569 /* For non-divisible grid we need pme_order iso pme_order-1 */
3570 /* In sum_qgrid_dd x overlap is copied in place: take padding into account.
3571 * y is always copied through a buffer: we don't need padding in z,
3572 * but we do need the overlap in x because of the communication order.
3574 init_overlap_comm(&pme->overlap[0], pme->pme_order,
3578 pme->nnodes_major, pme->nodeid_major,
3580 (div_round_up(pme->nky, pme->nnodes_minor)+pme->pme_order)*(pme->nkz+pme->pme_order-1));
3582 /* Along overlap dim 1 we can send in multiple pulses in sum_fftgrid_dd.
3583 * We do this with an offset buffer of equal size, so we need to allocate
3584 * extra for the offset. That's what the (+1)*pme->nkz is for.
3586 init_overlap_comm(&pme->overlap[1], pme->pme_order,
3590 pme->nnodes_minor, pme->nodeid_minor,
3592 (div_round_up(pme->nkx, pme->nnodes_major)+pme->pme_order+1)*pme->nkz);
3594 /* Double-check for a limitation of the (current) sum_fftgrid_dd code.
3595 * Note that gmx_pme_check_restrictions checked for this already.
3597 if (pme->bUseThreads && pme->overlap[0].noverlap_nodes > 1)
3599 gmx_incons("More than one communication pulse required for grid overlap communication along the major dimension while using threads");
3602 snew(pme->bsp_mod[XX], pme->nkx);
3603 snew(pme->bsp_mod[YY], pme->nky);
3604 snew(pme->bsp_mod[ZZ], pme->nkz);
3606 /* The required size of the interpolation grid, including overlap.
3607 * The allocated size (pmegrid_n?) might be slightly larger.
3609 pme->pmegrid_nx = pme->overlap[0].s2g1[pme->nodeid_major] -
3610 pme->overlap[0].s2g0[pme->nodeid_major];
3611 pme->pmegrid_ny = pme->overlap[1].s2g1[pme->nodeid_minor] -
3612 pme->overlap[1].s2g0[pme->nodeid_minor];
3613 pme->pmegrid_nz_base = pme->nkz;
3614 pme->pmegrid_nz = pme->pmegrid_nz_base + pme->pme_order - 1;
3615 set_grid_alignment(&pme->pmegrid_nz, pme->pme_order);
3617 pme->pmegrid_start_ix = pme->overlap[0].s2g0[pme->nodeid_major];
3618 pme->pmegrid_start_iy = pme->overlap[1].s2g0[pme->nodeid_minor];
3619 pme->pmegrid_start_iz = 0;
3621 make_gridindex5_to_localindex(pme->nkx,
3622 pme->pmegrid_start_ix,
3623 pme->pmegrid_nx - (pme->pme_order-1),
3624 &pme->nnx, &pme->fshx);
3625 make_gridindex5_to_localindex(pme->nky,
3626 pme->pmegrid_start_iy,
3627 pme->pmegrid_ny - (pme->pme_order-1),
3628 &pme->nny, &pme->fshy);
3629 make_gridindex5_to_localindex(pme->nkz,
3630 pme->pmegrid_start_iz,
3631 pme->pmegrid_nz_base,
3632 &pme->nnz, &pme->fshz);
3634 pme->spline_work = make_pme_spline_work(pme->pme_order);
3636 ndata[0] = pme->nkx;
3637 ndata[1] = pme->nky;
3638 ndata[2] = pme->nkz;
3639 /* It doesn't matter if we allocate too many grids here,
3640 * we only allocate and use the ones we need.
3642 if (EVDW_PME(ir->vdwtype))
3644 pme->ngrids = ((ir->ljpme_combination_rule == eljpmeLB) ? DO_Q_AND_LJ_LB : DO_Q_AND_LJ);
3650 snew(pme->fftgrid, pme->ngrids);
3651 snew(pme->cfftgrid, pme->ngrids);
3652 snew(pme->pfft_setup, pme->ngrids);
3654 for (i = 0; i < pme->ngrids; ++i)
3656 if ((i < DO_Q && EEL_PME(ir->coulombtype) && (i == 0 ||
3658 (i >= DO_Q && EVDW_PME(ir->vdwtype) && (i == 2 ||
3660 ir->ljpme_combination_rule == eljpmeLB)))
3662 pmegrids_init(&pme->pmegrid[i],
3663 pme->pmegrid_nx, pme->pmegrid_ny, pme->pmegrid_nz,
3664 pme->pmegrid_nz_base,
3668 pme->overlap[0].s2g1[pme->nodeid_major]-pme->overlap[0].s2g0[pme->nodeid_major+1],
3669 pme->overlap[1].s2g1[pme->nodeid_minor]-pme->overlap[1].s2g0[pme->nodeid_minor+1]);
3670 /* This routine will allocate the grid data to fit the FFTs */
3671 gmx_parallel_3dfft_init(&pme->pfft_setup[i], ndata,
3672 &pme->fftgrid[i], &pme->cfftgrid[i],
3674 bReproducible, pme->nthread);
3681 /* Use plain SPME B-spline interpolation */
3682 make_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
3686 /* Use the P3M grid-optimized influence function */
3687 make_p3m_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
3690 /* Use atc[0] for spreading */
3691 init_atomcomm(pme, &pme->atc[0], nnodes_major > 1 ? 0 : 1, TRUE);
3692 if (pme->ndecompdim >= 2)
3694 init_atomcomm(pme, &pme->atc[1], 1, FALSE);
3697 if (pme->nnodes == 1)
3699 pme->atc[0].n = homenr;
3700 pme_realloc_atomcomm_things(&pme->atc[0]);
3703 pme->lb_buf1 = NULL;
3704 pme->lb_buf2 = NULL;
3705 pme->lb_buf_nalloc = 0;
3710 /* Use fft5d, order after FFT is y major, z, x minor */
3712 snew(pme->work, pme->nthread);
3713 for (thread = 0; thread < pme->nthread; thread++)
3715 realloc_work(&pme->work[thread], pme->nkx);
3724 static void reuse_pmegrids(const pmegrids_t *old, pmegrids_t *new)
3728 for (d = 0; d < DIM; d++)
3730 if (new->grid.n[d] > old->grid.n[d])
3736 sfree_aligned(new->grid.grid);
3737 new->grid.grid = old->grid.grid;
3739 if (new->grid_th != NULL && new->nthread == old->nthread)
3741 sfree_aligned(new->grid_all);
3742 for (t = 0; t < new->nthread; t++)
3744 new->grid_th[t].grid = old->grid_th[t].grid;
3749 int gmx_pme_reinit(gmx_pme_t * pmedata,
3752 const t_inputrec * ir,
3760 irc.nkx = grid_size[XX];
3761 irc.nky = grid_size[YY];
3762 irc.nkz = grid_size[ZZ];
3764 if (pme_src->nnodes == 1)
3766 homenr = pme_src->atc[0].n;
3773 ret = gmx_pme_init(pmedata, cr, pme_src->nnodes_major, pme_src->nnodes_minor,
3774 &irc, homenr, pme_src->bFEP_q, pme_src->bFEP_lj, FALSE, pme_src->nthread);
3778 /* We can easily reuse the allocated pme grids in pme_src */
3779 reuse_pmegrids(&pme_src->pmegrid[PME_GRID_QA], &(*pmedata)->pmegrid[PME_GRID_QA]);
3780 /* We would like to reuse the fft grids, but that's harder */
3787 static void copy_local_grid(gmx_pme_t pme, pmegrids_t *pmegrids,
3788 int grid_index, int thread, real *fftgrid)
3790 ivec local_fft_ndata, local_fft_offset, local_fft_size;
3794 int offx, offy, offz, x, y, z, i0, i0t;
3799 gmx_parallel_3dfft_real_limits(pme->pfft_setup[grid_index],
3803 fft_my = local_fft_size[YY];
3804 fft_mz = local_fft_size[ZZ];
3806 pmegrid = &pmegrids->grid_th[thread];
3808 nsx = pmegrid->s[XX];
3809 nsy = pmegrid->s[YY];
3810 nsz = pmegrid->s[ZZ];
3812 for (d = 0; d < DIM; d++)
3814 nf[d] = min(pmegrid->n[d] - (pmegrid->order - 1),
3815 local_fft_ndata[d] - pmegrid->offset[d]);
3818 offx = pmegrid->offset[XX];
3819 offy = pmegrid->offset[YY];
3820 offz = pmegrid->offset[ZZ];
3822 /* Directly copy the non-overlapping parts of the local grids.
3823 * This also initializes the full grid.
3825 grid_th = pmegrid->grid;
3826 for (x = 0; x < nf[XX]; x++)
3828 for (y = 0; y < nf[YY]; y++)
3830 i0 = ((offx + x)*fft_my + (offy + y))*fft_mz + offz;
3831 i0t = (x*nsy + y)*nsz;
3832 for (z = 0; z < nf[ZZ]; z++)
3834 fftgrid[i0+z] = grid_th[i0t+z];
3841 reduce_threadgrid_overlap(gmx_pme_t pme,
3842 const pmegrids_t *pmegrids, int thread,
3843 real *fftgrid, real *commbuf_x, real *commbuf_y,
3846 ivec local_fft_ndata, local_fft_offset, local_fft_size;
3847 int fft_nx, fft_ny, fft_nz;
3852 int offx, offy, offz, x, y, z, i0, i0t;
3853 int sx, sy, sz, fx, fy, fz, tx1, ty1, tz1, ox, oy, oz;
3854 gmx_bool bClearBufX, bClearBufY, bClearBufXY, bClearBuf;
3855 gmx_bool bCommX, bCommY;
3858 const pmegrid_t *pmegrid, *pmegrid_g, *pmegrid_f;
3859 const real *grid_th;
3860 real *commbuf = NULL;
3862 gmx_parallel_3dfft_real_limits(pme->pfft_setup[grid_index],
3866 fft_nx = local_fft_ndata[XX];
3867 fft_ny = local_fft_ndata[YY];
3868 fft_nz = local_fft_ndata[ZZ];
3870 fft_my = local_fft_size[YY];
3871 fft_mz = local_fft_size[ZZ];
3873 /* This routine is called when all thread have finished spreading.
3874 * Here each thread sums grid contributions calculated by other threads
3875 * to the thread local grid volume.
3876 * To minimize the number of grid copying operations,
3877 * this routines sums immediately from the pmegrid to the fftgrid.
3880 /* Determine which part of the full node grid we should operate on,
3881 * this is our thread local part of the full grid.
3883 pmegrid = &pmegrids->grid_th[thread];
3885 for (d = 0; d < DIM; d++)
3887 /* Determine up to where our thread needs to copy from the
3888 * thread-local charge spreading grid to the rank-local FFT grid.
3889 * This is up to our spreading grid end minus order-1 and
3890 * not beyond the local FFT grid.
3893 min(pmegrid->offset[d]+pmegrid->n[d]-(pmegrid->order-1),
3894 local_fft_ndata[d]);
3897 offx = pmegrid->offset[XX];
3898 offy = pmegrid->offset[YY];
3899 offz = pmegrid->offset[ZZ];
3906 /* Now loop over all the thread data blocks that contribute
3907 * to the grid region we (our thread) are operating on.
3909 /* Note that fft_nx/y is equal to the number of grid points
3910 * between the first point of our node grid and the one of the next node.
3912 for (sx = 0; sx >= -pmegrids->nthread_comm[XX]; sx--)
3914 fx = pmegrid->ci[XX] + sx;
3919 fx += pmegrids->nc[XX];
3921 bCommX = (pme->nnodes_major > 1);
3923 pmegrid_g = &pmegrids->grid_th[fx*pmegrids->nc[YY]*pmegrids->nc[ZZ]];
3924 ox += pmegrid_g->offset[XX];
3925 /* Determine the end of our part of the source grid */
3928 /* Use our thread local source grid and target grid part */
3929 tx1 = min(ox + pmegrid_g->n[XX], localcopy_end[XX]);
3933 /* Use our thread local source grid and the spreading range */
3934 tx1 = min(ox + pmegrid_g->n[XX], pme->pme_order);
3937 for (sy = 0; sy >= -pmegrids->nthread_comm[YY]; sy--)
3939 fy = pmegrid->ci[YY] + sy;
3944 fy += pmegrids->nc[YY];
3946 bCommY = (pme->nnodes_minor > 1);
3948 pmegrid_g = &pmegrids->grid_th[fy*pmegrids->nc[ZZ]];
3949 oy += pmegrid_g->offset[YY];
3950 /* Determine the end of our part of the source grid */
3953 /* Use our thread local source grid and target grid part */
3954 ty1 = min(oy + pmegrid_g->n[YY], localcopy_end[YY]);
3958 /* Use our thread local source grid and the spreading range */
3959 ty1 = min(oy + pmegrid_g->n[YY], pme->pme_order);
3962 for (sz = 0; sz >= -pmegrids->nthread_comm[ZZ]; sz--)
3964 fz = pmegrid->ci[ZZ] + sz;
3968 fz += pmegrids->nc[ZZ];
3971 pmegrid_g = &pmegrids->grid_th[fz];
3972 oz += pmegrid_g->offset[ZZ];
3973 tz1 = min(oz + pmegrid_g->n[ZZ], localcopy_end[ZZ]);
3975 if (sx == 0 && sy == 0 && sz == 0)
3977 /* We have already added our local contribution
3978 * before calling this routine, so skip it here.
3983 thread_f = (fx*pmegrids->nc[YY] + fy)*pmegrids->nc[ZZ] + fz;
3985 pmegrid_f = &pmegrids->grid_th[thread_f];
3987 grid_th = pmegrid_f->grid;
3989 nsx = pmegrid_f->s[XX];
3990 nsy = pmegrid_f->s[YY];
3991 nsz = pmegrid_f->s[ZZ];
3993 #ifdef DEBUG_PME_REDUCE
3994 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",
3995 pme->nodeid, thread, thread_f,
3996 pme->pmegrid_start_ix,
3997 pme->pmegrid_start_iy,
3998 pme->pmegrid_start_iz,
4000 offx-ox, tx1-ox, offx, tx1,
4001 offy-oy, ty1-oy, offy, ty1,
4002 offz-oz, tz1-oz, offz, tz1);
4005 if (!(bCommX || bCommY))
4007 /* Copy from the thread local grid to the node grid */
4008 for (x = offx; x < tx1; x++)
4010 for (y = offy; y < ty1; y++)
4012 i0 = (x*fft_my + y)*fft_mz;
4013 i0t = ((x - ox)*nsy + (y - oy))*nsz - oz;
4014 for (z = offz; z < tz1; z++)
4016 fftgrid[i0+z] += grid_th[i0t+z];
4023 /* The order of this conditional decides
4024 * where the corner volume gets stored with x+y decomp.
4028 commbuf = commbuf_y;
4029 /* The y-size of the communication buffer is set by
4030 * the overlap of the grid part of our local slab
4031 * with the part starting at the next slab.
4034 pme->overlap[1].s2g1[pme->nodeid_minor] -
4035 pme->overlap[1].s2g0[pme->nodeid_minor+1];
4038 /* We index commbuf modulo the local grid size */
4039 commbuf += buf_my*fft_nx*fft_nz;
4041 bClearBuf = bClearBufXY;
4042 bClearBufXY = FALSE;
4046 bClearBuf = bClearBufY;
4052 commbuf = commbuf_x;
4054 bClearBuf = bClearBufX;
4058 /* Copy to the communication buffer */
4059 for (x = offx; x < tx1; x++)
4061 for (y = offy; y < ty1; y++)
4063 i0 = (x*buf_my + y)*fft_nz;
4064 i0t = ((x - ox)*nsy + (y - oy))*nsz - oz;
4068 /* First access of commbuf, initialize it */
4069 for (z = offz; z < tz1; z++)
4071 commbuf[i0+z] = grid_th[i0t+z];
4076 for (z = offz; z < tz1; z++)
4078 commbuf[i0+z] += grid_th[i0t+z];
4090 static void sum_fftgrid_dd(gmx_pme_t pme, real *fftgrid, int grid_index)
4092 ivec local_fft_ndata, local_fft_offset, local_fft_size;
4093 pme_overlap_t *overlap;
4094 int send_index0, send_nindex;
4099 int send_size_y, recv_size_y;
4100 int ipulse, send_id, recv_id, datasize, gridsize, size_yx;
4101 real *sendptr, *recvptr;
4102 int x, y, z, indg, indb;
4104 /* Note that this routine is only used for forward communication.
4105 * Since the force gathering, unlike the coefficient spreading,
4106 * can be trivially parallelized over the particles,
4107 * the backwards process is much simpler and can use the "old"
4108 * communication setup.
4111 gmx_parallel_3dfft_real_limits(pme->pfft_setup[grid_index],
4116 if (pme->nnodes_minor > 1)
4118 /* Major dimension */
4119 overlap = &pme->overlap[1];
4121 if (pme->nnodes_major > 1)
4123 size_yx = pme->overlap[0].comm_data[0].send_nindex;
4129 datasize = (local_fft_ndata[XX] + size_yx)*local_fft_ndata[ZZ];
4131 send_size_y = overlap->send_size;
4133 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
4135 send_id = overlap->send_id[ipulse];
4136 recv_id = overlap->recv_id[ipulse];
4138 overlap->comm_data[ipulse].send_index0 -
4139 overlap->comm_data[0].send_index0;
4140 send_nindex = overlap->comm_data[ipulse].send_nindex;
4141 /* We don't use recv_index0, as we always receive starting at 0 */
4142 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
4143 recv_size_y = overlap->comm_data[ipulse].recv_size;
4145 sendptr = overlap->sendbuf + send_index0*local_fft_ndata[ZZ];
4146 recvptr = overlap->recvbuf;
4150 fprintf(debug, "PME fftgrid comm y %2d x %2d x %2d\n",
4151 local_fft_ndata[XX], send_nindex, local_fft_ndata[ZZ]);
4155 MPI_Sendrecv(sendptr, send_size_y*datasize, GMX_MPI_REAL,
4157 recvptr, recv_size_y*datasize, GMX_MPI_REAL,
4159 overlap->mpi_comm, &stat);
4162 for (x = 0; x < local_fft_ndata[XX]; x++)
4164 for (y = 0; y < recv_nindex; y++)
4166 indg = (x*local_fft_size[YY] + y)*local_fft_size[ZZ];
4167 indb = (x*recv_size_y + y)*local_fft_ndata[ZZ];
4168 for (z = 0; z < local_fft_ndata[ZZ]; z++)
4170 fftgrid[indg+z] += recvptr[indb+z];
4175 if (pme->nnodes_major > 1)
4177 /* Copy from the received buffer to the send buffer for dim 0 */
4178 sendptr = pme->overlap[0].sendbuf;
4179 for (x = 0; x < size_yx; x++)
4181 for (y = 0; y < recv_nindex; y++)
4183 indg = (x*local_fft_ndata[YY] + y)*local_fft_ndata[ZZ];
4184 indb = ((local_fft_ndata[XX] + x)*recv_size_y + y)*local_fft_ndata[ZZ];
4185 for (z = 0; z < local_fft_ndata[ZZ]; z++)
4187 sendptr[indg+z] += recvptr[indb+z];
4195 /* We only support a single pulse here.
4196 * This is not a severe limitation, as this code is only used
4197 * with OpenMP and with OpenMP the (PME) domains can be larger.
4199 if (pme->nnodes_major > 1)
4201 /* Major dimension */
4202 overlap = &pme->overlap[0];
4204 datasize = local_fft_ndata[YY]*local_fft_ndata[ZZ];
4205 gridsize = local_fft_size[YY] *local_fft_size[ZZ];
4209 send_id = overlap->send_id[ipulse];
4210 recv_id = overlap->recv_id[ipulse];
4211 send_nindex = overlap->comm_data[ipulse].send_nindex;
4212 /* We don't use recv_index0, as we always receive starting at 0 */
4213 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
4215 sendptr = overlap->sendbuf;
4216 recvptr = overlap->recvbuf;
4220 fprintf(debug, "PME fftgrid comm x %2d x %2d x %2d\n",
4221 send_nindex, local_fft_ndata[YY], local_fft_ndata[ZZ]);
4225 MPI_Sendrecv(sendptr, send_nindex*datasize, GMX_MPI_REAL,
4227 recvptr, recv_nindex*datasize, GMX_MPI_REAL,
4229 overlap->mpi_comm, &stat);
4232 for (x = 0; x < recv_nindex; x++)
4234 for (y = 0; y < local_fft_ndata[YY]; y++)
4236 indg = (x*local_fft_size[YY] + y)*local_fft_size[ZZ];
4237 indb = (x*local_fft_ndata[YY] + y)*local_fft_ndata[ZZ];
4238 for (z = 0; z < local_fft_ndata[ZZ]; z++)
4240 fftgrid[indg+z] += recvptr[indb+z];
4248 static void spread_on_grid(gmx_pme_t pme,
4249 pme_atomcomm_t *atc, pmegrids_t *grids,
4250 gmx_bool bCalcSplines, gmx_bool bSpread,
4251 real *fftgrid, gmx_bool bDoSplines, int grid_index)
4253 int nthread, thread;
4254 #ifdef PME_TIME_THREADS
4255 gmx_cycles_t c1, c2, c3, ct1a, ct1b, ct1c;
4256 static double cs1 = 0, cs2 = 0, cs3 = 0;
4257 static double cs1a[6] = {0, 0, 0, 0, 0, 0};
4261 nthread = pme->nthread;
4262 assert(nthread > 0);
4264 #ifdef PME_TIME_THREADS
4265 c1 = omp_cyc_start();
4269 #pragma omp parallel for num_threads(nthread) schedule(static)
4270 for (thread = 0; thread < nthread; thread++)
4274 start = atc->n* thread /nthread;
4275 end = atc->n*(thread+1)/nthread;
4277 /* Compute fftgrid index for all atoms,
4278 * with help of some extra variables.
4280 calc_interpolation_idx(pme, atc, start, grid_index, end, thread);
4283 #ifdef PME_TIME_THREADS
4284 c1 = omp_cyc_end(c1);
4288 #ifdef PME_TIME_THREADS
4289 c2 = omp_cyc_start();
4291 #pragma omp parallel for num_threads(nthread) schedule(static)
4292 for (thread = 0; thread < nthread; thread++)
4294 splinedata_t *spline;
4295 pmegrid_t *grid = NULL;
4297 /* make local bsplines */
4298 if (grids == NULL || !pme->bUseThreads)
4300 spline = &atc->spline[0];
4306 grid = &grids->grid;
4311 spline = &atc->spline[thread];
4313 if (grids->nthread == 1)
4315 /* One thread, we operate on all coefficients */
4320 /* Get the indices our thread should operate on */
4321 make_thread_local_ind(atc, thread, spline);
4324 grid = &grids->grid_th[thread];
4329 make_bsplines(spline->theta, spline->dtheta, pme->pme_order,
4330 atc->fractx, spline->n, spline->ind, atc->coefficient, bDoSplines);
4335 /* put local atoms on grid. */
4336 #ifdef PME_TIME_SPREAD
4337 ct1a = omp_cyc_start();
4339 spread_coefficients_bsplines_thread(grid, atc, spline, pme->spline_work);
4341 if (pme->bUseThreads)
4343 copy_local_grid(pme, grids, grid_index, thread, fftgrid);
4345 #ifdef PME_TIME_SPREAD
4346 ct1a = omp_cyc_end(ct1a);
4347 cs1a[thread] += (double)ct1a;
4351 #ifdef PME_TIME_THREADS
4352 c2 = omp_cyc_end(c2);
4356 if (bSpread && pme->bUseThreads)
4358 #ifdef PME_TIME_THREADS
4359 c3 = omp_cyc_start();
4361 #pragma omp parallel for num_threads(grids->nthread) schedule(static)
4362 for (thread = 0; thread < grids->nthread; thread++)
4364 reduce_threadgrid_overlap(pme, grids, thread,
4366 pme->overlap[0].sendbuf,
4367 pme->overlap[1].sendbuf,
4370 #ifdef PME_TIME_THREADS
4371 c3 = omp_cyc_end(c3);
4375 if (pme->nnodes > 1)
4377 /* Communicate the overlapping part of the fftgrid.
4378 * For this communication call we need to check pme->bUseThreads
4379 * to have all ranks communicate here, regardless of pme->nthread.
4381 sum_fftgrid_dd(pme, fftgrid, grid_index);
4385 #ifdef PME_TIME_THREADS
4389 printf("idx %.2f spread %.2f red %.2f",
4390 cs1*1e-9, cs2*1e-9, cs3*1e-9);
4391 #ifdef PME_TIME_SPREAD
4392 for (thread = 0; thread < nthread; thread++)
4394 printf(" %.2f", cs1a[thread]*1e-9);
4403 static void dump_grid(FILE *fp,
4404 int sx, int sy, int sz, int nx, int ny, int nz,
4405 int my, int mz, const real *g)
4409 for (x = 0; x < nx; x++)
4411 for (y = 0; y < ny; y++)
4413 for (z = 0; z < nz; z++)
4415 fprintf(fp, "%2d %2d %2d %6.3f\n",
4416 sx+x, sy+y, sz+z, g[(x*my + y)*mz + z]);
4422 static void dump_local_fftgrid(gmx_pme_t pme, const real *fftgrid)
4424 ivec local_fft_ndata, local_fft_offset, local_fft_size;
4426 gmx_parallel_3dfft_real_limits(pme->pfft_setup[PME_GRID_QA],
4432 pme->pmegrid_start_ix,
4433 pme->pmegrid_start_iy,
4434 pme->pmegrid_start_iz,
4435 pme->pmegrid_nx-pme->pme_order+1,
4436 pme->pmegrid_ny-pme->pme_order+1,
4437 pme->pmegrid_nz-pme->pme_order+1,
4444 void gmx_pme_calc_energy(gmx_pme_t pme, int n, rvec *x, real *q, real *V)
4446 pme_atomcomm_t *atc;
4449 if (pme->nnodes > 1)
4451 gmx_incons("gmx_pme_calc_energy called in parallel");
4453 if (pme->bFEP_q > 1)
4455 gmx_incons("gmx_pme_calc_energy with free energy");
4458 atc = &pme->atc_energy;
4460 if (atc->spline == NULL)
4462 snew(atc->spline, atc->nthread);
4465 atc->bSpread = TRUE;
4466 atc->pme_order = pme->pme_order;
4468 pme_realloc_atomcomm_things(atc);
4470 atc->coefficient = q;
4472 /* We only use the A-charges grid */
4473 grid = &pme->pmegrid[PME_GRID_QA];
4475 /* Only calculate the spline coefficients, don't actually spread */
4476 spread_on_grid(pme, atc, NULL, TRUE, FALSE, pme->fftgrid[PME_GRID_QA], FALSE, PME_GRID_QA);
4478 *V = gather_energy_bsplines(pme, grid->grid.grid, atc);
4482 static void reset_pmeonly_counters(gmx_wallcycle_t wcycle,
4483 gmx_walltime_accounting_t walltime_accounting,
4484 t_nrnb *nrnb, t_inputrec *ir,
4487 /* Reset all the counters related to performance over the run */
4488 wallcycle_stop(wcycle, ewcRUN);
4489 wallcycle_reset_all(wcycle);
4491 if (ir->nsteps >= 0)
4493 /* ir->nsteps is not used here, but we update it for consistency */
4494 ir->nsteps -= step - ir->init_step;
4496 ir->init_step = step;
4497 wallcycle_start(wcycle, ewcRUN);
4498 walltime_accounting_start(walltime_accounting);
4502 static void gmx_pmeonly_switch(int *npmedata, gmx_pme_t **pmedata,
4504 t_commrec *cr, t_inputrec *ir,
4508 gmx_pme_t pme = NULL;
4511 while (ind < *npmedata)
4513 pme = (*pmedata)[ind];
4514 if (pme->nkx == grid_size[XX] &&
4515 pme->nky == grid_size[YY] &&
4516 pme->nkz == grid_size[ZZ])
4527 srenew(*pmedata, *npmedata);
4529 /* Generate a new PME data structure, copying part of the old pointers */
4530 gmx_pme_reinit(&((*pmedata)[ind]), cr, pme, ir, grid_size);
4532 *pme_ret = (*pmedata)[ind];
4535 int gmx_pmeonly(gmx_pme_t pme,
4536 t_commrec *cr, t_nrnb *mynrnb,
4537 gmx_wallcycle_t wcycle,
4538 gmx_walltime_accounting_t walltime_accounting,
4539 real ewaldcoeff_q, real ewaldcoeff_lj,
4544 gmx_pme_pp_t pme_pp;
4548 rvec *x_pp = NULL, *f_pp = NULL;
4549 real *chargeA = NULL, *chargeB = NULL;
4550 real *c6A = NULL, *c6B = NULL;
4551 real *sigmaA = NULL, *sigmaB = NULL;
4554 int maxshift_x = 0, maxshift_y = 0;
4555 real energy_q, energy_lj, dvdlambda_q, dvdlambda_lj;
4556 matrix vir_q, vir_lj;
4561 gmx_int64_t step, step_rel;
4564 /* This data will only use with PME tuning, i.e. switching PME grids */
4566 snew(pmedata, npmedata);
4569 pme_pp = gmx_pme_pp_init(cr);
4574 do /****** this is a quasi-loop over time steps! */
4576 /* The reason for having a loop here is PME grid tuning/switching */
4579 /* Domain decomposition */
4580 ret = gmx_pme_recv_coeffs_coords(pme_pp,
4586 &maxshift_x, &maxshift_y,
4587 &pme->bFEP_q, &pme->bFEP_lj,
4588 &lambda_q, &lambda_lj,
4592 grid_switch, &ewaldcoeff_q, &ewaldcoeff_lj);
4594 if (ret == pmerecvqxSWITCHGRID)
4596 /* Switch the PME grid to grid_switch */
4597 gmx_pmeonly_switch(&npmedata, &pmedata, grid_switch, cr, ir, &pme);
4600 if (ret == pmerecvqxRESETCOUNTERS)
4602 /* Reset the cycle and flop counters */
4603 reset_pmeonly_counters(wcycle, walltime_accounting, mynrnb, ir, step);
4606 while (ret == pmerecvqxSWITCHGRID || ret == pmerecvqxRESETCOUNTERS);
4608 if (ret == pmerecvqxFINISH)
4610 /* We should stop: break out of the loop */
4614 step_rel = step - ir->init_step;
4618 wallcycle_start(wcycle, ewcRUN);
4619 walltime_accounting_start(walltime_accounting);
4622 wallcycle_start(wcycle, ewcPMEMESH);
4629 gmx_pme_do(pme, 0, natoms, x_pp, f_pp,
4630 chargeA, chargeB, c6A, c6B, sigmaA, sigmaB, box,
4631 cr, maxshift_x, maxshift_y, mynrnb, wcycle,
4632 vir_q, ewaldcoeff_q, vir_lj, ewaldcoeff_lj,
4633 &energy_q, &energy_lj, lambda_q, lambda_lj, &dvdlambda_q, &dvdlambda_lj,
4634 pme_flags | GMX_PME_DO_ALL_F | (bEnerVir ? GMX_PME_CALC_ENER_VIR : 0));
4636 cycles = wallcycle_stop(wcycle, ewcPMEMESH);
4638 gmx_pme_send_force_vir_ener(pme_pp,
4639 f_pp, vir_q, energy_q, vir_lj, energy_lj,
4640 dvdlambda_q, dvdlambda_lj, cycles);
4643 } /***** end of quasi-loop, we stop with the break above */
4646 walltime_accounting_end(walltime_accounting);
4652 calc_initial_lb_coeffs(gmx_pme_t pme, real *local_c6, real *local_sigma)
4656 for (i = 0; i < pme->atc[0].n; ++i)
4660 sigma4 = local_sigma[i];
4661 sigma4 = sigma4*sigma4;
4662 sigma4 = sigma4*sigma4;
4663 pme->atc[0].coefficient[i] = local_c6[i] / sigma4;
4668 calc_next_lb_coeffs(gmx_pme_t pme, real *local_sigma)
4672 for (i = 0; i < pme->atc[0].n; ++i)
4674 pme->atc[0].coefficient[i] *= local_sigma[i];
4679 do_redist_pos_coeffs(gmx_pme_t pme, t_commrec *cr, int start, int homenr,
4680 gmx_bool bFirst, rvec x[], real *data)
4683 pme_atomcomm_t *atc;
4686 for (d = pme->ndecompdim - 1; d >= 0; d--)
4692 if (d == pme->ndecompdim - 1)
4700 n_d = pme->atc[d + 1].n;
4702 param_d = atc->coefficient;
4706 if (atc->npd > atc->pd_nalloc)
4708 atc->pd_nalloc = over_alloc_dd(atc->npd);
4709 srenew(atc->pd, atc->pd_nalloc);
4711 pme_calc_pidx_wrapper(n_d, pme->recipbox, x_d, atc);
4713 /* Redistribute x (only once) and qA/c6A or qB/c6B */
4714 if (DOMAINDECOMP(cr))
4716 dd_pmeredist_pos_coeffs(pme, n_d, bFirst, x_d, param_d, atc);
4721 int gmx_pme_do(gmx_pme_t pme,
4722 int start, int homenr,
4724 real *chargeA, real *chargeB,
4725 real *c6A, real *c6B,
4726 real *sigmaA, real *sigmaB,
4727 matrix box, t_commrec *cr,
4728 int maxshift_x, int maxshift_y,
4729 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
4730 matrix vir_q, real ewaldcoeff_q,
4731 matrix vir_lj, real ewaldcoeff_lj,
4732 real *energy_q, real *energy_lj,
4733 real lambda_q, real lambda_lj,
4734 real *dvdlambda_q, real *dvdlambda_lj,
4737 int d, i, j, k, ntot, npme, grid_index, max_grid_index;
4740 pme_atomcomm_t *atc = NULL;
4741 pmegrids_t *pmegrid = NULL;
4745 real *coefficient = NULL;
4750 gmx_parallel_3dfft_t pfft_setup;
4752 t_complex * cfftgrid;
4754 gmx_bool bFirst, bDoSplines;
4756 int fep_states_lj = pme->bFEP_lj ? 2 : 1;
4757 const gmx_bool bCalcEnerVir = flags & GMX_PME_CALC_ENER_VIR;
4758 const gmx_bool bCalcF = flags & GMX_PME_CALC_F;
4760 assert(pme->nnodes > 0);
4761 assert(pme->nnodes == 1 || pme->ndecompdim > 0);
4763 if (pme->nnodes > 1)
4767 if (atc->npd > atc->pd_nalloc)
4769 atc->pd_nalloc = over_alloc_dd(atc->npd);
4770 srenew(atc->pd, atc->pd_nalloc);
4772 for (d = pme->ndecompdim-1; d >= 0; d--)
4775 atc->maxshift = (atc->dimind == 0 ? maxshift_x : maxshift_y);
4781 /* This could be necessary for TPI */
4782 pme->atc[0].n = homenr;
4783 if (DOMAINDECOMP(cr))
4785 pme_realloc_atomcomm_things(atc);
4791 m_inv_ur0(box, pme->recipbox);
4794 /* For simplicity, we construct the splines for all particles if
4795 * more than one PME calculations is needed. Some optimization
4796 * could be done by keeping track of which atoms have splines
4797 * constructed, and construct new splines on each pass for atoms
4798 * that don't yet have them.
4801 bDoSplines = pme->bFEP || ((flags & GMX_PME_DO_COULOMB) && (flags & GMX_PME_DO_LJ));
4803 /* We need a maximum of four separate PME calculations:
4804 * grid_index=0: Coulomb PME with charges from state A
4805 * grid_index=1: Coulomb PME with charges from state B
4806 * grid_index=2: LJ PME with C6 from state A
4807 * grid_index=3: LJ PME with C6 from state B
4808 * For Lorentz-Berthelot combination rules, a separate loop is used to
4809 * calculate all the terms
4812 /* If we are doing LJ-PME with LB, we only do Q here */
4813 max_grid_index = (pme->ljpme_combination_rule == eljpmeLB) ? DO_Q : DO_Q_AND_LJ;
4815 for (grid_index = 0; grid_index < max_grid_index; ++grid_index)
4817 /* Check if we should do calculations at this grid_index
4818 * If grid_index is odd we should be doing FEP
4819 * If grid_index < 2 we should be doing electrostatic PME
4820 * If grid_index >= 2 we should be doing LJ-PME
4822 if ((grid_index < DO_Q && (!(flags & GMX_PME_DO_COULOMB) ||
4823 (grid_index == 1 && !pme->bFEP_q))) ||
4824 (grid_index >= DO_Q && (!(flags & GMX_PME_DO_LJ) ||
4825 (grid_index == 3 && !pme->bFEP_lj))))
4829 /* Unpack structure */
4830 pmegrid = &pme->pmegrid[grid_index];
4831 fftgrid = pme->fftgrid[grid_index];
4832 cfftgrid = pme->cfftgrid[grid_index];
4833 pfft_setup = pme->pfft_setup[grid_index];
4836 case 0: coefficient = chargeA + start; break;
4837 case 1: coefficient = chargeB + start; break;
4838 case 2: coefficient = c6A + start; break;
4839 case 3: coefficient = c6B + start; break;
4842 grid = pmegrid->grid.grid;
4846 fprintf(debug, "PME: number of ranks = %d, rank = %d\n",
4847 cr->nnodes, cr->nodeid);
4848 fprintf(debug, "Grid = %p\n", (void*)grid);
4851 gmx_fatal(FARGS, "No grid!");
4856 if (pme->nnodes == 1)
4858 atc->coefficient = coefficient;
4862 wallcycle_start(wcycle, ewcPME_REDISTXF);
4863 do_redist_pos_coeffs(pme, cr, start, homenr, bFirst, x, coefficient);
4866 wallcycle_stop(wcycle, ewcPME_REDISTXF);
4871 fprintf(debug, "Rank= %6d, pme local particles=%6d\n",
4872 cr->nodeid, atc->n);
4875 if (flags & GMX_PME_SPREAD)
4877 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
4879 /* Spread the coefficients on a grid */
4880 spread_on_grid(pme, &pme->atc[0], pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
4884 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
4886 inc_nrnb(nrnb, eNR_SPREADBSP,
4887 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
4889 if (!pme->bUseThreads)
4891 wrap_periodic_pmegrid(pme, grid);
4893 /* sum contributions to local grid from other nodes */
4895 if (pme->nnodes > 1)
4897 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
4902 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
4905 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
4908 dump_local_fftgrid(pme,fftgrid);
4913 /* Here we start a large thread parallel region */
4914 #pragma omp parallel num_threads(pme->nthread) private(thread)
4916 thread = gmx_omp_get_thread_num();
4917 if (flags & GMX_PME_SOLVE)
4924 wallcycle_start(wcycle, ewcPME_FFT);
4926 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
4930 wallcycle_stop(wcycle, ewcPME_FFT);
4934 /* solve in k-space for our local cells */
4937 wallcycle_start(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
4939 if (grid_index < DO_Q)
4942 solve_pme_yzx(pme, cfftgrid, ewaldcoeff_q,
4943 box[XX][XX]*box[YY][YY]*box[ZZ][ZZ],
4945 pme->nthread, thread);
4950 solve_pme_lj_yzx(pme, &cfftgrid, FALSE, ewaldcoeff_lj,
4951 box[XX][XX]*box[YY][YY]*box[ZZ][ZZ],
4953 pme->nthread, thread);
4958 wallcycle_stop(wcycle, (grid_index < DO_Q ? ewcPME_SOLVE : ewcLJPME));
4960 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
4970 wallcycle_start(wcycle, ewcPME_FFT);
4972 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
4976 wallcycle_stop(wcycle, ewcPME_FFT);
4980 if (pme->nodeid == 0)
4982 ntot = pme->nkx*pme->nky*pme->nkz;
4983 npme = ntot*log((real)ntot)/log(2.0);
4984 inc_nrnb(nrnb, eNR_FFT, 2*npme);
4987 /* Note: this wallcycle region is closed below
4988 outside an OpenMP region, so take care if
4989 refactoring code here. */
4990 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
4993 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
4996 /* End of thread parallel section.
4997 * With MPI we have to synchronize here before gmx_sum_qgrid_dd.
5002 /* distribute local grid to all nodes */
5004 if (pme->nnodes > 1)
5006 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
5011 unwrap_periodic_pmegrid(pme, grid);
5013 /* interpolate forces for our local atoms */
5017 /* If we are running without parallelization,
5018 * atc->f is the actual force array, not a buffer,
5019 * therefore we should not clear it.
5021 lambda = grid_index < DO_Q ? lambda_q : lambda_lj;
5022 bClearF = (bFirst && PAR(cr));
5023 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
5024 for (thread = 0; thread < pme->nthread; thread++)
5026 gather_f_bsplines(pme, grid, bClearF, atc,
5027 &atc->spline[thread],
5028 pme->bFEP ? (grid_index % 2 == 0 ? 1.0-lambda : lambda) : 1.0);
5033 inc_nrnb(nrnb, eNR_GATHERFBSP,
5034 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
5035 /* Note: this wallcycle region is opened above inside an OpenMP
5036 region, so take care if refactoring code here. */
5037 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
5042 /* This should only be called on the master thread
5043 * and after the threads have synchronized.
5047 get_pme_ener_vir_q(pme, pme->nthread, &energy_AB[grid_index], vir_AB[grid_index]);
5051 get_pme_ener_vir_lj(pme, pme->nthread, &energy_AB[grid_index], vir_AB[grid_index]);
5055 } /* of grid_index-loop */
5057 /* For Lorentz-Berthelot combination rules in LJ-PME, we need to calculate
5060 if ((flags & GMX_PME_DO_LJ) && pme->ljpme_combination_rule == eljpmeLB)
5062 /* Loop over A- and B-state if we are doing FEP */
5063 for (fep_state = 0; fep_state < fep_states_lj; ++fep_state)
5065 real *local_c6 = NULL, *local_sigma = NULL, *RedistC6 = NULL, *RedistSigma = NULL;
5066 if (pme->nnodes == 1)
5068 if (pme->lb_buf1 == NULL)
5070 pme->lb_buf_nalloc = pme->atc[0].n;
5071 snew(pme->lb_buf1, pme->lb_buf_nalloc);
5073 pme->atc[0].coefficient = pme->lb_buf1;
5078 local_sigma = sigmaA;
5082 local_sigma = sigmaB;
5085 gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
5095 RedistSigma = sigmaA;
5099 RedistSigma = sigmaB;
5102 gmx_incons("Trying to access wrong FEP-state in LJ-PME routine");
5104 wallcycle_start(wcycle, ewcPME_REDISTXF);
5106 do_redist_pos_coeffs(pme, cr, start, homenr, bFirst, x, RedistC6);
5107 if (pme->lb_buf_nalloc < atc->n)
5109 pme->lb_buf_nalloc = atc->nalloc;
5110 srenew(pme->lb_buf1, pme->lb_buf_nalloc);
5111 srenew(pme->lb_buf2, pme->lb_buf_nalloc);
5113 local_c6 = pme->lb_buf1;
5114 for (i = 0; i < atc->n; ++i)
5116 local_c6[i] = atc->coefficient[i];
5120 do_redist_pos_coeffs(pme, cr, start, homenr, FALSE, x, RedistSigma);
5121 local_sigma = pme->lb_buf2;
5122 for (i = 0; i < atc->n; ++i)
5124 local_sigma[i] = atc->coefficient[i];
5128 wallcycle_stop(wcycle, ewcPME_REDISTXF);
5130 calc_initial_lb_coeffs(pme, local_c6, local_sigma);
5132 /*Seven terms in LJ-PME with LB, grid_index < 2 reserved for electrostatics*/
5133 for (grid_index = 2; grid_index < 9; ++grid_index)
5135 /* Unpack structure */
5136 pmegrid = &pme->pmegrid[grid_index];
5137 fftgrid = pme->fftgrid[grid_index];
5138 cfftgrid = pme->cfftgrid[grid_index];
5139 pfft_setup = pme->pfft_setup[grid_index];
5140 calc_next_lb_coeffs(pme, local_sigma);
5141 grid = pmegrid->grid.grid;
5144 if (flags & GMX_PME_SPREAD)
5146 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
5147 /* Spread the c6 on a grid */
5148 spread_on_grid(pme, &pme->atc[0], pmegrid, bFirst, TRUE, fftgrid, bDoSplines, grid_index);
5152 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
5155 inc_nrnb(nrnb, eNR_SPREADBSP,
5156 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
5157 if (pme->nthread == 1)
5159 wrap_periodic_pmegrid(pme, grid);
5160 /* sum contributions to local grid from other nodes */
5162 if (pme->nnodes > 1)
5164 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_FORWARD);
5168 copy_pmegrid_to_fftgrid(pme, grid, fftgrid, grid_index);
5170 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
5172 /*Here we start a large thread parallel region*/
5173 #pragma omp parallel num_threads(pme->nthread) private(thread)
5175 thread = gmx_omp_get_thread_num();
5176 if (flags & GMX_PME_SOLVE)
5181 wallcycle_start(wcycle, ewcPME_FFT);
5184 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
5188 wallcycle_stop(wcycle, ewcPME_FFT);
5195 if (flags & GMX_PME_SOLVE)
5197 /* solve in k-space for our local cells */
5198 #pragma omp parallel num_threads(pme->nthread) private(thread)
5201 thread = gmx_omp_get_thread_num();
5204 wallcycle_start(wcycle, ewcLJPME);
5208 solve_pme_lj_yzx(pme, &pme->cfftgrid[2], TRUE, ewaldcoeff_lj,
5209 box[XX][XX]*box[YY][YY]*box[ZZ][ZZ],
5211 pme->nthread, thread);
5214 wallcycle_stop(wcycle, ewcLJPME);
5216 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
5223 /* This should only be called on the master thread and
5224 * after the threads have synchronized.
5226 get_pme_ener_vir_lj(pme, pme->nthread, &energy_AB[2+fep_state], vir_AB[2+fep_state]);
5231 bFirst = !(flags & GMX_PME_DO_COULOMB);
5232 calc_initial_lb_coeffs(pme, local_c6, local_sigma);
5233 for (grid_index = 8; grid_index >= 2; --grid_index)
5235 /* Unpack structure */
5236 pmegrid = &pme->pmegrid[grid_index];
5237 fftgrid = pme->fftgrid[grid_index];
5238 cfftgrid = pme->cfftgrid[grid_index];
5239 pfft_setup = pme->pfft_setup[grid_index];
5240 grid = pmegrid->grid.grid;
5241 calc_next_lb_coeffs(pme, local_sigma);
5243 #pragma omp parallel num_threads(pme->nthread) private(thread)
5245 thread = gmx_omp_get_thread_num();
5250 wallcycle_start(wcycle, ewcPME_FFT);
5253 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
5257 wallcycle_stop(wcycle, ewcPME_FFT);
5261 if (pme->nodeid == 0)
5263 ntot = pme->nkx*pme->nky*pme->nkz;
5264 npme = ntot*log((real)ntot)/log(2.0);
5265 inc_nrnb(nrnb, eNR_FFT, 2*npme);
5267 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
5270 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, grid_index, pme->nthread, thread);
5272 } /*#pragma omp parallel*/
5274 /* distribute local grid to all nodes */
5276 if (pme->nnodes > 1)
5278 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_GRID_BACKWARD);
5283 unwrap_periodic_pmegrid(pme, grid);
5285 /* interpolate forces for our local atoms */
5287 bClearF = (bFirst && PAR(cr));
5288 scale = pme->bFEP ? (fep_state < 1 ? 1.0-lambda_lj : lambda_lj) : 1.0;
5289 scale *= lb_scale_factor[grid_index-2];
5290 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
5291 for (thread = 0; thread < pme->nthread; thread++)
5293 gather_f_bsplines(pme, grid, bClearF, &pme->atc[0],
5294 &pme->atc[0].spline[thread],
5299 inc_nrnb(nrnb, eNR_GATHERFBSP,
5300 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
5301 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
5304 } /* for (grid_index = 8; grid_index >= 2; --grid_index) */
5306 } /* for (fep_state = 0; fep_state < fep_states_lj; ++fep_state) */
5307 } /* if ((flags & GMX_PME_DO_LJ) && pme->ljpme_combination_rule == eljpmeLB) */
5309 if (bCalcF && pme->nnodes > 1)
5311 wallcycle_start(wcycle, ewcPME_REDISTXF);
5312 for (d = 0; d < pme->ndecompdim; d++)
5315 if (d == pme->ndecompdim - 1)
5322 n_d = pme->atc[d+1].n;
5323 f_d = pme->atc[d+1].f;
5325 if (DOMAINDECOMP(cr))
5327 dd_pmeredist_f(pme, atc, n_d, f_d,
5328 d == pme->ndecompdim-1 && pme->bPPnode);
5332 wallcycle_stop(wcycle, ewcPME_REDISTXF);
5338 if (flags & GMX_PME_DO_COULOMB)
5342 *energy_q = energy_AB[0];
5343 m_add(vir_q, vir_AB[0], vir_q);
5347 *energy_q = (1.0-lambda_q)*energy_AB[0] + lambda_q*energy_AB[1];
5348 *dvdlambda_q += energy_AB[1] - energy_AB[0];
5349 for (i = 0; i < DIM; i++)
5351 for (j = 0; j < DIM; j++)
5353 vir_q[i][j] += (1.0-lambda_q)*vir_AB[0][i][j] +
5354 lambda_q*vir_AB[1][i][j];
5360 fprintf(debug, "Electrostatic PME mesh energy: %g\n", *energy_q);
5368 if (flags & GMX_PME_DO_LJ)
5372 *energy_lj = energy_AB[2];
5373 m_add(vir_lj, vir_AB[2], vir_lj);
5377 *energy_lj = (1.0-lambda_lj)*energy_AB[2] + lambda_lj*energy_AB[3];
5378 *dvdlambda_lj += energy_AB[3] - energy_AB[2];
5379 for (i = 0; i < DIM; i++)
5381 for (j = 0; j < DIM; j++)
5383 vir_lj[i][j] += (1.0-lambda_lj)*vir_AB[2][i][j] + lambda_lj*vir_AB[3][i][j];
5389 fprintf(debug, "Lennard-Jones PME mesh energy: %g\n", *energy_lj);