1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
36 /* IMPORTANT FOR DEVELOPERS:
38 * Triclinic pme stuff isn't entirely trivial, and we've experienced
39 * some bugs during development (many of them due to me). To avoid
40 * this in the future, please check the following things if you make
41 * changes in this file:
43 * 1. You should obtain identical (at least to the PME precision)
44 * energies, forces, and virial for
45 * a rectangular box and a triclinic one where the z (or y) axis is
46 * tilted a whole box side. For instance you could use these boxes:
48 * rectangular triclinic
53 * 2. You should check the energy conservation in a triclinic box.
55 * It might seem an overkill, but better safe than sorry.
63 #include "gromacs/fft/parallel_3dfft.h"
64 #include "gromacs/utility/gmxmpi.h"
73 #include "gmxcomplex.h"
77 #include "gmx_fatal.h"
82 #include "gmx_wallcycle.h"
84 #include "gmx_cyclecounter.h"
88 /* Single precision, with SSE2 or higher available */
89 #if defined(GMX_X86_SSE2) && !defined(GMX_DOUBLE)
91 #include "gmx_x86_simd_single.h"
94 /* Some old AMD processors could have problems with unaligned loads+stores */
96 #define PME_SSE_UNALIGNED
101 /* #define PRT_FORCE */
102 /* conditions for on the fly time-measurement */
103 /* #define TAKETIME (step > 1 && timesteps < 10) */
104 #define TAKETIME FALSE
106 /* #define PME_TIME_THREADS */
109 #define mpi_type MPI_DOUBLE
111 #define mpi_type MPI_FLOAT
114 /* GMX_CACHE_SEP should be a multiple of 16 to preserve alignment */
115 #define GMX_CACHE_SEP 64
117 /* We only define a maximum to be able to use local arrays without allocation.
118 * An order larger than 12 should never be needed, even for test cases.
119 * If needed it can be changed here.
121 #define PME_ORDER_MAX 12
123 /* Internal datastructures */
129 int recv_size; /* Receive buffer width, used with OpenMP */
140 int *send_id, *recv_id;
141 int send_size; /* Send buffer width, used with OpenMP */
142 pme_grid_comm_t *comm_data;
148 int *n; /* Cumulative counts of the number of particles per thread */
149 int nalloc; /* Allocation size of i */
150 int *i; /* Particle indices ordered on thread index (n) */
164 int dimind; /* The index of the dimension, 0=x, 1=y */
171 int *node_dest; /* The nodes to send x and q to with DD */
172 int *node_src; /* The nodes to receive x and q from with DD */
173 int *buf_index; /* Index for commnode into the buffers */
180 int *count; /* The number of atoms to send to each node */
182 int *rcount; /* The number of atoms to receive */
189 gmx_bool bSpread; /* These coordinates are used for spreading */
192 rvec *fractx; /* Fractional coordinate relative to the
193 * lower cell boundary
196 int *thread_idx; /* Which thread should spread which charge */
197 thread_plist_t *thread_plist;
198 splinedata_t *spline;
205 ivec ci; /* The spatial location of this grid */
206 ivec n; /* The used size of *grid, including order-1 */
207 ivec offset; /* The grid offset from the full node grid */
208 int order; /* PME spreading order */
209 ivec s; /* The allocated size of *grid, s >= n */
210 real *grid; /* The grid local thread, size n */
214 pmegrid_t grid; /* The full node grid (non thread-local) */
215 int nthread; /* The number of threads operating on this grid */
216 ivec nc; /* The local spatial decomposition over the threads */
217 pmegrid_t *grid_th; /* Array of grids for each thread */
218 real *grid_all; /* Allocated array for the grids in *grid_th */
219 int **g2t; /* The grid to thread index */
220 ivec nthread_comm; /* The number of threads to communicate with */
226 /* Masks for SSE aligned spreading and gathering */
227 __m128 mask_SSE0[6], mask_SSE1[6];
229 int dummy; /* C89 requires that struct has at least one member */
234 /* work data for solve_pme */
250 typedef struct gmx_pme {
251 int ndecompdim; /* The number of decomposition dimensions */
252 int nodeid; /* Our nodeid in mpi->mpi_comm */
255 int nnodes; /* The number of nodes doing PME */
260 MPI_Comm mpi_comm_d[2]; /* Indexed on dimension, 0=x, 1=y */
262 MPI_Datatype rvec_mpi; /* the pme vector's MPI type */
265 gmx_bool bUseThreads; /* Does any of the PME ranks have nthread>1 ? */
266 int nthread; /* The number of threads doing PME on our rank */
268 gmx_bool bPPnode; /* Node also does particle-particle forces */
269 gmx_bool bFEP; /* Compute Free energy contribution */
270 int nkx, nky, nkz; /* Grid dimensions */
271 gmx_bool bP3M; /* Do P3M: optimize the influence function */
275 pmegrids_t pmegridA; /* Grids on which we do spreading/interpolation, includes overlap */
277 /* The PME charge spreading grid sizes/strides, includes pme_order-1 */
278 int pmegrid_nx, pmegrid_ny, pmegrid_nz;
279 /* pmegrid_nz might be larger than strictly necessary to ensure
280 * memory alignment, pmegrid_nz_base gives the real base size.
283 /* The local PME grid starting indices */
284 int pmegrid_start_ix, pmegrid_start_iy, pmegrid_start_iz;
286 /* Work data for spreading and gathering */
287 pme_spline_work_t *spline_work;
289 real *fftgridA; /* Grids for FFT. With 1D FFT decomposition this can be a pointer */
290 real *fftgridB; /* inside the interpolation grid, but separate for 2D PME decomp. */
291 int fftgrid_nx, fftgrid_ny, fftgrid_nz;
293 t_complex *cfftgridA; /* Grids for complex FFT data */
294 t_complex *cfftgridB;
295 int cfftgrid_nx, cfftgrid_ny, cfftgrid_nz;
297 gmx_parallel_3dfft_t pfft_setupA;
298 gmx_parallel_3dfft_t pfft_setupB;
300 int *nnx, *nny, *nnz;
301 real *fshx, *fshy, *fshz;
303 pme_atomcomm_t atc[2]; /* Indexed on decomposition index */
307 pme_overlap_t overlap[2]; /* Indexed on dimension, 0=x, 1=y */
309 pme_atomcomm_t atc_energy; /* Only for gmx_pme_calc_energy */
311 rvec *bufv; /* Communication buffer */
312 real *bufr; /* Communication buffer */
313 int buf_nalloc; /* The communication buffer size */
315 /* thread local work data for solve_pme */
318 /* Work data for PME_redist */
319 gmx_bool redist_init;
327 int redist_buf_nalloc;
329 /* Work data for sum_qgrid */
330 real * sum_qgrid_tmp;
331 real * sum_qgrid_dd_tmp;
335 static void calc_interpolation_idx(gmx_pme_t pme, pme_atomcomm_t *atc,
336 int start, int end, int thread)
339 int *idxptr, tix, tiy, tiz;
340 real *xptr, *fptr, tx, ty, tz;
341 real rxx, ryx, ryy, rzx, rzy, rzz;
343 int start_ix, start_iy, start_iz;
344 int *g2tx, *g2ty, *g2tz;
346 int *thread_idx = NULL;
347 thread_plist_t *tpl = NULL;
355 start_ix = pme->pmegrid_start_ix;
356 start_iy = pme->pmegrid_start_iy;
357 start_iz = pme->pmegrid_start_iz;
359 rxx = pme->recipbox[XX][XX];
360 ryx = pme->recipbox[YY][XX];
361 ryy = pme->recipbox[YY][YY];
362 rzx = pme->recipbox[ZZ][XX];
363 rzy = pme->recipbox[ZZ][YY];
364 rzz = pme->recipbox[ZZ][ZZ];
366 g2tx = pme->pmegridA.g2t[XX];
367 g2ty = pme->pmegridA.g2t[YY];
368 g2tz = pme->pmegridA.g2t[ZZ];
370 bThreads = (atc->nthread > 1);
373 thread_idx = atc->thread_idx;
375 tpl = &atc->thread_plist[thread];
377 for (i = 0; i < atc->nthread; i++)
383 for (i = start; i < end; i++)
386 idxptr = atc->idx[i];
387 fptr = atc->fractx[i];
389 /* Fractional coordinates along box vectors, add 2.0 to make 100% sure we are positive for triclinic boxes */
390 tx = nx * ( xptr[XX] * rxx + xptr[YY] * ryx + xptr[ZZ] * rzx + 2.0 );
391 ty = ny * ( xptr[YY] * ryy + xptr[ZZ] * rzy + 2.0 );
392 tz = nz * ( xptr[ZZ] * rzz + 2.0 );
398 /* Because decomposition only occurs in x and y,
399 * we never have a fraction correction in z.
401 fptr[XX] = tx - tix + pme->fshx[tix];
402 fptr[YY] = ty - tiy + pme->fshy[tiy];
405 idxptr[XX] = pme->nnx[tix];
406 idxptr[YY] = pme->nny[tiy];
407 idxptr[ZZ] = pme->nnz[tiz];
410 range_check(idxptr[XX], 0, pme->pmegrid_nx);
411 range_check(idxptr[YY], 0, pme->pmegrid_ny);
412 range_check(idxptr[ZZ], 0, pme->pmegrid_nz);
417 thread_i = g2tx[idxptr[XX]] + g2ty[idxptr[YY]] + g2tz[idxptr[ZZ]];
418 thread_idx[i] = thread_i;
425 /* Make a list of particle indices sorted on thread */
427 /* Get the cumulative count */
428 for (i = 1; i < atc->nthread; i++)
430 tpl_n[i] += tpl_n[i-1];
432 /* The current implementation distributes particles equally
433 * over the threads, so we could actually allocate for that
434 * in pme_realloc_atomcomm_things.
436 if (tpl_n[atc->nthread-1] > tpl->nalloc)
438 tpl->nalloc = over_alloc_large(tpl_n[atc->nthread-1]);
439 srenew(tpl->i, tpl->nalloc);
441 /* Set tpl_n to the cumulative start */
442 for (i = atc->nthread-1; i >= 1; i--)
444 tpl_n[i] = tpl_n[i-1];
448 /* Fill our thread local array with indices sorted on thread */
449 for (i = start; i < end; i++)
451 tpl->i[tpl_n[atc->thread_idx[i]]++] = i;
453 /* Now tpl_n contains the cummulative count again */
457 static void make_thread_local_ind(pme_atomcomm_t *atc,
458 int thread, splinedata_t *spline)
460 int n, t, i, start, end;
463 /* Combine the indices made by each thread into one index */
467 for (t = 0; t < atc->nthread; t++)
469 tpl = &atc->thread_plist[t];
470 /* Copy our part (start - end) from the list of thread t */
473 start = tpl->n[thread-1];
475 end = tpl->n[thread];
476 for (i = start; i < end; i++)
478 spline->ind[n++] = tpl->i[i];
486 static void pme_calc_pidx(int start, int end,
487 matrix recipbox, rvec x[],
488 pme_atomcomm_t *atc, int *count)
493 real rxx, ryx, rzx, ryy, rzy;
496 /* Calculate PME task index (pidx) for each grid index.
497 * Here we always assign equally sized slabs to each node
498 * for load balancing reasons (the PME grid spacing is not used).
504 /* Reset the count */
505 for (i = 0; i < nslab; i++)
510 if (atc->dimind == 0)
512 rxx = recipbox[XX][XX];
513 ryx = recipbox[YY][XX];
514 rzx = recipbox[ZZ][XX];
515 /* Calculate the node index in x-dimension */
516 for (i = start; i < end; i++)
519 /* Fractional coordinates along box vectors */
520 s = nslab*(xptr[XX]*rxx + xptr[YY]*ryx + xptr[ZZ]*rzx);
521 si = (int)(s + 2*nslab) % nslab;
528 ryy = recipbox[YY][YY];
529 rzy = recipbox[ZZ][YY];
530 /* Calculate the node index in y-dimension */
531 for (i = start; i < end; i++)
534 /* Fractional coordinates along box vectors */
535 s = nslab*(xptr[YY]*ryy + xptr[ZZ]*rzy);
536 si = (int)(s + 2*nslab) % nslab;
543 static void pme_calc_pidx_wrapper(int natoms, matrix recipbox, rvec x[],
546 int nthread, thread, slab;
548 nthread = atc->nthread;
550 #pragma omp parallel for num_threads(nthread) schedule(static)
551 for (thread = 0; thread < nthread; thread++)
553 pme_calc_pidx(natoms* thread /nthread,
554 natoms*(thread+1)/nthread,
555 recipbox, x, atc, atc->count_thread[thread]);
557 /* Non-parallel reduction, since nslab is small */
559 for (thread = 1; thread < nthread; thread++)
561 for (slab = 0; slab < atc->nslab; slab++)
563 atc->count_thread[0][slab] += atc->count_thread[thread][slab];
568 static void realloc_splinevec(splinevec th, real **ptr_z, int nalloc)
570 const int padding = 4;
573 srenew(th[XX], nalloc);
574 srenew(th[YY], nalloc);
575 /* In z we add padding, this is only required for the aligned SSE code */
576 srenew(*ptr_z, nalloc+2*padding);
577 th[ZZ] = *ptr_z + padding;
579 for (i = 0; i < padding; i++)
582 (*ptr_z)[padding+nalloc+i] = 0;
586 static void pme_realloc_splinedata(splinedata_t *spline, pme_atomcomm_t *atc)
590 srenew(spline->ind, atc->nalloc);
591 /* Initialize the index to identity so it works without threads */
592 for (i = 0; i < atc->nalloc; i++)
597 realloc_splinevec(spline->theta, &spline->ptr_theta_z,
598 atc->pme_order*atc->nalloc);
599 realloc_splinevec(spline->dtheta, &spline->ptr_dtheta_z,
600 atc->pme_order*atc->nalloc);
603 static void pme_realloc_atomcomm_things(pme_atomcomm_t *atc)
605 int nalloc_old, i, j, nalloc_tpl;
607 /* We have to avoid a NULL pointer for atc->x to avoid
608 * possible fatal errors in MPI routines.
610 if (atc->n > atc->nalloc || atc->nalloc == 0)
612 nalloc_old = atc->nalloc;
613 atc->nalloc = over_alloc_dd(max(atc->n, 1));
617 srenew(atc->x, atc->nalloc);
618 srenew(atc->q, atc->nalloc);
619 srenew(atc->f, atc->nalloc);
620 for (i = nalloc_old; i < atc->nalloc; i++)
622 clear_rvec(atc->f[i]);
627 srenew(atc->fractx, atc->nalloc);
628 srenew(atc->idx, atc->nalloc);
630 if (atc->nthread > 1)
632 srenew(atc->thread_idx, atc->nalloc);
635 for (i = 0; i < atc->nthread; i++)
637 pme_realloc_splinedata(&atc->spline[i], atc);
643 static void pmeredist_pd(gmx_pme_t pme, gmx_bool forw,
644 int n, gmx_bool bXF, rvec *x_f, real *charge,
646 /* Redistribute particle data for PME calculation */
647 /* domain decomposition by x coordinate */
652 if (FALSE == pme->redist_init)
654 snew(pme->scounts, atc->nslab);
655 snew(pme->rcounts, atc->nslab);
656 snew(pme->sdispls, atc->nslab);
657 snew(pme->rdispls, atc->nslab);
658 snew(pme->sidx, atc->nslab);
659 pme->redist_init = TRUE;
661 if (n > pme->redist_buf_nalloc)
663 pme->redist_buf_nalloc = over_alloc_dd(n);
664 srenew(pme->redist_buf, pme->redist_buf_nalloc*DIM);
672 /* forward, redistribution from pp to pme */
674 /* Calculate send counts and exchange them with other nodes */
675 for (i = 0; (i < atc->nslab); i++)
679 for (i = 0; (i < n); i++)
681 pme->scounts[pme->idxa[i]]++;
683 MPI_Alltoall( pme->scounts, 1, MPI_INT, pme->rcounts, 1, MPI_INT, atc->mpi_comm);
685 /* Calculate send and receive displacements and index into send
690 for (i = 1; i < atc->nslab; i++)
692 pme->sdispls[i] = pme->sdispls[i-1]+pme->scounts[i-1];
693 pme->rdispls[i] = pme->rdispls[i-1]+pme->rcounts[i-1];
694 pme->sidx[i] = pme->sdispls[i];
696 /* Total # of particles to be received */
697 atc->n = pme->rdispls[atc->nslab-1] + pme->rcounts[atc->nslab-1];
699 pme_realloc_atomcomm_things(atc);
701 /* Copy particle coordinates into send buffer and exchange*/
702 for (i = 0; (i < n); i++)
704 ii = DIM*pme->sidx[pme->idxa[i]];
705 pme->sidx[pme->idxa[i]]++;
706 pme->redist_buf[ii+XX] = x_f[i][XX];
707 pme->redist_buf[ii+YY] = x_f[i][YY];
708 pme->redist_buf[ii+ZZ] = x_f[i][ZZ];
710 MPI_Alltoallv(pme->redist_buf, pme->scounts, pme->sdispls,
711 pme->rvec_mpi, atc->x, pme->rcounts, pme->rdispls,
712 pme->rvec_mpi, atc->mpi_comm);
716 /* Copy charge into send buffer and exchange*/
717 for (i = 0; i < atc->nslab; i++)
719 pme->sidx[i] = pme->sdispls[i];
721 for (i = 0; (i < n); i++)
723 ii = pme->sidx[pme->idxa[i]];
724 pme->sidx[pme->idxa[i]]++;
725 pme->redist_buf[ii] = charge[i];
727 MPI_Alltoallv(pme->redist_buf, pme->scounts, pme->sdispls, mpi_type,
728 atc->q, pme->rcounts, pme->rdispls, mpi_type,
731 else /* backward, redistribution from pme to pp */
733 MPI_Alltoallv(atc->f, pme->rcounts, pme->rdispls, pme->rvec_mpi,
734 pme->redist_buf, pme->scounts, pme->sdispls,
735 pme->rvec_mpi, atc->mpi_comm);
737 /* Copy data from receive buffer */
738 for (i = 0; i < atc->nslab; i++)
740 pme->sidx[i] = pme->sdispls[i];
742 for (i = 0; (i < n); i++)
744 ii = DIM*pme->sidx[pme->idxa[i]];
745 x_f[i][XX] += pme->redist_buf[ii+XX];
746 x_f[i][YY] += pme->redist_buf[ii+YY];
747 x_f[i][ZZ] += pme->redist_buf[ii+ZZ];
748 pme->sidx[pme->idxa[i]]++;
754 static void pme_dd_sendrecv(pme_atomcomm_t *atc,
755 gmx_bool bBackward, int shift,
756 void *buf_s, int nbyte_s,
757 void *buf_r, int nbyte_r)
763 if (bBackward == FALSE)
765 dest = atc->node_dest[shift];
766 src = atc->node_src[shift];
770 dest = atc->node_src[shift];
771 src = atc->node_dest[shift];
774 if (nbyte_s > 0 && nbyte_r > 0)
776 MPI_Sendrecv(buf_s, nbyte_s, MPI_BYTE,
778 buf_r, nbyte_r, MPI_BYTE,
780 atc->mpi_comm, &stat);
782 else if (nbyte_s > 0)
784 MPI_Send(buf_s, nbyte_s, MPI_BYTE,
788 else if (nbyte_r > 0)
790 MPI_Recv(buf_r, nbyte_r, MPI_BYTE,
792 atc->mpi_comm, &stat);
797 static void dd_pmeredist_x_q(gmx_pme_t pme,
798 int n, gmx_bool bX, rvec *x, real *charge,
801 int *commnode, *buf_index;
802 int nnodes_comm, i, nsend, local_pos, buf_pos, node, scount, rcount;
804 commnode = atc->node_dest;
805 buf_index = atc->buf_index;
807 nnodes_comm = min(2*atc->maxshift, atc->nslab-1);
810 for (i = 0; i < nnodes_comm; i++)
812 buf_index[commnode[i]] = nsend;
813 nsend += atc->count[commnode[i]];
817 if (atc->count[atc->nodeid] + nsend != n)
819 gmx_fatal(FARGS, "%d particles communicated to PME node %d are more than 2/3 times the cut-off out of the domain decomposition cell of their charge group in dimension %c.\n"
820 "This usually means that your system is not well equilibrated.",
821 n - (atc->count[atc->nodeid] + nsend),
822 pme->nodeid, 'x'+atc->dimind);
825 if (nsend > pme->buf_nalloc)
827 pme->buf_nalloc = over_alloc_dd(nsend);
828 srenew(pme->bufv, pme->buf_nalloc);
829 srenew(pme->bufr, pme->buf_nalloc);
832 atc->n = atc->count[atc->nodeid];
833 for (i = 0; i < nnodes_comm; i++)
835 scount = atc->count[commnode[i]];
836 /* Communicate the count */
839 fprintf(debug, "dimind %d PME node %d send to node %d: %d\n",
840 atc->dimind, atc->nodeid, commnode[i], scount);
842 pme_dd_sendrecv(atc, FALSE, i,
843 &scount, sizeof(int),
844 &atc->rcount[i], sizeof(int));
845 atc->n += atc->rcount[i];
848 pme_realloc_atomcomm_things(atc);
852 for (i = 0; i < n; i++)
855 if (node == atc->nodeid)
857 /* Copy direct to the receive buffer */
860 copy_rvec(x[i], atc->x[local_pos]);
862 atc->q[local_pos] = charge[i];
867 /* Copy to the send buffer */
870 copy_rvec(x[i], pme->bufv[buf_index[node]]);
872 pme->bufr[buf_index[node]] = charge[i];
878 for (i = 0; i < nnodes_comm; i++)
880 scount = atc->count[commnode[i]];
881 rcount = atc->rcount[i];
882 if (scount > 0 || rcount > 0)
886 /* Communicate the coordinates */
887 pme_dd_sendrecv(atc, FALSE, i,
888 pme->bufv[buf_pos], scount*sizeof(rvec),
889 atc->x[local_pos], rcount*sizeof(rvec));
891 /* Communicate the charges */
892 pme_dd_sendrecv(atc, FALSE, i,
893 pme->bufr+buf_pos, scount*sizeof(real),
894 atc->q+local_pos, rcount*sizeof(real));
896 local_pos += atc->rcount[i];
901 static void dd_pmeredist_f(gmx_pme_t pme, pme_atomcomm_t *atc,
905 int *commnode, *buf_index;
906 int nnodes_comm, local_pos, buf_pos, i, scount, rcount, node;
908 commnode = atc->node_dest;
909 buf_index = atc->buf_index;
911 nnodes_comm = min(2*atc->maxshift, atc->nslab-1);
913 local_pos = atc->count[atc->nodeid];
915 for (i = 0; i < nnodes_comm; i++)
917 scount = atc->rcount[i];
918 rcount = atc->count[commnode[i]];
919 if (scount > 0 || rcount > 0)
921 /* Communicate the forces */
922 pme_dd_sendrecv(atc, TRUE, i,
923 atc->f[local_pos], scount*sizeof(rvec),
924 pme->bufv[buf_pos], rcount*sizeof(rvec));
927 buf_index[commnode[i]] = buf_pos;
934 for (i = 0; i < n; i++)
937 if (node == atc->nodeid)
939 /* Add from the local force array */
940 rvec_inc(f[i], atc->f[local_pos]);
945 /* Add from the receive buffer */
946 rvec_inc(f[i], pme->bufv[buf_index[node]]);
953 for (i = 0; i < n; i++)
956 if (node == atc->nodeid)
958 /* Copy from the local force array */
959 copy_rvec(atc->f[local_pos], f[i]);
964 /* Copy from the receive buffer */
965 copy_rvec(pme->bufv[buf_index[node]], f[i]);
974 gmx_sum_qgrid_dd(gmx_pme_t pme, real *grid, int direction)
976 pme_overlap_t *overlap;
977 int send_index0, send_nindex;
978 int recv_index0, recv_nindex;
980 int i, j, k, ix, iy, iz, icnt;
981 int ipulse, send_id, recv_id, datasize;
983 real *sendptr, *recvptr;
985 /* Start with minor-rank communication. This is a bit of a pain since it is not contiguous */
986 overlap = &pme->overlap[1];
988 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
990 /* Since we have already (un)wrapped the overlap in the z-dimension,
991 * we only have to communicate 0 to nkz (not pmegrid_nz).
993 if (direction == GMX_SUM_QGRID_FORWARD)
995 send_id = overlap->send_id[ipulse];
996 recv_id = overlap->recv_id[ipulse];
997 send_index0 = overlap->comm_data[ipulse].send_index0;
998 send_nindex = overlap->comm_data[ipulse].send_nindex;
999 recv_index0 = overlap->comm_data[ipulse].recv_index0;
1000 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
1004 send_id = overlap->recv_id[ipulse];
1005 recv_id = overlap->send_id[ipulse];
1006 send_index0 = overlap->comm_data[ipulse].recv_index0;
1007 send_nindex = overlap->comm_data[ipulse].recv_nindex;
1008 recv_index0 = overlap->comm_data[ipulse].send_index0;
1009 recv_nindex = overlap->comm_data[ipulse].send_nindex;
1012 /* Copy data to contiguous send buffer */
1015 fprintf(debug, "PME send node %d %d -> %d grid start %d Communicating %d to %d\n",
1016 pme->nodeid, overlap->nodeid, send_id,
1017 pme->pmegrid_start_iy,
1018 send_index0-pme->pmegrid_start_iy,
1019 send_index0-pme->pmegrid_start_iy+send_nindex);
1022 for (i = 0; i < pme->pmegrid_nx; i++)
1025 for (j = 0; j < send_nindex; j++)
1027 iy = j + send_index0 - pme->pmegrid_start_iy;
1028 for (k = 0; k < pme->nkz; k++)
1031 overlap->sendbuf[icnt++] = grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz];
1036 datasize = pme->pmegrid_nx * pme->nkz;
1038 MPI_Sendrecv(overlap->sendbuf, send_nindex*datasize, GMX_MPI_REAL,
1040 overlap->recvbuf, recv_nindex*datasize, GMX_MPI_REAL,
1042 overlap->mpi_comm, &stat);
1044 /* Get data from contiguous recv buffer */
1047 fprintf(debug, "PME recv node %d %d <- %d grid start %d Communicating %d to %d\n",
1048 pme->nodeid, overlap->nodeid, recv_id,
1049 pme->pmegrid_start_iy,
1050 recv_index0-pme->pmegrid_start_iy,
1051 recv_index0-pme->pmegrid_start_iy+recv_nindex);
1054 for (i = 0; i < pme->pmegrid_nx; i++)
1057 for (j = 0; j < recv_nindex; j++)
1059 iy = j + recv_index0 - pme->pmegrid_start_iy;
1060 for (k = 0; k < pme->nkz; k++)
1063 if (direction == GMX_SUM_QGRID_FORWARD)
1065 grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz] += overlap->recvbuf[icnt++];
1069 grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz] = overlap->recvbuf[icnt++];
1076 /* Major dimension is easier, no copying required,
1077 * but we might have to sum to separate array.
1078 * Since we don't copy, we have to communicate up to pmegrid_nz,
1079 * not nkz as for the minor direction.
1081 overlap = &pme->overlap[0];
1083 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
1085 if (direction == GMX_SUM_QGRID_FORWARD)
1087 send_id = overlap->send_id[ipulse];
1088 recv_id = overlap->recv_id[ipulse];
1089 send_index0 = overlap->comm_data[ipulse].send_index0;
1090 send_nindex = overlap->comm_data[ipulse].send_nindex;
1091 recv_index0 = overlap->comm_data[ipulse].recv_index0;
1092 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
1093 recvptr = overlap->recvbuf;
1097 send_id = overlap->recv_id[ipulse];
1098 recv_id = overlap->send_id[ipulse];
1099 send_index0 = overlap->comm_data[ipulse].recv_index0;
1100 send_nindex = overlap->comm_data[ipulse].recv_nindex;
1101 recv_index0 = overlap->comm_data[ipulse].send_index0;
1102 recv_nindex = overlap->comm_data[ipulse].send_nindex;
1103 recvptr = grid + (recv_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1106 sendptr = grid + (send_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1107 datasize = pme->pmegrid_ny * pme->pmegrid_nz;
1111 fprintf(debug, "PME send node %d %d -> %d grid start %d Communicating %d to %d\n",
1112 pme->nodeid, overlap->nodeid, send_id,
1113 pme->pmegrid_start_ix,
1114 send_index0-pme->pmegrid_start_ix,
1115 send_index0-pme->pmegrid_start_ix+send_nindex);
1116 fprintf(debug, "PME recv node %d %d <- %d grid start %d Communicating %d to %d\n",
1117 pme->nodeid, overlap->nodeid, recv_id,
1118 pme->pmegrid_start_ix,
1119 recv_index0-pme->pmegrid_start_ix,
1120 recv_index0-pme->pmegrid_start_ix+recv_nindex);
1123 MPI_Sendrecv(sendptr, send_nindex*datasize, GMX_MPI_REAL,
1125 recvptr, recv_nindex*datasize, GMX_MPI_REAL,
1127 overlap->mpi_comm, &stat);
1129 /* ADD data from contiguous recv buffer */
1130 if (direction == GMX_SUM_QGRID_FORWARD)
1132 p = grid + (recv_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1133 for (i = 0; i < recv_nindex*datasize; i++)
1135 p[i] += overlap->recvbuf[i];
1144 copy_pmegrid_to_fftgrid(gmx_pme_t pme, real *pmegrid, real *fftgrid)
1146 ivec local_fft_ndata, local_fft_offset, local_fft_size;
1147 ivec local_pme_size;
1151 /* Dimensions should be identical for A/B grid, so we just use A here */
1152 gmx_parallel_3dfft_real_limits(pme->pfft_setupA,
1157 local_pme_size[0] = pme->pmegrid_nx;
1158 local_pme_size[1] = pme->pmegrid_ny;
1159 local_pme_size[2] = pme->pmegrid_nz;
1161 /* The fftgrid is always 'justified' to the lower-left corner of the PME grid,
1162 the offset is identical, and the PME grid always has more data (due to overlap)
1167 char fn[STRLEN], format[STRLEN];
1169 sprintf(fn, "pmegrid%d.pdb", pme->nodeid);
1170 fp = ffopen(fn, "w");
1171 sprintf(fn, "pmegrid%d.txt", pme->nodeid);
1172 fp2 = ffopen(fn, "w");
1173 sprintf(format, "%s%s\n", pdbformat, "%6.2f%6.2f");
1176 for (ix = 0; ix < local_fft_ndata[XX]; ix++)
1178 for (iy = 0; iy < local_fft_ndata[YY]; iy++)
1180 for (iz = 0; iz < local_fft_ndata[ZZ]; iz++)
1182 pmeidx = ix*(local_pme_size[YY]*local_pme_size[ZZ])+iy*(local_pme_size[ZZ])+iz;
1183 fftidx = ix*(local_fft_size[YY]*local_fft_size[ZZ])+iy*(local_fft_size[ZZ])+iz;
1184 fftgrid[fftidx] = pmegrid[pmeidx];
1186 val = 100*pmegrid[pmeidx];
1187 if (pmegrid[pmeidx] != 0)
1189 fprintf(fp, format, "ATOM", pmeidx, "CA", "GLY", ' ', pmeidx, ' ',
1190 5.0*ix, 5.0*iy, 5.0*iz, 1.0, val);
1192 if (pmegrid[pmeidx] != 0)
1194 fprintf(fp2, "%-12s %5d %5d %5d %12.5e\n",
1196 pme->pmegrid_start_ix + ix,
1197 pme->pmegrid_start_iy + iy,
1198 pme->pmegrid_start_iz + iz,
1214 static gmx_cycles_t omp_cyc_start()
1216 return gmx_cycles_read();
1219 static gmx_cycles_t omp_cyc_end(gmx_cycles_t c)
1221 return gmx_cycles_read() - c;
1226 copy_fftgrid_to_pmegrid(gmx_pme_t pme, const real *fftgrid, real *pmegrid,
1227 int nthread, int thread)
1229 ivec local_fft_ndata, local_fft_offset, local_fft_size;
1230 ivec local_pme_size;
1231 int ixy0, ixy1, ixy, ix, iy, iz;
1233 #ifdef PME_TIME_THREADS
1235 static double cs1 = 0;
1239 #ifdef PME_TIME_THREADS
1240 c1 = omp_cyc_start();
1242 /* Dimensions should be identical for A/B grid, so we just use A here */
1243 gmx_parallel_3dfft_real_limits(pme->pfft_setupA,
1248 local_pme_size[0] = pme->pmegrid_nx;
1249 local_pme_size[1] = pme->pmegrid_ny;
1250 local_pme_size[2] = pme->pmegrid_nz;
1252 /* The fftgrid is always 'justified' to the lower-left corner of the PME grid,
1253 the offset is identical, and the PME grid always has more data (due to overlap)
1255 ixy0 = ((thread )*local_fft_ndata[XX]*local_fft_ndata[YY])/nthread;
1256 ixy1 = ((thread+1)*local_fft_ndata[XX]*local_fft_ndata[YY])/nthread;
1258 for (ixy = ixy0; ixy < ixy1; ixy++)
1260 ix = ixy/local_fft_ndata[YY];
1261 iy = ixy - ix*local_fft_ndata[YY];
1263 pmeidx = (ix*local_pme_size[YY] + iy)*local_pme_size[ZZ];
1264 fftidx = (ix*local_fft_size[YY] + iy)*local_fft_size[ZZ];
1265 for (iz = 0; iz < local_fft_ndata[ZZ]; iz++)
1267 pmegrid[pmeidx+iz] = fftgrid[fftidx+iz];
1271 #ifdef PME_TIME_THREADS
1272 c1 = omp_cyc_end(c1);
1277 printf("copy %.2f\n", cs1*1e-9);
1286 wrap_periodic_pmegrid(gmx_pme_t pme, real *pmegrid)
1288 int nx, ny, nz, pnx, pny, pnz, ny_x, overlap, ix, iy, iz;
1294 pnx = pme->pmegrid_nx;
1295 pny = pme->pmegrid_ny;
1296 pnz = pme->pmegrid_nz;
1298 overlap = pme->pme_order - 1;
1300 /* Add periodic overlap in z */
1301 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1303 for (iy = 0; iy < pme->pmegrid_ny; iy++)
1305 for (iz = 0; iz < overlap; iz++)
1307 pmegrid[(ix*pny+iy)*pnz+iz] +=
1308 pmegrid[(ix*pny+iy)*pnz+nz+iz];
1313 if (pme->nnodes_minor == 1)
1315 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1317 for (iy = 0; iy < overlap; iy++)
1319 for (iz = 0; iz < nz; iz++)
1321 pmegrid[(ix*pny+iy)*pnz+iz] +=
1322 pmegrid[(ix*pny+ny+iy)*pnz+iz];
1328 if (pme->nnodes_major == 1)
1330 ny_x = (pme->nnodes_minor == 1 ? ny : pme->pmegrid_ny);
1332 for (ix = 0; ix < overlap; ix++)
1334 for (iy = 0; iy < ny_x; iy++)
1336 for (iz = 0; iz < nz; iz++)
1338 pmegrid[(ix*pny+iy)*pnz+iz] +=
1339 pmegrid[((nx+ix)*pny+iy)*pnz+iz];
1348 unwrap_periodic_pmegrid(gmx_pme_t pme, real *pmegrid)
1350 int nx, ny, nz, pnx, pny, pnz, ny_x, overlap, ix;
1356 pnx = pme->pmegrid_nx;
1357 pny = pme->pmegrid_ny;
1358 pnz = pme->pmegrid_nz;
1360 overlap = pme->pme_order - 1;
1362 if (pme->nnodes_major == 1)
1364 ny_x = (pme->nnodes_minor == 1 ? ny : pme->pmegrid_ny);
1366 for (ix = 0; ix < overlap; ix++)
1370 for (iy = 0; iy < ny_x; iy++)
1372 for (iz = 0; iz < nz; iz++)
1374 pmegrid[((nx+ix)*pny+iy)*pnz+iz] =
1375 pmegrid[(ix*pny+iy)*pnz+iz];
1381 if (pme->nnodes_minor == 1)
1383 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1384 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1388 for (iy = 0; iy < overlap; iy++)
1390 for (iz = 0; iz < nz; iz++)
1392 pmegrid[(ix*pny+ny+iy)*pnz+iz] =
1393 pmegrid[(ix*pny+iy)*pnz+iz];
1399 /* Copy periodic overlap in z */
1400 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1401 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1405 for (iy = 0; iy < pme->pmegrid_ny; iy++)
1407 for (iz = 0; iz < overlap; iz++)
1409 pmegrid[(ix*pny+iy)*pnz+nz+iz] =
1410 pmegrid[(ix*pny+iy)*pnz+iz];
1416 static void clear_grid(int nx, int ny, int nz, real *grid,
1418 int fx, int fy, int fz,
1422 int fsx, fsy, fsz, gx, gy, gz, g0x, g0y, x, y, z;
1425 nc = 2 + (order - 2)/FLBS;
1426 ncz = 2 + (order - 2)/FLBSZ;
1428 for (fsx = fx; fsx < fx+nc; fsx++)
1430 for (fsy = fy; fsy < fy+nc; fsy++)
1432 for (fsz = fz; fsz < fz+ncz; fsz++)
1434 flind = (fsx*fs[YY] + fsy)*fs[ZZ] + fsz;
1435 if (flag[flind] == 0)
1440 g0x = (gx*ny + gy)*nz + gz;
1441 for (x = 0; x < FLBS; x++)
1444 for (y = 0; y < FLBS; y++)
1446 for (z = 0; z < FLBSZ; z++)
1462 /* This has to be a macro to enable full compiler optimization with xlC (and probably others too) */
1463 #define DO_BSPLINE(order) \
1464 for (ithx = 0; (ithx < order); ithx++) \
1466 index_x = (i0+ithx)*pny*pnz; \
1467 valx = qn*thx[ithx]; \
1469 for (ithy = 0; (ithy < order); ithy++) \
1471 valxy = valx*thy[ithy]; \
1472 index_xy = index_x+(j0+ithy)*pnz; \
1474 for (ithz = 0; (ithz < order); ithz++) \
1476 index_xyz = index_xy+(k0+ithz); \
1477 grid[index_xyz] += valxy*thz[ithz]; \
1483 static void spread_q_bsplines_thread(pmegrid_t *pmegrid,
1484 pme_atomcomm_t *atc, splinedata_t *spline,
1485 pme_spline_work_t *work)
1488 /* spread charges from home atoms to local grid */
1491 int b, i, nn, n, ithx, ithy, ithz, i0, j0, k0;
1493 int order, norder, index_x, index_xy, index_xyz;
1494 real valx, valxy, qn;
1495 real *thx, *thy, *thz;
1496 int localsize, bndsize;
1497 int pnx, pny, pnz, ndatatot;
1498 int offx, offy, offz;
1500 pnx = pmegrid->s[XX];
1501 pny = pmegrid->s[YY];
1502 pnz = pmegrid->s[ZZ];
1504 offx = pmegrid->offset[XX];
1505 offy = pmegrid->offset[YY];
1506 offz = pmegrid->offset[ZZ];
1508 ndatatot = pnx*pny*pnz;
1509 grid = pmegrid->grid;
1510 for (i = 0; i < ndatatot; i++)
1515 order = pmegrid->order;
1517 for (nn = 0; nn < spline->n; nn++)
1519 n = spline->ind[nn];
1524 idxptr = atc->idx[n];
1527 i0 = idxptr[XX] - offx;
1528 j0 = idxptr[YY] - offy;
1529 k0 = idxptr[ZZ] - offz;
1531 thx = spline->theta[XX] + norder;
1532 thy = spline->theta[YY] + norder;
1533 thz = spline->theta[ZZ] + norder;
1539 #ifdef PME_SSE_UNALIGNED
1540 #define PME_SPREAD_SSE_ORDER4
1542 #define PME_SPREAD_SSE_ALIGNED
1545 #include "pme_sse_single.h"
1552 #define PME_SPREAD_SSE_ALIGNED
1554 #include "pme_sse_single.h"
1567 static void set_grid_alignment(int *pmegrid_nz, int pme_order)
1571 #ifndef PME_SSE_UNALIGNED
1576 /* Round nz up to a multiple of 4 to ensure alignment */
1577 *pmegrid_nz = ((*pmegrid_nz + 3) & ~3);
1582 static void set_gridsize_alignment(int *gridsize, int pme_order)
1585 #ifndef PME_SSE_UNALIGNED
1588 /* Add extra elements to ensured aligned operations do not go
1589 * beyond the allocated grid size.
1590 * Note that for pme_order=5, the pme grid z-size alignment
1591 * ensures that we will not go beyond the grid size.
1599 static void pmegrid_init(pmegrid_t *grid,
1600 int cx, int cy, int cz,
1601 int x0, int y0, int z0,
1602 int x1, int y1, int z1,
1603 gmx_bool set_alignment,
1612 grid->offset[XX] = x0;
1613 grid->offset[YY] = y0;
1614 grid->offset[ZZ] = z0;
1615 grid->n[XX] = x1 - x0 + pme_order - 1;
1616 grid->n[YY] = y1 - y0 + pme_order - 1;
1617 grid->n[ZZ] = z1 - z0 + pme_order - 1;
1618 copy_ivec(grid->n, grid->s);
1621 set_grid_alignment(&nz, pme_order);
1626 else if (nz != grid->s[ZZ])
1628 gmx_incons("pmegrid_init call with an unaligned z size");
1631 grid->order = pme_order;
1634 gridsize = grid->s[XX]*grid->s[YY]*grid->s[ZZ];
1635 set_gridsize_alignment(&gridsize, pme_order);
1636 snew_aligned(grid->grid, gridsize, 16);
1644 static int div_round_up(int enumerator, int denominator)
1646 return (enumerator + denominator - 1)/denominator;
1649 static void make_subgrid_division(const ivec n, int ovl, int nthread,
1652 int gsize_opt, gsize;
1657 for (nsx = 1; nsx <= nthread; nsx++)
1659 if (nthread % nsx == 0)
1661 for (nsy = 1; nsy <= nthread; nsy++)
1663 if (nsx*nsy <= nthread && nthread % (nsx*nsy) == 0)
1665 nsz = nthread/(nsx*nsy);
1667 /* Determine the number of grid points per thread */
1669 (div_round_up(n[XX], nsx) + ovl)*
1670 (div_round_up(n[YY], nsy) + ovl)*
1671 (div_round_up(n[ZZ], nsz) + ovl);
1673 /* Minimize the number of grids points per thread
1674 * and, secondarily, the number of cuts in minor dimensions.
1676 if (gsize_opt == -1 ||
1677 gsize < gsize_opt ||
1678 (gsize == gsize_opt &&
1679 (nsz < nsub[ZZ] || (nsz == nsub[ZZ] && nsy < nsub[YY]))))
1691 env = getenv("GMX_PME_THREAD_DIVISION");
1694 sscanf(env, "%d %d %d", &nsub[XX], &nsub[YY], &nsub[ZZ]);
1697 if (nsub[XX]*nsub[YY]*nsub[ZZ] != nthread)
1699 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);
1703 static void pmegrids_init(pmegrids_t *grids,
1704 int nx, int ny, int nz, int nz_base,
1706 gmx_bool bUseThreads,
1711 ivec n, n_base, g0, g1;
1712 int t, x, y, z, d, i, tfac;
1713 int max_comm_lines = -1;
1715 n[XX] = nx - (pme_order - 1);
1716 n[YY] = ny - (pme_order - 1);
1717 n[ZZ] = nz - (pme_order - 1);
1719 copy_ivec(n, n_base);
1720 n_base[ZZ] = nz_base;
1722 pmegrid_init(&grids->grid, 0, 0, 0, 0, 0, 0, n[XX], n[YY], n[ZZ], FALSE, pme_order,
1725 grids->nthread = nthread;
1727 make_subgrid_division(n_base, pme_order-1, grids->nthread, grids->nc);
1734 for (d = 0; d < DIM; d++)
1736 nst[d] = div_round_up(n[d], grids->nc[d]) + pme_order - 1;
1738 set_grid_alignment(&nst[ZZ], pme_order);
1742 fprintf(debug, "pmegrid thread local division: %d x %d x %d\n",
1743 grids->nc[XX], grids->nc[YY], grids->nc[ZZ]);
1744 fprintf(debug, "pmegrid %d %d %d max thread pmegrid %d %d %d\n",
1746 nst[XX], nst[YY], nst[ZZ]);
1749 snew(grids->grid_th, grids->nthread);
1751 gridsize = nst[XX]*nst[YY]*nst[ZZ];
1752 set_gridsize_alignment(&gridsize, pme_order);
1753 snew_aligned(grids->grid_all,
1754 grids->nthread*gridsize+(grids->nthread+1)*GMX_CACHE_SEP,
1757 for (x = 0; x < grids->nc[XX]; x++)
1759 for (y = 0; y < grids->nc[YY]; y++)
1761 for (z = 0; z < grids->nc[ZZ]; z++)
1763 pmegrid_init(&grids->grid_th[t],
1765 (n[XX]*(x ))/grids->nc[XX],
1766 (n[YY]*(y ))/grids->nc[YY],
1767 (n[ZZ]*(z ))/grids->nc[ZZ],
1768 (n[XX]*(x+1))/grids->nc[XX],
1769 (n[YY]*(y+1))/grids->nc[YY],
1770 (n[ZZ]*(z+1))/grids->nc[ZZ],
1773 grids->grid_all+GMX_CACHE_SEP+t*(gridsize+GMX_CACHE_SEP));
1781 grids->grid_th = NULL;
1784 snew(grids->g2t, DIM);
1786 for (d = DIM-1; d >= 0; d--)
1788 snew(grids->g2t[d], n[d]);
1790 for (i = 0; i < n[d]; i++)
1792 /* The second check should match the parameters
1793 * of the pmegrid_init call above.
1795 while (t + 1 < grids->nc[d] && i >= (n[d]*(t+1))/grids->nc[d])
1799 grids->g2t[d][i] = t*tfac;
1802 tfac *= grids->nc[d];
1806 case XX: max_comm_lines = overlap_x; break;
1807 case YY: max_comm_lines = overlap_y; break;
1808 case ZZ: max_comm_lines = pme_order - 1; break;
1810 grids->nthread_comm[d] = 0;
1811 while ((n[d]*grids->nthread_comm[d])/grids->nc[d] < max_comm_lines &&
1812 grids->nthread_comm[d] < grids->nc[d])
1814 grids->nthread_comm[d]++;
1818 fprintf(debug, "pmegrid thread grid communication range in %c: %d\n",
1819 'x'+d, grids->nthread_comm[d]);
1821 /* It should be possible to make grids->nthread_comm[d]==grids->nc[d]
1822 * work, but this is not a problematic restriction.
1824 if (grids->nc[d] > 1 && grids->nthread_comm[d] > grids->nc[d])
1826 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);
1832 static void pmegrids_destroy(pmegrids_t *grids)
1836 if (grids->grid.grid != NULL)
1838 sfree(grids->grid.grid);
1840 if (grids->nthread > 0)
1842 for (t = 0; t < grids->nthread; t++)
1844 sfree(grids->grid_th[t].grid);
1846 sfree(grids->grid_th);
1852 static void realloc_work(pme_work_t *work, int nkx)
1854 if (nkx > work->nalloc)
1857 srenew(work->mhx, work->nalloc);
1858 srenew(work->mhy, work->nalloc);
1859 srenew(work->mhz, work->nalloc);
1860 srenew(work->m2, work->nalloc);
1861 /* Allocate an aligned pointer for SSE operations, including 3 extra
1862 * elements at the end since SSE operates on 4 elements at a time.
1864 sfree_aligned(work->denom);
1865 sfree_aligned(work->tmp1);
1866 sfree_aligned(work->eterm);
1867 snew_aligned(work->denom, work->nalloc+3, 16);
1868 snew_aligned(work->tmp1, work->nalloc+3, 16);
1869 snew_aligned(work->eterm, work->nalloc+3, 16);
1870 srenew(work->m2inv, work->nalloc);
1875 static void free_work(pme_work_t *work)
1881 sfree_aligned(work->denom);
1882 sfree_aligned(work->tmp1);
1883 sfree_aligned(work->eterm);
1889 /* Calculate exponentials through SSE in float precision */
1890 inline static void calc_exponentials(int start, int end, real f, real *d_aligned, real *r_aligned, real *e_aligned)
1893 const __m128 two = _mm_set_ps(2.0f, 2.0f, 2.0f, 2.0f);
1896 __m128 tmp_d1, d_inv, tmp_r, tmp_e;
1898 f_sse = _mm_load1_ps(&f);
1899 for (kx = 0; kx < end; kx += 4)
1901 tmp_d1 = _mm_load_ps(d_aligned+kx);
1902 lu = _mm_rcp_ps(tmp_d1);
1903 d_inv = _mm_mul_ps(lu, _mm_sub_ps(two, _mm_mul_ps(lu, tmp_d1)));
1904 tmp_r = _mm_load_ps(r_aligned+kx);
1905 tmp_r = gmx_mm_exp_ps(tmp_r);
1906 tmp_e = _mm_mul_ps(f_sse, d_inv);
1907 tmp_e = _mm_mul_ps(tmp_e, tmp_r);
1908 _mm_store_ps(e_aligned+kx, tmp_e);
1913 inline static void calc_exponentials(int start, int end, real f, real *d, real *r, real *e)
1916 for (kx = start; kx < end; kx++)
1920 for (kx = start; kx < end; kx++)
1924 for (kx = start; kx < end; kx++)
1926 e[kx] = f*r[kx]*d[kx];
1932 static int solve_pme_yzx(gmx_pme_t pme, t_complex *grid,
1933 real ewaldcoeff, real vol,
1935 int nthread, int thread)
1937 /* do recip sum over local cells in grid */
1938 /* y major, z middle, x minor or continuous */
1940 int kx, ky, kz, maxkx, maxky, maxkz;
1941 int nx, ny, nz, iyz0, iyz1, iyz, iy, iz, kxstart, kxend;
1943 real factor = M_PI*M_PI/(ewaldcoeff*ewaldcoeff);
1944 real ets2, struct2, vfactor, ets2vf;
1945 real d1, d2, energy = 0;
1947 real virxx = 0, virxy = 0, virxz = 0, viryy = 0, viryz = 0, virzz = 0;
1948 real rxx, ryx, ryy, rzx, rzy, rzz;
1950 real *mhx, *mhy, *mhz, *m2, *denom, *tmp1, *eterm, *m2inv;
1951 real mhxk, mhyk, mhzk, m2k;
1954 ivec local_ndata, local_offset, local_size;
1957 elfac = ONE_4PI_EPS0/pme->epsilon_r;
1963 /* Dimensions should be identical for A/B grid, so we just use A here */
1964 gmx_parallel_3dfft_complex_limits(pme->pfft_setupA,
1970 rxx = pme->recipbox[XX][XX];
1971 ryx = pme->recipbox[YY][XX];
1972 ryy = pme->recipbox[YY][YY];
1973 rzx = pme->recipbox[ZZ][XX];
1974 rzy = pme->recipbox[ZZ][YY];
1975 rzz = pme->recipbox[ZZ][ZZ];
1981 work = &pme->work[thread];
1986 denom = work->denom;
1988 eterm = work->eterm;
1989 m2inv = work->m2inv;
1991 iyz0 = local_ndata[YY]*local_ndata[ZZ]* thread /nthread;
1992 iyz1 = local_ndata[YY]*local_ndata[ZZ]*(thread+1)/nthread;
1994 for (iyz = iyz0; iyz < iyz1; iyz++)
1996 iy = iyz/local_ndata[ZZ];
1997 iz = iyz - iy*local_ndata[ZZ];
1999 ky = iy + local_offset[YY];
2010 by = M_PI*vol*pme->bsp_mod[YY][ky];
2012 kz = iz + local_offset[ZZ];
2016 bz = pme->bsp_mod[ZZ][kz];
2018 /* 0.5 correction for corner points */
2020 if (kz == 0 || kz == (nz+1)/2)
2025 p0 = grid + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2027 /* We should skip the k-space point (0,0,0) */
2028 if (local_offset[XX] > 0 || ky > 0 || kz > 0)
2030 kxstart = local_offset[XX];
2034 kxstart = local_offset[XX] + 1;
2037 kxend = local_offset[XX] + local_ndata[XX];
2041 /* More expensive inner loop, especially because of the storage
2042 * of the mh elements in array's.
2043 * Because x is the minor grid index, all mh elements
2044 * depend on kx for triclinic unit cells.
2047 /* Two explicit loops to avoid a conditional inside the loop */
2048 for (kx = kxstart; kx < maxkx; kx++)
2053 mhyk = mx * ryx + my * ryy;
2054 mhzk = mx * rzx + my * rzy + mz * rzz;
2055 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2060 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2061 tmp1[kx] = -factor*m2k;
2064 for (kx = maxkx; kx < kxend; kx++)
2069 mhyk = mx * ryx + my * ryy;
2070 mhzk = mx * rzx + my * rzy + mz * rzz;
2071 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2076 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2077 tmp1[kx] = -factor*m2k;
2080 for (kx = kxstart; kx < kxend; kx++)
2082 m2inv[kx] = 1.0/m2[kx];
2085 calc_exponentials(kxstart, kxend, elfac, denom, tmp1, eterm);
2087 for (kx = kxstart; kx < kxend; kx++, p0++)
2092 p0->re = d1*eterm[kx];
2093 p0->im = d2*eterm[kx];
2095 struct2 = 2.0*(d1*d1+d2*d2);
2097 tmp1[kx] = eterm[kx]*struct2;
2100 for (kx = kxstart; kx < kxend; kx++)
2102 ets2 = corner_fac*tmp1[kx];
2103 vfactor = (factor*m2[kx] + 1.0)*2.0*m2inv[kx];
2106 ets2vf = ets2*vfactor;
2107 virxx += ets2vf*mhx[kx]*mhx[kx] - ets2;
2108 virxy += ets2vf*mhx[kx]*mhy[kx];
2109 virxz += ets2vf*mhx[kx]*mhz[kx];
2110 viryy += ets2vf*mhy[kx]*mhy[kx] - ets2;
2111 viryz += ets2vf*mhy[kx]*mhz[kx];
2112 virzz += ets2vf*mhz[kx]*mhz[kx] - ets2;
2117 /* We don't need to calculate the energy and the virial.
2118 * In this case the triclinic overhead is small.
2121 /* Two explicit loops to avoid a conditional inside the loop */
2123 for (kx = kxstart; kx < maxkx; kx++)
2128 mhyk = mx * ryx + my * ryy;
2129 mhzk = mx * rzx + my * rzy + mz * rzz;
2130 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2131 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2132 tmp1[kx] = -factor*m2k;
2135 for (kx = maxkx; kx < kxend; kx++)
2140 mhyk = mx * ryx + my * ryy;
2141 mhzk = mx * rzx + my * rzy + mz * rzz;
2142 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2143 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2144 tmp1[kx] = -factor*m2k;
2147 calc_exponentials(kxstart, kxend, elfac, denom, tmp1, eterm);
2149 for (kx = kxstart; kx < kxend; kx++, p0++)
2154 p0->re = d1*eterm[kx];
2155 p0->im = d2*eterm[kx];
2162 /* Update virial with local values.
2163 * The virial is symmetric by definition.
2164 * this virial seems ok for isotropic scaling, but I'm
2165 * experiencing problems on semiisotropic membranes.
2166 * IS THAT COMMENT STILL VALID??? (DvdS, 2001/02/07).
2168 work->vir[XX][XX] = 0.25*virxx;
2169 work->vir[YY][YY] = 0.25*viryy;
2170 work->vir[ZZ][ZZ] = 0.25*virzz;
2171 work->vir[XX][YY] = work->vir[YY][XX] = 0.25*virxy;
2172 work->vir[XX][ZZ] = work->vir[ZZ][XX] = 0.25*virxz;
2173 work->vir[YY][ZZ] = work->vir[ZZ][YY] = 0.25*viryz;
2175 /* This energy should be corrected for a charged system */
2176 work->energy = 0.5*energy;
2179 /* Return the loop count */
2180 return local_ndata[YY]*local_ndata[XX];
2183 static void get_pme_ener_vir(const gmx_pme_t pme, int nthread,
2184 real *mesh_energy, matrix vir)
2186 /* This function sums output over threads
2187 * and should therefore only be called after thread synchronization.
2191 *mesh_energy = pme->work[0].energy;
2192 copy_mat(pme->work[0].vir, vir);
2194 for (thread = 1; thread < nthread; thread++)
2196 *mesh_energy += pme->work[thread].energy;
2197 m_add(vir, pme->work[thread].vir, vir);
2201 #define DO_FSPLINE(order) \
2202 for (ithx = 0; (ithx < order); ithx++) \
2204 index_x = (i0+ithx)*pny*pnz; \
2208 for (ithy = 0; (ithy < order); ithy++) \
2210 index_xy = index_x+(j0+ithy)*pnz; \
2215 for (ithz = 0; (ithz < order); ithz++) \
2217 gval = grid[index_xy+(k0+ithz)]; \
2218 fxy1 += thz[ithz]*gval; \
2219 fz1 += dthz[ithz]*gval; \
2228 static void gather_f_bsplines(gmx_pme_t pme, real *grid,
2229 gmx_bool bClearF, pme_atomcomm_t *atc,
2230 splinedata_t *spline,
2233 /* sum forces for local particles */
2234 int nn, n, ithx, ithy, ithz, i0, j0, k0;
2235 int index_x, index_xy;
2236 int nx, ny, nz, pnx, pny, pnz;
2238 real tx, ty, dx, dy, qn;
2239 real fx, fy, fz, gval;
2241 real *thx, *thy, *thz, *dthx, *dthy, *dthz;
2243 real rxx, ryx, ryy, rzx, rzy, rzz;
2246 pme_spline_work_t *work;
2248 work = pme->spline_work;
2250 order = pme->pme_order;
2251 thx = spline->theta[XX];
2252 thy = spline->theta[YY];
2253 thz = spline->theta[ZZ];
2254 dthx = spline->dtheta[XX];
2255 dthy = spline->dtheta[YY];
2256 dthz = spline->dtheta[ZZ];
2260 pnx = pme->pmegrid_nx;
2261 pny = pme->pmegrid_ny;
2262 pnz = pme->pmegrid_nz;
2264 rxx = pme->recipbox[XX][XX];
2265 ryx = pme->recipbox[YY][XX];
2266 ryy = pme->recipbox[YY][YY];
2267 rzx = pme->recipbox[ZZ][XX];
2268 rzy = pme->recipbox[ZZ][YY];
2269 rzz = pme->recipbox[ZZ][ZZ];
2271 for (nn = 0; nn < spline->n; nn++)
2273 n = spline->ind[nn];
2274 qn = scale*atc->q[n];
2287 idxptr = atc->idx[n];
2294 /* Pointer arithmetic alert, next six statements */
2295 thx = spline->theta[XX] + norder;
2296 thy = spline->theta[YY] + norder;
2297 thz = spline->theta[ZZ] + norder;
2298 dthx = spline->dtheta[XX] + norder;
2299 dthy = spline->dtheta[YY] + norder;
2300 dthz = spline->dtheta[ZZ] + norder;
2306 #ifdef PME_SSE_UNALIGNED
2307 #define PME_GATHER_F_SSE_ORDER4
2309 #define PME_GATHER_F_SSE_ALIGNED
2312 #include "pme_sse_single.h"
2319 #define PME_GATHER_F_SSE_ALIGNED
2321 #include "pme_sse_single.h"
2331 atc->f[n][XX] += -qn*( fx*nx*rxx );
2332 atc->f[n][YY] += -qn*( fx*nx*ryx + fy*ny*ryy );
2333 atc->f[n][ZZ] += -qn*( fx*nx*rzx + fy*ny*rzy + fz*nz*rzz );
2336 /* Since the energy and not forces are interpolated
2337 * the net force might not be exactly zero.
2338 * This can be solved by also interpolating F, but
2339 * that comes at a cost.
2340 * A better hack is to remove the net force every
2341 * step, but that must be done at a higher level
2342 * since this routine doesn't see all atoms if running
2343 * in parallel. Don't know how important it is? EL 990726
2348 static real gather_energy_bsplines(gmx_pme_t pme, real *grid,
2349 pme_atomcomm_t *atc)
2351 splinedata_t *spline;
2352 int n, ithx, ithy, ithz, i0, j0, k0;
2353 int index_x, index_xy;
2355 real energy, pot, tx, ty, qn, gval;
2356 real *thx, *thy, *thz;
2360 spline = &atc->spline[0];
2362 order = pme->pme_order;
2365 for (n = 0; (n < atc->n); n++)
2371 idxptr = atc->idx[n];
2378 /* Pointer arithmetic alert, next three statements */
2379 thx = spline->theta[XX] + norder;
2380 thy = spline->theta[YY] + norder;
2381 thz = spline->theta[ZZ] + norder;
2384 for (ithx = 0; (ithx < order); ithx++)
2386 index_x = (i0+ithx)*pme->pmegrid_ny*pme->pmegrid_nz;
2389 for (ithy = 0; (ithy < order); ithy++)
2391 index_xy = index_x+(j0+ithy)*pme->pmegrid_nz;
2394 for (ithz = 0; (ithz < order); ithz++)
2396 gval = grid[index_xy+(k0+ithz)];
2397 pot += tx*ty*thz[ithz]*gval;
2410 /* Macro to force loop unrolling by fixing order.
2411 * This gives a significant performance gain.
2413 #define CALC_SPLINE(order) \
2417 real data[PME_ORDER_MAX]; \
2418 real ddata[PME_ORDER_MAX]; \
2420 for (j = 0; (j < DIM); j++) \
2424 /* dr is relative offset from lower cell limit */ \
2425 data[order-1] = 0; \
2429 for (k = 3; (k < order); k++) \
2431 div = 1.0/(k - 1.0); \
2432 data[k-1] = div*dr*data[k-2]; \
2433 for (l = 1; (l < (k-1)); l++) \
2435 data[k-l-1] = div*((dr+l)*data[k-l-2]+(k-l-dr)* \
2438 data[0] = div*(1-dr)*data[0]; \
2440 /* differentiate */ \
2441 ddata[0] = -data[0]; \
2442 for (k = 1; (k < order); k++) \
2444 ddata[k] = data[k-1] - data[k]; \
2447 div = 1.0/(order - 1); \
2448 data[order-1] = div*dr*data[order-2]; \
2449 for (l = 1; (l < (order-1)); l++) \
2451 data[order-l-1] = div*((dr+l)*data[order-l-2]+ \
2452 (order-l-dr)*data[order-l-1]); \
2454 data[0] = div*(1 - dr)*data[0]; \
2456 for (k = 0; k < order; k++) \
2458 theta[j][i*order+k] = data[k]; \
2459 dtheta[j][i*order+k] = ddata[k]; \
2464 void make_bsplines(splinevec theta, splinevec dtheta, int order,
2465 rvec fractx[], int nr, int ind[], real charge[],
2466 gmx_bool bFreeEnergy)
2468 /* construct splines for local atoms */
2472 for (i = 0; i < nr; i++)
2474 /* With free energy we do not use the charge check.
2475 * In most cases this will be more efficient than calling make_bsplines
2476 * twice, since usually more than half the particles have charges.
2479 if (bFreeEnergy || charge[ii] != 0.0)
2484 case 4: CALC_SPLINE(4); break;
2485 case 5: CALC_SPLINE(5); break;
2486 default: CALC_SPLINE(order); break;
2493 void make_dft_mod(real *mod, real *data, int ndata)
2498 for (i = 0; i < ndata; i++)
2501 for (j = 0; j < ndata; j++)
2503 arg = (2.0*M_PI*i*j)/ndata;
2504 sc += data[j]*cos(arg);
2505 ss += data[j]*sin(arg);
2507 mod[i] = sc*sc+ss*ss;
2509 for (i = 0; i < ndata; i++)
2513 mod[i] = (mod[i-1]+mod[i+1])*0.5;
2519 static void make_bspline_moduli(splinevec bsp_mod,
2520 int nx, int ny, int nz, int order)
2522 int nmax = max(nx, max(ny, nz));
2523 real *data, *ddata, *bsp_data;
2529 snew(bsp_data, nmax);
2535 for (k = 3; k < order; k++)
2539 for (l = 1; l < (k-1); l++)
2541 data[k-l-1] = div*(l*data[k-l-2]+(k-l)*data[k-l-1]);
2543 data[0] = div*data[0];
2546 ddata[0] = -data[0];
2547 for (k = 1; k < order; k++)
2549 ddata[k] = data[k-1]-data[k];
2551 div = 1.0/(order-1);
2553 for (l = 1; l < (order-1); l++)
2555 data[order-l-1] = div*(l*data[order-l-2]+(order-l)*data[order-l-1]);
2557 data[0] = div*data[0];
2559 for (i = 0; i < nmax; i++)
2563 for (i = 1; i <= order; i++)
2565 bsp_data[i] = data[i-1];
2568 make_dft_mod(bsp_mod[XX], bsp_data, nx);
2569 make_dft_mod(bsp_mod[YY], bsp_data, ny);
2570 make_dft_mod(bsp_mod[ZZ], bsp_data, nz);
2578 /* Return the P3M optimal influence function */
2579 static double do_p3m_influence(double z, int order)
2586 /* The formula and most constants can be found in:
2587 * Ballenegger et al., JCTC 8, 936 (2012)
2592 return 1.0 - 2.0*z2/3.0;
2595 return 1.0 - z2 + 2.0*z4/15.0;
2598 return 1.0 - 4.0*z2/3.0 + 2.0*z4/5.0 + 4.0*z2*z4/315.0;
2601 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;
2604 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;
2607 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;
2609 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;
2616 /* Calculate the P3M B-spline moduli for one dimension */
2617 static void make_p3m_bspline_moduli_dim(real *bsp_mod, int n, int order)
2619 double zarg, zai, sinzai, infl;
2624 gmx_fatal(FARGS, "The current P3M code only supports orders up to 8");
2631 for (i = -maxk; i < 0; i++)
2635 infl = do_p3m_influence(sinzai, order);
2636 bsp_mod[n+i] = infl*infl*pow(sinzai/zai, -2.0*order);
2639 for (i = 1; i < maxk; i++)
2643 infl = do_p3m_influence(sinzai, order);
2644 bsp_mod[i] = infl*infl*pow(sinzai/zai, -2.0*order);
2648 /* Calculate the P3M B-spline moduli */
2649 static void make_p3m_bspline_moduli(splinevec bsp_mod,
2650 int nx, int ny, int nz, int order)
2652 make_p3m_bspline_moduli_dim(bsp_mod[XX], nx, order);
2653 make_p3m_bspline_moduli_dim(bsp_mod[YY], ny, order);
2654 make_p3m_bspline_moduli_dim(bsp_mod[ZZ], nz, order);
2658 static void setup_coordinate_communication(pme_atomcomm_t *atc)
2666 for (i = 1; i <= nslab/2; i++)
2668 fw = (atc->nodeid + i) % nslab;
2669 bw = (atc->nodeid - i + nslab) % nslab;
2672 atc->node_dest[n] = fw;
2673 atc->node_src[n] = bw;
2678 atc->node_dest[n] = bw;
2679 atc->node_src[n] = fw;
2685 int gmx_pme_destroy(FILE *log, gmx_pme_t *pmedata)
2691 fprintf(log, "Destroying PME data structures.\n");
2694 sfree((*pmedata)->nnx);
2695 sfree((*pmedata)->nny);
2696 sfree((*pmedata)->nnz);
2698 pmegrids_destroy(&(*pmedata)->pmegridA);
2700 sfree((*pmedata)->fftgridA);
2701 sfree((*pmedata)->cfftgridA);
2702 gmx_parallel_3dfft_destroy((*pmedata)->pfft_setupA);
2704 if ((*pmedata)->pmegridB.grid.grid != NULL)
2706 pmegrids_destroy(&(*pmedata)->pmegridB);
2707 sfree((*pmedata)->fftgridB);
2708 sfree((*pmedata)->cfftgridB);
2709 gmx_parallel_3dfft_destroy((*pmedata)->pfft_setupB);
2711 for (thread = 0; thread < (*pmedata)->nthread; thread++)
2713 free_work(&(*pmedata)->work[thread]);
2715 sfree((*pmedata)->work);
2723 static int mult_up(int n, int f)
2725 return ((n + f - 1)/f)*f;
2729 static double pme_load_imbalance(gmx_pme_t pme)
2734 nma = pme->nnodes_major;
2735 nmi = pme->nnodes_minor;
2737 n1 = mult_up(pme->nkx, nma)*mult_up(pme->nky, nmi)*pme->nkz;
2738 n2 = mult_up(pme->nkx, nma)*mult_up(pme->nkz, nmi)*pme->nky;
2739 n3 = mult_up(pme->nky, nma)*mult_up(pme->nkz, nmi)*pme->nkx;
2741 /* pme_solve is roughly double the cost of an fft */
2743 return (n1 + n2 + 3*n3)/(double)(6*pme->nkx*pme->nky*pme->nkz);
2746 static void init_atomcomm(gmx_pme_t pme, pme_atomcomm_t *atc, t_commrec *cr,
2747 int dimind, gmx_bool bSpread)
2749 int nk, k, s, thread;
2751 atc->dimind = dimind;
2756 if (pme->nnodes > 1)
2758 atc->mpi_comm = pme->mpi_comm_d[dimind];
2759 MPI_Comm_size(atc->mpi_comm, &atc->nslab);
2760 MPI_Comm_rank(atc->mpi_comm, &atc->nodeid);
2764 fprintf(debug, "For PME atom communication in dimind %d: nslab %d rank %d\n", atc->dimind, atc->nslab, atc->nodeid);
2768 atc->bSpread = bSpread;
2769 atc->pme_order = pme->pme_order;
2773 /* These three allocations are not required for particle decomp. */
2774 snew(atc->node_dest, atc->nslab);
2775 snew(atc->node_src, atc->nslab);
2776 setup_coordinate_communication(atc);
2778 snew(atc->count_thread, pme->nthread);
2779 for (thread = 0; thread < pme->nthread; thread++)
2781 snew(atc->count_thread[thread], atc->nslab);
2783 atc->count = atc->count_thread[0];
2784 snew(atc->rcount, atc->nslab);
2785 snew(atc->buf_index, atc->nslab);
2788 atc->nthread = pme->nthread;
2789 if (atc->nthread > 1)
2791 snew(atc->thread_plist, atc->nthread);
2793 snew(atc->spline, atc->nthread);
2794 for (thread = 0; thread < atc->nthread; thread++)
2796 if (atc->nthread > 1)
2798 snew(atc->thread_plist[thread].n, atc->nthread+2*GMX_CACHE_SEP);
2799 atc->thread_plist[thread].n += GMX_CACHE_SEP;
2801 snew(atc->spline[thread].thread_one, pme->nthread);
2802 atc->spline[thread].thread_one[thread] = 1;
2807 init_overlap_comm(pme_overlap_t * ol,
2817 int lbnd, rbnd, maxlr, b, i;
2820 pme_grid_comm_t *pgc;
2822 int fft_start, fft_end, send_index1, recv_index1;
2826 ol->mpi_comm = comm;
2829 ol->nnodes = nnodes;
2830 ol->nodeid = nodeid;
2832 /* Linear translation of the PME grid won't affect reciprocal space
2833 * calculations, so to optimize we only interpolate "upwards",
2834 * which also means we only have to consider overlap in one direction.
2835 * I.e., particles on this node might also be spread to grid indices
2836 * that belong to higher nodes (modulo nnodes)
2839 snew(ol->s2g0, ol->nnodes+1);
2840 snew(ol->s2g1, ol->nnodes);
2843 fprintf(debug, "PME slab boundaries:");
2845 for (i = 0; i < nnodes; i++)
2847 /* s2g0 the local interpolation grid start.
2848 * s2g1 the local interpolation grid end.
2849 * Because grid overlap communication only goes forward,
2850 * the grid the slabs for fft's should be rounded down.
2852 ol->s2g0[i] = ( i *ndata + 0 )/nnodes;
2853 ol->s2g1[i] = ((i+1)*ndata + nnodes-1)/nnodes + norder - 1;
2857 fprintf(debug, " %3d %3d", ol->s2g0[i], ol->s2g1[i]);
2860 ol->s2g0[nnodes] = ndata;
2863 fprintf(debug, "\n");
2866 /* Determine with how many nodes we need to communicate the grid overlap */
2872 for (i = 0; i < nnodes; i++)
2874 if ((i+b < nnodes && ol->s2g1[i] > ol->s2g0[i+b]) ||
2875 (i+b >= nnodes && ol->s2g1[i] > ol->s2g0[i+b-nnodes] + ndata))
2881 while (bCont && b < nnodes);
2882 ol->noverlap_nodes = b - 1;
2884 snew(ol->send_id, ol->noverlap_nodes);
2885 snew(ol->recv_id, ol->noverlap_nodes);
2886 for (b = 0; b < ol->noverlap_nodes; b++)
2888 ol->send_id[b] = (ol->nodeid + (b + 1)) % ol->nnodes;
2889 ol->recv_id[b] = (ol->nodeid - (b + 1) + ol->nnodes) % ol->nnodes;
2891 snew(ol->comm_data, ol->noverlap_nodes);
2894 for (b = 0; b < ol->noverlap_nodes; b++)
2896 pgc = &ol->comm_data[b];
2898 fft_start = ol->s2g0[ol->send_id[b]];
2899 fft_end = ol->s2g0[ol->send_id[b]+1];
2900 if (ol->send_id[b] < nodeid)
2905 send_index1 = ol->s2g1[nodeid];
2906 send_index1 = min(send_index1, fft_end);
2907 pgc->send_index0 = fft_start;
2908 pgc->send_nindex = max(0, send_index1 - pgc->send_index0);
2909 ol->send_size += pgc->send_nindex;
2911 /* We always start receiving to the first index of our slab */
2912 fft_start = ol->s2g0[ol->nodeid];
2913 fft_end = ol->s2g0[ol->nodeid+1];
2914 recv_index1 = ol->s2g1[ol->recv_id[b]];
2915 if (ol->recv_id[b] > nodeid)
2917 recv_index1 -= ndata;
2919 recv_index1 = min(recv_index1, fft_end);
2920 pgc->recv_index0 = fft_start;
2921 pgc->recv_nindex = max(0, recv_index1 - pgc->recv_index0);
2925 /* Communicate the buffer sizes to receive */
2926 for (b = 0; b < ol->noverlap_nodes; b++)
2928 MPI_Sendrecv(&ol->send_size, 1, MPI_INT, ol->send_id[b], b,
2929 &ol->comm_data[b].recv_size, 1, MPI_INT, ol->recv_id[b], b,
2930 ol->mpi_comm, &stat);
2934 /* For non-divisible grid we need pme_order iso pme_order-1 */
2935 snew(ol->sendbuf, norder*commplainsize);
2936 snew(ol->recvbuf, norder*commplainsize);
2940 make_gridindex5_to_localindex(int n, int local_start, int local_range,
2941 int **global_to_local,
2942 real **fraction_shift)
2950 for (i = 0; (i < 5*n); i++)
2952 /* Determine the global to local grid index */
2953 gtl[i] = (i - local_start + n) % n;
2954 /* For coordinates that fall within the local grid the fraction
2955 * is correct, we don't need to shift it.
2958 if (local_range < n)
2960 /* Due to rounding issues i could be 1 beyond the lower or
2961 * upper boundary of the local grid. Correct the index for this.
2962 * If we shift the index, we need to shift the fraction by
2963 * the same amount in the other direction to not affect
2965 * Note that due to this shifting the weights at the end of
2966 * the spline might change, but that will only involve values
2967 * between zero and values close to the precision of a real,
2968 * which is anyhow the accuracy of the whole mesh calculation.
2970 /* With local_range=0 we should not change i=local_start */
2971 if (i % n != local_start)
2978 else if (gtl[i] == local_range)
2980 gtl[i] = local_range - 1;
2987 *global_to_local = gtl;
2988 *fraction_shift = fsh;
2991 static pme_spline_work_t *make_pme_spline_work(int order)
2993 pme_spline_work_t *work;
3000 snew_aligned(work, 1, 16);
3002 zero_SSE = _mm_setzero_ps();
3004 /* Generate bit masks to mask out the unused grid entries,
3005 * as we only operate on order of the 8 grid entries that are
3006 * load into 2 SSE float registers.
3008 for (of = 0; of < 8-(order-1); of++)
3010 for (i = 0; i < 8; i++)
3012 tmp[i] = (i >= of && i < of+order ? 1 : 0);
3014 work->mask_SSE0[of] = _mm_loadu_ps(tmp);
3015 work->mask_SSE1[of] = _mm_loadu_ps(tmp+4);
3016 work->mask_SSE0[of] = _mm_cmpgt_ps(work->mask_SSE0[of], zero_SSE);
3017 work->mask_SSE1[of] = _mm_cmpgt_ps(work->mask_SSE1[of], zero_SSE);
3026 void gmx_pme_check_restrictions(int pme_order,
3027 int nkx, int nky, int nkz,
3030 gmx_bool bUseThreads,
3032 gmx_bool *bValidSettings)
3034 if (pme_order > PME_ORDER_MAX)
3038 *bValidSettings = FALSE;
3041 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.",
3042 pme_order, PME_ORDER_MAX);
3045 if (nkx <= pme_order*(nnodes_major > 1 ? 2 : 1) ||
3046 nky <= pme_order*(nnodes_minor > 1 ? 2 : 1) ||
3051 *bValidSettings = FALSE;
3054 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",
3058 /* Check for a limitation of the (current) sum_fftgrid_dd code.
3059 * We only allow multiple communication pulses in dim 1, not in dim 0.
3061 if (bUseThreads && (nkx < nnodes_major*pme_order &&
3062 nkx != nnodes_major*(pme_order - 1)))
3066 *bValidSettings = FALSE;
3069 gmx_fatal(FARGS, "The number of PME grid lines per node along x is %g. But when using OpenMP threads, the number of grid lines per node along x should be >= pme_order (%d) or = pmeorder-1. To resolve this issue, use less nodes along x (and possibly more along y and/or z) by specifying -dd manually.",
3070 nkx/(double)nnodes_major, pme_order);
3073 if (bValidSettings != NULL)
3075 *bValidSettings = TRUE;
3081 int gmx_pme_init(gmx_pme_t * pmedata,
3087 gmx_bool bFreeEnergy,
3088 gmx_bool bReproducible,
3091 gmx_pme_t pme = NULL;
3093 int use_threads, sum_use_threads;
3098 fprintf(debug, "Creating PME data structures.\n");
3102 pme->redist_init = FALSE;
3103 pme->sum_qgrid_tmp = NULL;
3104 pme->sum_qgrid_dd_tmp = NULL;
3105 pme->buf_nalloc = 0;
3106 pme->redist_buf_nalloc = 0;
3109 pme->bPPnode = TRUE;
3111 pme->nnodes_major = nnodes_major;
3112 pme->nnodes_minor = nnodes_minor;
3115 if (nnodes_major*nnodes_minor > 1)
3117 pme->mpi_comm = cr->mpi_comm_mygroup;
3119 MPI_Comm_rank(pme->mpi_comm, &pme->nodeid);
3120 MPI_Comm_size(pme->mpi_comm, &pme->nnodes);
3121 if (pme->nnodes != nnodes_major*nnodes_minor)
3123 gmx_incons("PME node count mismatch");
3128 pme->mpi_comm = MPI_COMM_NULL;
3132 if (pme->nnodes == 1)
3135 pme->mpi_comm_d[0] = MPI_COMM_NULL;
3136 pme->mpi_comm_d[1] = MPI_COMM_NULL;
3138 pme->ndecompdim = 0;
3139 pme->nodeid_major = 0;
3140 pme->nodeid_minor = 0;
3142 pme->mpi_comm_d[0] = pme->mpi_comm_d[1] = MPI_COMM_NULL;
3147 if (nnodes_minor == 1)
3150 pme->mpi_comm_d[0] = pme->mpi_comm;
3151 pme->mpi_comm_d[1] = MPI_COMM_NULL;
3153 pme->ndecompdim = 1;
3154 pme->nodeid_major = pme->nodeid;
3155 pme->nodeid_minor = 0;
3158 else if (nnodes_major == 1)
3161 pme->mpi_comm_d[0] = MPI_COMM_NULL;
3162 pme->mpi_comm_d[1] = pme->mpi_comm;
3164 pme->ndecompdim = 1;
3165 pme->nodeid_major = 0;
3166 pme->nodeid_minor = pme->nodeid;
3170 if (pme->nnodes % nnodes_major != 0)
3172 gmx_incons("For 2D PME decomposition, #PME nodes must be divisible by the number of nodes in the major dimension");
3174 pme->ndecompdim = 2;
3177 MPI_Comm_split(pme->mpi_comm, pme->nodeid % nnodes_minor,
3178 pme->nodeid, &pme->mpi_comm_d[0]); /* My communicator along major dimension */
3179 MPI_Comm_split(pme->mpi_comm, pme->nodeid/nnodes_minor,
3180 pme->nodeid, &pme->mpi_comm_d[1]); /* My communicator along minor dimension */
3182 MPI_Comm_rank(pme->mpi_comm_d[0], &pme->nodeid_major);
3183 MPI_Comm_size(pme->mpi_comm_d[0], &pme->nnodes_major);
3184 MPI_Comm_rank(pme->mpi_comm_d[1], &pme->nodeid_minor);
3185 MPI_Comm_size(pme->mpi_comm_d[1], &pme->nnodes_minor);
3188 pme->bPPnode = (cr->duty & DUTY_PP);
3191 pme->nthread = nthread;
3193 /* Check if any of the PME MPI ranks uses threads */
3194 use_threads = (pme->nthread > 1 ? 1 : 0);
3196 if (pme->nnodes > 1)
3198 MPI_Allreduce(&use_threads, &sum_use_threads, 1, MPI_INT,
3199 MPI_SUM, pme->mpi_comm);
3204 sum_use_threads = use_threads;
3206 pme->bUseThreads = (sum_use_threads > 0);
3208 if (ir->ePBC == epbcSCREW)
3210 gmx_fatal(FARGS, "pme does not (yet) work with pbc = screw");
3213 pme->bFEP = ((ir->efep != efepNO) && bFreeEnergy);
3217 pme->bP3M = (ir->coulombtype == eelP3M_AD || getenv("GMX_PME_P3M") != NULL);
3218 pme->pme_order = ir->pme_order;
3219 pme->epsilon_r = ir->epsilon_r;
3221 /* If we violate restrictions, generate a fatal error here */
3222 gmx_pme_check_restrictions(pme->pme_order,
3223 pme->nkx, pme->nky, pme->nkz,
3230 if (pme->nnodes > 1)
3235 MPI_Type_contiguous(DIM, mpi_type, &(pme->rvec_mpi));
3236 MPI_Type_commit(&(pme->rvec_mpi));
3239 /* Note that the charge spreading and force gathering, which usually
3240 * takes about the same amount of time as FFT+solve_pme,
3241 * is always fully load balanced
3242 * (unless the charge distribution is inhomogeneous).
3245 imbal = pme_load_imbalance(pme);
3246 if (imbal >= 1.2 && pme->nodeid_major == 0 && pme->nodeid_minor == 0)
3250 "NOTE: The load imbalance in PME FFT and solve is %d%%.\n"
3251 " For optimal PME load balancing\n"
3252 " PME grid_x (%d) and grid_y (%d) should be divisible by #PME_nodes_x (%d)\n"
3253 " and PME grid_y (%d) and grid_z (%d) should be divisible by #PME_nodes_y (%d)\n"
3255 (int)((imbal-1)*100 + 0.5),
3256 pme->nkx, pme->nky, pme->nnodes_major,
3257 pme->nky, pme->nkz, pme->nnodes_minor);
3261 /* For non-divisible grid we need pme_order iso pme_order-1 */
3262 /* In sum_qgrid_dd x overlap is copied in place: take padding into account.
3263 * y is always copied through a buffer: we don't need padding in z,
3264 * but we do need the overlap in x because of the communication order.
3266 init_overlap_comm(&pme->overlap[0], pme->pme_order,
3270 pme->nnodes_major, pme->nodeid_major,
3272 (div_round_up(pme->nky, pme->nnodes_minor)+pme->pme_order)*(pme->nkz+pme->pme_order-1));
3274 /* Along overlap dim 1 we can send in multiple pulses in sum_fftgrid_dd.
3275 * We do this with an offset buffer of equal size, so we need to allocate
3276 * extra for the offset. That's what the (+1)*pme->nkz is for.
3278 init_overlap_comm(&pme->overlap[1], pme->pme_order,
3282 pme->nnodes_minor, pme->nodeid_minor,
3284 (div_round_up(pme->nkx, pme->nnodes_major)+pme->pme_order+1)*pme->nkz);
3286 /* Double-check for a limitation of the (current) sum_fftgrid_dd code.
3287 * Note that gmx_pme_check_restrictions checked for this already.
3289 if (pme->bUseThreads && pme->overlap[0].noverlap_nodes > 1)
3291 gmx_incons("More than one communication pulse required for grid overlap communication along the major dimension while using threads");
3294 snew(pme->bsp_mod[XX], pme->nkx);
3295 snew(pme->bsp_mod[YY], pme->nky);
3296 snew(pme->bsp_mod[ZZ], pme->nkz);
3298 /* The required size of the interpolation grid, including overlap.
3299 * The allocated size (pmegrid_n?) might be slightly larger.
3301 pme->pmegrid_nx = pme->overlap[0].s2g1[pme->nodeid_major] -
3302 pme->overlap[0].s2g0[pme->nodeid_major];
3303 pme->pmegrid_ny = pme->overlap[1].s2g1[pme->nodeid_minor] -
3304 pme->overlap[1].s2g0[pme->nodeid_minor];
3305 pme->pmegrid_nz_base = pme->nkz;
3306 pme->pmegrid_nz = pme->pmegrid_nz_base + pme->pme_order - 1;
3307 set_grid_alignment(&pme->pmegrid_nz, pme->pme_order);
3309 pme->pmegrid_start_ix = pme->overlap[0].s2g0[pme->nodeid_major];
3310 pme->pmegrid_start_iy = pme->overlap[1].s2g0[pme->nodeid_minor];
3311 pme->pmegrid_start_iz = 0;
3313 make_gridindex5_to_localindex(pme->nkx,
3314 pme->pmegrid_start_ix,
3315 pme->pmegrid_nx - (pme->pme_order-1),
3316 &pme->nnx, &pme->fshx);
3317 make_gridindex5_to_localindex(pme->nky,
3318 pme->pmegrid_start_iy,
3319 pme->pmegrid_ny - (pme->pme_order-1),
3320 &pme->nny, &pme->fshy);
3321 make_gridindex5_to_localindex(pme->nkz,
3322 pme->pmegrid_start_iz,
3323 pme->pmegrid_nz_base,
3324 &pme->nnz, &pme->fshz);
3326 pmegrids_init(&pme->pmegridA,
3327 pme->pmegrid_nx, pme->pmegrid_ny, pme->pmegrid_nz,
3328 pme->pmegrid_nz_base,
3332 pme->overlap[0].s2g1[pme->nodeid_major]-pme->overlap[0].s2g0[pme->nodeid_major+1],
3333 pme->overlap[1].s2g1[pme->nodeid_minor]-pme->overlap[1].s2g0[pme->nodeid_minor+1]);
3335 pme->spline_work = make_pme_spline_work(pme->pme_order);
3337 ndata[0] = pme->nkx;
3338 ndata[1] = pme->nky;
3339 ndata[2] = pme->nkz;
3341 /* This routine will allocate the grid data to fit the FFTs */
3342 gmx_parallel_3dfft_init(&pme->pfft_setupA, ndata,
3343 &pme->fftgridA, &pme->cfftgridA,
3345 pme->overlap[0].s2g0, pme->overlap[1].s2g0,
3346 bReproducible, pme->nthread);
3350 pmegrids_init(&pme->pmegridB,
3351 pme->pmegrid_nx, pme->pmegrid_ny, pme->pmegrid_nz,
3352 pme->pmegrid_nz_base,
3356 pme->nkx % pme->nnodes_major != 0,
3357 pme->nky % pme->nnodes_minor != 0);
3359 gmx_parallel_3dfft_init(&pme->pfft_setupB, ndata,
3360 &pme->fftgridB, &pme->cfftgridB,
3362 pme->overlap[0].s2g0, pme->overlap[1].s2g0,
3363 bReproducible, pme->nthread);
3367 pme->pmegridB.grid.grid = NULL;
3368 pme->fftgridB = NULL;
3369 pme->cfftgridB = NULL;
3374 /* Use plain SPME B-spline interpolation */
3375 make_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
3379 /* Use the P3M grid-optimized influence function */
3380 make_p3m_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
3383 /* Use atc[0] for spreading */
3384 init_atomcomm(pme, &pme->atc[0], cr, nnodes_major > 1 ? 0 : 1, TRUE);
3385 if (pme->ndecompdim >= 2)
3387 init_atomcomm(pme, &pme->atc[1], cr, 1, FALSE);
3390 if (pme->nnodes == 1)
3392 pme->atc[0].n = homenr;
3393 pme_realloc_atomcomm_things(&pme->atc[0]);
3399 /* Use fft5d, order after FFT is y major, z, x minor */
3401 snew(pme->work, pme->nthread);
3402 for (thread = 0; thread < pme->nthread; thread++)
3404 realloc_work(&pme->work[thread], pme->nkx);
3413 static void reuse_pmegrids(const pmegrids_t *old, pmegrids_t *new)
3417 for (d = 0; d < DIM; d++)
3419 if (new->grid.n[d] > old->grid.n[d])
3425 sfree_aligned(new->grid.grid);
3426 new->grid.grid = old->grid.grid;
3428 if (new->grid_th != NULL && new->nthread == old->nthread)
3430 sfree_aligned(new->grid_all);
3431 for (t = 0; t < new->nthread; t++)
3433 new->grid_th[t].grid = old->grid_th[t].grid;
3438 int gmx_pme_reinit(gmx_pme_t * pmedata,
3441 const t_inputrec * ir,
3449 irc.nkx = grid_size[XX];
3450 irc.nky = grid_size[YY];
3451 irc.nkz = grid_size[ZZ];
3453 if (pme_src->nnodes == 1)
3455 homenr = pme_src->atc[0].n;
3462 ret = gmx_pme_init(pmedata, cr, pme_src->nnodes_major, pme_src->nnodes_minor,
3463 &irc, homenr, pme_src->bFEP, FALSE, pme_src->nthread);
3467 /* We can easily reuse the allocated pme grids in pme_src */
3468 reuse_pmegrids(&pme_src->pmegridA, &(*pmedata)->pmegridA);
3469 /* We would like to reuse the fft grids, but that's harder */
3476 static void copy_local_grid(gmx_pme_t pme,
3477 pmegrids_t *pmegrids, int thread, real *fftgrid)
3479 ivec local_fft_ndata, local_fft_offset, local_fft_size;
3483 int offx, offy, offz, x, y, z, i0, i0t;
3488 gmx_parallel_3dfft_real_limits(pme->pfft_setupA,
3492 fft_my = local_fft_size[YY];
3493 fft_mz = local_fft_size[ZZ];
3495 pmegrid = &pmegrids->grid_th[thread];
3497 nsx = pmegrid->s[XX];
3498 nsy = pmegrid->s[YY];
3499 nsz = pmegrid->s[ZZ];
3501 for (d = 0; d < DIM; d++)
3503 nf[d] = min(pmegrid->n[d] - (pmegrid->order - 1),
3504 local_fft_ndata[d] - pmegrid->offset[d]);
3507 offx = pmegrid->offset[XX];
3508 offy = pmegrid->offset[YY];
3509 offz = pmegrid->offset[ZZ];
3511 /* Directly copy the non-overlapping parts of the local grids.
3512 * This also initializes the full grid.
3514 grid_th = pmegrid->grid;
3515 for (x = 0; x < nf[XX]; x++)
3517 for (y = 0; y < nf[YY]; y++)
3519 i0 = ((offx + x)*fft_my + (offy + y))*fft_mz + offz;
3520 i0t = (x*nsy + y)*nsz;
3521 for (z = 0; z < nf[ZZ]; z++)
3523 fftgrid[i0+z] = grid_th[i0t+z];
3530 reduce_threadgrid_overlap(gmx_pme_t pme,
3531 const pmegrids_t *pmegrids, int thread,
3532 real *fftgrid, real *commbuf_x, real *commbuf_y)
3534 ivec local_fft_ndata, local_fft_offset, local_fft_size;
3535 int fft_nx, fft_ny, fft_nz;
3540 int offx, offy, offz, x, y, z, i0, i0t;
3541 int sx, sy, sz, fx, fy, fz, tx1, ty1, tz1, ox, oy, oz;
3542 gmx_bool bClearBufX, bClearBufY, bClearBufXY, bClearBuf;
3543 gmx_bool bCommX, bCommY;
3546 const pmegrid_t *pmegrid, *pmegrid_g, *pmegrid_f;
3547 const real *grid_th;
3548 real *commbuf = NULL;
3550 gmx_parallel_3dfft_real_limits(pme->pfft_setupA,
3554 fft_nx = local_fft_ndata[XX];
3555 fft_ny = local_fft_ndata[YY];
3556 fft_nz = local_fft_ndata[ZZ];
3558 fft_my = local_fft_size[YY];
3559 fft_mz = local_fft_size[ZZ];
3561 /* This routine is called when all thread have finished spreading.
3562 * Here each thread sums grid contributions calculated by other threads
3563 * to the thread local grid volume.
3564 * To minimize the number of grid copying operations,
3565 * this routines sums immediately from the pmegrid to the fftgrid.
3568 /* Determine which part of the full node grid we should operate on,
3569 * this is our thread local part of the full grid.
3571 pmegrid = &pmegrids->grid_th[thread];
3573 for (d = 0; d < DIM; d++)
3575 ne[d] = min(pmegrid->offset[d]+pmegrid->n[d]-(pmegrid->order-1),
3576 local_fft_ndata[d]);
3579 offx = pmegrid->offset[XX];
3580 offy = pmegrid->offset[YY];
3581 offz = pmegrid->offset[ZZ];
3588 /* Now loop over all the thread data blocks that contribute
3589 * to the grid region we (our thread) are operating on.
3591 /* Note that ffy_nx/y is equal to the number of grid points
3592 * between the first point of our node grid and the one of the next node.
3594 for (sx = 0; sx >= -pmegrids->nthread_comm[XX]; sx--)
3596 fx = pmegrid->ci[XX] + sx;
3601 fx += pmegrids->nc[XX];
3603 bCommX = (pme->nnodes_major > 1);
3605 pmegrid_g = &pmegrids->grid_th[fx*pmegrids->nc[YY]*pmegrids->nc[ZZ]];
3606 ox += pmegrid_g->offset[XX];
3609 tx1 = min(ox + pmegrid_g->n[XX], ne[XX]);
3613 tx1 = min(ox + pmegrid_g->n[XX], pme->pme_order);
3616 for (sy = 0; sy >= -pmegrids->nthread_comm[YY]; sy--)
3618 fy = pmegrid->ci[YY] + sy;
3623 fy += pmegrids->nc[YY];
3625 bCommY = (pme->nnodes_minor > 1);
3627 pmegrid_g = &pmegrids->grid_th[fy*pmegrids->nc[ZZ]];
3628 oy += pmegrid_g->offset[YY];
3631 ty1 = min(oy + pmegrid_g->n[YY], ne[YY]);
3635 ty1 = min(oy + pmegrid_g->n[YY], pme->pme_order);
3638 for (sz = 0; sz >= -pmegrids->nthread_comm[ZZ]; sz--)
3640 fz = pmegrid->ci[ZZ] + sz;
3644 fz += pmegrids->nc[ZZ];
3647 pmegrid_g = &pmegrids->grid_th[fz];
3648 oz += pmegrid_g->offset[ZZ];
3649 tz1 = min(oz + pmegrid_g->n[ZZ], ne[ZZ]);
3651 if (sx == 0 && sy == 0 && sz == 0)
3653 /* We have already added our local contribution
3654 * before calling this routine, so skip it here.
3659 thread_f = (fx*pmegrids->nc[YY] + fy)*pmegrids->nc[ZZ] + fz;
3661 pmegrid_f = &pmegrids->grid_th[thread_f];
3663 grid_th = pmegrid_f->grid;
3665 nsx = pmegrid_f->s[XX];
3666 nsy = pmegrid_f->s[YY];
3667 nsz = pmegrid_f->s[ZZ];
3669 #ifdef DEBUG_PME_REDUCE
3670 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",
3671 pme->nodeid, thread, thread_f,
3672 pme->pmegrid_start_ix,
3673 pme->pmegrid_start_iy,
3674 pme->pmegrid_start_iz,
3676 offx-ox, tx1-ox, offx, tx1,
3677 offy-oy, ty1-oy, offy, ty1,
3678 offz-oz, tz1-oz, offz, tz1);
3681 if (!(bCommX || bCommY))
3683 /* Copy from the thread local grid to the node grid */
3684 for (x = offx; x < tx1; x++)
3686 for (y = offy; y < ty1; y++)
3688 i0 = (x*fft_my + y)*fft_mz;
3689 i0t = ((x - ox)*nsy + (y - oy))*nsz - oz;
3690 for (z = offz; z < tz1; z++)
3692 fftgrid[i0+z] += grid_th[i0t+z];
3699 /* The order of this conditional decides
3700 * where the corner volume gets stored with x+y decomp.
3704 commbuf = commbuf_y;
3705 buf_my = ty1 - offy;
3708 /* We index commbuf modulo the local grid size */
3709 commbuf += buf_my*fft_nx*fft_nz;
3711 bClearBuf = bClearBufXY;
3712 bClearBufXY = FALSE;
3716 bClearBuf = bClearBufY;
3722 commbuf = commbuf_x;
3724 bClearBuf = bClearBufX;
3728 /* Copy to the communication buffer */
3729 for (x = offx; x < tx1; x++)
3731 for (y = offy; y < ty1; y++)
3733 i0 = (x*buf_my + y)*fft_nz;
3734 i0t = ((x - ox)*nsy + (y - oy))*nsz - oz;
3738 /* First access of commbuf, initialize it */
3739 for (z = offz; z < tz1; z++)
3741 commbuf[i0+z] = grid_th[i0t+z];
3746 for (z = offz; z < tz1; z++)
3748 commbuf[i0+z] += grid_th[i0t+z];
3760 static void sum_fftgrid_dd(gmx_pme_t pme, real *fftgrid)
3762 ivec local_fft_ndata, local_fft_offset, local_fft_size;
3763 pme_overlap_t *overlap;
3764 int send_index0, send_nindex;
3769 int send_size_y, recv_size_y;
3770 int ipulse, send_id, recv_id, datasize, gridsize, size_yx;
3771 real *sendptr, *recvptr;
3772 int x, y, z, indg, indb;
3774 /* Note that this routine is only used for forward communication.
3775 * Since the force gathering, unlike the charge spreading,
3776 * can be trivially parallelized over the particles,
3777 * the backwards process is much simpler and can use the "old"
3778 * communication setup.
3781 gmx_parallel_3dfft_real_limits(pme->pfft_setupA,
3786 if (pme->nnodes_minor > 1)
3788 /* Major dimension */
3789 overlap = &pme->overlap[1];
3791 if (pme->nnodes_major > 1)
3793 size_yx = pme->overlap[0].comm_data[0].send_nindex;
3799 datasize = (local_fft_ndata[XX] + size_yx)*local_fft_ndata[ZZ];
3801 send_size_y = overlap->send_size;
3803 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
3805 send_id = overlap->send_id[ipulse];
3806 recv_id = overlap->recv_id[ipulse];
3808 overlap->comm_data[ipulse].send_index0 -
3809 overlap->comm_data[0].send_index0;
3810 send_nindex = overlap->comm_data[ipulse].send_nindex;
3811 /* We don't use recv_index0, as we always receive starting at 0 */
3812 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
3813 recv_size_y = overlap->comm_data[ipulse].recv_size;
3815 sendptr = overlap->sendbuf + send_index0*local_fft_ndata[ZZ];
3816 recvptr = overlap->recvbuf;
3819 MPI_Sendrecv(sendptr, send_size_y*datasize, GMX_MPI_REAL,
3821 recvptr, recv_size_y*datasize, GMX_MPI_REAL,
3823 overlap->mpi_comm, &stat);
3826 for (x = 0; x < local_fft_ndata[XX]; x++)
3828 for (y = 0; y < recv_nindex; y++)
3830 indg = (x*local_fft_size[YY] + y)*local_fft_size[ZZ];
3831 indb = (x*recv_size_y + y)*local_fft_ndata[ZZ];
3832 for (z = 0; z < local_fft_ndata[ZZ]; z++)
3834 fftgrid[indg+z] += recvptr[indb+z];
3839 if (pme->nnodes_major > 1)
3841 /* Copy from the received buffer to the send buffer for dim 0 */
3842 sendptr = pme->overlap[0].sendbuf;
3843 for (x = 0; x < size_yx; x++)
3845 for (y = 0; y < recv_nindex; y++)
3847 indg = (x*local_fft_ndata[YY] + y)*local_fft_ndata[ZZ];
3848 indb = ((local_fft_ndata[XX] + x)*recv_size_y + y)*local_fft_ndata[ZZ];
3849 for (z = 0; z < local_fft_ndata[ZZ]; z++)
3851 sendptr[indg+z] += recvptr[indb+z];
3859 /* We only support a single pulse here.
3860 * This is not a severe limitation, as this code is only used
3861 * with OpenMP and with OpenMP the (PME) domains can be larger.
3863 if (pme->nnodes_major > 1)
3865 /* Major dimension */
3866 overlap = &pme->overlap[0];
3868 datasize = local_fft_ndata[YY]*local_fft_ndata[ZZ];
3869 gridsize = local_fft_size[YY] *local_fft_size[ZZ];
3873 send_id = overlap->send_id[ipulse];
3874 recv_id = overlap->recv_id[ipulse];
3875 send_nindex = overlap->comm_data[ipulse].send_nindex;
3876 /* We don't use recv_index0, as we always receive starting at 0 */
3877 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
3879 sendptr = overlap->sendbuf;
3880 recvptr = overlap->recvbuf;
3884 fprintf(debug, "PME fftgrid comm %2d x %2d x %2d\n",
3885 send_nindex, local_fft_ndata[YY], local_fft_ndata[ZZ]);
3889 MPI_Sendrecv(sendptr, send_nindex*datasize, GMX_MPI_REAL,
3891 recvptr, recv_nindex*datasize, GMX_MPI_REAL,
3893 overlap->mpi_comm, &stat);
3896 for (x = 0; x < recv_nindex; x++)
3898 for (y = 0; y < local_fft_ndata[YY]; y++)
3900 indg = (x*local_fft_size[YY] + y)*local_fft_size[ZZ];
3901 indb = (x*local_fft_ndata[YY] + y)*local_fft_ndata[ZZ];
3902 for (z = 0; z < local_fft_ndata[ZZ]; z++)
3904 fftgrid[indg+z] += recvptr[indb+z];
3912 static void spread_on_grid(gmx_pme_t pme,
3913 pme_atomcomm_t *atc, pmegrids_t *grids,
3914 gmx_bool bCalcSplines, gmx_bool bSpread,
3917 int nthread, thread;
3918 #ifdef PME_TIME_THREADS
3919 gmx_cycles_t c1, c2, c3, ct1a, ct1b, ct1c;
3920 static double cs1 = 0, cs2 = 0, cs3 = 0;
3921 static double cs1a[6] = {0, 0, 0, 0, 0, 0};
3925 nthread = pme->nthread;
3926 assert(nthread > 0);
3928 #ifdef PME_TIME_THREADS
3929 c1 = omp_cyc_start();
3933 #pragma omp parallel for num_threads(nthread) schedule(static)
3934 for (thread = 0; thread < nthread; thread++)
3938 start = atc->n* thread /nthread;
3939 end = atc->n*(thread+1)/nthread;
3941 /* Compute fftgrid index for all atoms,
3942 * with help of some extra variables.
3944 calc_interpolation_idx(pme, atc, start, end, thread);
3947 #ifdef PME_TIME_THREADS
3948 c1 = omp_cyc_end(c1);
3952 #ifdef PME_TIME_THREADS
3953 c2 = omp_cyc_start();
3955 #pragma omp parallel for num_threads(nthread) schedule(static)
3956 for (thread = 0; thread < nthread; thread++)
3958 splinedata_t *spline;
3959 pmegrid_t *grid = NULL;
3961 /* make local bsplines */
3962 if (grids == NULL || !pme->bUseThreads)
3964 spline = &atc->spline[0];
3970 grid = &grids->grid;
3975 spline = &atc->spline[thread];
3977 if (grids->nthread == 1)
3979 /* One thread, we operate on all charges */
3984 /* Get the indices our thread should operate on */
3985 make_thread_local_ind(atc, thread, spline);
3988 grid = &grids->grid_th[thread];
3993 make_bsplines(spline->theta, spline->dtheta, pme->pme_order,
3994 atc->fractx, spline->n, spline->ind, atc->q, pme->bFEP);
3999 /* put local atoms on grid. */
4000 #ifdef PME_TIME_SPREAD
4001 ct1a = omp_cyc_start();
4003 spread_q_bsplines_thread(grid, atc, spline, pme->spline_work);
4005 if (pme->bUseThreads)
4007 copy_local_grid(pme, grids, thread, fftgrid);
4009 #ifdef PME_TIME_SPREAD
4010 ct1a = omp_cyc_end(ct1a);
4011 cs1a[thread] += (double)ct1a;
4015 #ifdef PME_TIME_THREADS
4016 c2 = omp_cyc_end(c2);
4020 if (bSpread && pme->bUseThreads)
4022 #ifdef PME_TIME_THREADS
4023 c3 = omp_cyc_start();
4025 #pragma omp parallel for num_threads(grids->nthread) schedule(static)
4026 for (thread = 0; thread < grids->nthread; thread++)
4028 reduce_threadgrid_overlap(pme, grids, thread,
4030 pme->overlap[0].sendbuf,
4031 pme->overlap[1].sendbuf);
4033 #ifdef PME_TIME_THREADS
4034 c3 = omp_cyc_end(c3);
4038 if (pme->nnodes > 1)
4040 /* Communicate the overlapping part of the fftgrid.
4041 * For this communication call we need to check pme->bUseThreads
4042 * to have all ranks communicate here, regardless of pme->nthread.
4044 sum_fftgrid_dd(pme, fftgrid);
4048 #ifdef PME_TIME_THREADS
4052 printf("idx %.2f spread %.2f red %.2f",
4053 cs1*1e-9, cs2*1e-9, cs3*1e-9);
4054 #ifdef PME_TIME_SPREAD
4055 for (thread = 0; thread < nthread; thread++)
4057 printf(" %.2f", cs1a[thread]*1e-9);
4066 static void dump_grid(FILE *fp,
4067 int sx, int sy, int sz, int nx, int ny, int nz,
4068 int my, int mz, const real *g)
4072 for (x = 0; x < nx; x++)
4074 for (y = 0; y < ny; y++)
4076 for (z = 0; z < nz; z++)
4078 fprintf(fp, "%2d %2d %2d %6.3f\n",
4079 sx+x, sy+y, sz+z, g[(x*my + y)*mz + z]);
4085 static void dump_local_fftgrid(gmx_pme_t pme, const real *fftgrid)
4087 ivec local_fft_ndata, local_fft_offset, local_fft_size;
4089 gmx_parallel_3dfft_real_limits(pme->pfft_setupA,
4095 pme->pmegrid_start_ix,
4096 pme->pmegrid_start_iy,
4097 pme->pmegrid_start_iz,
4098 pme->pmegrid_nx-pme->pme_order+1,
4099 pme->pmegrid_ny-pme->pme_order+1,
4100 pme->pmegrid_nz-pme->pme_order+1,
4107 void gmx_pme_calc_energy(gmx_pme_t pme, int n, rvec *x, real *q, real *V)
4109 pme_atomcomm_t *atc;
4112 if (pme->nnodes > 1)
4114 gmx_incons("gmx_pme_calc_energy called in parallel");
4118 gmx_incons("gmx_pme_calc_energy with free energy");
4121 atc = &pme->atc_energy;
4123 if (atc->spline == NULL)
4125 snew(atc->spline, atc->nthread);
4128 atc->bSpread = TRUE;
4129 atc->pme_order = pme->pme_order;
4131 pme_realloc_atomcomm_things(atc);
4135 /* We only use the A-charges grid */
4136 grid = &pme->pmegridA;
4138 /* Only calculate the spline coefficients, don't actually spread */
4139 spread_on_grid(pme, atc, NULL, TRUE, FALSE, pme->fftgridA);
4141 *V = gather_energy_bsplines(pme, grid->grid.grid, atc);
4145 static void reset_pmeonly_counters(t_commrec *cr, gmx_wallcycle_t wcycle,
4146 t_nrnb *nrnb, t_inputrec *ir,
4147 gmx_large_int_t step)
4149 /* Reset all the counters related to performance over the run */
4150 wallcycle_stop(wcycle, ewcRUN);
4151 wallcycle_reset_all(wcycle);
4153 if (ir->nsteps >= 0)
4155 /* ir->nsteps is not used here, but we update it for consistency */
4156 ir->nsteps -= step - ir->init_step;
4158 ir->init_step = step;
4159 wallcycle_start(wcycle, ewcRUN);
4163 static void gmx_pmeonly_switch(int *npmedata, gmx_pme_t **pmedata,
4165 t_commrec *cr, t_inputrec *ir,
4169 gmx_pme_t pme = NULL;
4172 while (ind < *npmedata)
4174 pme = (*pmedata)[ind];
4175 if (pme->nkx == grid_size[XX] &&
4176 pme->nky == grid_size[YY] &&
4177 pme->nkz == grid_size[ZZ])
4188 srenew(*pmedata, *npmedata);
4190 /* Generate a new PME data structure, copying part of the old pointers */
4191 gmx_pme_reinit(&((*pmedata)[ind]), cr, pme, ir, grid_size);
4193 *pme_ret = (*pmedata)[ind];
4197 int gmx_pmeonly(gmx_pme_t pme,
4198 t_commrec *cr, t_nrnb *nrnb,
4199 gmx_wallcycle_t wcycle,
4200 real ewaldcoeff, gmx_bool bGatherOnly,
4205 gmx_pme_pp_t pme_pp;
4209 rvec *x_pp = NULL, *f_pp = NULL;
4210 real *chargeA = NULL, *chargeB = NULL;
4212 int maxshift_x = 0, maxshift_y = 0;
4213 real energy, dvdlambda;
4218 gmx_large_int_t step, step_rel;
4221 /* This data will only use with PME tuning, i.e. switching PME grids */
4223 snew(pmedata, npmedata);
4226 pme_pp = gmx_pme_pp_init(cr);
4231 do /****** this is a quasi-loop over time steps! */
4233 /* The reason for having a loop here is PME grid tuning/switching */
4236 /* Domain decomposition */
4237 ret = gmx_pme_recv_q_x(pme_pp,
4239 &chargeA, &chargeB, box, &x_pp, &f_pp,
4240 &maxshift_x, &maxshift_y,
4241 &pme->bFEP, &lambda,
4244 grid_switch, &ewaldcoeff);
4246 if (ret == pmerecvqxSWITCHGRID)
4248 /* Switch the PME grid to grid_switch */
4249 gmx_pmeonly_switch(&npmedata, &pmedata, grid_switch, cr, ir, &pme);
4252 if (ret == pmerecvqxRESETCOUNTERS)
4254 /* Reset the cycle and flop counters */
4255 reset_pmeonly_counters(cr, wcycle, nrnb, ir, step);
4258 while (ret == pmerecvqxSWITCHGRID || ret == pmerecvqxRESETCOUNTERS);
4260 if (ret == pmerecvqxFINISH)
4262 /* We should stop: break out of the loop */
4266 step_rel = step - ir->init_step;
4270 wallcycle_start(wcycle, ewcRUN);
4273 wallcycle_start(wcycle, ewcPMEMESH);
4277 gmx_pme_do(pme, 0, natoms, x_pp, f_pp, chargeA, chargeB, box,
4278 cr, maxshift_x, maxshift_y, nrnb, wcycle, vir, ewaldcoeff,
4279 &energy, lambda, &dvdlambda,
4280 GMX_PME_DO_ALL_F | (bEnerVir ? GMX_PME_CALC_ENER_VIR : 0));
4282 cycles = wallcycle_stop(wcycle, ewcPMEMESH);
4284 gmx_pme_send_force_vir_ener(pme_pp,
4285 f_pp, vir, energy, dvdlambda,
4289 } /***** end of quasi-loop, we stop with the break above */
4295 int gmx_pme_do(gmx_pme_t pme,
4296 int start, int homenr,
4298 real *chargeA, real *chargeB,
4299 matrix box, t_commrec *cr,
4300 int maxshift_x, int maxshift_y,
4301 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
4302 matrix vir, real ewaldcoeff,
4303 real *energy, real lambda,
4304 real *dvdlambda, int flags)
4306 int q, d, i, j, ntot, npme;
4309 pme_atomcomm_t *atc = NULL;
4310 pmegrids_t *pmegrid = NULL;
4314 real *charge = NULL, *q_d;
4318 gmx_parallel_3dfft_t pfft_setup;
4320 t_complex * cfftgrid;
4322 const gmx_bool bCalcEnerVir = flags & GMX_PME_CALC_ENER_VIR;
4323 const gmx_bool bCalcF = flags & GMX_PME_CALC_F;
4325 assert(pme->nnodes > 0);
4326 assert(pme->nnodes == 1 || pme->ndecompdim > 0);
4328 if (pme->nnodes > 1)
4332 if (atc->npd > atc->pd_nalloc)
4334 atc->pd_nalloc = over_alloc_dd(atc->npd);
4335 srenew(atc->pd, atc->pd_nalloc);
4337 atc->maxshift = (atc->dimind == 0 ? maxshift_x : maxshift_y);
4341 /* This could be necessary for TPI */
4342 pme->atc[0].n = homenr;
4345 for (q = 0; q < (pme->bFEP ? 2 : 1); q++)
4349 pmegrid = &pme->pmegridA;
4350 fftgrid = pme->fftgridA;
4351 cfftgrid = pme->cfftgridA;
4352 pfft_setup = pme->pfft_setupA;
4353 charge = chargeA+start;
4357 pmegrid = &pme->pmegridB;
4358 fftgrid = pme->fftgridB;
4359 cfftgrid = pme->cfftgridB;
4360 pfft_setup = pme->pfft_setupB;
4361 charge = chargeB+start;
4363 grid = pmegrid->grid.grid;
4364 /* Unpack structure */
4367 fprintf(debug, "PME: nnodes = %d, nodeid = %d\n",
4368 cr->nnodes, cr->nodeid);
4369 fprintf(debug, "Grid = %p\n", (void*)grid);
4372 gmx_fatal(FARGS, "No grid!");
4377 m_inv_ur0(box, pme->recipbox);
4379 if (pme->nnodes == 1)
4382 if (DOMAINDECOMP(cr))
4385 pme_realloc_atomcomm_things(atc);
4393 wallcycle_start(wcycle, ewcPME_REDISTXF);
4394 for (d = pme->ndecompdim-1; d >= 0; d--)
4396 if (d == pme->ndecompdim-1)
4404 n_d = pme->atc[d+1].n;
4410 if (atc->npd > atc->pd_nalloc)
4412 atc->pd_nalloc = over_alloc_dd(atc->npd);
4413 srenew(atc->pd, atc->pd_nalloc);
4415 atc->maxshift = (atc->dimind == 0 ? maxshift_x : maxshift_y);
4416 pme_calc_pidx_wrapper(n_d, pme->recipbox, x_d, atc);
4419 /* Redistribute x (only once) and qA or qB */
4420 if (DOMAINDECOMP(cr))
4422 dd_pmeredist_x_q(pme, n_d, q == 0, x_d, q_d, atc);
4426 pmeredist_pd(pme, TRUE, n_d, q == 0, x_d, q_d, atc);
4431 wallcycle_stop(wcycle, ewcPME_REDISTXF);
4436 fprintf(debug, "Node= %6d, pme local particles=%6d\n",
4437 cr->nodeid, atc->n);
4440 if (flags & GMX_PME_SPREAD_Q)
4442 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
4444 /* Spread the charges on a grid */
4445 spread_on_grid(pme, &pme->atc[0], pmegrid, q == 0, TRUE, fftgrid);
4449 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
4451 inc_nrnb(nrnb, eNR_SPREADQBSP,
4452 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
4454 if (!pme->bUseThreads)
4456 wrap_periodic_pmegrid(pme, grid);
4458 /* sum contributions to local grid from other nodes */
4460 if (pme->nnodes > 1)
4462 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_QGRID_FORWARD);
4467 copy_pmegrid_to_fftgrid(pme, grid, fftgrid);
4470 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
4473 dump_local_fftgrid(pme,fftgrid);
4478 /* Here we start a large thread parallel region */
4479 #pragma omp parallel num_threads(pme->nthread) private(thread)
4481 thread = gmx_omp_get_thread_num();
4482 if (flags & GMX_PME_SOLVE)
4489 wallcycle_start(wcycle, ewcPME_FFT);
4491 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
4492 fftgrid, cfftgrid, thread, wcycle);
4495 wallcycle_stop(wcycle, ewcPME_FFT);
4499 /* solve in k-space for our local cells */
4502 wallcycle_start(wcycle, ewcPME_SOLVE);
4505 solve_pme_yzx(pme, cfftgrid, ewaldcoeff,
4506 box[XX][XX]*box[YY][YY]*box[ZZ][ZZ],
4508 pme->nthread, thread);
4511 wallcycle_stop(wcycle, ewcPME_SOLVE);
4513 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
4523 wallcycle_start(wcycle, ewcPME_FFT);
4525 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
4526 cfftgrid, fftgrid, thread, wcycle);
4529 wallcycle_stop(wcycle, ewcPME_FFT);
4533 if (pme->nodeid == 0)
4535 ntot = pme->nkx*pme->nky*pme->nkz;
4536 npme = ntot*log((real)ntot)/log(2.0);
4537 inc_nrnb(nrnb, eNR_FFT, 2*npme);
4540 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
4543 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, pme->nthread, thread);
4546 /* End of thread parallel section.
4547 * With MPI we have to synchronize here before gmx_sum_qgrid_dd.
4552 /* distribute local grid to all nodes */
4554 if (pme->nnodes > 1)
4556 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_QGRID_BACKWARD);
4561 unwrap_periodic_pmegrid(pme, grid);
4563 /* interpolate forces for our local atoms */
4567 /* If we are running without parallelization,
4568 * atc->f is the actual force array, not a buffer,
4569 * therefore we should not clear it.
4571 bClearF = (q == 0 && PAR(cr));
4572 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
4573 for (thread = 0; thread < pme->nthread; thread++)
4575 gather_f_bsplines(pme, grid, bClearF, atc,
4576 &atc->spline[thread],
4577 pme->bFEP ? (q == 0 ? 1.0-lambda : lambda) : 1.0);
4582 inc_nrnb(nrnb, eNR_GATHERFBSP,
4583 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
4584 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
4589 /* This should only be called on the master thread
4590 * and after the threads have synchronized.
4592 get_pme_ener_vir(pme, pme->nthread, &energy_AB[q], vir_AB[q]);
4596 if (bCalcF && pme->nnodes > 1)
4598 wallcycle_start(wcycle, ewcPME_REDISTXF);
4599 for (d = 0; d < pme->ndecompdim; d++)
4602 if (d == pme->ndecompdim - 1)
4609 n_d = pme->atc[d+1].n;
4610 f_d = pme->atc[d+1].f;
4612 if (DOMAINDECOMP(cr))
4614 dd_pmeredist_f(pme, atc, n_d, f_d,
4615 d == pme->ndecompdim-1 && pme->bPPnode);
4619 pmeredist_pd(pme, FALSE, n_d, TRUE, f_d, NULL, atc);
4623 wallcycle_stop(wcycle, ewcPME_REDISTXF);
4631 *energy = energy_AB[0];
4632 m_add(vir, vir_AB[0], vir);
4636 *energy = (1.0-lambda)*energy_AB[0] + lambda*energy_AB[1];
4637 *dvdlambda += energy_AB[1] - energy_AB[0];
4638 for (i = 0; i < DIM; i++)
4640 for (j = 0; j < DIM; j++)
4642 vir[i][j] += (1.0-lambda)*vir_AB[0][i][j] +
4643 lambda*vir_AB[1][i][j];
4655 fprintf(debug, "PME mesh energy: %g\n", *energy);