1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * Gromacs Copyright (c) 1991-2005
5 * David van der Spoel, Erik Lindahl, University of Groningen.
7 * This program is free software; you can redistribute it and/or
8 * modify it under the terms of the GNU General Public License
9 * as published by the Free Software Foundation; either version 2
10 * of the License, or (at your option) any later version.
12 * To help us fund GROMACS development, we humbly ask that you cite
13 * the research papers on the package. Check out http://www.gromacs.org
16 * Gnomes, ROck Monsters And Chili Sauce
34 #include "gmx_parallel_3dfft.h"
36 #include "gmxcomplex.h"
37 #include "gmx_fatal.h"
41 struct gmx_parallel_3dfft {
46 static int *copy_int_array(int n,int *src)
50 dest = (int*)malloc(n*sizeof(int));
57 static int *make_slab2grid(int nnodes,int ngrid)
61 s2g = (int*)malloc((nnodes+1)*sizeof(int));
62 for(i=0; i<nnodes+1; i++) {
63 /* We always round up */
64 s2g[i] = (i*ngrid + nnodes - 1)/nnodes;
71 gmx_parallel_3dfft_init (gmx_parallel_3dfft_t * pfft_setup,
74 t_complex ** complex_data,
76 int * slab2index_major,
77 int * slab2index_minor,
78 gmx_bool bReproducible,
81 int rN=ndata[2],M=ndata[1],K=ndata[0];
82 int flags = FFT5D_REALCOMPLEX | FFT5D_ORDER_YZ; /* FFT5D_DEBUG */
83 MPI_Comm rcomm[]={comm[1],comm[0]};
84 int Nb,Mb,Kb; /* dimension for backtransform (in starting order) */
85 t_complex *buf1, *buf2; /*intermediate buffers - used internally.*/
88 if (bReproducible) flags |= FFT5D_NOMEASURE;
90 if (!(flags&FFT5D_ORDER_YZ)) {
93 Nb=K;Mb=rN;Kb=M; /* currently always true because ORDER_YZ always set */
96 (*pfft_setup)->p1 = fft5d_plan_3d(rN,M,K,rcomm, flags, (t_complex**)real_data, complex_data, &buf1, &buf2, nthreads);
98 (*pfft_setup)->p2 = fft5d_plan_3d(Nb,Mb,Kb,rcomm,
99 (flags|FFT5D_BACKWARD|FFT5D_NOMALLOC)^FFT5D_ORDER_YZ, complex_data, (t_complex**)real_data, &buf1, &buf2, nthreads);
101 return (*pfft_setup)->p1 != 0 && (*pfft_setup)->p2 !=0;
106 fft5d_limits(fft5d_plan p,
111 int N1,M0,K0,K1,*coor;
112 fft5d_local_size(p,&N1,&M0,&K0,&K1,&coor); /* M0=MG/P[0], K1=KG/P[1], NG,MG,KG global sizes */
115 local_offset[1]=p->oM[0]; /*=p->coor[0]*p->MG/p->P[0]; */
116 local_offset[0]=p->oK[0]; /*=p->coor[1]*p->KG/p->P[1]; */
118 local_ndata[2]=p->rC[0];
119 local_ndata[1]=p->pM[0];
120 local_ndata[0]=p->pK[0];
122 if ((!(p->flags&FFT5D_BACKWARD)) && (p->flags&FFT5D_REALCOMPLEX)) {
123 local_size[2]=p->C[0]*2;
125 local_size[2]=p->C[0];
127 local_size[1]=p->pM[0];
128 local_size[0]=p->pK[0];
133 gmx_parallel_3dfft_real_limits(gmx_parallel_3dfft_t pfft_setup,
137 return fft5d_limits(pfft_setup->p1,local_ndata,local_offset,local_size);
140 static void reorder_ivec_yzx(ivec v)
151 gmx_parallel_3dfft_complex_limits(gmx_parallel_3dfft_t pfft_setup,
159 /* For now everything is in-order, but prepare to save communication by avoiding transposes */
160 complex_order[0] = 0;
161 complex_order[1] = 1;
162 complex_order[2] = 2;
164 ret = fft5d_limits(pfft_setup->p2,local_ndata,local_offset,local_size);
166 reorder_ivec_yzx(local_ndata);
167 reorder_ivec_yzx(local_offset);
168 reorder_ivec_yzx(local_size);
175 gmx_parallel_3dfft_execute(gmx_parallel_3dfft_t pfft_setup,
176 enum gmx_fft_direction dir,
180 gmx_wallcycle_t wcycle) {
181 if ((!(pfft_setup->p1->flags&FFT5D_REALCOMPLEX)) ^ (dir==GMX_FFT_FORWARD ||dir==GMX_FFT_BACKWARD))
183 gmx_fatal(FARGS,"Invalid transform. Plan and execution don't match regarding reel/complex");
185 if (dir==GMX_FFT_FORWARD || dir==GMX_FFT_REAL_TO_COMPLEX) {
186 fft5d_execute(pfft_setup->p1,thread,wcycle);
188 fft5d_execute(pfft_setup->p2,thread,wcycle);
194 gmx_parallel_3dfft_destroy(gmx_parallel_3dfft_t pfft_setup) {
195 fft5d_destroy(pfft_setup->p2);
196 fft5d_destroy(pfft_setup->p1);