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/utility/gmxmpi.h"
72 #include "gmxcomplex.h"
76 #include "gmx_fatal.h"
82 #include "gmx_wallcycle.h"
83 #include "gmx_parallel_3dfft.h"
85 #include "gmx_cyclecounter.h"
89 /* Single precision, with SSE2 or higher available */
90 #if defined(GMX_X86_SSE2) && !defined(GMX_DOUBLE)
92 #include "gmx_x86_simd_single.h"
95 /* Some old AMD processors could have problems with unaligned loads+stores */
97 #define PME_SSE_UNALIGNED
102 /* #define PRT_FORCE */
103 /* conditions for on the fly time-measurement */
104 /* #define TAKETIME (step > 1 && timesteps < 10) */
105 #define TAKETIME FALSE
107 /* #define PME_TIME_THREADS */
110 #define mpi_type MPI_DOUBLE
112 #define mpi_type MPI_FLOAT
115 /* GMX_CACHE_SEP should be a multiple of 16 to preserve alignment */
116 #define GMX_CACHE_SEP 64
118 /* We only define a maximum to be able to use local arrays without allocation.
119 * An order larger than 12 should never be needed, even for test cases.
120 * If needed it can be changed here.
122 #define PME_ORDER_MAX 12
124 /* Internal datastructures */
130 int recv_size; /* Receive buffer width, used with OpenMP */
141 int *send_id, *recv_id;
142 int send_size; /* Send buffer width, used with OpenMP */
143 pme_grid_comm_t *comm_data;
149 int *n; /* Cumulative counts of the number of particles per thread */
150 int nalloc; /* Allocation size of i */
151 int *i; /* Particle indices ordered on thread index (n) */
165 int dimind; /* The index of the dimension, 0=x, 1=y */
172 int *node_dest; /* The nodes to send x and q to with DD */
173 int *node_src; /* The nodes to receive x and q from with DD */
174 int *buf_index; /* Index for commnode into the buffers */
181 int *count; /* The number of atoms to send to each node */
183 int *rcount; /* The number of atoms to receive */
190 gmx_bool bSpread; /* These coordinates are used for spreading */
193 rvec *fractx; /* Fractional coordinate relative to the
194 * lower cell boundary
197 int *thread_idx; /* Which thread should spread which charge */
198 thread_plist_t *thread_plist;
199 splinedata_t *spline;
206 ivec ci; /* The spatial location of this grid */
207 ivec n; /* The used size of *grid, including order-1 */
208 ivec offset; /* The grid offset from the full node grid */
209 int order; /* PME spreading order */
210 ivec s; /* The allocated size of *grid, s >= n */
211 real *grid; /* The grid local thread, size n */
215 pmegrid_t grid; /* The full node grid (non thread-local) */
216 int nthread; /* The number of threads operating on this grid */
217 ivec nc; /* The local spatial decomposition over the threads */
218 pmegrid_t *grid_th; /* Array of grids for each thread */
219 real *grid_all; /* Allocated array for the grids in *grid_th */
220 int **g2t; /* The grid to thread index */
221 ivec nthread_comm; /* The number of threads to communicate with */
227 /* Masks for SSE aligned spreading and gathering */
228 __m128 mask_SSE0[6], mask_SSE1[6];
230 int dummy; /* C89 requires that struct has at least one member */
235 /* work data for solve_pme */
251 typedef struct gmx_pme {
252 int ndecompdim; /* The number of decomposition dimensions */
253 int nodeid; /* Our nodeid in mpi->mpi_comm */
256 int nnodes; /* The number of nodes doing PME */
261 MPI_Comm mpi_comm_d[2]; /* Indexed on dimension, 0=x, 1=y */
263 MPI_Datatype rvec_mpi; /* the pme vector's MPI type */
266 gmx_bool bUseThreads; /* Does any of the PME ranks have nthread>1 ? */
267 int nthread; /* The number of threads doing PME on our rank */
269 gmx_bool bPPnode; /* Node also does particle-particle forces */
270 gmx_bool bFEP; /* Compute Free energy contribution */
271 int nkx, nky, nkz; /* Grid dimensions */
272 gmx_bool bP3M; /* Do P3M: optimize the influence function */
276 pmegrids_t pmegridA; /* Grids on which we do spreading/interpolation, includes overlap */
278 /* The PME charge spreading grid sizes/strides, includes pme_order-1 */
279 int pmegrid_nx, pmegrid_ny, pmegrid_nz;
280 /* pmegrid_nz might be larger than strictly necessary to ensure
281 * memory alignment, pmegrid_nz_base gives the real base size.
284 /* The local PME grid starting indices */
285 int pmegrid_start_ix, pmegrid_start_iy, pmegrid_start_iz;
287 /* Work data for spreading and gathering */
288 pme_spline_work_t *spline_work;
290 real *fftgridA; /* Grids for FFT. With 1D FFT decomposition this can be a pointer */
291 real *fftgridB; /* inside the interpolation grid, but separate for 2D PME decomp. */
292 int fftgrid_nx, fftgrid_ny, fftgrid_nz;
294 t_complex *cfftgridA; /* Grids for complex FFT data */
295 t_complex *cfftgridB;
296 int cfftgrid_nx, cfftgrid_ny, cfftgrid_nz;
298 gmx_parallel_3dfft_t pfft_setupA;
299 gmx_parallel_3dfft_t pfft_setupB;
301 int *nnx, *nny, *nnz;
302 real *fshx, *fshy, *fshz;
304 pme_atomcomm_t atc[2]; /* Indexed on decomposition index */
308 pme_overlap_t overlap[2]; /* Indexed on dimension, 0=x, 1=y */
310 pme_atomcomm_t atc_energy; /* Only for gmx_pme_calc_energy */
312 rvec *bufv; /* Communication buffer */
313 real *bufr; /* Communication buffer */
314 int buf_nalloc; /* The communication buffer size */
316 /* thread local work data for solve_pme */
319 /* Work data for PME_redist */
320 gmx_bool redist_init;
328 int redist_buf_nalloc;
330 /* Work data for sum_qgrid */
331 real * sum_qgrid_tmp;
332 real * sum_qgrid_dd_tmp;
336 static void calc_interpolation_idx(gmx_pme_t pme, pme_atomcomm_t *atc,
337 int start, int end, int thread)
340 int *idxptr, tix, tiy, tiz;
341 real *xptr, *fptr, tx, ty, tz;
342 real rxx, ryx, ryy, rzx, rzy, rzz;
344 int start_ix, start_iy, start_iz;
345 int *g2tx, *g2ty, *g2tz;
347 int *thread_idx = NULL;
348 thread_plist_t *tpl = NULL;
356 start_ix = pme->pmegrid_start_ix;
357 start_iy = pme->pmegrid_start_iy;
358 start_iz = pme->pmegrid_start_iz;
360 rxx = pme->recipbox[XX][XX];
361 ryx = pme->recipbox[YY][XX];
362 ryy = pme->recipbox[YY][YY];
363 rzx = pme->recipbox[ZZ][XX];
364 rzy = pme->recipbox[ZZ][YY];
365 rzz = pme->recipbox[ZZ][ZZ];
367 g2tx = pme->pmegridA.g2t[XX];
368 g2ty = pme->pmegridA.g2t[YY];
369 g2tz = pme->pmegridA.g2t[ZZ];
371 bThreads = (atc->nthread > 1);
374 thread_idx = atc->thread_idx;
376 tpl = &atc->thread_plist[thread];
378 for (i = 0; i < atc->nthread; i++)
384 for (i = start; i < end; i++)
387 idxptr = atc->idx[i];
388 fptr = atc->fractx[i];
390 /* Fractional coordinates along box vectors, add 2.0 to make 100% sure we are positive for triclinic boxes */
391 tx = nx * ( xptr[XX] * rxx + xptr[YY] * ryx + xptr[ZZ] * rzx + 2.0 );
392 ty = ny * ( xptr[YY] * ryy + xptr[ZZ] * rzy + 2.0 );
393 tz = nz * ( xptr[ZZ] * rzz + 2.0 );
399 /* Because decomposition only occurs in x and y,
400 * we never have a fraction correction in z.
402 fptr[XX] = tx - tix + pme->fshx[tix];
403 fptr[YY] = ty - tiy + pme->fshy[tiy];
406 idxptr[XX] = pme->nnx[tix];
407 idxptr[YY] = pme->nny[tiy];
408 idxptr[ZZ] = pme->nnz[tiz];
411 range_check(idxptr[XX], 0, pme->pmegrid_nx);
412 range_check(idxptr[YY], 0, pme->pmegrid_ny);
413 range_check(idxptr[ZZ], 0, pme->pmegrid_nz);
418 thread_i = g2tx[idxptr[XX]] + g2ty[idxptr[YY]] + g2tz[idxptr[ZZ]];
419 thread_idx[i] = thread_i;
426 /* Make a list of particle indices sorted on thread */
428 /* Get the cumulative count */
429 for (i = 1; i < atc->nthread; i++)
431 tpl_n[i] += tpl_n[i-1];
433 /* The current implementation distributes particles equally
434 * over the threads, so we could actually allocate for that
435 * in pme_realloc_atomcomm_things.
437 if (tpl_n[atc->nthread-1] > tpl->nalloc)
439 tpl->nalloc = over_alloc_large(tpl_n[atc->nthread-1]);
440 srenew(tpl->i, tpl->nalloc);
442 /* Set tpl_n to the cumulative start */
443 for (i = atc->nthread-1; i >= 1; i--)
445 tpl_n[i] = tpl_n[i-1];
449 /* Fill our thread local array with indices sorted on thread */
450 for (i = start; i < end; i++)
452 tpl->i[tpl_n[atc->thread_idx[i]]++] = i;
454 /* Now tpl_n contains the cummulative count again */
458 static void make_thread_local_ind(pme_atomcomm_t *atc,
459 int thread, splinedata_t *spline)
461 int n, t, i, start, end;
464 /* Combine the indices made by each thread into one index */
468 for (t = 0; t < atc->nthread; t++)
470 tpl = &atc->thread_plist[t];
471 /* Copy our part (start - end) from the list of thread t */
474 start = tpl->n[thread-1];
476 end = tpl->n[thread];
477 for (i = start; i < end; i++)
479 spline->ind[n++] = tpl->i[i];
487 static void pme_calc_pidx(int start, int end,
488 matrix recipbox, rvec x[],
489 pme_atomcomm_t *atc, int *count)
494 real rxx, ryx, rzx, ryy, rzy;
497 /* Calculate PME task index (pidx) for each grid index.
498 * Here we always assign equally sized slabs to each node
499 * for load balancing reasons (the PME grid spacing is not used).
505 /* Reset the count */
506 for (i = 0; i < nslab; i++)
511 if (atc->dimind == 0)
513 rxx = recipbox[XX][XX];
514 ryx = recipbox[YY][XX];
515 rzx = recipbox[ZZ][XX];
516 /* Calculate the node index in x-dimension */
517 for (i = start; i < end; i++)
520 /* Fractional coordinates along box vectors */
521 s = nslab*(xptr[XX]*rxx + xptr[YY]*ryx + xptr[ZZ]*rzx);
522 si = (int)(s + 2*nslab) % nslab;
529 ryy = recipbox[YY][YY];
530 rzy = recipbox[ZZ][YY];
531 /* Calculate the node index in y-dimension */
532 for (i = start; i < end; i++)
535 /* Fractional coordinates along box vectors */
536 s = nslab*(xptr[YY]*ryy + xptr[ZZ]*rzy);
537 si = (int)(s + 2*nslab) % nslab;
544 static void pme_calc_pidx_wrapper(int natoms, matrix recipbox, rvec x[],
547 int nthread, thread, slab;
549 nthread = atc->nthread;
551 #pragma omp parallel for num_threads(nthread) schedule(static)
552 for (thread = 0; thread < nthread; thread++)
554 pme_calc_pidx(natoms* thread /nthread,
555 natoms*(thread+1)/nthread,
556 recipbox, x, atc, atc->count_thread[thread]);
558 /* Non-parallel reduction, since nslab is small */
560 for (thread = 1; thread < nthread; thread++)
562 for (slab = 0; slab < atc->nslab; slab++)
564 atc->count_thread[0][slab] += atc->count_thread[thread][slab];
569 static void realloc_splinevec(splinevec th, real **ptr_z, int nalloc)
571 const int padding = 4;
574 srenew(th[XX], nalloc);
575 srenew(th[YY], nalloc);
576 /* In z we add padding, this is only required for the aligned SSE code */
577 srenew(*ptr_z, nalloc+2*padding);
578 th[ZZ] = *ptr_z + padding;
580 for (i = 0; i < padding; i++)
583 (*ptr_z)[padding+nalloc+i] = 0;
587 static void pme_realloc_splinedata(splinedata_t *spline, pme_atomcomm_t *atc)
591 srenew(spline->ind, atc->nalloc);
592 /* Initialize the index to identity so it works without threads */
593 for (i = 0; i < atc->nalloc; i++)
598 realloc_splinevec(spline->theta, &spline->ptr_theta_z,
599 atc->pme_order*atc->nalloc);
600 realloc_splinevec(spline->dtheta, &spline->ptr_dtheta_z,
601 atc->pme_order*atc->nalloc);
604 static void pme_realloc_atomcomm_things(pme_atomcomm_t *atc)
606 int nalloc_old, i, j, nalloc_tpl;
608 /* We have to avoid a NULL pointer for atc->x to avoid
609 * possible fatal errors in MPI routines.
611 if (atc->n > atc->nalloc || atc->nalloc == 0)
613 nalloc_old = atc->nalloc;
614 atc->nalloc = over_alloc_dd(max(atc->n, 1));
618 srenew(atc->x, atc->nalloc);
619 srenew(atc->q, atc->nalloc);
620 srenew(atc->f, atc->nalloc);
621 for (i = nalloc_old; i < atc->nalloc; i++)
623 clear_rvec(atc->f[i]);
628 srenew(atc->fractx, atc->nalloc);
629 srenew(atc->idx, atc->nalloc);
631 if (atc->nthread > 1)
633 srenew(atc->thread_idx, atc->nalloc);
636 for (i = 0; i < atc->nthread; i++)
638 pme_realloc_splinedata(&atc->spline[i], atc);
644 static void pmeredist_pd(gmx_pme_t pme, gmx_bool forw,
645 int n, gmx_bool bXF, rvec *x_f, real *charge,
647 /* Redistribute particle data for PME calculation */
648 /* domain decomposition by x coordinate */
653 if (FALSE == pme->redist_init)
655 snew(pme->scounts, atc->nslab);
656 snew(pme->rcounts, atc->nslab);
657 snew(pme->sdispls, atc->nslab);
658 snew(pme->rdispls, atc->nslab);
659 snew(pme->sidx, atc->nslab);
660 pme->redist_init = TRUE;
662 if (n > pme->redist_buf_nalloc)
664 pme->redist_buf_nalloc = over_alloc_dd(n);
665 srenew(pme->redist_buf, pme->redist_buf_nalloc*DIM);
673 /* forward, redistribution from pp to pme */
675 /* Calculate send counts and exchange them with other nodes */
676 for (i = 0; (i < atc->nslab); i++)
680 for (i = 0; (i < n); i++)
682 pme->scounts[pme->idxa[i]]++;
684 MPI_Alltoall( pme->scounts, 1, MPI_INT, pme->rcounts, 1, MPI_INT, atc->mpi_comm);
686 /* Calculate send and receive displacements and index into send
691 for (i = 1; i < atc->nslab; i++)
693 pme->sdispls[i] = pme->sdispls[i-1]+pme->scounts[i-1];
694 pme->rdispls[i] = pme->rdispls[i-1]+pme->rcounts[i-1];
695 pme->sidx[i] = pme->sdispls[i];
697 /* Total # of particles to be received */
698 atc->n = pme->rdispls[atc->nslab-1] + pme->rcounts[atc->nslab-1];
700 pme_realloc_atomcomm_things(atc);
702 /* Copy particle coordinates into send buffer and exchange*/
703 for (i = 0; (i < n); i++)
705 ii = DIM*pme->sidx[pme->idxa[i]];
706 pme->sidx[pme->idxa[i]]++;
707 pme->redist_buf[ii+XX] = x_f[i][XX];
708 pme->redist_buf[ii+YY] = x_f[i][YY];
709 pme->redist_buf[ii+ZZ] = x_f[i][ZZ];
711 MPI_Alltoallv(pme->redist_buf, pme->scounts, pme->sdispls,
712 pme->rvec_mpi, atc->x, pme->rcounts, pme->rdispls,
713 pme->rvec_mpi, atc->mpi_comm);
717 /* Copy charge into send buffer and exchange*/
718 for (i = 0; i < atc->nslab; i++)
720 pme->sidx[i] = pme->sdispls[i];
722 for (i = 0; (i < n); i++)
724 ii = pme->sidx[pme->idxa[i]];
725 pme->sidx[pme->idxa[i]]++;
726 pme->redist_buf[ii] = charge[i];
728 MPI_Alltoallv(pme->redist_buf, pme->scounts, pme->sdispls, mpi_type,
729 atc->q, pme->rcounts, pme->rdispls, mpi_type,
732 else /* backward, redistribution from pme to pp */
734 MPI_Alltoallv(atc->f, pme->rcounts, pme->rdispls, pme->rvec_mpi,
735 pme->redist_buf, pme->scounts, pme->sdispls,
736 pme->rvec_mpi, atc->mpi_comm);
738 /* Copy data from receive buffer */
739 for (i = 0; i < atc->nslab; i++)
741 pme->sidx[i] = pme->sdispls[i];
743 for (i = 0; (i < n); i++)
745 ii = DIM*pme->sidx[pme->idxa[i]];
746 x_f[i][XX] += pme->redist_buf[ii+XX];
747 x_f[i][YY] += pme->redist_buf[ii+YY];
748 x_f[i][ZZ] += pme->redist_buf[ii+ZZ];
749 pme->sidx[pme->idxa[i]]++;
755 static void pme_dd_sendrecv(pme_atomcomm_t *atc,
756 gmx_bool bBackward, int shift,
757 void *buf_s, int nbyte_s,
758 void *buf_r, int nbyte_r)
764 if (bBackward == FALSE)
766 dest = atc->node_dest[shift];
767 src = atc->node_src[shift];
771 dest = atc->node_src[shift];
772 src = atc->node_dest[shift];
775 if (nbyte_s > 0 && nbyte_r > 0)
777 MPI_Sendrecv(buf_s, nbyte_s, MPI_BYTE,
779 buf_r, nbyte_r, MPI_BYTE,
781 atc->mpi_comm, &stat);
783 else if (nbyte_s > 0)
785 MPI_Send(buf_s, nbyte_s, MPI_BYTE,
789 else if (nbyte_r > 0)
791 MPI_Recv(buf_r, nbyte_r, MPI_BYTE,
793 atc->mpi_comm, &stat);
798 static void dd_pmeredist_x_q(gmx_pme_t pme,
799 int n, gmx_bool bX, rvec *x, real *charge,
802 int *commnode, *buf_index;
803 int nnodes_comm, i, nsend, local_pos, buf_pos, node, scount, rcount;
805 commnode = atc->node_dest;
806 buf_index = atc->buf_index;
808 nnodes_comm = min(2*atc->maxshift, atc->nslab-1);
811 for (i = 0; i < nnodes_comm; i++)
813 buf_index[commnode[i]] = nsend;
814 nsend += atc->count[commnode[i]];
818 if (atc->count[atc->nodeid] + nsend != n)
820 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"
821 "This usually means that your system is not well equilibrated.",
822 n - (atc->count[atc->nodeid] + nsend),
823 pme->nodeid, 'x'+atc->dimind);
826 if (nsend > pme->buf_nalloc)
828 pme->buf_nalloc = over_alloc_dd(nsend);
829 srenew(pme->bufv, pme->buf_nalloc);
830 srenew(pme->bufr, pme->buf_nalloc);
833 atc->n = atc->count[atc->nodeid];
834 for (i = 0; i < nnodes_comm; i++)
836 scount = atc->count[commnode[i]];
837 /* Communicate the count */
840 fprintf(debug, "dimind %d PME node %d send to node %d: %d\n",
841 atc->dimind, atc->nodeid, commnode[i], scount);
843 pme_dd_sendrecv(atc, FALSE, i,
844 &scount, sizeof(int),
845 &atc->rcount[i], sizeof(int));
846 atc->n += atc->rcount[i];
849 pme_realloc_atomcomm_things(atc);
853 for (i = 0; i < n; i++)
856 if (node == atc->nodeid)
858 /* Copy direct to the receive buffer */
861 copy_rvec(x[i], atc->x[local_pos]);
863 atc->q[local_pos] = charge[i];
868 /* Copy to the send buffer */
871 copy_rvec(x[i], pme->bufv[buf_index[node]]);
873 pme->bufr[buf_index[node]] = charge[i];
879 for (i = 0; i < nnodes_comm; i++)
881 scount = atc->count[commnode[i]];
882 rcount = atc->rcount[i];
883 if (scount > 0 || rcount > 0)
887 /* Communicate the coordinates */
888 pme_dd_sendrecv(atc, FALSE, i,
889 pme->bufv[buf_pos], scount*sizeof(rvec),
890 atc->x[local_pos], rcount*sizeof(rvec));
892 /* Communicate the charges */
893 pme_dd_sendrecv(atc, FALSE, i,
894 pme->bufr+buf_pos, scount*sizeof(real),
895 atc->q+local_pos, rcount*sizeof(real));
897 local_pos += atc->rcount[i];
902 static void dd_pmeredist_f(gmx_pme_t pme, pme_atomcomm_t *atc,
906 int *commnode, *buf_index;
907 int nnodes_comm, local_pos, buf_pos, i, scount, rcount, node;
909 commnode = atc->node_dest;
910 buf_index = atc->buf_index;
912 nnodes_comm = min(2*atc->maxshift, atc->nslab-1);
914 local_pos = atc->count[atc->nodeid];
916 for (i = 0; i < nnodes_comm; i++)
918 scount = atc->rcount[i];
919 rcount = atc->count[commnode[i]];
920 if (scount > 0 || rcount > 0)
922 /* Communicate the forces */
923 pme_dd_sendrecv(atc, TRUE, i,
924 atc->f[local_pos], scount*sizeof(rvec),
925 pme->bufv[buf_pos], rcount*sizeof(rvec));
928 buf_index[commnode[i]] = buf_pos;
935 for (i = 0; i < n; i++)
938 if (node == atc->nodeid)
940 /* Add from the local force array */
941 rvec_inc(f[i], atc->f[local_pos]);
946 /* Add from the receive buffer */
947 rvec_inc(f[i], pme->bufv[buf_index[node]]);
954 for (i = 0; i < n; i++)
957 if (node == atc->nodeid)
959 /* Copy from the local force array */
960 copy_rvec(atc->f[local_pos], f[i]);
965 /* Copy from the receive buffer */
966 copy_rvec(pme->bufv[buf_index[node]], f[i]);
975 gmx_sum_qgrid_dd(gmx_pme_t pme, real *grid, int direction)
977 pme_overlap_t *overlap;
978 int send_index0, send_nindex;
979 int recv_index0, recv_nindex;
981 int i, j, k, ix, iy, iz, icnt;
982 int ipulse, send_id, recv_id, datasize;
984 real *sendptr, *recvptr;
986 /* Start with minor-rank communication. This is a bit of a pain since it is not contiguous */
987 overlap = &pme->overlap[1];
989 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
991 /* Since we have already (un)wrapped the overlap in the z-dimension,
992 * we only have to communicate 0 to nkz (not pmegrid_nz).
994 if (direction == GMX_SUM_QGRID_FORWARD)
996 send_id = overlap->send_id[ipulse];
997 recv_id = overlap->recv_id[ipulse];
998 send_index0 = overlap->comm_data[ipulse].send_index0;
999 send_nindex = overlap->comm_data[ipulse].send_nindex;
1000 recv_index0 = overlap->comm_data[ipulse].recv_index0;
1001 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
1005 send_id = overlap->recv_id[ipulse];
1006 recv_id = overlap->send_id[ipulse];
1007 send_index0 = overlap->comm_data[ipulse].recv_index0;
1008 send_nindex = overlap->comm_data[ipulse].recv_nindex;
1009 recv_index0 = overlap->comm_data[ipulse].send_index0;
1010 recv_nindex = overlap->comm_data[ipulse].send_nindex;
1013 /* Copy data to contiguous send buffer */
1016 fprintf(debug, "PME send node %d %d -> %d grid start %d Communicating %d to %d\n",
1017 pme->nodeid, overlap->nodeid, send_id,
1018 pme->pmegrid_start_iy,
1019 send_index0-pme->pmegrid_start_iy,
1020 send_index0-pme->pmegrid_start_iy+send_nindex);
1023 for (i = 0; i < pme->pmegrid_nx; i++)
1026 for (j = 0; j < send_nindex; j++)
1028 iy = j + send_index0 - pme->pmegrid_start_iy;
1029 for (k = 0; k < pme->nkz; k++)
1032 overlap->sendbuf[icnt++] = grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz];
1037 datasize = pme->pmegrid_nx * pme->nkz;
1039 MPI_Sendrecv(overlap->sendbuf, send_nindex*datasize, GMX_MPI_REAL,
1041 overlap->recvbuf, recv_nindex*datasize, GMX_MPI_REAL,
1043 overlap->mpi_comm, &stat);
1045 /* Get data from contiguous recv buffer */
1048 fprintf(debug, "PME recv node %d %d <- %d grid start %d Communicating %d to %d\n",
1049 pme->nodeid, overlap->nodeid, recv_id,
1050 pme->pmegrid_start_iy,
1051 recv_index0-pme->pmegrid_start_iy,
1052 recv_index0-pme->pmegrid_start_iy+recv_nindex);
1055 for (i = 0; i < pme->pmegrid_nx; i++)
1058 for (j = 0; j < recv_nindex; j++)
1060 iy = j + recv_index0 - pme->pmegrid_start_iy;
1061 for (k = 0; k < pme->nkz; k++)
1064 if (direction == GMX_SUM_QGRID_FORWARD)
1066 grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz] += overlap->recvbuf[icnt++];
1070 grid[ix*(pme->pmegrid_ny*pme->pmegrid_nz)+iy*(pme->pmegrid_nz)+iz] = overlap->recvbuf[icnt++];
1077 /* Major dimension is easier, no copying required,
1078 * but we might have to sum to separate array.
1079 * Since we don't copy, we have to communicate up to pmegrid_nz,
1080 * not nkz as for the minor direction.
1082 overlap = &pme->overlap[0];
1084 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
1086 if (direction == GMX_SUM_QGRID_FORWARD)
1088 send_id = overlap->send_id[ipulse];
1089 recv_id = overlap->recv_id[ipulse];
1090 send_index0 = overlap->comm_data[ipulse].send_index0;
1091 send_nindex = overlap->comm_data[ipulse].send_nindex;
1092 recv_index0 = overlap->comm_data[ipulse].recv_index0;
1093 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
1094 recvptr = overlap->recvbuf;
1098 send_id = overlap->recv_id[ipulse];
1099 recv_id = overlap->send_id[ipulse];
1100 send_index0 = overlap->comm_data[ipulse].recv_index0;
1101 send_nindex = overlap->comm_data[ipulse].recv_nindex;
1102 recv_index0 = overlap->comm_data[ipulse].send_index0;
1103 recv_nindex = overlap->comm_data[ipulse].send_nindex;
1104 recvptr = grid + (recv_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1107 sendptr = grid + (send_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1108 datasize = pme->pmegrid_ny * pme->pmegrid_nz;
1112 fprintf(debug, "PME send node %d %d -> %d grid start %d Communicating %d to %d\n",
1113 pme->nodeid, overlap->nodeid, send_id,
1114 pme->pmegrid_start_ix,
1115 send_index0-pme->pmegrid_start_ix,
1116 send_index0-pme->pmegrid_start_ix+send_nindex);
1117 fprintf(debug, "PME recv node %d %d <- %d grid start %d Communicating %d to %d\n",
1118 pme->nodeid, overlap->nodeid, recv_id,
1119 pme->pmegrid_start_ix,
1120 recv_index0-pme->pmegrid_start_ix,
1121 recv_index0-pme->pmegrid_start_ix+recv_nindex);
1124 MPI_Sendrecv(sendptr, send_nindex*datasize, GMX_MPI_REAL,
1126 recvptr, recv_nindex*datasize, GMX_MPI_REAL,
1128 overlap->mpi_comm, &stat);
1130 /* ADD data from contiguous recv buffer */
1131 if (direction == GMX_SUM_QGRID_FORWARD)
1133 p = grid + (recv_index0-pme->pmegrid_start_ix)*(pme->pmegrid_ny*pme->pmegrid_nz);
1134 for (i = 0; i < recv_nindex*datasize; i++)
1136 p[i] += overlap->recvbuf[i];
1145 copy_pmegrid_to_fftgrid(gmx_pme_t pme, real *pmegrid, real *fftgrid)
1147 ivec local_fft_ndata, local_fft_offset, local_fft_size;
1148 ivec local_pme_size;
1152 /* Dimensions should be identical for A/B grid, so we just use A here */
1153 gmx_parallel_3dfft_real_limits(pme->pfft_setupA,
1158 local_pme_size[0] = pme->pmegrid_nx;
1159 local_pme_size[1] = pme->pmegrid_ny;
1160 local_pme_size[2] = pme->pmegrid_nz;
1162 /* The fftgrid is always 'justified' to the lower-left corner of the PME grid,
1163 the offset is identical, and the PME grid always has more data (due to overlap)
1168 char fn[STRLEN], format[STRLEN];
1170 sprintf(fn, "pmegrid%d.pdb", pme->nodeid);
1171 fp = ffopen(fn, "w");
1172 sprintf(fn, "pmegrid%d.txt", pme->nodeid);
1173 fp2 = ffopen(fn, "w");
1174 sprintf(format, "%s%s\n", pdbformat, "%6.2f%6.2f");
1177 for (ix = 0; ix < local_fft_ndata[XX]; ix++)
1179 for (iy = 0; iy < local_fft_ndata[YY]; iy++)
1181 for (iz = 0; iz < local_fft_ndata[ZZ]; iz++)
1183 pmeidx = ix*(local_pme_size[YY]*local_pme_size[ZZ])+iy*(local_pme_size[ZZ])+iz;
1184 fftidx = ix*(local_fft_size[YY]*local_fft_size[ZZ])+iy*(local_fft_size[ZZ])+iz;
1185 fftgrid[fftidx] = pmegrid[pmeidx];
1187 val = 100*pmegrid[pmeidx];
1188 if (pmegrid[pmeidx] != 0)
1190 fprintf(fp, format, "ATOM", pmeidx, "CA", "GLY", ' ', pmeidx, ' ',
1191 5.0*ix, 5.0*iy, 5.0*iz, 1.0, val);
1193 if (pmegrid[pmeidx] != 0)
1195 fprintf(fp2, "%-12s %5d %5d %5d %12.5e\n",
1197 pme->pmegrid_start_ix + ix,
1198 pme->pmegrid_start_iy + iy,
1199 pme->pmegrid_start_iz + iz,
1215 static gmx_cycles_t omp_cyc_start()
1217 return gmx_cycles_read();
1220 static gmx_cycles_t omp_cyc_end(gmx_cycles_t c)
1222 return gmx_cycles_read() - c;
1227 copy_fftgrid_to_pmegrid(gmx_pme_t pme, const real *fftgrid, real *pmegrid,
1228 int nthread, int thread)
1230 ivec local_fft_ndata, local_fft_offset, local_fft_size;
1231 ivec local_pme_size;
1232 int ixy0, ixy1, ixy, ix, iy, iz;
1234 #ifdef PME_TIME_THREADS
1236 static double cs1 = 0;
1240 #ifdef PME_TIME_THREADS
1241 c1 = omp_cyc_start();
1243 /* Dimensions should be identical for A/B grid, so we just use A here */
1244 gmx_parallel_3dfft_real_limits(pme->pfft_setupA,
1249 local_pme_size[0] = pme->pmegrid_nx;
1250 local_pme_size[1] = pme->pmegrid_ny;
1251 local_pme_size[2] = pme->pmegrid_nz;
1253 /* The fftgrid is always 'justified' to the lower-left corner of the PME grid,
1254 the offset is identical, and the PME grid always has more data (due to overlap)
1256 ixy0 = ((thread )*local_fft_ndata[XX]*local_fft_ndata[YY])/nthread;
1257 ixy1 = ((thread+1)*local_fft_ndata[XX]*local_fft_ndata[YY])/nthread;
1259 for (ixy = ixy0; ixy < ixy1; ixy++)
1261 ix = ixy/local_fft_ndata[YY];
1262 iy = ixy - ix*local_fft_ndata[YY];
1264 pmeidx = (ix*local_pme_size[YY] + iy)*local_pme_size[ZZ];
1265 fftidx = (ix*local_fft_size[YY] + iy)*local_fft_size[ZZ];
1266 for (iz = 0; iz < local_fft_ndata[ZZ]; iz++)
1268 pmegrid[pmeidx+iz] = fftgrid[fftidx+iz];
1272 #ifdef PME_TIME_THREADS
1273 c1 = omp_cyc_end(c1);
1278 printf("copy %.2f\n", cs1*1e-9);
1287 wrap_periodic_pmegrid(gmx_pme_t pme, real *pmegrid)
1289 int nx, ny, nz, pnx, pny, pnz, ny_x, overlap, ix, iy, iz;
1295 pnx = pme->pmegrid_nx;
1296 pny = pme->pmegrid_ny;
1297 pnz = pme->pmegrid_nz;
1299 overlap = pme->pme_order - 1;
1301 /* Add periodic overlap in z */
1302 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1304 for (iy = 0; iy < pme->pmegrid_ny; iy++)
1306 for (iz = 0; iz < overlap; iz++)
1308 pmegrid[(ix*pny+iy)*pnz+iz] +=
1309 pmegrid[(ix*pny+iy)*pnz+nz+iz];
1314 if (pme->nnodes_minor == 1)
1316 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1318 for (iy = 0; iy < overlap; iy++)
1320 for (iz = 0; iz < nz; iz++)
1322 pmegrid[(ix*pny+iy)*pnz+iz] +=
1323 pmegrid[(ix*pny+ny+iy)*pnz+iz];
1329 if (pme->nnodes_major == 1)
1331 ny_x = (pme->nnodes_minor == 1 ? ny : pme->pmegrid_ny);
1333 for (ix = 0; ix < overlap; ix++)
1335 for (iy = 0; iy < ny_x; iy++)
1337 for (iz = 0; iz < nz; iz++)
1339 pmegrid[(ix*pny+iy)*pnz+iz] +=
1340 pmegrid[((nx+ix)*pny+iy)*pnz+iz];
1349 unwrap_periodic_pmegrid(gmx_pme_t pme, real *pmegrid)
1351 int nx, ny, nz, pnx, pny, pnz, ny_x, overlap, ix;
1357 pnx = pme->pmegrid_nx;
1358 pny = pme->pmegrid_ny;
1359 pnz = pme->pmegrid_nz;
1361 overlap = pme->pme_order - 1;
1363 if (pme->nnodes_major == 1)
1365 ny_x = (pme->nnodes_minor == 1 ? ny : pme->pmegrid_ny);
1367 for (ix = 0; ix < overlap; ix++)
1371 for (iy = 0; iy < ny_x; iy++)
1373 for (iz = 0; iz < nz; iz++)
1375 pmegrid[((nx+ix)*pny+iy)*pnz+iz] =
1376 pmegrid[(ix*pny+iy)*pnz+iz];
1382 if (pme->nnodes_minor == 1)
1384 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1385 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1389 for (iy = 0; iy < overlap; iy++)
1391 for (iz = 0; iz < nz; iz++)
1393 pmegrid[(ix*pny+ny+iy)*pnz+iz] =
1394 pmegrid[(ix*pny+iy)*pnz+iz];
1400 /* Copy periodic overlap in z */
1401 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
1402 for (ix = 0; ix < pme->pmegrid_nx; ix++)
1406 for (iy = 0; iy < pme->pmegrid_ny; iy++)
1408 for (iz = 0; iz < overlap; iz++)
1410 pmegrid[(ix*pny+iy)*pnz+nz+iz] =
1411 pmegrid[(ix*pny+iy)*pnz+iz];
1417 static void clear_grid(int nx, int ny, int nz, real *grid,
1419 int fx, int fy, int fz,
1423 int fsx, fsy, fsz, gx, gy, gz, g0x, g0y, x, y, z;
1426 nc = 2 + (order - 2)/FLBS;
1427 ncz = 2 + (order - 2)/FLBSZ;
1429 for (fsx = fx; fsx < fx+nc; fsx++)
1431 for (fsy = fy; fsy < fy+nc; fsy++)
1433 for (fsz = fz; fsz < fz+ncz; fsz++)
1435 flind = (fsx*fs[YY] + fsy)*fs[ZZ] + fsz;
1436 if (flag[flind] == 0)
1441 g0x = (gx*ny + gy)*nz + gz;
1442 for (x = 0; x < FLBS; x++)
1445 for (y = 0; y < FLBS; y++)
1447 for (z = 0; z < FLBSZ; z++)
1463 /* This has to be a macro to enable full compiler optimization with xlC (and probably others too) */
1464 #define DO_BSPLINE(order) \
1465 for (ithx = 0; (ithx < order); ithx++) \
1467 index_x = (i0+ithx)*pny*pnz; \
1468 valx = qn*thx[ithx]; \
1470 for (ithy = 0; (ithy < order); ithy++) \
1472 valxy = valx*thy[ithy]; \
1473 index_xy = index_x+(j0+ithy)*pnz; \
1475 for (ithz = 0; (ithz < order); ithz++) \
1477 index_xyz = index_xy+(k0+ithz); \
1478 grid[index_xyz] += valxy*thz[ithz]; \
1484 static void spread_q_bsplines_thread(pmegrid_t *pmegrid,
1485 pme_atomcomm_t *atc, splinedata_t *spline,
1486 pme_spline_work_t *work)
1489 /* spread charges from home atoms to local grid */
1492 int b, i, nn, n, ithx, ithy, ithz, i0, j0, k0;
1494 int order, norder, index_x, index_xy, index_xyz;
1495 real valx, valxy, qn;
1496 real *thx, *thy, *thz;
1497 int localsize, bndsize;
1498 int pnx, pny, pnz, ndatatot;
1499 int offx, offy, offz;
1501 pnx = pmegrid->s[XX];
1502 pny = pmegrid->s[YY];
1503 pnz = pmegrid->s[ZZ];
1505 offx = pmegrid->offset[XX];
1506 offy = pmegrid->offset[YY];
1507 offz = pmegrid->offset[ZZ];
1509 ndatatot = pnx*pny*pnz;
1510 grid = pmegrid->grid;
1511 for (i = 0; i < ndatatot; i++)
1516 order = pmegrid->order;
1518 for (nn = 0; nn < spline->n; nn++)
1520 n = spline->ind[nn];
1525 idxptr = atc->idx[n];
1528 i0 = idxptr[XX] - offx;
1529 j0 = idxptr[YY] - offy;
1530 k0 = idxptr[ZZ] - offz;
1532 thx = spline->theta[XX] + norder;
1533 thy = spline->theta[YY] + norder;
1534 thz = spline->theta[ZZ] + norder;
1540 #ifdef PME_SSE_UNALIGNED
1541 #define PME_SPREAD_SSE_ORDER4
1543 #define PME_SPREAD_SSE_ALIGNED
1546 #include "pme_sse_single.h"
1553 #define PME_SPREAD_SSE_ALIGNED
1555 #include "pme_sse_single.h"
1568 static void set_grid_alignment(int *pmegrid_nz, int pme_order)
1572 #ifndef PME_SSE_UNALIGNED
1577 /* Round nz up to a multiple of 4 to ensure alignment */
1578 *pmegrid_nz = ((*pmegrid_nz + 3) & ~3);
1583 static void set_gridsize_alignment(int *gridsize, int pme_order)
1586 #ifndef PME_SSE_UNALIGNED
1589 /* Add extra elements to ensured aligned operations do not go
1590 * beyond the allocated grid size.
1591 * Note that for pme_order=5, the pme grid z-size alignment
1592 * ensures that we will not go beyond the grid size.
1600 static void pmegrid_init(pmegrid_t *grid,
1601 int cx, int cy, int cz,
1602 int x0, int y0, int z0,
1603 int x1, int y1, int z1,
1604 gmx_bool set_alignment,
1613 grid->offset[XX] = x0;
1614 grid->offset[YY] = y0;
1615 grid->offset[ZZ] = z0;
1616 grid->n[XX] = x1 - x0 + pme_order - 1;
1617 grid->n[YY] = y1 - y0 + pme_order - 1;
1618 grid->n[ZZ] = z1 - z0 + pme_order - 1;
1619 copy_ivec(grid->n, grid->s);
1622 set_grid_alignment(&nz, pme_order);
1627 else if (nz != grid->s[ZZ])
1629 gmx_incons("pmegrid_init call with an unaligned z size");
1632 grid->order = pme_order;
1635 gridsize = grid->s[XX]*grid->s[YY]*grid->s[ZZ];
1636 set_gridsize_alignment(&gridsize, pme_order);
1637 snew_aligned(grid->grid, gridsize, 16);
1645 static int div_round_up(int enumerator, int denominator)
1647 return (enumerator + denominator - 1)/denominator;
1650 static void make_subgrid_division(const ivec n, int ovl, int nthread,
1653 int gsize_opt, gsize;
1658 for (nsx = 1; nsx <= nthread; nsx++)
1660 if (nthread % nsx == 0)
1662 for (nsy = 1; nsy <= nthread; nsy++)
1664 if (nsx*nsy <= nthread && nthread % (nsx*nsy) == 0)
1666 nsz = nthread/(nsx*nsy);
1668 /* Determine the number of grid points per thread */
1670 (div_round_up(n[XX], nsx) + ovl)*
1671 (div_round_up(n[YY], nsy) + ovl)*
1672 (div_round_up(n[ZZ], nsz) + ovl);
1674 /* Minimize the number of grids points per thread
1675 * and, secondarily, the number of cuts in minor dimensions.
1677 if (gsize_opt == -1 ||
1678 gsize < gsize_opt ||
1679 (gsize == gsize_opt &&
1680 (nsz < nsub[ZZ] || (nsz == nsub[ZZ] && nsy < nsub[YY]))))
1692 env = getenv("GMX_PME_THREAD_DIVISION");
1695 sscanf(env, "%d %d %d", &nsub[XX], &nsub[YY], &nsub[ZZ]);
1698 if (nsub[XX]*nsub[YY]*nsub[ZZ] != nthread)
1700 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);
1704 static void pmegrids_init(pmegrids_t *grids,
1705 int nx, int ny, int nz, int nz_base,
1707 gmx_bool bUseThreads,
1712 ivec n, n_base, g0, g1;
1713 int t, x, y, z, d, i, tfac;
1714 int max_comm_lines = -1;
1716 n[XX] = nx - (pme_order - 1);
1717 n[YY] = ny - (pme_order - 1);
1718 n[ZZ] = nz - (pme_order - 1);
1720 copy_ivec(n, n_base);
1721 n_base[ZZ] = nz_base;
1723 pmegrid_init(&grids->grid, 0, 0, 0, 0, 0, 0, n[XX], n[YY], n[ZZ], FALSE, pme_order,
1726 grids->nthread = nthread;
1728 make_subgrid_division(n_base, pme_order-1, grids->nthread, grids->nc);
1735 for (d = 0; d < DIM; d++)
1737 nst[d] = div_round_up(n[d], grids->nc[d]) + pme_order - 1;
1739 set_grid_alignment(&nst[ZZ], pme_order);
1743 fprintf(debug, "pmegrid thread local division: %d x %d x %d\n",
1744 grids->nc[XX], grids->nc[YY], grids->nc[ZZ]);
1745 fprintf(debug, "pmegrid %d %d %d max thread pmegrid %d %d %d\n",
1747 nst[XX], nst[YY], nst[ZZ]);
1750 snew(grids->grid_th, grids->nthread);
1752 gridsize = nst[XX]*nst[YY]*nst[ZZ];
1753 set_gridsize_alignment(&gridsize, pme_order);
1754 snew_aligned(grids->grid_all,
1755 grids->nthread*gridsize+(grids->nthread+1)*GMX_CACHE_SEP,
1758 for (x = 0; x < grids->nc[XX]; x++)
1760 for (y = 0; y < grids->nc[YY]; y++)
1762 for (z = 0; z < grids->nc[ZZ]; z++)
1764 pmegrid_init(&grids->grid_th[t],
1766 (n[XX]*(x ))/grids->nc[XX],
1767 (n[YY]*(y ))/grids->nc[YY],
1768 (n[ZZ]*(z ))/grids->nc[ZZ],
1769 (n[XX]*(x+1))/grids->nc[XX],
1770 (n[YY]*(y+1))/grids->nc[YY],
1771 (n[ZZ]*(z+1))/grids->nc[ZZ],
1774 grids->grid_all+GMX_CACHE_SEP+t*(gridsize+GMX_CACHE_SEP));
1782 grids->grid_th = NULL;
1785 snew(grids->g2t, DIM);
1787 for (d = DIM-1; d >= 0; d--)
1789 snew(grids->g2t[d], n[d]);
1791 for (i = 0; i < n[d]; i++)
1793 /* The second check should match the parameters
1794 * of the pmegrid_init call above.
1796 while (t + 1 < grids->nc[d] && i >= (n[d]*(t+1))/grids->nc[d])
1800 grids->g2t[d][i] = t*tfac;
1803 tfac *= grids->nc[d];
1807 case XX: max_comm_lines = overlap_x; break;
1808 case YY: max_comm_lines = overlap_y; break;
1809 case ZZ: max_comm_lines = pme_order - 1; break;
1811 grids->nthread_comm[d] = 0;
1812 while ((n[d]*grids->nthread_comm[d])/grids->nc[d] < max_comm_lines &&
1813 grids->nthread_comm[d] < grids->nc[d])
1815 grids->nthread_comm[d]++;
1819 fprintf(debug, "pmegrid thread grid communication range in %c: %d\n",
1820 'x'+d, grids->nthread_comm[d]);
1822 /* It should be possible to make grids->nthread_comm[d]==grids->nc[d]
1823 * work, but this is not a problematic restriction.
1825 if (grids->nc[d] > 1 && grids->nthread_comm[d] > grids->nc[d])
1827 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);
1833 static void pmegrids_destroy(pmegrids_t *grids)
1837 if (grids->grid.grid != NULL)
1839 sfree(grids->grid.grid);
1841 if (grids->nthread > 0)
1843 for (t = 0; t < grids->nthread; t++)
1845 sfree(grids->grid_th[t].grid);
1847 sfree(grids->grid_th);
1853 static void realloc_work(pme_work_t *work, int nkx)
1855 if (nkx > work->nalloc)
1858 srenew(work->mhx, work->nalloc);
1859 srenew(work->mhy, work->nalloc);
1860 srenew(work->mhz, work->nalloc);
1861 srenew(work->m2, work->nalloc);
1862 /* Allocate an aligned pointer for SSE operations, including 3 extra
1863 * elements at the end since SSE operates on 4 elements at a time.
1865 sfree_aligned(work->denom);
1866 sfree_aligned(work->tmp1);
1867 sfree_aligned(work->eterm);
1868 snew_aligned(work->denom, work->nalloc+3, 16);
1869 snew_aligned(work->tmp1, work->nalloc+3, 16);
1870 snew_aligned(work->eterm, work->nalloc+3, 16);
1871 srenew(work->m2inv, work->nalloc);
1876 static void free_work(pme_work_t *work)
1882 sfree_aligned(work->denom);
1883 sfree_aligned(work->tmp1);
1884 sfree_aligned(work->eterm);
1890 /* Calculate exponentials through SSE in float precision */
1891 inline static void calc_exponentials(int start, int end, real f, real *d_aligned, real *r_aligned, real *e_aligned)
1894 const __m128 two = _mm_set_ps(2.0f, 2.0f, 2.0f, 2.0f);
1897 __m128 tmp_d1, d_inv, tmp_r, tmp_e;
1899 f_sse = _mm_load1_ps(&f);
1900 for (kx = 0; kx < end; kx += 4)
1902 tmp_d1 = _mm_load_ps(d_aligned+kx);
1903 lu = _mm_rcp_ps(tmp_d1);
1904 d_inv = _mm_mul_ps(lu, _mm_sub_ps(two, _mm_mul_ps(lu, tmp_d1)));
1905 tmp_r = _mm_load_ps(r_aligned+kx);
1906 tmp_r = gmx_mm_exp_ps(tmp_r);
1907 tmp_e = _mm_mul_ps(f_sse, d_inv);
1908 tmp_e = _mm_mul_ps(tmp_e, tmp_r);
1909 _mm_store_ps(e_aligned+kx, tmp_e);
1914 inline static void calc_exponentials(int start, int end, real f, real *d, real *r, real *e)
1917 for (kx = start; kx < end; kx++)
1921 for (kx = start; kx < end; kx++)
1925 for (kx = start; kx < end; kx++)
1927 e[kx] = f*r[kx]*d[kx];
1933 static int solve_pme_yzx(gmx_pme_t pme, t_complex *grid,
1934 real ewaldcoeff, real vol,
1936 int nthread, int thread)
1938 /* do recip sum over local cells in grid */
1939 /* y major, z middle, x minor or continuous */
1941 int kx, ky, kz, maxkx, maxky, maxkz;
1942 int nx, ny, nz, iyz0, iyz1, iyz, iy, iz, kxstart, kxend;
1944 real factor = M_PI*M_PI/(ewaldcoeff*ewaldcoeff);
1945 real ets2, struct2, vfactor, ets2vf;
1946 real d1, d2, energy = 0;
1948 real virxx = 0, virxy = 0, virxz = 0, viryy = 0, viryz = 0, virzz = 0;
1949 real rxx, ryx, ryy, rzx, rzy, rzz;
1951 real *mhx, *mhy, *mhz, *m2, *denom, *tmp1, *eterm, *m2inv;
1952 real mhxk, mhyk, mhzk, m2k;
1955 ivec local_ndata, local_offset, local_size;
1958 elfac = ONE_4PI_EPS0/pme->epsilon_r;
1964 /* Dimensions should be identical for A/B grid, so we just use A here */
1965 gmx_parallel_3dfft_complex_limits(pme->pfft_setupA,
1971 rxx = pme->recipbox[XX][XX];
1972 ryx = pme->recipbox[YY][XX];
1973 ryy = pme->recipbox[YY][YY];
1974 rzx = pme->recipbox[ZZ][XX];
1975 rzy = pme->recipbox[ZZ][YY];
1976 rzz = pme->recipbox[ZZ][ZZ];
1982 work = &pme->work[thread];
1987 denom = work->denom;
1989 eterm = work->eterm;
1990 m2inv = work->m2inv;
1992 iyz0 = local_ndata[YY]*local_ndata[ZZ]* thread /nthread;
1993 iyz1 = local_ndata[YY]*local_ndata[ZZ]*(thread+1)/nthread;
1995 for (iyz = iyz0; iyz < iyz1; iyz++)
1997 iy = iyz/local_ndata[ZZ];
1998 iz = iyz - iy*local_ndata[ZZ];
2000 ky = iy + local_offset[YY];
2011 by = M_PI*vol*pme->bsp_mod[YY][ky];
2013 kz = iz + local_offset[ZZ];
2017 bz = pme->bsp_mod[ZZ][kz];
2019 /* 0.5 correction for corner points */
2021 if (kz == 0 || kz == (nz+1)/2)
2026 p0 = grid + iy*local_size[ZZ]*local_size[XX] + iz*local_size[XX];
2028 /* We should skip the k-space point (0,0,0) */
2029 if (local_offset[XX] > 0 || ky > 0 || kz > 0)
2031 kxstart = local_offset[XX];
2035 kxstart = local_offset[XX] + 1;
2038 kxend = local_offset[XX] + local_ndata[XX];
2042 /* More expensive inner loop, especially because of the storage
2043 * of the mh elements in array's.
2044 * Because x is the minor grid index, all mh elements
2045 * depend on kx for triclinic unit cells.
2048 /* Two explicit loops to avoid a conditional inside the loop */
2049 for (kx = kxstart; kx < maxkx; kx++)
2054 mhyk = mx * ryx + my * ryy;
2055 mhzk = mx * rzx + my * rzy + mz * rzz;
2056 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2061 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2062 tmp1[kx] = -factor*m2k;
2065 for (kx = maxkx; kx < kxend; kx++)
2070 mhyk = mx * ryx + my * ryy;
2071 mhzk = mx * rzx + my * rzy + mz * rzz;
2072 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2077 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2078 tmp1[kx] = -factor*m2k;
2081 for (kx = kxstart; kx < kxend; kx++)
2083 m2inv[kx] = 1.0/m2[kx];
2086 calc_exponentials(kxstart, kxend, elfac, denom, tmp1, eterm);
2088 for (kx = kxstart; kx < kxend; kx++, p0++)
2093 p0->re = d1*eterm[kx];
2094 p0->im = d2*eterm[kx];
2096 struct2 = 2.0*(d1*d1+d2*d2);
2098 tmp1[kx] = eterm[kx]*struct2;
2101 for (kx = kxstart; kx < kxend; kx++)
2103 ets2 = corner_fac*tmp1[kx];
2104 vfactor = (factor*m2[kx] + 1.0)*2.0*m2inv[kx];
2107 ets2vf = ets2*vfactor;
2108 virxx += ets2vf*mhx[kx]*mhx[kx] - ets2;
2109 virxy += ets2vf*mhx[kx]*mhy[kx];
2110 virxz += ets2vf*mhx[kx]*mhz[kx];
2111 viryy += ets2vf*mhy[kx]*mhy[kx] - ets2;
2112 viryz += ets2vf*mhy[kx]*mhz[kx];
2113 virzz += ets2vf*mhz[kx]*mhz[kx] - ets2;
2118 /* We don't need to calculate the energy and the virial.
2119 * In this case the triclinic overhead is small.
2122 /* Two explicit loops to avoid a conditional inside the loop */
2124 for (kx = kxstart; kx < maxkx; kx++)
2129 mhyk = mx * ryx + my * ryy;
2130 mhzk = mx * rzx + my * rzy + mz * rzz;
2131 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2132 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2133 tmp1[kx] = -factor*m2k;
2136 for (kx = maxkx; kx < kxend; kx++)
2141 mhyk = mx * ryx + my * ryy;
2142 mhzk = mx * rzx + my * rzy + mz * rzz;
2143 m2k = mhxk*mhxk + mhyk*mhyk + mhzk*mhzk;
2144 denom[kx] = m2k*bz*by*pme->bsp_mod[XX][kx];
2145 tmp1[kx] = -factor*m2k;
2148 calc_exponentials(kxstart, kxend, elfac, denom, tmp1, eterm);
2150 for (kx = kxstart; kx < kxend; kx++, p0++)
2155 p0->re = d1*eterm[kx];
2156 p0->im = d2*eterm[kx];
2163 /* Update virial with local values.
2164 * The virial is symmetric by definition.
2165 * this virial seems ok for isotropic scaling, but I'm
2166 * experiencing problems on semiisotropic membranes.
2167 * IS THAT COMMENT STILL VALID??? (DvdS, 2001/02/07).
2169 work->vir[XX][XX] = 0.25*virxx;
2170 work->vir[YY][YY] = 0.25*viryy;
2171 work->vir[ZZ][ZZ] = 0.25*virzz;
2172 work->vir[XX][YY] = work->vir[YY][XX] = 0.25*virxy;
2173 work->vir[XX][ZZ] = work->vir[ZZ][XX] = 0.25*virxz;
2174 work->vir[YY][ZZ] = work->vir[ZZ][YY] = 0.25*viryz;
2176 /* This energy should be corrected for a charged system */
2177 work->energy = 0.5*energy;
2180 /* Return the loop count */
2181 return local_ndata[YY]*local_ndata[XX];
2184 static void get_pme_ener_vir(const gmx_pme_t pme, int nthread,
2185 real *mesh_energy, matrix vir)
2187 /* This function sums output over threads
2188 * and should therefore only be called after thread synchronization.
2192 *mesh_energy = pme->work[0].energy;
2193 copy_mat(pme->work[0].vir, vir);
2195 for (thread = 1; thread < nthread; thread++)
2197 *mesh_energy += pme->work[thread].energy;
2198 m_add(vir, pme->work[thread].vir, vir);
2202 #define DO_FSPLINE(order) \
2203 for (ithx = 0; (ithx < order); ithx++) \
2205 index_x = (i0+ithx)*pny*pnz; \
2209 for (ithy = 0; (ithy < order); ithy++) \
2211 index_xy = index_x+(j0+ithy)*pnz; \
2216 for (ithz = 0; (ithz < order); ithz++) \
2218 gval = grid[index_xy+(k0+ithz)]; \
2219 fxy1 += thz[ithz]*gval; \
2220 fz1 += dthz[ithz]*gval; \
2229 static void gather_f_bsplines(gmx_pme_t pme, real *grid,
2230 gmx_bool bClearF, pme_atomcomm_t *atc,
2231 splinedata_t *spline,
2234 /* sum forces for local particles */
2235 int nn, n, ithx, ithy, ithz, i0, j0, k0;
2236 int index_x, index_xy;
2237 int nx, ny, nz, pnx, pny, pnz;
2239 real tx, ty, dx, dy, qn;
2240 real fx, fy, fz, gval;
2242 real *thx, *thy, *thz, *dthx, *dthy, *dthz;
2244 real rxx, ryx, ryy, rzx, rzy, rzz;
2247 pme_spline_work_t *work;
2249 work = pme->spline_work;
2251 order = pme->pme_order;
2252 thx = spline->theta[XX];
2253 thy = spline->theta[YY];
2254 thz = spline->theta[ZZ];
2255 dthx = spline->dtheta[XX];
2256 dthy = spline->dtheta[YY];
2257 dthz = spline->dtheta[ZZ];
2261 pnx = pme->pmegrid_nx;
2262 pny = pme->pmegrid_ny;
2263 pnz = pme->pmegrid_nz;
2265 rxx = pme->recipbox[XX][XX];
2266 ryx = pme->recipbox[YY][XX];
2267 ryy = pme->recipbox[YY][YY];
2268 rzx = pme->recipbox[ZZ][XX];
2269 rzy = pme->recipbox[ZZ][YY];
2270 rzz = pme->recipbox[ZZ][ZZ];
2272 for (nn = 0; nn < spline->n; nn++)
2274 n = spline->ind[nn];
2275 qn = scale*atc->q[n];
2288 idxptr = atc->idx[n];
2295 /* Pointer arithmetic alert, next six statements */
2296 thx = spline->theta[XX] + norder;
2297 thy = spline->theta[YY] + norder;
2298 thz = spline->theta[ZZ] + norder;
2299 dthx = spline->dtheta[XX] + norder;
2300 dthy = spline->dtheta[YY] + norder;
2301 dthz = spline->dtheta[ZZ] + norder;
2307 #ifdef PME_SSE_UNALIGNED
2308 #define PME_GATHER_F_SSE_ORDER4
2310 #define PME_GATHER_F_SSE_ALIGNED
2313 #include "pme_sse_single.h"
2320 #define PME_GATHER_F_SSE_ALIGNED
2322 #include "pme_sse_single.h"
2332 atc->f[n][XX] += -qn*( fx*nx*rxx );
2333 atc->f[n][YY] += -qn*( fx*nx*ryx + fy*ny*ryy );
2334 atc->f[n][ZZ] += -qn*( fx*nx*rzx + fy*ny*rzy + fz*nz*rzz );
2337 /* Since the energy and not forces are interpolated
2338 * the net force might not be exactly zero.
2339 * This can be solved by also interpolating F, but
2340 * that comes at a cost.
2341 * A better hack is to remove the net force every
2342 * step, but that must be done at a higher level
2343 * since this routine doesn't see all atoms if running
2344 * in parallel. Don't know how important it is? EL 990726
2349 static real gather_energy_bsplines(gmx_pme_t pme, real *grid,
2350 pme_atomcomm_t *atc)
2352 splinedata_t *spline;
2353 int n, ithx, ithy, ithz, i0, j0, k0;
2354 int index_x, index_xy;
2356 real energy, pot, tx, ty, qn, gval;
2357 real *thx, *thy, *thz;
2361 spline = &atc->spline[0];
2363 order = pme->pme_order;
2366 for (n = 0; (n < atc->n); n++)
2372 idxptr = atc->idx[n];
2379 /* Pointer arithmetic alert, next three statements */
2380 thx = spline->theta[XX] + norder;
2381 thy = spline->theta[YY] + norder;
2382 thz = spline->theta[ZZ] + norder;
2385 for (ithx = 0; (ithx < order); ithx++)
2387 index_x = (i0+ithx)*pme->pmegrid_ny*pme->pmegrid_nz;
2390 for (ithy = 0; (ithy < order); ithy++)
2392 index_xy = index_x+(j0+ithy)*pme->pmegrid_nz;
2395 for (ithz = 0; (ithz < order); ithz++)
2397 gval = grid[index_xy+(k0+ithz)];
2398 pot += tx*ty*thz[ithz]*gval;
2411 /* Macro to force loop unrolling by fixing order.
2412 * This gives a significant performance gain.
2414 #define CALC_SPLINE(order) \
2418 real data[PME_ORDER_MAX]; \
2419 real ddata[PME_ORDER_MAX]; \
2421 for (j = 0; (j < DIM); j++) \
2425 /* dr is relative offset from lower cell limit */ \
2426 data[order-1] = 0; \
2430 for (k = 3; (k < order); k++) \
2432 div = 1.0/(k - 1.0); \
2433 data[k-1] = div*dr*data[k-2]; \
2434 for (l = 1; (l < (k-1)); l++) \
2436 data[k-l-1] = div*((dr+l)*data[k-l-2]+(k-l-dr)* \
2439 data[0] = div*(1-dr)*data[0]; \
2441 /* differentiate */ \
2442 ddata[0] = -data[0]; \
2443 for (k = 1; (k < order); k++) \
2445 ddata[k] = data[k-1] - data[k]; \
2448 div = 1.0/(order - 1); \
2449 data[order-1] = div*dr*data[order-2]; \
2450 for (l = 1; (l < (order-1)); l++) \
2452 data[order-l-1] = div*((dr+l)*data[order-l-2]+ \
2453 (order-l-dr)*data[order-l-1]); \
2455 data[0] = div*(1 - dr)*data[0]; \
2457 for (k = 0; k < order; k++) \
2459 theta[j][i*order+k] = data[k]; \
2460 dtheta[j][i*order+k] = ddata[k]; \
2465 void make_bsplines(splinevec theta, splinevec dtheta, int order,
2466 rvec fractx[], int nr, int ind[], real charge[],
2467 gmx_bool bFreeEnergy)
2469 /* construct splines for local atoms */
2473 for (i = 0; i < nr; i++)
2475 /* With free energy we do not use the charge check.
2476 * In most cases this will be more efficient than calling make_bsplines
2477 * twice, since usually more than half the particles have charges.
2480 if (bFreeEnergy || charge[ii] != 0.0)
2485 case 4: CALC_SPLINE(4); break;
2486 case 5: CALC_SPLINE(5); break;
2487 default: CALC_SPLINE(order); break;
2494 void make_dft_mod(real *mod, real *data, int ndata)
2499 for (i = 0; i < ndata; i++)
2502 for (j = 0; j < ndata; j++)
2504 arg = (2.0*M_PI*i*j)/ndata;
2505 sc += data[j]*cos(arg);
2506 ss += data[j]*sin(arg);
2508 mod[i] = sc*sc+ss*ss;
2510 for (i = 0; i < ndata; i++)
2514 mod[i] = (mod[i-1]+mod[i+1])*0.5;
2520 static void make_bspline_moduli(splinevec bsp_mod,
2521 int nx, int ny, int nz, int order)
2523 int nmax = max(nx, max(ny, nz));
2524 real *data, *ddata, *bsp_data;
2530 snew(bsp_data, nmax);
2536 for (k = 3; k < order; k++)
2540 for (l = 1; l < (k-1); l++)
2542 data[k-l-1] = div*(l*data[k-l-2]+(k-l)*data[k-l-1]);
2544 data[0] = div*data[0];
2547 ddata[0] = -data[0];
2548 for (k = 1; k < order; k++)
2550 ddata[k] = data[k-1]-data[k];
2552 div = 1.0/(order-1);
2554 for (l = 1; l < (order-1); l++)
2556 data[order-l-1] = div*(l*data[order-l-2]+(order-l)*data[order-l-1]);
2558 data[0] = div*data[0];
2560 for (i = 0; i < nmax; i++)
2564 for (i = 1; i <= order; i++)
2566 bsp_data[i] = data[i-1];
2569 make_dft_mod(bsp_mod[XX], bsp_data, nx);
2570 make_dft_mod(bsp_mod[YY], bsp_data, ny);
2571 make_dft_mod(bsp_mod[ZZ], bsp_data, nz);
2579 /* Return the P3M optimal influence function */
2580 static double do_p3m_influence(double z, int order)
2587 /* The formula and most constants can be found in:
2588 * Ballenegger et al., JCTC 8, 936 (2012)
2593 return 1.0 - 2.0*z2/3.0;
2596 return 1.0 - z2 + 2.0*z4/15.0;
2599 return 1.0 - 4.0*z2/3.0 + 2.0*z4/5.0 + 4.0*z2*z4/315.0;
2602 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;
2605 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;
2608 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;
2610 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;
2617 /* Calculate the P3M B-spline moduli for one dimension */
2618 static void make_p3m_bspline_moduli_dim(real *bsp_mod, int n, int order)
2620 double zarg, zai, sinzai, infl;
2625 gmx_fatal(FARGS, "The current P3M code only supports orders up to 8");
2632 for (i = -maxk; i < 0; i++)
2636 infl = do_p3m_influence(sinzai, order);
2637 bsp_mod[n+i] = infl*infl*pow(sinzai/zai, -2.0*order);
2640 for (i = 1; i < maxk; i++)
2644 infl = do_p3m_influence(sinzai, order);
2645 bsp_mod[i] = infl*infl*pow(sinzai/zai, -2.0*order);
2649 /* Calculate the P3M B-spline moduli */
2650 static void make_p3m_bspline_moduli(splinevec bsp_mod,
2651 int nx, int ny, int nz, int order)
2653 make_p3m_bspline_moduli_dim(bsp_mod[XX], nx, order);
2654 make_p3m_bspline_moduli_dim(bsp_mod[YY], ny, order);
2655 make_p3m_bspline_moduli_dim(bsp_mod[ZZ], nz, order);
2659 static void setup_coordinate_communication(pme_atomcomm_t *atc)
2667 for (i = 1; i <= nslab/2; i++)
2669 fw = (atc->nodeid + i) % nslab;
2670 bw = (atc->nodeid - i + nslab) % nslab;
2673 atc->node_dest[n] = fw;
2674 atc->node_src[n] = bw;
2679 atc->node_dest[n] = bw;
2680 atc->node_src[n] = fw;
2686 int gmx_pme_destroy(FILE *log, gmx_pme_t *pmedata)
2692 fprintf(log, "Destroying PME data structures.\n");
2695 sfree((*pmedata)->nnx);
2696 sfree((*pmedata)->nny);
2697 sfree((*pmedata)->nnz);
2699 pmegrids_destroy(&(*pmedata)->pmegridA);
2701 sfree((*pmedata)->fftgridA);
2702 sfree((*pmedata)->cfftgridA);
2703 gmx_parallel_3dfft_destroy((*pmedata)->pfft_setupA);
2705 if ((*pmedata)->pmegridB.grid.grid != NULL)
2707 pmegrids_destroy(&(*pmedata)->pmegridB);
2708 sfree((*pmedata)->fftgridB);
2709 sfree((*pmedata)->cfftgridB);
2710 gmx_parallel_3dfft_destroy((*pmedata)->pfft_setupB);
2712 for (thread = 0; thread < (*pmedata)->nthread; thread++)
2714 free_work(&(*pmedata)->work[thread]);
2716 sfree((*pmedata)->work);
2724 static int mult_up(int n, int f)
2726 return ((n + f - 1)/f)*f;
2730 static double pme_load_imbalance(gmx_pme_t pme)
2735 nma = pme->nnodes_major;
2736 nmi = pme->nnodes_minor;
2738 n1 = mult_up(pme->nkx, nma)*mult_up(pme->nky, nmi)*pme->nkz;
2739 n2 = mult_up(pme->nkx, nma)*mult_up(pme->nkz, nmi)*pme->nky;
2740 n3 = mult_up(pme->nky, nma)*mult_up(pme->nkz, nmi)*pme->nkx;
2742 /* pme_solve is roughly double the cost of an fft */
2744 return (n1 + n2 + 3*n3)/(double)(6*pme->nkx*pme->nky*pme->nkz);
2747 static void init_atomcomm(gmx_pme_t pme, pme_atomcomm_t *atc, t_commrec *cr,
2748 int dimind, gmx_bool bSpread)
2750 int nk, k, s, thread;
2752 atc->dimind = dimind;
2757 if (pme->nnodes > 1)
2759 atc->mpi_comm = pme->mpi_comm_d[dimind];
2760 MPI_Comm_size(atc->mpi_comm, &atc->nslab);
2761 MPI_Comm_rank(atc->mpi_comm, &atc->nodeid);
2765 fprintf(debug, "For PME atom communication in dimind %d: nslab %d rank %d\n", atc->dimind, atc->nslab, atc->nodeid);
2769 atc->bSpread = bSpread;
2770 atc->pme_order = pme->pme_order;
2774 /* These three allocations are not required for particle decomp. */
2775 snew(atc->node_dest, atc->nslab);
2776 snew(atc->node_src, atc->nslab);
2777 setup_coordinate_communication(atc);
2779 snew(atc->count_thread, pme->nthread);
2780 for (thread = 0; thread < pme->nthread; thread++)
2782 snew(atc->count_thread[thread], atc->nslab);
2784 atc->count = atc->count_thread[0];
2785 snew(atc->rcount, atc->nslab);
2786 snew(atc->buf_index, atc->nslab);
2789 atc->nthread = pme->nthread;
2790 if (atc->nthread > 1)
2792 snew(atc->thread_plist, atc->nthread);
2794 snew(atc->spline, atc->nthread);
2795 for (thread = 0; thread < atc->nthread; thread++)
2797 if (atc->nthread > 1)
2799 snew(atc->thread_plist[thread].n, atc->nthread+2*GMX_CACHE_SEP);
2800 atc->thread_plist[thread].n += GMX_CACHE_SEP;
2802 snew(atc->spline[thread].thread_one, pme->nthread);
2803 atc->spline[thread].thread_one[thread] = 1;
2808 init_overlap_comm(pme_overlap_t * ol,
2818 int lbnd, rbnd, maxlr, b, i;
2821 pme_grid_comm_t *pgc;
2823 int fft_start, fft_end, send_index1, recv_index1;
2827 ol->mpi_comm = comm;
2830 ol->nnodes = nnodes;
2831 ol->nodeid = nodeid;
2833 /* Linear translation of the PME grid won't affect reciprocal space
2834 * calculations, so to optimize we only interpolate "upwards",
2835 * which also means we only have to consider overlap in one direction.
2836 * I.e., particles on this node might also be spread to grid indices
2837 * that belong to higher nodes (modulo nnodes)
2840 snew(ol->s2g0, ol->nnodes+1);
2841 snew(ol->s2g1, ol->nnodes);
2844 fprintf(debug, "PME slab boundaries:");
2846 for (i = 0; i < nnodes; i++)
2848 /* s2g0 the local interpolation grid start.
2849 * s2g1 the local interpolation grid end.
2850 * Because grid overlap communication only goes forward,
2851 * the grid the slabs for fft's should be rounded down.
2853 ol->s2g0[i] = ( i *ndata + 0 )/nnodes;
2854 ol->s2g1[i] = ((i+1)*ndata + nnodes-1)/nnodes + norder - 1;
2858 fprintf(debug, " %3d %3d", ol->s2g0[i], ol->s2g1[i]);
2861 ol->s2g0[nnodes] = ndata;
2864 fprintf(debug, "\n");
2867 /* Determine with how many nodes we need to communicate the grid overlap */
2873 for (i = 0; i < nnodes; i++)
2875 if ((i+b < nnodes && ol->s2g1[i] > ol->s2g0[i+b]) ||
2876 (i+b >= nnodes && ol->s2g1[i] > ol->s2g0[i+b-nnodes] + ndata))
2882 while (bCont && b < nnodes);
2883 ol->noverlap_nodes = b - 1;
2885 snew(ol->send_id, ol->noverlap_nodes);
2886 snew(ol->recv_id, ol->noverlap_nodes);
2887 for (b = 0; b < ol->noverlap_nodes; b++)
2889 ol->send_id[b] = (ol->nodeid + (b + 1)) % ol->nnodes;
2890 ol->recv_id[b] = (ol->nodeid - (b + 1) + ol->nnodes) % ol->nnodes;
2892 snew(ol->comm_data, ol->noverlap_nodes);
2895 for (b = 0; b < ol->noverlap_nodes; b++)
2897 pgc = &ol->comm_data[b];
2899 fft_start = ol->s2g0[ol->send_id[b]];
2900 fft_end = ol->s2g0[ol->send_id[b]+1];
2901 if (ol->send_id[b] < nodeid)
2906 send_index1 = ol->s2g1[nodeid];
2907 send_index1 = min(send_index1, fft_end);
2908 pgc->send_index0 = fft_start;
2909 pgc->send_nindex = max(0, send_index1 - pgc->send_index0);
2910 ol->send_size += pgc->send_nindex;
2912 /* We always start receiving to the first index of our slab */
2913 fft_start = ol->s2g0[ol->nodeid];
2914 fft_end = ol->s2g0[ol->nodeid+1];
2915 recv_index1 = ol->s2g1[ol->recv_id[b]];
2916 if (ol->recv_id[b] > nodeid)
2918 recv_index1 -= ndata;
2920 recv_index1 = min(recv_index1, fft_end);
2921 pgc->recv_index0 = fft_start;
2922 pgc->recv_nindex = max(0, recv_index1 - pgc->recv_index0);
2926 /* Communicate the buffer sizes to receive */
2927 for (b = 0; b < ol->noverlap_nodes; b++)
2929 MPI_Sendrecv(&ol->send_size, 1, MPI_INT, ol->send_id[b], b,
2930 &ol->comm_data[b].recv_size, 1, MPI_INT, ol->recv_id[b], b,
2931 ol->mpi_comm, &stat);
2935 /* For non-divisible grid we need pme_order iso pme_order-1 */
2936 snew(ol->sendbuf, norder*commplainsize);
2937 snew(ol->recvbuf, norder*commplainsize);
2941 make_gridindex5_to_localindex(int n, int local_start, int local_range,
2942 int **global_to_local,
2943 real **fraction_shift)
2951 for (i = 0; (i < 5*n); i++)
2953 /* Determine the global to local grid index */
2954 gtl[i] = (i - local_start + n) % n;
2955 /* For coordinates that fall within the local grid the fraction
2956 * is correct, we don't need to shift it.
2959 if (local_range < n)
2961 /* Due to rounding issues i could be 1 beyond the lower or
2962 * upper boundary of the local grid. Correct the index for this.
2963 * If we shift the index, we need to shift the fraction by
2964 * the same amount in the other direction to not affect
2966 * Note that due to this shifting the weights at the end of
2967 * the spline might change, but that will only involve values
2968 * between zero and values close to the precision of a real,
2969 * which is anyhow the accuracy of the whole mesh calculation.
2971 /* With local_range=0 we should not change i=local_start */
2972 if (i % n != local_start)
2979 else if (gtl[i] == local_range)
2981 gtl[i] = local_range - 1;
2988 *global_to_local = gtl;
2989 *fraction_shift = fsh;
2992 static pme_spline_work_t *make_pme_spline_work(int order)
2994 pme_spline_work_t *work;
3001 snew_aligned(work, 1, 16);
3003 zero_SSE = _mm_setzero_ps();
3005 /* Generate bit masks to mask out the unused grid entries,
3006 * as we only operate on order of the 8 grid entries that are
3007 * load into 2 SSE float registers.
3009 for (of = 0; of < 8-(order-1); of++)
3011 for (i = 0; i < 8; i++)
3013 tmp[i] = (i >= of && i < of+order ? 1 : 0);
3015 work->mask_SSE0[of] = _mm_loadu_ps(tmp);
3016 work->mask_SSE1[of] = _mm_loadu_ps(tmp+4);
3017 work->mask_SSE0[of] = _mm_cmpgt_ps(work->mask_SSE0[of], zero_SSE);
3018 work->mask_SSE1[of] = _mm_cmpgt_ps(work->mask_SSE1[of], zero_SSE);
3027 void gmx_pme_check_restrictions(int pme_order,
3028 int nkx, int nky, int nkz,
3031 gmx_bool bUseThreads,
3033 gmx_bool *bValidSettings)
3035 if (pme_order > PME_ORDER_MAX)
3039 *bValidSettings = FALSE;
3042 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.",
3043 pme_order, PME_ORDER_MAX);
3046 if (nkx <= pme_order*(nnodes_major > 1 ? 2 : 1) ||
3047 nky <= pme_order*(nnodes_minor > 1 ? 2 : 1) ||
3052 *bValidSettings = FALSE;
3055 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",
3059 /* Check for a limitation of the (current) sum_fftgrid_dd code.
3060 * We only allow multiple communication pulses in dim 1, not in dim 0.
3062 if (bUseThreads && (nkx < nnodes_major*pme_order &&
3063 nkx != nnodes_major*(pme_order - 1)))
3067 *bValidSettings = FALSE;
3070 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.",
3071 nkx/(double)nnodes_major, pme_order);
3074 if (bValidSettings != NULL)
3076 *bValidSettings = TRUE;
3082 int gmx_pme_init(gmx_pme_t * pmedata,
3088 gmx_bool bFreeEnergy,
3089 gmx_bool bReproducible,
3092 gmx_pme_t pme = NULL;
3094 int use_threads, sum_use_threads;
3099 fprintf(debug, "Creating PME data structures.\n");
3103 pme->redist_init = FALSE;
3104 pme->sum_qgrid_tmp = NULL;
3105 pme->sum_qgrid_dd_tmp = NULL;
3106 pme->buf_nalloc = 0;
3107 pme->redist_buf_nalloc = 0;
3110 pme->bPPnode = TRUE;
3112 pme->nnodes_major = nnodes_major;
3113 pme->nnodes_minor = nnodes_minor;
3116 if (nnodes_major*nnodes_minor > 1)
3118 pme->mpi_comm = cr->mpi_comm_mygroup;
3120 MPI_Comm_rank(pme->mpi_comm, &pme->nodeid);
3121 MPI_Comm_size(pme->mpi_comm, &pme->nnodes);
3122 if (pme->nnodes != nnodes_major*nnodes_minor)
3124 gmx_incons("PME node count mismatch");
3129 pme->mpi_comm = MPI_COMM_NULL;
3133 if (pme->nnodes == 1)
3136 pme->mpi_comm_d[0] = MPI_COMM_NULL;
3137 pme->mpi_comm_d[1] = MPI_COMM_NULL;
3139 pme->ndecompdim = 0;
3140 pme->nodeid_major = 0;
3141 pme->nodeid_minor = 0;
3143 pme->mpi_comm_d[0] = pme->mpi_comm_d[1] = MPI_COMM_NULL;
3148 if (nnodes_minor == 1)
3151 pme->mpi_comm_d[0] = pme->mpi_comm;
3152 pme->mpi_comm_d[1] = MPI_COMM_NULL;
3154 pme->ndecompdim = 1;
3155 pme->nodeid_major = pme->nodeid;
3156 pme->nodeid_minor = 0;
3159 else if (nnodes_major == 1)
3162 pme->mpi_comm_d[0] = MPI_COMM_NULL;
3163 pme->mpi_comm_d[1] = pme->mpi_comm;
3165 pme->ndecompdim = 1;
3166 pme->nodeid_major = 0;
3167 pme->nodeid_minor = pme->nodeid;
3171 if (pme->nnodes % nnodes_major != 0)
3173 gmx_incons("For 2D PME decomposition, #PME nodes must be divisible by the number of nodes in the major dimension");
3175 pme->ndecompdim = 2;
3178 MPI_Comm_split(pme->mpi_comm, pme->nodeid % nnodes_minor,
3179 pme->nodeid, &pme->mpi_comm_d[0]); /* My communicator along major dimension */
3180 MPI_Comm_split(pme->mpi_comm, pme->nodeid/nnodes_minor,
3181 pme->nodeid, &pme->mpi_comm_d[1]); /* My communicator along minor dimension */
3183 MPI_Comm_rank(pme->mpi_comm_d[0], &pme->nodeid_major);
3184 MPI_Comm_size(pme->mpi_comm_d[0], &pme->nnodes_major);
3185 MPI_Comm_rank(pme->mpi_comm_d[1], &pme->nodeid_minor);
3186 MPI_Comm_size(pme->mpi_comm_d[1], &pme->nnodes_minor);
3189 pme->bPPnode = (cr->duty & DUTY_PP);
3192 pme->nthread = nthread;
3194 /* Check if any of the PME MPI ranks uses threads */
3195 use_threads = (pme->nthread > 1 ? 1 : 0);
3197 if (pme->nnodes > 1)
3199 MPI_Allreduce(&use_threads, &sum_use_threads, 1, MPI_INT,
3200 MPI_SUM, pme->mpi_comm);
3205 sum_use_threads = use_threads;
3207 pme->bUseThreads = (sum_use_threads > 0);
3209 if (ir->ePBC == epbcSCREW)
3211 gmx_fatal(FARGS, "pme does not (yet) work with pbc = screw");
3214 pme->bFEP = ((ir->efep != efepNO) && bFreeEnergy);
3218 pme->bP3M = (ir->coulombtype == eelP3M_AD || getenv("GMX_PME_P3M") != NULL);
3219 pme->pme_order = ir->pme_order;
3220 pme->epsilon_r = ir->epsilon_r;
3222 /* If we violate restrictions, generate a fatal error here */
3223 gmx_pme_check_restrictions(pme->pme_order,
3224 pme->nkx, pme->nky, pme->nkz,
3231 if (pme->nnodes > 1)
3236 MPI_Type_contiguous(DIM, mpi_type, &(pme->rvec_mpi));
3237 MPI_Type_commit(&(pme->rvec_mpi));
3240 /* Note that the charge spreading and force gathering, which usually
3241 * takes about the same amount of time as FFT+solve_pme,
3242 * is always fully load balanced
3243 * (unless the charge distribution is inhomogeneous).
3246 imbal = pme_load_imbalance(pme);
3247 if (imbal >= 1.2 && pme->nodeid_major == 0 && pme->nodeid_minor == 0)
3251 "NOTE: The load imbalance in PME FFT and solve is %d%%.\n"
3252 " For optimal PME load balancing\n"
3253 " PME grid_x (%d) and grid_y (%d) should be divisible by #PME_nodes_x (%d)\n"
3254 " and PME grid_y (%d) and grid_z (%d) should be divisible by #PME_nodes_y (%d)\n"
3256 (int)((imbal-1)*100 + 0.5),
3257 pme->nkx, pme->nky, pme->nnodes_major,
3258 pme->nky, pme->nkz, pme->nnodes_minor);
3262 /* For non-divisible grid we need pme_order iso pme_order-1 */
3263 /* In sum_qgrid_dd x overlap is copied in place: take padding into account.
3264 * y is always copied through a buffer: we don't need padding in z,
3265 * but we do need the overlap in x because of the communication order.
3267 init_overlap_comm(&pme->overlap[0], pme->pme_order,
3271 pme->nnodes_major, pme->nodeid_major,
3273 (div_round_up(pme->nky, pme->nnodes_minor)+pme->pme_order)*(pme->nkz+pme->pme_order-1));
3275 /* Along overlap dim 1 we can send in multiple pulses in sum_fftgrid_dd.
3276 * We do this with an offset buffer of equal size, so we need to allocate
3277 * extra for the offset. That's what the (+1)*pme->nkz is for.
3279 init_overlap_comm(&pme->overlap[1], pme->pme_order,
3283 pme->nnodes_minor, pme->nodeid_minor,
3285 (div_round_up(pme->nkx, pme->nnodes_major)+pme->pme_order+1)*pme->nkz);
3287 /* Double-check for a limitation of the (current) sum_fftgrid_dd code.
3288 * Note that gmx_pme_check_restrictions checked for this already.
3290 if (pme->bUseThreads && pme->overlap[0].noverlap_nodes > 1)
3292 gmx_incons("More than one communication pulse required for grid overlap communication along the major dimension while using threads");
3295 snew(pme->bsp_mod[XX], pme->nkx);
3296 snew(pme->bsp_mod[YY], pme->nky);
3297 snew(pme->bsp_mod[ZZ], pme->nkz);
3299 /* The required size of the interpolation grid, including overlap.
3300 * The allocated size (pmegrid_n?) might be slightly larger.
3302 pme->pmegrid_nx = pme->overlap[0].s2g1[pme->nodeid_major] -
3303 pme->overlap[0].s2g0[pme->nodeid_major];
3304 pme->pmegrid_ny = pme->overlap[1].s2g1[pme->nodeid_minor] -
3305 pme->overlap[1].s2g0[pme->nodeid_minor];
3306 pme->pmegrid_nz_base = pme->nkz;
3307 pme->pmegrid_nz = pme->pmegrid_nz_base + pme->pme_order - 1;
3308 set_grid_alignment(&pme->pmegrid_nz, pme->pme_order);
3310 pme->pmegrid_start_ix = pme->overlap[0].s2g0[pme->nodeid_major];
3311 pme->pmegrid_start_iy = pme->overlap[1].s2g0[pme->nodeid_minor];
3312 pme->pmegrid_start_iz = 0;
3314 make_gridindex5_to_localindex(pme->nkx,
3315 pme->pmegrid_start_ix,
3316 pme->pmegrid_nx - (pme->pme_order-1),
3317 &pme->nnx, &pme->fshx);
3318 make_gridindex5_to_localindex(pme->nky,
3319 pme->pmegrid_start_iy,
3320 pme->pmegrid_ny - (pme->pme_order-1),
3321 &pme->nny, &pme->fshy);
3322 make_gridindex5_to_localindex(pme->nkz,
3323 pme->pmegrid_start_iz,
3324 pme->pmegrid_nz_base,
3325 &pme->nnz, &pme->fshz);
3327 pmegrids_init(&pme->pmegridA,
3328 pme->pmegrid_nx, pme->pmegrid_ny, pme->pmegrid_nz,
3329 pme->pmegrid_nz_base,
3333 pme->overlap[0].s2g1[pme->nodeid_major]-pme->overlap[0].s2g0[pme->nodeid_major+1],
3334 pme->overlap[1].s2g1[pme->nodeid_minor]-pme->overlap[1].s2g0[pme->nodeid_minor+1]);
3336 pme->spline_work = make_pme_spline_work(pme->pme_order);
3338 ndata[0] = pme->nkx;
3339 ndata[1] = pme->nky;
3340 ndata[2] = pme->nkz;
3342 /* This routine will allocate the grid data to fit the FFTs */
3343 gmx_parallel_3dfft_init(&pme->pfft_setupA, ndata,
3344 &pme->fftgridA, &pme->cfftgridA,
3346 pme->overlap[0].s2g0, pme->overlap[1].s2g0,
3347 bReproducible, pme->nthread);
3351 pmegrids_init(&pme->pmegridB,
3352 pme->pmegrid_nx, pme->pmegrid_ny, pme->pmegrid_nz,
3353 pme->pmegrid_nz_base,
3357 pme->nkx % pme->nnodes_major != 0,
3358 pme->nky % pme->nnodes_minor != 0);
3360 gmx_parallel_3dfft_init(&pme->pfft_setupB, ndata,
3361 &pme->fftgridB, &pme->cfftgridB,
3363 pme->overlap[0].s2g0, pme->overlap[1].s2g0,
3364 bReproducible, pme->nthread);
3368 pme->pmegridB.grid.grid = NULL;
3369 pme->fftgridB = NULL;
3370 pme->cfftgridB = NULL;
3375 /* Use plain SPME B-spline interpolation */
3376 make_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
3380 /* Use the P3M grid-optimized influence function */
3381 make_p3m_bspline_moduli(pme->bsp_mod, pme->nkx, pme->nky, pme->nkz, pme->pme_order);
3384 /* Use atc[0] for spreading */
3385 init_atomcomm(pme, &pme->atc[0], cr, nnodes_major > 1 ? 0 : 1, TRUE);
3386 if (pme->ndecompdim >= 2)
3388 init_atomcomm(pme, &pme->atc[1], cr, 1, FALSE);
3391 if (pme->nnodes == 1)
3393 pme->atc[0].n = homenr;
3394 pme_realloc_atomcomm_things(&pme->atc[0]);
3400 /* Use fft5d, order after FFT is y major, z, x minor */
3402 snew(pme->work, pme->nthread);
3403 for (thread = 0; thread < pme->nthread; thread++)
3405 realloc_work(&pme->work[thread], pme->nkx);
3414 static void reuse_pmegrids(const pmegrids_t *old, pmegrids_t *new)
3418 for (d = 0; d < DIM; d++)
3420 if (new->grid.n[d] > old->grid.n[d])
3426 sfree_aligned(new->grid.grid);
3427 new->grid.grid = old->grid.grid;
3429 if (new->grid_th != NULL && new->nthread == old->nthread)
3431 sfree_aligned(new->grid_all);
3432 for (t = 0; t < new->nthread; t++)
3434 new->grid_th[t].grid = old->grid_th[t].grid;
3439 int gmx_pme_reinit(gmx_pme_t * pmedata,
3442 const t_inputrec * ir,
3450 irc.nkx = grid_size[XX];
3451 irc.nky = grid_size[YY];
3452 irc.nkz = grid_size[ZZ];
3454 if (pme_src->nnodes == 1)
3456 homenr = pme_src->atc[0].n;
3463 ret = gmx_pme_init(pmedata, cr, pme_src->nnodes_major, pme_src->nnodes_minor,
3464 &irc, homenr, pme_src->bFEP, FALSE, pme_src->nthread);
3468 /* We can easily reuse the allocated pme grids in pme_src */
3469 reuse_pmegrids(&pme_src->pmegridA, &(*pmedata)->pmegridA);
3470 /* We would like to reuse the fft grids, but that's harder */
3477 static void copy_local_grid(gmx_pme_t pme,
3478 pmegrids_t *pmegrids, int thread, real *fftgrid)
3480 ivec local_fft_ndata, local_fft_offset, local_fft_size;
3484 int offx, offy, offz, x, y, z, i0, i0t;
3489 gmx_parallel_3dfft_real_limits(pme->pfft_setupA,
3493 fft_my = local_fft_size[YY];
3494 fft_mz = local_fft_size[ZZ];
3496 pmegrid = &pmegrids->grid_th[thread];
3498 nsx = pmegrid->s[XX];
3499 nsy = pmegrid->s[YY];
3500 nsz = pmegrid->s[ZZ];
3502 for (d = 0; d < DIM; d++)
3504 nf[d] = min(pmegrid->n[d] - (pmegrid->order - 1),
3505 local_fft_ndata[d] - pmegrid->offset[d]);
3508 offx = pmegrid->offset[XX];
3509 offy = pmegrid->offset[YY];
3510 offz = pmegrid->offset[ZZ];
3512 /* Directly copy the non-overlapping parts of the local grids.
3513 * This also initializes the full grid.
3515 grid_th = pmegrid->grid;
3516 for (x = 0; x < nf[XX]; x++)
3518 for (y = 0; y < nf[YY]; y++)
3520 i0 = ((offx + x)*fft_my + (offy + y))*fft_mz + offz;
3521 i0t = (x*nsy + y)*nsz;
3522 for (z = 0; z < nf[ZZ]; z++)
3524 fftgrid[i0+z] = grid_th[i0t+z];
3531 reduce_threadgrid_overlap(gmx_pme_t pme,
3532 const pmegrids_t *pmegrids, int thread,
3533 real *fftgrid, real *commbuf_x, real *commbuf_y)
3535 ivec local_fft_ndata, local_fft_offset, local_fft_size;
3536 int fft_nx, fft_ny, fft_nz;
3541 int offx, offy, offz, x, y, z, i0, i0t;
3542 int sx, sy, sz, fx, fy, fz, tx1, ty1, tz1, ox, oy, oz;
3543 gmx_bool bClearBufX, bClearBufY, bClearBufXY, bClearBuf;
3544 gmx_bool bCommX, bCommY;
3547 const pmegrid_t *pmegrid, *pmegrid_g, *pmegrid_f;
3548 const real *grid_th;
3549 real *commbuf = NULL;
3551 gmx_parallel_3dfft_real_limits(pme->pfft_setupA,
3555 fft_nx = local_fft_ndata[XX];
3556 fft_ny = local_fft_ndata[YY];
3557 fft_nz = local_fft_ndata[ZZ];
3559 fft_my = local_fft_size[YY];
3560 fft_mz = local_fft_size[ZZ];
3562 /* This routine is called when all thread have finished spreading.
3563 * Here each thread sums grid contributions calculated by other threads
3564 * to the thread local grid volume.
3565 * To minimize the number of grid copying operations,
3566 * this routines sums immediately from the pmegrid to the fftgrid.
3569 /* Determine which part of the full node grid we should operate on,
3570 * this is our thread local part of the full grid.
3572 pmegrid = &pmegrids->grid_th[thread];
3574 for (d = 0; d < DIM; d++)
3576 ne[d] = min(pmegrid->offset[d]+pmegrid->n[d]-(pmegrid->order-1),
3577 local_fft_ndata[d]);
3580 offx = pmegrid->offset[XX];
3581 offy = pmegrid->offset[YY];
3582 offz = pmegrid->offset[ZZ];
3589 /* Now loop over all the thread data blocks that contribute
3590 * to the grid region we (our thread) are operating on.
3592 /* Note that ffy_nx/y is equal to the number of grid points
3593 * between the first point of our node grid and the one of the next node.
3595 for (sx = 0; sx >= -pmegrids->nthread_comm[XX]; sx--)
3597 fx = pmegrid->ci[XX] + sx;
3602 fx += pmegrids->nc[XX];
3604 bCommX = (pme->nnodes_major > 1);
3606 pmegrid_g = &pmegrids->grid_th[fx*pmegrids->nc[YY]*pmegrids->nc[ZZ]];
3607 ox += pmegrid_g->offset[XX];
3610 tx1 = min(ox + pmegrid_g->n[XX], ne[XX]);
3614 tx1 = min(ox + pmegrid_g->n[XX], pme->pme_order);
3617 for (sy = 0; sy >= -pmegrids->nthread_comm[YY]; sy--)
3619 fy = pmegrid->ci[YY] + sy;
3624 fy += pmegrids->nc[YY];
3626 bCommY = (pme->nnodes_minor > 1);
3628 pmegrid_g = &pmegrids->grid_th[fy*pmegrids->nc[ZZ]];
3629 oy += pmegrid_g->offset[YY];
3632 ty1 = min(oy + pmegrid_g->n[YY], ne[YY]);
3636 ty1 = min(oy + pmegrid_g->n[YY], pme->pme_order);
3639 for (sz = 0; sz >= -pmegrids->nthread_comm[ZZ]; sz--)
3641 fz = pmegrid->ci[ZZ] + sz;
3645 fz += pmegrids->nc[ZZ];
3648 pmegrid_g = &pmegrids->grid_th[fz];
3649 oz += pmegrid_g->offset[ZZ];
3650 tz1 = min(oz + pmegrid_g->n[ZZ], ne[ZZ]);
3652 if (sx == 0 && sy == 0 && sz == 0)
3654 /* We have already added our local contribution
3655 * before calling this routine, so skip it here.
3660 thread_f = (fx*pmegrids->nc[YY] + fy)*pmegrids->nc[ZZ] + fz;
3662 pmegrid_f = &pmegrids->grid_th[thread_f];
3664 grid_th = pmegrid_f->grid;
3666 nsx = pmegrid_f->s[XX];
3667 nsy = pmegrid_f->s[YY];
3668 nsz = pmegrid_f->s[ZZ];
3670 #ifdef DEBUG_PME_REDUCE
3671 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",
3672 pme->nodeid, thread, thread_f,
3673 pme->pmegrid_start_ix,
3674 pme->pmegrid_start_iy,
3675 pme->pmegrid_start_iz,
3677 offx-ox, tx1-ox, offx, tx1,
3678 offy-oy, ty1-oy, offy, ty1,
3679 offz-oz, tz1-oz, offz, tz1);
3682 if (!(bCommX || bCommY))
3684 /* Copy from the thread local grid to the node grid */
3685 for (x = offx; x < tx1; x++)
3687 for (y = offy; y < ty1; y++)
3689 i0 = (x*fft_my + y)*fft_mz;
3690 i0t = ((x - ox)*nsy + (y - oy))*nsz - oz;
3691 for (z = offz; z < tz1; z++)
3693 fftgrid[i0+z] += grid_th[i0t+z];
3700 /* The order of this conditional decides
3701 * where the corner volume gets stored with x+y decomp.
3705 commbuf = commbuf_y;
3706 buf_my = ty1 - offy;
3709 /* We index commbuf modulo the local grid size */
3710 commbuf += buf_my*fft_nx*fft_nz;
3712 bClearBuf = bClearBufXY;
3713 bClearBufXY = FALSE;
3717 bClearBuf = bClearBufY;
3723 commbuf = commbuf_x;
3725 bClearBuf = bClearBufX;
3729 /* Copy to the communication buffer */
3730 for (x = offx; x < tx1; x++)
3732 for (y = offy; y < ty1; y++)
3734 i0 = (x*buf_my + y)*fft_nz;
3735 i0t = ((x - ox)*nsy + (y - oy))*nsz - oz;
3739 /* First access of commbuf, initialize it */
3740 for (z = offz; z < tz1; z++)
3742 commbuf[i0+z] = grid_th[i0t+z];
3747 for (z = offz; z < tz1; z++)
3749 commbuf[i0+z] += grid_th[i0t+z];
3761 static void sum_fftgrid_dd(gmx_pme_t pme, real *fftgrid)
3763 ivec local_fft_ndata, local_fft_offset, local_fft_size;
3764 pme_overlap_t *overlap;
3765 int send_index0, send_nindex;
3770 int send_size_y, recv_size_y;
3771 int ipulse, send_id, recv_id, datasize, gridsize, size_yx;
3772 real *sendptr, *recvptr;
3773 int x, y, z, indg, indb;
3775 /* Note that this routine is only used for forward communication.
3776 * Since the force gathering, unlike the charge spreading,
3777 * can be trivially parallelized over the particles,
3778 * the backwards process is much simpler and can use the "old"
3779 * communication setup.
3782 gmx_parallel_3dfft_real_limits(pme->pfft_setupA,
3787 if (pme->nnodes_minor > 1)
3789 /* Major dimension */
3790 overlap = &pme->overlap[1];
3792 if (pme->nnodes_major > 1)
3794 size_yx = pme->overlap[0].comm_data[0].send_nindex;
3800 datasize = (local_fft_ndata[XX] + size_yx)*local_fft_ndata[ZZ];
3802 send_size_y = overlap->send_size;
3804 for (ipulse = 0; ipulse < overlap->noverlap_nodes; ipulse++)
3806 send_id = overlap->send_id[ipulse];
3807 recv_id = overlap->recv_id[ipulse];
3809 overlap->comm_data[ipulse].send_index0 -
3810 overlap->comm_data[0].send_index0;
3811 send_nindex = overlap->comm_data[ipulse].send_nindex;
3812 /* We don't use recv_index0, as we always receive starting at 0 */
3813 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
3814 recv_size_y = overlap->comm_data[ipulse].recv_size;
3816 sendptr = overlap->sendbuf + send_index0*local_fft_ndata[ZZ];
3817 recvptr = overlap->recvbuf;
3820 MPI_Sendrecv(sendptr, send_size_y*datasize, GMX_MPI_REAL,
3822 recvptr, recv_size_y*datasize, GMX_MPI_REAL,
3824 overlap->mpi_comm, &stat);
3827 for (x = 0; x < local_fft_ndata[XX]; x++)
3829 for (y = 0; y < recv_nindex; y++)
3831 indg = (x*local_fft_size[YY] + y)*local_fft_size[ZZ];
3832 indb = (x*recv_size_y + y)*local_fft_ndata[ZZ];
3833 for (z = 0; z < local_fft_ndata[ZZ]; z++)
3835 fftgrid[indg+z] += recvptr[indb+z];
3840 if (pme->nnodes_major > 1)
3842 /* Copy from the received buffer to the send buffer for dim 0 */
3843 sendptr = pme->overlap[0].sendbuf;
3844 for (x = 0; x < size_yx; x++)
3846 for (y = 0; y < recv_nindex; y++)
3848 indg = (x*local_fft_ndata[YY] + y)*local_fft_ndata[ZZ];
3849 indb = ((local_fft_ndata[XX] + x)*recv_size_y + y)*local_fft_ndata[ZZ];
3850 for (z = 0; z < local_fft_ndata[ZZ]; z++)
3852 sendptr[indg+z] += recvptr[indb+z];
3860 /* We only support a single pulse here.
3861 * This is not a severe limitation, as this code is only used
3862 * with OpenMP and with OpenMP the (PME) domains can be larger.
3864 if (pme->nnodes_major > 1)
3866 /* Major dimension */
3867 overlap = &pme->overlap[0];
3869 datasize = local_fft_ndata[YY]*local_fft_ndata[ZZ];
3870 gridsize = local_fft_size[YY] *local_fft_size[ZZ];
3874 send_id = overlap->send_id[ipulse];
3875 recv_id = overlap->recv_id[ipulse];
3876 send_nindex = overlap->comm_data[ipulse].send_nindex;
3877 /* We don't use recv_index0, as we always receive starting at 0 */
3878 recv_nindex = overlap->comm_data[ipulse].recv_nindex;
3880 sendptr = overlap->sendbuf;
3881 recvptr = overlap->recvbuf;
3885 fprintf(debug, "PME fftgrid comm %2d x %2d x %2d\n",
3886 send_nindex, local_fft_ndata[YY], local_fft_ndata[ZZ]);
3890 MPI_Sendrecv(sendptr, send_nindex*datasize, GMX_MPI_REAL,
3892 recvptr, recv_nindex*datasize, GMX_MPI_REAL,
3894 overlap->mpi_comm, &stat);
3897 for (x = 0; x < recv_nindex; x++)
3899 for (y = 0; y < local_fft_ndata[YY]; y++)
3901 indg = (x*local_fft_size[YY] + y)*local_fft_size[ZZ];
3902 indb = (x*local_fft_ndata[YY] + y)*local_fft_ndata[ZZ];
3903 for (z = 0; z < local_fft_ndata[ZZ]; z++)
3905 fftgrid[indg+z] += recvptr[indb+z];
3913 static void spread_on_grid(gmx_pme_t pme,
3914 pme_atomcomm_t *atc, pmegrids_t *grids,
3915 gmx_bool bCalcSplines, gmx_bool bSpread,
3918 int nthread, thread;
3919 #ifdef PME_TIME_THREADS
3920 gmx_cycles_t c1, c2, c3, ct1a, ct1b, ct1c;
3921 static double cs1 = 0, cs2 = 0, cs3 = 0;
3922 static double cs1a[6] = {0, 0, 0, 0, 0, 0};
3926 nthread = pme->nthread;
3927 assert(nthread > 0);
3929 #ifdef PME_TIME_THREADS
3930 c1 = omp_cyc_start();
3934 #pragma omp parallel for num_threads(nthread) schedule(static)
3935 for (thread = 0; thread < nthread; thread++)
3939 start = atc->n* thread /nthread;
3940 end = atc->n*(thread+1)/nthread;
3942 /* Compute fftgrid index for all atoms,
3943 * with help of some extra variables.
3945 calc_interpolation_idx(pme, atc, start, end, thread);
3948 #ifdef PME_TIME_THREADS
3949 c1 = omp_cyc_end(c1);
3953 #ifdef PME_TIME_THREADS
3954 c2 = omp_cyc_start();
3956 #pragma omp parallel for num_threads(nthread) schedule(static)
3957 for (thread = 0; thread < nthread; thread++)
3959 splinedata_t *spline;
3960 pmegrid_t *grid = NULL;
3962 /* make local bsplines */
3963 if (grids == NULL || !pme->bUseThreads)
3965 spline = &atc->spline[0];
3971 grid = &grids->grid;
3976 spline = &atc->spline[thread];
3978 if (grids->nthread == 1)
3980 /* One thread, we operate on all charges */
3985 /* Get the indices our thread should operate on */
3986 make_thread_local_ind(atc, thread, spline);
3989 grid = &grids->grid_th[thread];
3994 make_bsplines(spline->theta, spline->dtheta, pme->pme_order,
3995 atc->fractx, spline->n, spline->ind, atc->q, pme->bFEP);
4000 /* put local atoms on grid. */
4001 #ifdef PME_TIME_SPREAD
4002 ct1a = omp_cyc_start();
4004 spread_q_bsplines_thread(grid, atc, spline, pme->spline_work);
4006 if (pme->bUseThreads)
4008 copy_local_grid(pme, grids, thread, fftgrid);
4010 #ifdef PME_TIME_SPREAD
4011 ct1a = omp_cyc_end(ct1a);
4012 cs1a[thread] += (double)ct1a;
4016 #ifdef PME_TIME_THREADS
4017 c2 = omp_cyc_end(c2);
4021 if (bSpread && pme->bUseThreads)
4023 #ifdef PME_TIME_THREADS
4024 c3 = omp_cyc_start();
4026 #pragma omp parallel for num_threads(grids->nthread) schedule(static)
4027 for (thread = 0; thread < grids->nthread; thread++)
4029 reduce_threadgrid_overlap(pme, grids, thread,
4031 pme->overlap[0].sendbuf,
4032 pme->overlap[1].sendbuf);
4034 #ifdef PME_TIME_THREADS
4035 c3 = omp_cyc_end(c3);
4039 if (pme->nnodes > 1)
4041 /* Communicate the overlapping part of the fftgrid.
4042 * For this communication call we need to check pme->bUseThreads
4043 * to have all ranks communicate here, regardless of pme->nthread.
4045 sum_fftgrid_dd(pme, fftgrid);
4049 #ifdef PME_TIME_THREADS
4053 printf("idx %.2f spread %.2f red %.2f",
4054 cs1*1e-9, cs2*1e-9, cs3*1e-9);
4055 #ifdef PME_TIME_SPREAD
4056 for (thread = 0; thread < nthread; thread++)
4058 printf(" %.2f", cs1a[thread]*1e-9);
4067 static void dump_grid(FILE *fp,
4068 int sx, int sy, int sz, int nx, int ny, int nz,
4069 int my, int mz, const real *g)
4073 for (x = 0; x < nx; x++)
4075 for (y = 0; y < ny; y++)
4077 for (z = 0; z < nz; z++)
4079 fprintf(fp, "%2d %2d %2d %6.3f\n",
4080 sx+x, sy+y, sz+z, g[(x*my + y)*mz + z]);
4086 static void dump_local_fftgrid(gmx_pme_t pme, const real *fftgrid)
4088 ivec local_fft_ndata, local_fft_offset, local_fft_size;
4090 gmx_parallel_3dfft_real_limits(pme->pfft_setupA,
4096 pme->pmegrid_start_ix,
4097 pme->pmegrid_start_iy,
4098 pme->pmegrid_start_iz,
4099 pme->pmegrid_nx-pme->pme_order+1,
4100 pme->pmegrid_ny-pme->pme_order+1,
4101 pme->pmegrid_nz-pme->pme_order+1,
4108 void gmx_pme_calc_energy(gmx_pme_t pme, int n, rvec *x, real *q, real *V)
4110 pme_atomcomm_t *atc;
4113 if (pme->nnodes > 1)
4115 gmx_incons("gmx_pme_calc_energy called in parallel");
4119 gmx_incons("gmx_pme_calc_energy with free energy");
4122 atc = &pme->atc_energy;
4124 if (atc->spline == NULL)
4126 snew(atc->spline, atc->nthread);
4129 atc->bSpread = TRUE;
4130 atc->pme_order = pme->pme_order;
4132 pme_realloc_atomcomm_things(atc);
4136 /* We only use the A-charges grid */
4137 grid = &pme->pmegridA;
4139 /* Only calculate the spline coefficients, don't actually spread */
4140 spread_on_grid(pme, atc, NULL, TRUE, FALSE, pme->fftgridA);
4142 *V = gather_energy_bsplines(pme, grid->grid.grid, atc);
4146 static void reset_pmeonly_counters(t_commrec *cr, gmx_wallcycle_t wcycle,
4147 t_nrnb *nrnb, t_inputrec *ir,
4148 gmx_large_int_t step)
4150 /* Reset all the counters related to performance over the run */
4151 wallcycle_stop(wcycle, ewcRUN);
4152 wallcycle_reset_all(wcycle);
4154 if (ir->nsteps >= 0)
4156 /* ir->nsteps is not used here, but we update it for consistency */
4157 ir->nsteps -= step - ir->init_step;
4159 ir->init_step = step;
4160 wallcycle_start(wcycle, ewcRUN);
4164 static void gmx_pmeonly_switch(int *npmedata, gmx_pme_t **pmedata,
4166 t_commrec *cr, t_inputrec *ir,
4170 gmx_pme_t pme = NULL;
4173 while (ind < *npmedata)
4175 pme = (*pmedata)[ind];
4176 if (pme->nkx == grid_size[XX] &&
4177 pme->nky == grid_size[YY] &&
4178 pme->nkz == grid_size[ZZ])
4189 srenew(*pmedata, *npmedata);
4191 /* Generate a new PME data structure, copying part of the old pointers */
4192 gmx_pme_reinit(&((*pmedata)[ind]), cr, pme, ir, grid_size);
4194 *pme_ret = (*pmedata)[ind];
4198 int gmx_pmeonly(gmx_pme_t pme,
4199 t_commrec *cr, t_nrnb *nrnb,
4200 gmx_wallcycle_t wcycle,
4201 real ewaldcoeff, gmx_bool bGatherOnly,
4206 gmx_pme_pp_t pme_pp;
4210 rvec *x_pp = NULL, *f_pp = NULL;
4211 real *chargeA = NULL, *chargeB = NULL;
4213 int maxshift_x = 0, maxshift_y = 0;
4214 real energy, dvdlambda;
4219 gmx_large_int_t step, step_rel;
4222 /* This data will only use with PME tuning, i.e. switching PME grids */
4224 snew(pmedata, npmedata);
4227 pme_pp = gmx_pme_pp_init(cr);
4232 do /****** this is a quasi-loop over time steps! */
4234 /* The reason for having a loop here is PME grid tuning/switching */
4237 /* Domain decomposition */
4238 ret = gmx_pme_recv_q_x(pme_pp,
4240 &chargeA, &chargeB, box, &x_pp, &f_pp,
4241 &maxshift_x, &maxshift_y,
4242 &pme->bFEP, &lambda,
4245 grid_switch, &ewaldcoeff);
4247 if (ret == pmerecvqxSWITCHGRID)
4249 /* Switch the PME grid to grid_switch */
4250 gmx_pmeonly_switch(&npmedata, &pmedata, grid_switch, cr, ir, &pme);
4253 if (ret == pmerecvqxRESETCOUNTERS)
4255 /* Reset the cycle and flop counters */
4256 reset_pmeonly_counters(cr, wcycle, nrnb, ir, step);
4259 while (ret == pmerecvqxSWITCHGRID || ret == pmerecvqxRESETCOUNTERS);
4261 if (ret == pmerecvqxFINISH)
4263 /* We should stop: break out of the loop */
4267 step_rel = step - ir->init_step;
4271 wallcycle_start(wcycle, ewcRUN);
4274 wallcycle_start(wcycle, ewcPMEMESH);
4278 gmx_pme_do(pme, 0, natoms, x_pp, f_pp, chargeA, chargeB, box,
4279 cr, maxshift_x, maxshift_y, nrnb, wcycle, vir, ewaldcoeff,
4280 &energy, lambda, &dvdlambda,
4281 GMX_PME_DO_ALL_F | (bEnerVir ? GMX_PME_CALC_ENER_VIR : 0));
4283 cycles = wallcycle_stop(wcycle, ewcPMEMESH);
4285 gmx_pme_send_force_vir_ener(pme_pp,
4286 f_pp, vir, energy, dvdlambda,
4290 } /***** end of quasi-loop, we stop with the break above */
4296 int gmx_pme_do(gmx_pme_t pme,
4297 int start, int homenr,
4299 real *chargeA, real *chargeB,
4300 matrix box, t_commrec *cr,
4301 int maxshift_x, int maxshift_y,
4302 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
4303 matrix vir, real ewaldcoeff,
4304 real *energy, real lambda,
4305 real *dvdlambda, int flags)
4307 int q, d, i, j, ntot, npme;
4310 pme_atomcomm_t *atc = NULL;
4311 pmegrids_t *pmegrid = NULL;
4315 real *charge = NULL, *q_d;
4319 gmx_parallel_3dfft_t pfft_setup;
4321 t_complex * cfftgrid;
4323 const gmx_bool bCalcEnerVir = flags & GMX_PME_CALC_ENER_VIR;
4324 const gmx_bool bCalcF = flags & GMX_PME_CALC_F;
4326 assert(pme->nnodes > 0);
4327 assert(pme->nnodes == 1 || pme->ndecompdim > 0);
4329 if (pme->nnodes > 1)
4333 if (atc->npd > atc->pd_nalloc)
4335 atc->pd_nalloc = over_alloc_dd(atc->npd);
4336 srenew(atc->pd, atc->pd_nalloc);
4338 atc->maxshift = (atc->dimind == 0 ? maxshift_x : maxshift_y);
4342 /* This could be necessary for TPI */
4343 pme->atc[0].n = homenr;
4346 for (q = 0; q < (pme->bFEP ? 2 : 1); q++)
4350 pmegrid = &pme->pmegridA;
4351 fftgrid = pme->fftgridA;
4352 cfftgrid = pme->cfftgridA;
4353 pfft_setup = pme->pfft_setupA;
4354 charge = chargeA+start;
4358 pmegrid = &pme->pmegridB;
4359 fftgrid = pme->fftgridB;
4360 cfftgrid = pme->cfftgridB;
4361 pfft_setup = pme->pfft_setupB;
4362 charge = chargeB+start;
4364 grid = pmegrid->grid.grid;
4365 /* Unpack structure */
4368 fprintf(debug, "PME: nnodes = %d, nodeid = %d\n",
4369 cr->nnodes, cr->nodeid);
4370 fprintf(debug, "Grid = %p\n", (void*)grid);
4373 gmx_fatal(FARGS, "No grid!");
4378 m_inv_ur0(box, pme->recipbox);
4380 if (pme->nnodes == 1)
4383 if (DOMAINDECOMP(cr))
4386 pme_realloc_atomcomm_things(atc);
4394 wallcycle_start(wcycle, ewcPME_REDISTXF);
4395 for (d = pme->ndecompdim-1; d >= 0; d--)
4397 if (d == pme->ndecompdim-1)
4405 n_d = pme->atc[d+1].n;
4411 if (atc->npd > atc->pd_nalloc)
4413 atc->pd_nalloc = over_alloc_dd(atc->npd);
4414 srenew(atc->pd, atc->pd_nalloc);
4416 atc->maxshift = (atc->dimind == 0 ? maxshift_x : maxshift_y);
4417 pme_calc_pidx_wrapper(n_d, pme->recipbox, x_d, atc);
4420 /* Redistribute x (only once) and qA or qB */
4421 if (DOMAINDECOMP(cr))
4423 dd_pmeredist_x_q(pme, n_d, q == 0, x_d, q_d, atc);
4427 pmeredist_pd(pme, TRUE, n_d, q == 0, x_d, q_d, atc);
4432 wallcycle_stop(wcycle, ewcPME_REDISTXF);
4437 fprintf(debug, "Node= %6d, pme local particles=%6d\n",
4438 cr->nodeid, atc->n);
4441 if (flags & GMX_PME_SPREAD_Q)
4443 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
4445 /* Spread the charges on a grid */
4446 spread_on_grid(pme, &pme->atc[0], pmegrid, q == 0, TRUE, fftgrid);
4450 inc_nrnb(nrnb, eNR_WEIGHTS, DIM*atc->n);
4452 inc_nrnb(nrnb, eNR_SPREADQBSP,
4453 pme->pme_order*pme->pme_order*pme->pme_order*atc->n);
4455 if (!pme->bUseThreads)
4457 wrap_periodic_pmegrid(pme, grid);
4459 /* sum contributions to local grid from other nodes */
4461 if (pme->nnodes > 1)
4463 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_QGRID_FORWARD);
4468 copy_pmegrid_to_fftgrid(pme, grid, fftgrid);
4471 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
4474 dump_local_fftgrid(pme,fftgrid);
4479 /* Here we start a large thread parallel region */
4480 #pragma omp parallel num_threads(pme->nthread) private(thread)
4482 thread = gmx_omp_get_thread_num();
4483 if (flags & GMX_PME_SOLVE)
4490 wallcycle_start(wcycle, ewcPME_FFT);
4492 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_REAL_TO_COMPLEX,
4493 fftgrid, cfftgrid, thread, wcycle);
4496 wallcycle_stop(wcycle, ewcPME_FFT);
4500 /* solve in k-space for our local cells */
4503 wallcycle_start(wcycle, ewcPME_SOLVE);
4506 solve_pme_yzx(pme, cfftgrid, ewaldcoeff,
4507 box[XX][XX]*box[YY][YY]*box[ZZ][ZZ],
4509 pme->nthread, thread);
4512 wallcycle_stop(wcycle, ewcPME_SOLVE);
4514 inc_nrnb(nrnb, eNR_SOLVEPME, loop_count);
4524 wallcycle_start(wcycle, ewcPME_FFT);
4526 gmx_parallel_3dfft_execute(pfft_setup, GMX_FFT_COMPLEX_TO_REAL,
4527 cfftgrid, fftgrid, thread, wcycle);
4530 wallcycle_stop(wcycle, ewcPME_FFT);
4534 if (pme->nodeid == 0)
4536 ntot = pme->nkx*pme->nky*pme->nkz;
4537 npme = ntot*log((real)ntot)/log(2.0);
4538 inc_nrnb(nrnb, eNR_FFT, 2*npme);
4541 wallcycle_start(wcycle, ewcPME_SPREADGATHER);
4544 copy_fftgrid_to_pmegrid(pme, fftgrid, grid, pme->nthread, thread);
4547 /* End of thread parallel section.
4548 * With MPI we have to synchronize here before gmx_sum_qgrid_dd.
4553 /* distribute local grid to all nodes */
4555 if (pme->nnodes > 1)
4557 gmx_sum_qgrid_dd(pme, grid, GMX_SUM_QGRID_BACKWARD);
4562 unwrap_periodic_pmegrid(pme, grid);
4564 /* interpolate forces for our local atoms */
4568 /* If we are running without parallelization,
4569 * atc->f is the actual force array, not a buffer,
4570 * therefore we should not clear it.
4572 bClearF = (q == 0 && PAR(cr));
4573 #pragma omp parallel for num_threads(pme->nthread) schedule(static)
4574 for (thread = 0; thread < pme->nthread; thread++)
4576 gather_f_bsplines(pme, grid, bClearF, atc,
4577 &atc->spline[thread],
4578 pme->bFEP ? (q == 0 ? 1.0-lambda : lambda) : 1.0);
4583 inc_nrnb(nrnb, eNR_GATHERFBSP,
4584 pme->pme_order*pme->pme_order*pme->pme_order*pme->atc[0].n);
4585 wallcycle_stop(wcycle, ewcPME_SPREADGATHER);
4590 /* This should only be called on the master thread
4591 * and after the threads have synchronized.
4593 get_pme_ener_vir(pme, pme->nthread, &energy_AB[q], vir_AB[q]);
4597 if (bCalcF && pme->nnodes > 1)
4599 wallcycle_start(wcycle, ewcPME_REDISTXF);
4600 for (d = 0; d < pme->ndecompdim; d++)
4603 if (d == pme->ndecompdim - 1)
4610 n_d = pme->atc[d+1].n;
4611 f_d = pme->atc[d+1].f;
4613 if (DOMAINDECOMP(cr))
4615 dd_pmeredist_f(pme, atc, n_d, f_d,
4616 d == pme->ndecompdim-1 && pme->bPPnode);
4620 pmeredist_pd(pme, FALSE, n_d, TRUE, f_d, NULL, atc);
4624 wallcycle_stop(wcycle, ewcPME_REDISTXF);
4632 *energy = energy_AB[0];
4633 m_add(vir, vir_AB[0], vir);
4637 *energy = (1.0-lambda)*energy_AB[0] + lambda*energy_AB[1];
4638 *dvdlambda += energy_AB[1] - energy_AB[0];
4639 for (i = 0; i < DIM; i++)
4641 for (j = 0; j < DIM; j++)
4643 vir[i][j] += (1.0-lambda)*vir_AB[0][i][j] +
4644 lambda*vir_AB[1][i][j];
4656 fprintf(debug, "PME mesh energy: %g\n", *energy);