xref: /petsc/src/mat/impls/fft/fftw/fftw.c (revision 0c9b18e4ba9674abf7df50cf4600152bbc2a9f57)
1b2b77a04SHong Zhang 
2b2b77a04SHong Zhang /*
3b2b77a04SHong Zhang     Provides an interface to the FFTW package.
4b2b77a04SHong Zhang     Testing examples can be found in ~src/mat/examples/tests
5b2b77a04SHong Zhang */
6b2b77a04SHong Zhang 
7c6db04a5SJed Brown #include <../src/mat/impls/fft/fft.h>   /*I "petscmat.h" I*/
8b2b77a04SHong Zhang EXTERN_C_BEGIN
9c6db04a5SJed Brown #include <fftw3-mpi.h>
10b2b77a04SHong Zhang EXTERN_C_END
11b2b77a04SHong Zhang 
12b2b77a04SHong Zhang typedef struct {
13b9ae062cSHong Zhang   ptrdiff_t   ndim_fftw,*dim_fftw;
14e5d7f247SAmlan Barua   PetscInt   partial_dim;
15b2b77a04SHong Zhang   fftw_plan   p_forward,p_backward;
16b2b77a04SHong Zhang   unsigned    p_flag; /* planner flags, FFTW_ESTIMATE,FFTW_MEASURE, FFTW_PATIENT, FFTW_EXHAUSTIVE */
17b2b77a04SHong Zhang   PetscScalar *finarray,*foutarray,*binarray,*boutarray; /* keep track of arrays becaue fftw plan should be
18b2b77a04SHong Zhang                                                             executed for the arrays with which the plan was created */
19b2b77a04SHong Zhang } Mat_FFTW;
20b2b77a04SHong Zhang 
21b2b77a04SHong Zhang extern PetscErrorCode MatMult_SeqFFTW(Mat,Vec,Vec);
22b2b77a04SHong Zhang extern PetscErrorCode MatMultTranspose_SeqFFTW(Mat,Vec,Vec);
23b2b77a04SHong Zhang extern PetscErrorCode MatMult_MPIFFTW(Mat,Vec,Vec);
24b2b77a04SHong Zhang extern PetscErrorCode MatMultTranspose_MPIFFTW(Mat,Vec,Vec);
25b2b77a04SHong Zhang extern PetscErrorCode MatDestroy_FFTW(Mat);
26b2b77a04SHong Zhang extern PetscErrorCode VecDestroy_MPIFFTW(Vec);
27b2b77a04SHong Zhang extern PetscErrorCode MatGetVecs_FFTW(Mat,Vec*,Vec*);
284f7415efSAmlan Barua extern PetscErrorCode MatGetVecsFFT_FFTW(Mat,Vec*,Vec*,Vec*);
29b2b77a04SHong Zhang 
30b2b77a04SHong Zhang #undef __FUNCT__
31b2b77a04SHong Zhang #define __FUNCT__ "MatMult_SeqFFTW"
32b2b77a04SHong Zhang PetscErrorCode MatMult_SeqFFTW(Mat A,Vec x,Vec y)
33b2b77a04SHong Zhang {
34b2b77a04SHong Zhang   PetscErrorCode ierr;
35b2b77a04SHong Zhang   Mat_FFT        *fft  = (Mat_FFT*)A->data;
36b2b77a04SHong Zhang   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
37b2b77a04SHong Zhang   PetscScalar    *x_array,*y_array;
38b2b77a04SHong Zhang   PetscInt       ndim=fft->ndim,*dim=fft->dim;
39b2b77a04SHong Zhang 
40b2b77a04SHong Zhang   PetscFunctionBegin;
41b2b77a04SHong Zhang   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
42b2b77a04SHong Zhang   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
43b2b77a04SHong Zhang   if (!fftw->p_forward){ /* create a plan, then excute it */
44b2b77a04SHong Zhang     switch (ndim){
45b2b77a04SHong Zhang     case 1:
4658a912c5SAmlan Barua #if defined(PETSC_USE_COMPLEX)
47b2b77a04SHong Zhang       fftw->p_forward = fftw_plan_dft_1d(dim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_FORWARD,fftw->p_flag);
4858a912c5SAmlan Barua #else
4958a912c5SAmlan Barua       fftw->p_forward = fftw_plan_dft_r2c_1d(dim[0],(double *)x_array,(fftw_complex*)y_array,fftw->p_flag);
5058a912c5SAmlan Barua #endif
51b2b77a04SHong Zhang       break;
52b2b77a04SHong Zhang     case 2:
5358a912c5SAmlan Barua #if defined(PETSC_USE_COMPLEX)
54b2b77a04SHong Zhang       fftw->p_forward = fftw_plan_dft_2d(dim[0],dim[1],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_FORWARD,fftw->p_flag);
5558a912c5SAmlan Barua #else
5658a912c5SAmlan Barua       fftw->p_forward = fftw_plan_dft_r2c_2d(dim[0],dim[1],(double *)x_array,(fftw_complex*)y_array,fftw->p_flag);
5758a912c5SAmlan Barua #endif
58b2b77a04SHong Zhang       break;
59b2b77a04SHong Zhang     case 3:
6058a912c5SAmlan Barua #if defined(PETSC_USE_COMPLEX)
61b2b77a04SHong Zhang       fftw->p_forward = fftw_plan_dft_3d(dim[0],dim[1],dim[2],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_FORWARD,fftw->p_flag);
6258a912c5SAmlan Barua #else
6358a912c5SAmlan Barua       fftw->p_forward = fftw_plan_dft_r2c_3d(dim[0],dim[1],dim[2],(double *)x_array,(fftw_complex*)y_array,fftw->p_flag);
6458a912c5SAmlan Barua #endif
65b2b77a04SHong Zhang       break;
66b2b77a04SHong Zhang     default:
6758a912c5SAmlan Barua #if defined(PETSC_USE_COMPLEX)
68b2b77a04SHong Zhang       fftw->p_forward = fftw_plan_dft(ndim,dim,(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_FORWARD,fftw->p_flag);
6958a912c5SAmlan Barua #else
7058a912c5SAmlan Barua       fftw->p_forward = fftw_plan_dft_r2c(ndim,dim,(double *)x_array,(fftw_complex*)y_array,fftw->p_flag);
7158a912c5SAmlan Barua #endif
72b2b77a04SHong Zhang       break;
73b2b77a04SHong Zhang     }
74b2b77a04SHong Zhang     fftw->finarray  = x_array;
75b2b77a04SHong Zhang     fftw->foutarray = y_array;
76b2b77a04SHong Zhang     /* Warning: if (fftw->p_flag!==FFTW_ESTIMATE) The data in the in/out arrays is overwritten!
77b2b77a04SHong Zhang                 planning should be done before x is initialized! See FFTW manual sec2.1 or sec4 */
78b2b77a04SHong Zhang     fftw_execute(fftw->p_forward);
79b2b77a04SHong Zhang   } else { /* use existing plan */
80b2b77a04SHong Zhang     if (fftw->finarray != x_array || fftw->foutarray != y_array){ /* use existing plan on new arrays */
81b2b77a04SHong Zhang       fftw_execute_dft(fftw->p_forward,(fftw_complex*)x_array,(fftw_complex*)y_array);
82b2b77a04SHong Zhang     } else {
83b2b77a04SHong Zhang       fftw_execute(fftw->p_forward);
84b2b77a04SHong Zhang     }
85b2b77a04SHong Zhang   }
86b2b77a04SHong Zhang   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
87b2b77a04SHong Zhang   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
88b2b77a04SHong Zhang   PetscFunctionReturn(0);
89b2b77a04SHong Zhang }
90b2b77a04SHong Zhang 
91b2b77a04SHong Zhang #undef __FUNCT__
92b2b77a04SHong Zhang #define __FUNCT__ "MatMultTranspose_SeqFFTW"
93b2b77a04SHong Zhang PetscErrorCode MatMultTranspose_SeqFFTW(Mat A,Vec x,Vec y)
94b2b77a04SHong Zhang {
95b2b77a04SHong Zhang   PetscErrorCode ierr;
96b2b77a04SHong Zhang   Mat_FFT        *fft = (Mat_FFT*)A->data;
97b2b77a04SHong Zhang   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
98b2b77a04SHong Zhang   PetscScalar    *x_array,*y_array;
99b2b77a04SHong Zhang   PetscInt       ndim=fft->ndim,*dim=fft->dim;
100b2b77a04SHong Zhang 
101b2b77a04SHong Zhang   PetscFunctionBegin;
102b2b77a04SHong Zhang   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
103b2b77a04SHong Zhang   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
104b2b77a04SHong Zhang   if (!fftw->p_backward){ /* create a plan, then excute it */
105b2b77a04SHong Zhang     switch (ndim){
106b2b77a04SHong Zhang     case 1:
10758a912c5SAmlan Barua #if defined(PETSC_USE_COMPLEX)
108b2b77a04SHong Zhang       fftw->p_backward = fftw_plan_dft_1d(dim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_BACKWARD,fftw->p_flag);
10958a912c5SAmlan Barua #else
110e81bb599SAmlan Barua       fftw->p_backward= fftw_plan_dft_c2r_1d(dim[0],(fftw_complex*)x_array,(double *)y_array,fftw->p_flag);
11158a912c5SAmlan Barua #endif
112b2b77a04SHong Zhang       break;
113b2b77a04SHong Zhang     case 2:
11458a912c5SAmlan Barua #if defined(PETSC_USE_COMPLEX)
115b2b77a04SHong Zhang       fftw->p_backward = fftw_plan_dft_2d(dim[0],dim[1],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_BACKWARD,fftw->p_flag);
11658a912c5SAmlan Barua #else
117e81bb599SAmlan Barua       fftw->p_backward= fftw_plan_dft_c2r_2d(dim[0],dim[1],(fftw_complex*)x_array,(double *)y_array,fftw->p_flag);
11858a912c5SAmlan Barua #endif
119b2b77a04SHong Zhang       break;
120b2b77a04SHong Zhang     case 3:
12158a912c5SAmlan Barua #if defined(PETSC_USE_COMPLEX)
122b2b77a04SHong Zhang       fftw->p_backward = fftw_plan_dft_3d(dim[0],dim[1],dim[2],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_BACKWARD,fftw->p_flag);
12358a912c5SAmlan Barua #else
124e81bb599SAmlan Barua       fftw->p_backward= fftw_plan_dft_c2r_3d(dim[0],dim[1],dim[2],(fftw_complex*)x_array,(double *)y_array,fftw->p_flag);
12558a912c5SAmlan Barua #endif
126b2b77a04SHong Zhang       break;
127b2b77a04SHong Zhang     default:
12858a912c5SAmlan Barua #if defined(PETSC_USE_COMPLEX)
129b2b77a04SHong Zhang       fftw->p_backward = fftw_plan_dft(ndim,dim,(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_BACKWARD,fftw->p_flag);
13058a912c5SAmlan Barua #else
131e81bb599SAmlan Barua       fftw->p_backward= fftw_plan_dft_c2r(ndim,dim,(fftw_complex*)x_array,(double *)y_array,fftw->p_flag);
13258a912c5SAmlan Barua #endif
133b2b77a04SHong Zhang       break;
134b2b77a04SHong Zhang     }
135b2b77a04SHong Zhang     fftw->binarray  = x_array;
136b2b77a04SHong Zhang     fftw->boutarray = y_array;
137b2b77a04SHong Zhang     fftw_execute(fftw->p_backward);CHKERRQ(ierr);
138b2b77a04SHong Zhang   } else { /* use existing plan */
139b2b77a04SHong Zhang     if (fftw->binarray != x_array || fftw->boutarray != y_array){ /* use existing plan on new arrays */
140b2b77a04SHong Zhang       fftw_execute_dft(fftw->p_backward,(fftw_complex*)x_array,(fftw_complex*)y_array);
141b2b77a04SHong Zhang     } else {
142b2b77a04SHong Zhang       fftw_execute(fftw->p_backward);CHKERRQ(ierr);
143b2b77a04SHong Zhang     }
144b2b77a04SHong Zhang   }
145b2b77a04SHong Zhang   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
146b2b77a04SHong Zhang   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
147b2b77a04SHong Zhang   PetscFunctionReturn(0);
148b2b77a04SHong Zhang }
149b2b77a04SHong Zhang 
150b2b77a04SHong Zhang #undef __FUNCT__
151b2b77a04SHong Zhang #define __FUNCT__ "MatMult_MPIFFTW"
152b2b77a04SHong Zhang PetscErrorCode MatMult_MPIFFTW(Mat A,Vec x,Vec y)
153b2b77a04SHong Zhang {
154b2b77a04SHong Zhang   PetscErrorCode ierr;
155b2b77a04SHong Zhang   Mat_FFT        *fft  = (Mat_FFT*)A->data;
156b2b77a04SHong Zhang   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
157b2b77a04SHong Zhang   PetscScalar    *x_array,*y_array;
158c92418ccSAmlan Barua   PetscInt       ndim=fft->ndim,*dim=fft->dim;
159b2b77a04SHong Zhang   MPI_Comm       comm=((PetscObject)A)->comm;
160b2b77a04SHong Zhang 
161b2b77a04SHong Zhang   PetscFunctionBegin;
162b2b77a04SHong Zhang   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
163b2b77a04SHong Zhang   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
164b2b77a04SHong Zhang   if (!fftw->p_forward){ /* create a plan, then excute it */
165b2b77a04SHong Zhang     switch (ndim){
166b2b77a04SHong Zhang     case 1:
1673c3a9c75SAmlan Barua #if defined(PETSC_USE_COMPLEX)
168b2b77a04SHong Zhang       fftw->p_forward = fftw_mpi_plan_dft_1d(dim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_FORWARD,fftw->p_flag);
169ae0a50aaSHong Zhang #else
170ae0a50aaSHong Zhang       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers yet");
1713c3a9c75SAmlan Barua #endif
172b2b77a04SHong Zhang       break;
173b2b77a04SHong Zhang     case 2:
1743c3a9c75SAmlan Barua #if defined(PETSC_USE_COMPLEX)
175b2b77a04SHong Zhang       fftw->p_forward = fftw_mpi_plan_dft_2d(dim[0],dim[1],(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_FORWARD,fftw->p_flag);
1763c3a9c75SAmlan Barua #else
1773c3a9c75SAmlan Barua       fftw->p_forward = fftw_mpi_plan_dft_r2c_2d(dim[0],dim[1],(double *)x_array,(fftw_complex*)y_array,comm,FFTW_ESTIMATE);
1783c3a9c75SAmlan Barua #endif
179b2b77a04SHong Zhang       break;
180b2b77a04SHong Zhang     case 3:
1813c3a9c75SAmlan Barua #if defined(PETSC_USE_COMPLEX)
182b2b77a04SHong Zhang       fftw->p_forward = fftw_mpi_plan_dft_3d(dim[0],dim[1],dim[2],(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_FORWARD,fftw->p_flag);
1833c3a9c75SAmlan Barua #else
18451d1eed7SAmlan Barua       fftw->p_forward = fftw_mpi_plan_dft_r2c_3d(dim[0],dim[1],dim[2],(double *)x_array,(fftw_complex*)y_array,comm,FFTW_ESTIMATE);
1853c3a9c75SAmlan Barua #endif
186b2b77a04SHong Zhang       break;
187b2b77a04SHong Zhang     default:
1883c3a9c75SAmlan Barua #if defined(PETSC_USE_COMPLEX)
189c92418ccSAmlan Barua       fftw->p_forward = fftw_mpi_plan_dft(fftw->ndim_fftw,fftw->dim_fftw,(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_FORWARD,fftw->p_flag);
1903c3a9c75SAmlan Barua #else
1913c3a9c75SAmlan Barua       fftw->p_forward = fftw_mpi_plan_dft_r2c(fftw->ndim_fftw,fftw->dim_fftw,(double *)x_array,(fftw_complex*)y_array,comm,FFTW_ESTIMATE);
1923c3a9c75SAmlan Barua #endif
193b2b77a04SHong Zhang       break;
194b2b77a04SHong Zhang     }
195b2b77a04SHong Zhang     fftw->finarray  = x_array;
196b2b77a04SHong Zhang     fftw->foutarray = y_array;
197b2b77a04SHong Zhang     /* Warning: if (fftw->p_flag!==FFTW_ESTIMATE) The data in the in/out arrays is overwritten!
198b2b77a04SHong Zhang                 planning should be done before x is initialized! See FFTW manual sec2.1 or sec4 */
199b2b77a04SHong Zhang     fftw_execute(fftw->p_forward);
200b2b77a04SHong Zhang   } else { /* use existing plan */
201b2b77a04SHong Zhang     if (fftw->finarray != x_array || fftw->foutarray != y_array){ /* use existing plan on new arrays */
202b2b77a04SHong Zhang       fftw_execute_dft(fftw->p_forward,(fftw_complex*)x_array,(fftw_complex*)y_array);
203b2b77a04SHong Zhang     } else {
204b2b77a04SHong Zhang       fftw_execute(fftw->p_forward);
205b2b77a04SHong Zhang     }
206b2b77a04SHong Zhang   }
207b2b77a04SHong Zhang   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
208b2b77a04SHong Zhang   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
209b2b77a04SHong Zhang   PetscFunctionReturn(0);
210b2b77a04SHong Zhang }
211b2b77a04SHong Zhang 
212b2b77a04SHong Zhang #undef __FUNCT__
213b2b77a04SHong Zhang #define __FUNCT__ "MatMultTranspose_MPIFFTW"
214b2b77a04SHong Zhang PetscErrorCode MatMultTranspose_MPIFFTW(Mat A,Vec x,Vec y)
215b2b77a04SHong Zhang {
216b2b77a04SHong Zhang   PetscErrorCode ierr;
217b2b77a04SHong Zhang   Mat_FFT        *fft  = (Mat_FFT*)A->data;
218b2b77a04SHong Zhang   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
219b2b77a04SHong Zhang   PetscScalar    *x_array,*y_array;
220c92418ccSAmlan Barua   PetscInt       ndim=fft->ndim,*dim=fft->dim;
221b2b77a04SHong Zhang   MPI_Comm       comm=((PetscObject)A)->comm;
222b2b77a04SHong Zhang 
223b2b77a04SHong Zhang   PetscFunctionBegin;
224b2b77a04SHong Zhang   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
225b2b77a04SHong Zhang   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
226b2b77a04SHong Zhang   if (!fftw->p_backward){ /* create a plan, then excute it */
227b2b77a04SHong Zhang     switch (ndim){
228b2b77a04SHong Zhang     case 1:
2293c3a9c75SAmlan Barua #if defined(PETSC_USE_COMPLEX)
230b2b77a04SHong Zhang       fftw->p_backward = fftw_mpi_plan_dft_1d(dim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_BACKWARD,fftw->p_flag);
231ae0a50aaSHong Zhang #else
232ae0a50aaSHong Zhang       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers yet");
2333c3a9c75SAmlan Barua #endif
234b2b77a04SHong Zhang       break;
235b2b77a04SHong Zhang     case 2:
2363c3a9c75SAmlan Barua #if defined(PETSC_USE_COMPLEX)
237b2b77a04SHong Zhang       fftw->p_backward = fftw_mpi_plan_dft_2d(dim[0],dim[1],(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_BACKWARD,fftw->p_flag);
2383c3a9c75SAmlan Barua #else
2393c3a9c75SAmlan Barua       fftw->p_backward = fftw_mpi_plan_dft_c2r_2d(dim[0],dim[1],(fftw_complex*)x_array,(double *)y_array,comm,FFTW_ESTIMATE);
2403c3a9c75SAmlan Barua #endif
241b2b77a04SHong Zhang       break;
242b2b77a04SHong Zhang     case 3:
2433c3a9c75SAmlan Barua #if defined(PETSC_USE_COMPLEX)
244b2b77a04SHong Zhang       fftw->p_backward = fftw_mpi_plan_dft_3d(dim[0],dim[1],dim[2],(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_BACKWARD,fftw->p_flag);
2453c3a9c75SAmlan Barua #else
2463c3a9c75SAmlan Barua       fftw->p_backward = fftw_mpi_plan_dft_c2r_3d(dim[0],dim[1],dim[2],(fftw_complex*)x_array,(double *)y_array,comm,FFTW_ESTIMATE);
2473c3a9c75SAmlan Barua #endif
248b2b77a04SHong Zhang       break;
249b2b77a04SHong Zhang     default:
2503c3a9c75SAmlan Barua #if defined(PETSC_USE_COMPLEX)
251c92418ccSAmlan Barua       fftw->p_backward = fftw_mpi_plan_dft(fftw->ndim_fftw,fftw->dim_fftw,(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_BACKWARD,fftw->p_flag);
2523c3a9c75SAmlan Barua #else
2533c3a9c75SAmlan Barua       fftw->p_backward = fftw_mpi_plan_dft_c2r(fftw->ndim_fftw,fftw->dim_fftw,(fftw_complex*)x_array,(double *)y_array,comm,FFTW_ESTIMATE);
2543c3a9c75SAmlan Barua #endif
255b2b77a04SHong Zhang       break;
256b2b77a04SHong Zhang     }
257b2b77a04SHong Zhang     fftw->binarray  = x_array;
258b2b77a04SHong Zhang     fftw->boutarray = y_array;
259b2b77a04SHong Zhang     fftw_execute(fftw->p_backward);CHKERRQ(ierr);
260b2b77a04SHong Zhang   } else { /* use existing plan */
261b2b77a04SHong Zhang     if (fftw->binarray != x_array || fftw->boutarray != y_array){ /* use existing plan on new arrays */
262b2b77a04SHong Zhang       fftw_execute_dft(fftw->p_backward,(fftw_complex*)x_array,(fftw_complex*)y_array);
263b2b77a04SHong Zhang     } else {
264b2b77a04SHong Zhang       fftw_execute(fftw->p_backward);CHKERRQ(ierr);
265b2b77a04SHong Zhang     }
266b2b77a04SHong Zhang   }
267b2b77a04SHong Zhang   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
268b2b77a04SHong Zhang   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
269b2b77a04SHong Zhang   PetscFunctionReturn(0);
270b2b77a04SHong Zhang }
271b2b77a04SHong Zhang 
272b2b77a04SHong Zhang #undef __FUNCT__
273b2b77a04SHong Zhang #define __FUNCT__ "MatDestroy_FFTW"
274b2b77a04SHong Zhang PetscErrorCode MatDestroy_FFTW(Mat A)
275b2b77a04SHong Zhang {
276b2b77a04SHong Zhang   Mat_FFT        *fft = (Mat_FFT*)A->data;
277b2b77a04SHong Zhang   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
278b2b77a04SHong Zhang   PetscErrorCode ierr;
279b2b77a04SHong Zhang 
280b2b77a04SHong Zhang   PetscFunctionBegin;
281b2b77a04SHong Zhang   fftw_destroy_plan(fftw->p_forward);
282b2b77a04SHong Zhang   fftw_destroy_plan(fftw->p_backward);
283b1301b2eSHong Zhang   ierr = PetscFree(fftw->dim_fftw);CHKERRQ(ierr);
284bf0cc555SLisandro Dalcin   ierr = PetscFree(fft->data);CHKERRQ(ierr);
285b2b77a04SHong Zhang   PetscFunctionReturn(0);
286b2b77a04SHong Zhang }
287b2b77a04SHong Zhang 
288c6db04a5SJed Brown #include <../src/vec/vec/impls/mpi/pvecimpl.h>   /*I  "petscvec.h"   I*/
289b2b77a04SHong Zhang #undef __FUNCT__
290b2b77a04SHong Zhang #define __FUNCT__ "VecDestroy_MPIFFTW"
291b2b77a04SHong Zhang PetscErrorCode VecDestroy_MPIFFTW(Vec v)
292b2b77a04SHong Zhang {
293b2b77a04SHong Zhang   PetscErrorCode  ierr;
294b2b77a04SHong Zhang   PetscScalar     *array;
295b2b77a04SHong Zhang 
296b2b77a04SHong Zhang   PetscFunctionBegin;
297b2b77a04SHong Zhang   ierr = VecGetArray(v,&array);CHKERRQ(ierr);
298b2b77a04SHong Zhang   fftw_free((fftw_complex*)array);CHKERRQ(ierr);
299b2b77a04SHong Zhang   ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
300b2b77a04SHong Zhang   ierr = VecDestroy_MPI(v);CHKERRQ(ierr);
301b2b77a04SHong Zhang   PetscFunctionReturn(0);
302b2b77a04SHong Zhang }
303b2b77a04SHong Zhang 
304b2b77a04SHong Zhang #undef __FUNCT__
3053c3a9c75SAmlan Barua #define __FUNCT__ "MatGetVecs1DC_FFTW"
306c92418ccSAmlan Barua /*
307c92418ccSAmlan Barua    MatGetVecs_FFTW1D - Get Vectors(s) compatible with matrix, i.e. with the
308c92418ccSAmlan Barua      parallel layout determined by FFTW-1D
309c92418ccSAmlan Barua 
310c92418ccSAmlan Barua */
3116a506ed8SAmlan Barua PetscErrorCode  MatGetVecs_FFTW1D(Mat A,Vec *fin,Vec *fout,Vec *bout)
312c92418ccSAmlan Barua {
313c92418ccSAmlan Barua   PetscErrorCode ierr;
314c92418ccSAmlan Barua   PetscMPIInt    size,rank;
315c92418ccSAmlan Barua   MPI_Comm       comm=((PetscObject)A)->comm;
316c92418ccSAmlan Barua   Mat_FFT        *fft = (Mat_FFT*)A->data;
317c92418ccSAmlan Barua //  Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
318c92418ccSAmlan Barua   PetscInt       N=fft->N;
319c92418ccSAmlan Barua   PetscInt       ndim=fft->ndim,*dim=fft->dim;
320c92418ccSAmlan Barua   ptrdiff_t      f_alloc_local,f_local_n0,f_local_0_start;
321c92418ccSAmlan Barua   ptrdiff_t      f_local_n1,f_local_1_end;
322c92418ccSAmlan Barua   ptrdiff_t      b_alloc_local,b_local_n0,b_local_0_start;
323c92418ccSAmlan Barua   ptrdiff_t      b_local_n1,b_local_1_end;
324ae0a50aaSHong Zhang   fftw_complex   *data_fin,*data_fout,*data_bout;
325c92418ccSAmlan Barua 
326c92418ccSAmlan Barua   PetscFunctionBegin;
327c92418ccSAmlan Barua #if !defined(PETSC_USE_COMPLEX)
328c92418ccSAmlan Barua   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers");
329c92418ccSAmlan Barua #endif
330c92418ccSAmlan Barua   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
331c92418ccSAmlan Barua   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
332c92418ccSAmlan Barua   if (size == 1){
333c92418ccSAmlan Barua     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Works only for parallel 1D");
334c92418ccSAmlan Barua   }
335c92418ccSAmlan Barua   else {
336c92418ccSAmlan Barua       if (ndim>1){
337c92418ccSAmlan Barua         SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Works only for parallel 1D");}
338c92418ccSAmlan Barua       else {
339c92418ccSAmlan Barua           f_alloc_local = fftw_mpi_local_size_1d(dim[0],comm,FFTW_FORWARD,FFTW_ESTIMATE,&f_local_n0,&f_local_0_start,&f_local_n1,&f_local_1_end);
340c92418ccSAmlan Barua           b_alloc_local = fftw_mpi_local_size_1d(dim[0],comm,FFTW_BACKWARD,FFTW_ESTIMATE,&b_local_n0,&b_local_0_start,&b_local_n1,&b_local_1_end);
341c92418ccSAmlan Barua           if (fin) {
342c92418ccSAmlan Barua             data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*f_alloc_local);
343c92418ccSAmlan Barua             ierr = VecCreateMPIWithArray(comm,f_local_n0,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
344c92418ccSAmlan Barua             (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
345c92418ccSAmlan Barua           }
346c92418ccSAmlan Barua           if (fout) {
347c92418ccSAmlan Barua             data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*f_alloc_local);
348c92418ccSAmlan Barua             ierr = VecCreateMPIWithArray(comm,f_local_n1,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
349c92418ccSAmlan Barua             (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
350c92418ccSAmlan Barua           }
351c92418ccSAmlan Barua           if (bout) {
352c92418ccSAmlan Barua             data_bout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*b_alloc_local);
353c92418ccSAmlan Barua             ierr = VecCreateMPIWithArray(comm,b_local_n1,N,(const PetscScalar*)data_bout,bout);CHKERRQ(ierr);
354c92418ccSAmlan Barua             (*bout)->ops->destroy   = VecDestroy_MPIFFTW;
355c92418ccSAmlan Barua           }
356c92418ccSAmlan Barua   }
357c92418ccSAmlan Barua   if (fin){
358c92418ccSAmlan Barua     ierr = PetscLayoutReference(A->cmap,&(*fin)->map);CHKERRQ(ierr);
359c92418ccSAmlan Barua   }
360c92418ccSAmlan Barua   if (fout){
361c92418ccSAmlan Barua     ierr = PetscLayoutReference(A->rmap,&(*fout)->map);CHKERRQ(ierr);
362c92418ccSAmlan Barua   }
363c92418ccSAmlan Barua   if (bout){
364c92418ccSAmlan Barua     ierr = PetscLayoutReference(A->rmap,&(*bout)->map);CHKERRQ(ierr);
365c92418ccSAmlan Barua   }
366c92418ccSAmlan Barua   PetscFunctionReturn(0);
367c92418ccSAmlan Barua }
368c92418ccSAmlan Barua 
369c92418ccSAmlan Barua 
370c92418ccSAmlan Barua }
3713c3a9c75SAmlan Barua 
3724f7415efSAmlan Barua 
3734f7415efSAmlan Barua #undef __FUNCT__
3744f7415efSAmlan Barua #define __FUNCT__ "MatGetVecsFFT"
3754f7415efSAmlan Barua PetscErrorCode MatGetVecsFFT(Mat A,Vec *x,Vec *y,Vec *z)
3764f7415efSAmlan Barua {
3774f7415efSAmlan Barua   PetscErrorCode ierr;
3784f7415efSAmlan Barua   PetscFunctionBegin;
3794f7415efSAmlan Barua   ierr = PetscTryMethod(A,"MatGetVecsFFT_C",(Mat,Vec*,Vec*,Vec*),(A,x,y,z));CHKERRQ(ierr);
3804f7415efSAmlan Barua   PetscFunctionReturn(0);
3814f7415efSAmlan Barua }
3824f7415efSAmlan Barua 
3834f7415efSAmlan Barua #undef __FUNCT__
3844f7415efSAmlan Barua #define __FUNCT__ "MatGetVecsFFT_FFTW"
3854f7415efSAmlan Barua /*
3864f7415efSAmlan Barua    MatGetVecs_FFTW - Get vector(s) compatible with the matrix, i.e. with the
3874f7415efSAmlan Barua      parallel layout determined by FFTW
3884f7415efSAmlan Barua 
3894f7415efSAmlan Barua    Collective on Mat
3904f7415efSAmlan Barua 
3914f7415efSAmlan Barua    Input Parameter:
3924f7415efSAmlan Barua .  mat - the matrix
3934f7415efSAmlan Barua 
3944f7415efSAmlan Barua    Output Parameter:
3954f7415efSAmlan Barua +   fin - (optional) input vector of forward FFTW
3964f7415efSAmlan Barua -   fout - (optional) output vector of forward FFTW
3974f7415efSAmlan Barua 
3984f7415efSAmlan Barua   Level: advanced
3994f7415efSAmlan Barua 
4004f7415efSAmlan Barua .seealso: MatCreateFFTW()
4014f7415efSAmlan Barua */
4024f7415efSAmlan Barua EXTERN_C_BEGIN
4034f7415efSAmlan Barua PetscErrorCode  MatGetVecsFFT_FFTW(Mat A,Vec *fin,Vec *fout,Vec *bout)
4044f7415efSAmlan Barua {
4054f7415efSAmlan Barua   PetscErrorCode ierr;
4064f7415efSAmlan Barua   PetscMPIInt    size,rank;
4074f7415efSAmlan Barua   MPI_Comm       comm=((PetscObject)A)->comm;
4084f7415efSAmlan Barua   Mat_FFT        *fft = (Mat_FFT*)A->data;
4094f7415efSAmlan Barua   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
4104f7415efSAmlan Barua   PetscInt       N=fft->N, N1, n1,vsize;
4114f7415efSAmlan Barua   PetscInt       ndim=fft->ndim,*dim=fft->dim,n=fft->n;
4124f7415efSAmlan Barua 
4134f7415efSAmlan Barua   PetscFunctionBegin;
4144f7415efSAmlan Barua   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
4154f7415efSAmlan Barua   PetscValidType(A,1);
4164f7415efSAmlan Barua 
4174f7415efSAmlan Barua   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr);
4184f7415efSAmlan Barua   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr);
4194f7415efSAmlan Barua   printf("size is %d\n",size);
4204f7415efSAmlan Barua   if (size == 1){ /* sequential case */
4214f7415efSAmlan Barua   printf("Routine is getting called\n");
4224f7415efSAmlan Barua #if defined(PETSC_USE_COMPLEX)
4234f7415efSAmlan Barua     if (fin) {ierr = VecCreateSeq(PETSC_COMM_SELF,N,fin);CHKERRQ(ierr);}
4244f7415efSAmlan Barua     if (fout){ierr = VecCreateSeq(PETSC_COMM_SELF,N,fout);CHKERRQ(ierr);}
4254f7415efSAmlan Barua     if (bout){ierr = VecCreateSeq(PETSC_COMM_SELF,N,bout);CHKERRQ(ierr);}
4264f7415efSAmlan Barua #else
4274f7415efSAmlan Barua     if (fin) {ierr = VecCreateSeq(PETSC_COMM_SELF,N*2*(dim[ndim-1]/2+1)/dim[ndim-1],fin);CHKERRQ(ierr);}
4284f7415efSAmlan Barua     if (fout){ierr = VecCreateSeq(PETSC_COMM_SELF,2*N*(dim[ndim-1]/2+1)/dim[ndim-1],fout);CHKERRQ(ierr);}
4294f7415efSAmlan Barua     if (bout){ierr = VecCreateSeq(PETSC_COMM_SELF,2*N*(dim[ndim-1]/2+1)/dim[ndim-1],bout);CHKERRQ(ierr);}
4304f7415efSAmlan Barua #endif
4314f7415efSAmlan Barua   } else {
4324f7415efSAmlan Barua     ptrdiff_t      alloc_local,local_n0,local_0_start;
4334f7415efSAmlan Barua     ptrdiff_t      local_n1,local_1_end;
434*0c9b18e4SAmlan Barua     fftw_complex   *data_fin,*data_fout,*data_bout;
4354f7415efSAmlan Barua     double         *data_finr,*data_boutr ;
4364f7415efSAmlan Barua     ptrdiff_t      local_1_start,temp;
4374f7415efSAmlan Barua     switch (ndim){
4384f7415efSAmlan Barua           case 1:
4394f7415efSAmlan Barua                  SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Not implemented yet");
4404f7415efSAmlan Barua           case 2:
4414f7415efSAmlan Barua #if !defined(PETSC_USE_COMPLEX)
4424f7415efSAmlan Barua                  alloc_local =  fftw_mpi_local_size_2d_transposed(dim[0],dim[1]/2+1,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
4434f7415efSAmlan Barua                  N1 = 2*dim[0]*(dim[1]/2+1); n1 = 2*local_n0*(dim[1]/2+1);
4444f7415efSAmlan Barua                  if (fin) {
4454f7415efSAmlan Barua                            data_finr=(double *)fftw_malloc(sizeof(double)*alloc_local*2);
4464f7415efSAmlan Barua                            ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,(PetscInt)n1,N1,(PetscScalar*)data_finr,fin);CHKERRQ(ierr);
4474f7415efSAmlan Barua                            ierr = VecGetSize(*fin,&vsize);CHKERRQ(ierr);
4484f7415efSAmlan Barua                            //printf("The code comes here with vector size %d\n",vsize);
4494f7415efSAmlan Barua                            (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
4504f7415efSAmlan Barua                           }
4514f7415efSAmlan Barua                  if (fout) {
4524f7415efSAmlan Barua                             data_fout=(fftw_complex *)fftw_malloc(sizeof(fftw_complex)*alloc_local);
4534f7415efSAmlan Barua                             ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n1,N1,(PetscScalar*)data_fout,fout);CHKERRQ(ierr);
4544f7415efSAmlan Barua                             (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
4554f7415efSAmlan Barua                            }
4564f7415efSAmlan Barua                  if (bout) {
4574f7415efSAmlan Barua                            data_boutr=(double *)fftw_malloc(sizeof(double)*alloc_local*2);
4584f7415efSAmlan Barua                            ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,(PetscInt)n1,N1,(PetscScalar*)data_boutr,bout);CHKERRQ(ierr);
4594f7415efSAmlan Barua                            ierr = VecGetSize(*bout,&vsize);CHKERRQ(ierr);
4604f7415efSAmlan Barua                            //printf("The code comes here with vector size %d\n",vsize);
4614f7415efSAmlan Barua                            (*bout)->ops->destroy   = VecDestroy_MPIFFTW;
4624f7415efSAmlan Barua                           }
4634f7415efSAmlan Barua 
4644f7415efSAmlan Barua                            //printf("Vector size from fftw.c is  given by %d, %d\n",n1,N1);
4654f7415efSAmlan Barua 
4664f7415efSAmlan Barua #else
4674f7415efSAmlan Barua       /* Get local size */
468*0c9b18e4SAmlan Barua                  printf("Code works for paralllel 2d complex DFT\n");
4694f7415efSAmlan Barua                  alloc_local = fftw_mpi_local_size_2d(dim[0],dim[1],comm,&local_n0,&local_0_start);
4704f7415efSAmlan Barua                  if (fin) {
4714f7415efSAmlan Barua                            data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
4724f7415efSAmlan Barua                            ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
4734f7415efSAmlan Barua                            (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
4744f7415efSAmlan Barua                           }
4754f7415efSAmlan Barua                  if (fout) {
4764f7415efSAmlan Barua                             data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
4774f7415efSAmlan Barua                             ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
4784f7415efSAmlan Barua                             (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
4794f7415efSAmlan Barua                            }
4804f7415efSAmlan Barua                  if (bout) {
4814f7415efSAmlan Barua                            data_bout  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
4824f7415efSAmlan Barua                            ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_bout,bout);CHKERRQ(ierr);
4834f7415efSAmlan Barua                            (*bout)->ops->destroy   = VecDestroy_MPIFFTW;
4844f7415efSAmlan Barua                           }
4854f7415efSAmlan Barua 
4864f7415efSAmlan Barua      //printf("Hope this does not come here");
4874f7415efSAmlan Barua #endif
4884f7415efSAmlan Barua       break;
4894f7415efSAmlan Barua 
4904f7415efSAmlan Barua           case 3:
4914f7415efSAmlan Barua #if !defined(PETSC_USE_COMPLEX)
4924f7415efSAmlan Barua                  alloc_local =  fftw_mpi_local_size_3d_transposed(dim[0],dim[1],dim[2]/2+1,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
4934f7415efSAmlan Barua                  N1 = 2*dim[0]*dim[1]*(dim[2]/2+1); n1 = 2*local_n0*dim[1]*(dim[2]/2+1);
4944f7415efSAmlan Barua                  if (fin) {
4954f7415efSAmlan Barua                          data_finr=(double *)fftw_malloc(sizeof(double)*alloc_local*2);
4964f7415efSAmlan Barua                          ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,(PetscInt)n1,N1,(PetscScalar*)data_finr,fin);CHKERRQ(ierr);
4974f7415efSAmlan Barua                          ierr = VecGetSize(*fin,&vsize);CHKERRQ(ierr);
4984f7415efSAmlan Barua         //printf("The code comes here with vector size %d\n",vsize);
4994f7415efSAmlan Barua                          (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
5004f7415efSAmlan Barua                          }
5014f7415efSAmlan Barua                  if (fout) {
5024f7415efSAmlan Barua                           data_fout=(fftw_complex *)fftw_malloc(sizeof(fftw_complex)*alloc_local);
5034f7415efSAmlan Barua                           ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n1,N1,(PetscScalar*)data_fout,fout);CHKERRQ(ierr);
5044f7415efSAmlan Barua                           (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
5054f7415efSAmlan Barua                           }
5064f7415efSAmlan Barua                  if (bout) {
5074f7415efSAmlan Barua                          data_boutr=(double *)fftw_malloc(sizeof(double)*alloc_local*2);
5084f7415efSAmlan Barua                          ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,(PetscInt)n1,N1,(PetscScalar*)data_boutr,bout);CHKERRQ(ierr);
5094f7415efSAmlan Barua                         // ierr = VecGetSize(*fin,&vsize);CHKERRQ(ierr);
5104f7415efSAmlan Barua         //printf("The code comes here with vector size %d\n",vsize);
5114f7415efSAmlan Barua                          (*bout)->ops->destroy   = VecDestroy_MPIFFTW;
5124f7415efSAmlan Barua                           }
5134f7415efSAmlan Barua 
5144f7415efSAmlan Barua       //printf("Vector size from fftw.c is  given by %d, %d\n",n1,N1);
5154f7415efSAmlan Barua 
5164f7415efSAmlan Barua #else
517*0c9b18e4SAmlan Barua                  alloc_local = fftw_mpi_local_size_3d(dim[0],dim[1],dim[2],comm,&local_n0,&local_0_start);
518*0c9b18e4SAmlan Barua //      printf("The quantity n is %d",n);
519*0c9b18e4SAmlan Barua                  if (fin) {
520*0c9b18e4SAmlan Barua                          data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
521*0c9b18e4SAmlan Barua                          ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
522*0c9b18e4SAmlan Barua                          (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
523*0c9b18e4SAmlan Barua                          }
524*0c9b18e4SAmlan Barua                  if (fout) {
525*0c9b18e4SAmlan Barua                           data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
526*0c9b18e4SAmlan Barua                           ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
527*0c9b18e4SAmlan Barua                           (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
528*0c9b18e4SAmlan Barua                          }
529*0c9b18e4SAmlan Barua                  if (bout) {
530*0c9b18e4SAmlan Barua                          data_bout  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
531*0c9b18e4SAmlan Barua                          ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_bout,bout);CHKERRQ(ierr);
532*0c9b18e4SAmlan Barua                          (*bout)->ops->destroy   = VecDestroy_MPIFFTW;
533*0c9b18e4SAmlan Barua                          }
534*0c9b18e4SAmlan Barua 
5354f7415efSAmlan Barua #endif
5364f7415efSAmlan Barua                  break;
5374f7415efSAmlan Barua           default:
5384f7415efSAmlan Barua #if !defined(PETSC_USE_COMPLEX)
5394f7415efSAmlan Barua                  temp = (fftw->dim_fftw)[fftw->ndim_fftw-1];
5404f7415efSAmlan Barua                  printf("The value of temp is %ld\n",temp);
5414f7415efSAmlan Barua                  (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp/2 + 1;
5424f7415efSAmlan Barua                  alloc_local = fftw_mpi_local_size_transposed(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
5434f7415efSAmlan Barua                  N1 = 2*N*(PetscInt)((fftw->dim_fftw)[fftw->ndim_fftw-1])/((PetscInt) temp);
5444f7415efSAmlan Barua                  (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp;
5454f7415efSAmlan Barua 
5464f7415efSAmlan Barua                  if (fin) {
5474f7415efSAmlan Barua                          data_finr=(double *)fftw_malloc(sizeof(double)*alloc_local*2);
5484f7415efSAmlan Barua                          ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,(PetscInt)n,N1,(PetscScalar*)data_finr,fin);CHKERRQ(ierr);
5494f7415efSAmlan Barua                          ierr = VecGetSize(*fin,&vsize);CHKERRQ(ierr);
5504f7415efSAmlan Barua         //printf("The code comes here with vector size %d\n",vsize);
5514f7415efSAmlan Barua                          (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
5524f7415efSAmlan Barua                         }
5534f7415efSAmlan Barua                  if (fout) {
5544f7415efSAmlan Barua                           data_fout=(fftw_complex *)fftw_malloc(sizeof(fftw_complex)*alloc_local);
5554f7415efSAmlan Barua                           ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n,N1,(PetscScalar*)data_fout,fout);CHKERRQ(ierr);
5564f7415efSAmlan Barua                           (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
5574f7415efSAmlan Barua                         }
5584f7415efSAmlan Barua                  if (bout) {
5594f7415efSAmlan Barua                          data_boutr=(double *)fftw_malloc(sizeof(double)*alloc_local*2);
5604f7415efSAmlan Barua                          ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,(PetscInt)n,N1,(PetscScalar*)data_boutr,bout);CHKERRQ(ierr);
5614f7415efSAmlan Barua                          //ierr = VecGetSize(*fin,&vsize);CHKERRQ(ierr);
5624f7415efSAmlan Barua         //printf("The code comes here with vector size %d\n",vsize);
5634f7415efSAmlan Barua                          (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
5644f7415efSAmlan Barua                         }
5654f7415efSAmlan Barua 
5664f7415efSAmlan Barua #else
567*0c9b18e4SAmlan Barua                 alloc_local = fftw_mpi_local_size(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start);
568*0c9b18e4SAmlan Barua                 if (fin) {
569*0c9b18e4SAmlan Barua                        data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
570*0c9b18e4SAmlan Barua                        ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
571*0c9b18e4SAmlan Barua                        (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
572*0c9b18e4SAmlan Barua                        }
573*0c9b18e4SAmlan Barua                 if (fout) {
574*0c9b18e4SAmlan Barua                          data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
575*0c9b18e4SAmlan Barua                          ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
576*0c9b18e4SAmlan Barua                          (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
577*0c9b18e4SAmlan Barua                        }
578*0c9b18e4SAmlan Barua                 if (bout) {
579*0c9b18e4SAmlan Barua                        data_bout  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
580*0c9b18e4SAmlan Barua                        ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_bout,bout);CHKERRQ(ierr);
581*0c9b18e4SAmlan Barua                        (*bout)->ops->destroy   = VecDestroy_MPIFFTW;
582*0c9b18e4SAmlan Barua                        }
5834f7415efSAmlan Barua 
584*0c9b18e4SAmlan Barua 
585*0c9b18e4SAmlan Barua 
5864f7415efSAmlan Barua #endif
5874f7415efSAmlan Barua             break;
5884f7415efSAmlan Barua           }
5894f7415efSAmlan Barua   }
5904f7415efSAmlan Barua   PetscFunctionReturn(0);
5914f7415efSAmlan Barua }
5924f7415efSAmlan Barua EXTERN_C_END
5934f7415efSAmlan Barua 
5944f7415efSAmlan Barua 
5954f7415efSAmlan Barua 
596c92418ccSAmlan Barua #undef __FUNCT__
597b2b77a04SHong Zhang #define __FUNCT__ "MatGetVecs_FFTW"
598b2b77a04SHong Zhang /*
599b2b77a04SHong Zhang    MatGetVecs_FFTW - Get vector(s) compatible with the matrix, i.e. with the
600b2b77a04SHong Zhang      parallel layout determined by FFTW
601b2b77a04SHong Zhang 
602b2b77a04SHong Zhang    Collective on Mat
603b2b77a04SHong Zhang 
604b2b77a04SHong Zhang    Input Parameter:
605b2b77a04SHong Zhang .  mat - the matrix
606b2b77a04SHong Zhang 
607b2b77a04SHong Zhang    Output Parameter:
608b2b77a04SHong Zhang +   fin - (optional) input vector of forward FFTW
609b2b77a04SHong Zhang -   fout - (optional) output vector of forward FFTW
610b2b77a04SHong Zhang 
611b2b77a04SHong Zhang   Level: advanced
612b2b77a04SHong Zhang 
613b2b77a04SHong Zhang .seealso: MatCreateFFTW()
614b2b77a04SHong Zhang */
615b2b77a04SHong Zhang PetscErrorCode  MatGetVecs_FFTW(Mat A,Vec *fin,Vec *fout)
616b2b77a04SHong Zhang {
617b2b77a04SHong Zhang   PetscErrorCode ierr;
618b2b77a04SHong Zhang   PetscMPIInt    size,rank;
619b2b77a04SHong Zhang   MPI_Comm       comm=((PetscObject)A)->comm;
620b2b77a04SHong Zhang   Mat_FFT        *fft = (Mat_FFT*)A->data;
621c92418ccSAmlan Barua   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
6223c3a9c75SAmlan Barua   PetscInt       N=fft->N, N1, n1,vsize;
623e81bb599SAmlan Barua   PetscInt       ndim=fft->ndim,*dim=fft->dim,n=fft->n;
624b2b77a04SHong Zhang 
625b2b77a04SHong Zhang   PetscFunctionBegin;
626b2b77a04SHong Zhang   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
627b2b77a04SHong Zhang   PetscValidType(A,1);
628b2b77a04SHong Zhang 
629b2b77a04SHong Zhang   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
630b2b77a04SHong Zhang   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
631b2b77a04SHong Zhang   if (size == 1){ /* sequential case */
632e5d7f247SAmlan Barua #if defined(PETSC_USE_COMPLEX)
633b2b77a04SHong Zhang     if (fin) {ierr = VecCreateSeq(PETSC_COMM_SELF,N,fin);CHKERRQ(ierr);}
634b2b77a04SHong Zhang     if (fout){ierr = VecCreateSeq(PETSC_COMM_SELF,N,fout);CHKERRQ(ierr);}
635e5d7f247SAmlan Barua #else
636e81bb599SAmlan Barua     if (fin) {ierr = VecCreateSeq(PETSC_COMM_SELF,N*2*(dim[ndim-1]/2+1)/dim[ndim-1],fin);CHKERRQ(ierr);}
637e81bb599SAmlan Barua     if (fout){ierr = VecCreateSeq(PETSC_COMM_SELF,2*N*(dim[ndim-1]/2+1)/dim[ndim-1],fout);CHKERRQ(ierr);}
638e81bb599SAmlan Barua #endif
639b2b77a04SHong Zhang   } else {        /* mpi case */
640b2b77a04SHong Zhang     ptrdiff_t      alloc_local,local_n0,local_0_start;
6416882372aSHong Zhang     ptrdiff_t      local_n1,local_1_end;
642b2b77a04SHong Zhang     fftw_complex   *data_fin,*data_fout;
64351d1eed7SAmlan Barua     double         *data_finr ;
644b3a17365SAmlan Barua     ptrdiff_t      local_1_start,temp;
645c92418ccSAmlan Barua //    PetscInt ctr;
646c92418ccSAmlan Barua //    ptrdiff_t      ndim1,*pdim;
647c92418ccSAmlan Barua //    ndim1=(ptrdiff_t) ndim;
648c92418ccSAmlan Barua //    pdim = (ptrdiff_t *)calloc(ndim,sizeof(ptrdiff_t));
64911902ff2SHong Zhang 
650c92418ccSAmlan Barua //    for(ctr=0;ctr<ndim;ctr++)
651f76f76beSAmlan Barua //        {k
652c92418ccSAmlan Barua //           pdim[ctr] = dim[ctr];
653c92418ccSAmlan Barua //       }
654b2b77a04SHong Zhang 
655f76f76beSAmlan Barua 
656f76f76beSAmlan Barua 
657b2b77a04SHong Zhang     switch (ndim){
658b2b77a04SHong Zhang     case 1:
6596882372aSHong Zhang       /* Get local size */
660e5d7f247SAmlan Barua       /* We need to write an error message here saying that one cannot call this routine when doing parallel 1D real FFTW */
661c92418ccSAmlan Barua //      SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Works only for parallel Multi-dimensional FFTW, Dimension>1. Check Documentation for MatGetVecs_FFTW1D routine");
6626882372aSHong Zhang       alloc_local = fftw_mpi_local_size_1d(dim[0],comm,FFTW_FORWARD,FFTW_ESTIMATE,&local_n0,&local_0_start,&local_n1,&local_1_end);
6636882372aSHong Zhang       if (fin) {
6646882372aSHong Zhang         data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
6656882372aSHong Zhang         ierr = VecCreateMPIWithArray(comm,local_n0,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
6666882372aSHong Zhang         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
6676882372aSHong Zhang       }
6686882372aSHong Zhang       if (fout) {
6696882372aSHong Zhang         data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
6706882372aSHong Zhang         ierr = VecCreateMPIWithArray(comm,local_n1,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
6716882372aSHong Zhang         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
6726882372aSHong Zhang       }
673b2b77a04SHong Zhang       break;
674b2b77a04SHong Zhang     case 2:
6753c3a9c75SAmlan Barua #if !defined(PETSC_USE_COMPLEX)
6763c3a9c75SAmlan Barua       alloc_local =  fftw_mpi_local_size_2d_transposed(dim[0],dim[1]/2+1,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
6773c3a9c75SAmlan Barua       N1 = 2*dim[0]*(dim[1]/2+1); n1 = 2*local_n0*(dim[1]/2+1);
6783c3a9c75SAmlan Barua       if (fin) {
6793c3a9c75SAmlan Barua         data_finr=(double *)fftw_malloc(sizeof(double)*alloc_local*2);
68054dd5118SAmlan Barua         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,(PetscInt)n1,N1,(PetscScalar*)data_finr,fin);CHKERRQ(ierr);
6813c3a9c75SAmlan Barua         ierr = VecGetSize(*fin,&vsize);CHKERRQ(ierr);
68209dd8a53SAmlan Barua         //printf("The code comes here with vector size %d\n",vsize);
6833c3a9c75SAmlan Barua         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
6843c3a9c75SAmlan Barua       }
6853c3a9c75SAmlan Barua       if (fout) {
68609dd8a53SAmlan Barua         data_fout=(fftw_complex *)fftw_malloc(sizeof(fftw_complex)*alloc_local);
68709dd8a53SAmlan Barua         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n1,N1,(PetscScalar*)data_fout,fout);CHKERRQ(ierr);
6883c3a9c75SAmlan Barua         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
6893c3a9c75SAmlan Barua       }
690f76f76beSAmlan Barua       printf("Vector size from fftw.c is  given by %d, %d\n",n1,N1);
6913c3a9c75SAmlan Barua 
6923c3a9c75SAmlan Barua #else
693b2b77a04SHong Zhang       /* Get local size */
69464657d84SAmlan Barua      //printf("Hope this does not come here");
695b2b77a04SHong Zhang       alloc_local = fftw_mpi_local_size_2d(dim[0],dim[1],comm,&local_n0,&local_0_start);
696b2b77a04SHong Zhang       if (fin) {
697b2b77a04SHong Zhang         data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
698b2b77a04SHong Zhang         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
699b2b77a04SHong Zhang         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
700b2b77a04SHong Zhang       }
701b2b77a04SHong Zhang       if (fout) {
702b2b77a04SHong Zhang         data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
703b2b77a04SHong Zhang         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
704b2b77a04SHong Zhang         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
705b2b77a04SHong Zhang       }
70664657d84SAmlan Barua      //printf("Hope this does not come here");
7073c3a9c75SAmlan Barua #endif
708b2b77a04SHong Zhang       break;
709b2b77a04SHong Zhang     case 3:
7106882372aSHong Zhang       /* Get local size */
7113c3a9c75SAmlan Barua #if !defined(PETSC_USE_COMPLEX)
71251d1eed7SAmlan Barua       alloc_local =  fftw_mpi_local_size_3d_transposed(dim[0],dim[1],dim[2]/2+1,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
71351d1eed7SAmlan Barua       N1 = 2*dim[0]*dim[1]*(dim[2]/2+1); n1 = 2*local_n0*dim[1]*(dim[2]/2+1);
71451d1eed7SAmlan Barua       if (fin) {
71551d1eed7SAmlan Barua         data_finr=(double *)fftw_malloc(sizeof(double)*alloc_local*2);
71651d1eed7SAmlan Barua         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,(PetscInt)n1,N1,(PetscScalar*)data_finr,fin);CHKERRQ(ierr);
71751d1eed7SAmlan Barua         ierr = VecGetSize(*fin,&vsize);CHKERRQ(ierr);
71851d1eed7SAmlan Barua         //printf("The code comes here with vector size %d\n",vsize);
71951d1eed7SAmlan Barua         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
72051d1eed7SAmlan Barua       }
72151d1eed7SAmlan Barua       if (fout) {
72251d1eed7SAmlan Barua         data_fout=(fftw_complex *)fftw_malloc(sizeof(fftw_complex)*alloc_local);
72351d1eed7SAmlan Barua         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n1,N1,(PetscScalar*)data_fout,fout);CHKERRQ(ierr);
72451d1eed7SAmlan Barua         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
72551d1eed7SAmlan Barua       }
72651d1eed7SAmlan Barua       printf("Vector size from fftw.c is  given by %d, %d\n",n1,N1);
72751d1eed7SAmlan Barua 
72851d1eed7SAmlan Barua 
7293c3a9c75SAmlan Barua #else
7306882372aSHong Zhang       alloc_local = fftw_mpi_local_size_3d(dim[0],dim[1],dim[2],comm,&local_n0,&local_0_start);
73111902ff2SHong Zhang //      printf("The quantity n is %d",n);
7326882372aSHong Zhang       if (fin) {
7336882372aSHong Zhang         data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
7346882372aSHong Zhang         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
7356882372aSHong Zhang         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
7366882372aSHong Zhang       }
7376882372aSHong Zhang       if (fout) {
7386882372aSHong Zhang         data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
7396882372aSHong Zhang         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
7406882372aSHong Zhang         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
7416882372aSHong Zhang       }
7423c3a9c75SAmlan Barua #endif
743b2b77a04SHong Zhang       break;
744b2b77a04SHong Zhang     default:
7456882372aSHong Zhang       /* Get local size */
7463c3a9c75SAmlan Barua #if !defined(PETSC_USE_COMPLEX)
747b3a17365SAmlan Barua       temp = (fftw->dim_fftw)[fftw->ndim_fftw-1];
748b3a17365SAmlan Barua       printf("The value of temp is %ld\n",temp);
749b3a17365SAmlan Barua       (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp/2 + 1;
750b3a17365SAmlan Barua       alloc_local = fftw_mpi_local_size_transposed(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
751b3a17365SAmlan Barua       N1 = 2*N*(PetscInt)((fftw->dim_fftw)[fftw->ndim_fftw-1])/((PetscInt) temp);
752b3a17365SAmlan Barua       (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp;
753b3a17365SAmlan Barua       if (fin) {
754b3a17365SAmlan Barua         data_finr=(double *)fftw_malloc(sizeof(double)*alloc_local*2);
755b3a17365SAmlan Barua         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,(PetscInt)n,N1,(PetscScalar*)data_finr,fin);CHKERRQ(ierr);
756b3a17365SAmlan Barua         ierr = VecGetSize(*fin,&vsize);CHKERRQ(ierr);
757b3a17365SAmlan Barua         //printf("The code comes here with vector size %d\n",vsize);
758b3a17365SAmlan Barua         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
759b3a17365SAmlan Barua       }
760b3a17365SAmlan Barua       if (fout) {
761b3a17365SAmlan Barua         data_fout=(fftw_complex *)fftw_malloc(sizeof(fftw_complex)*alloc_local);
762b3a17365SAmlan Barua         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n,N1,(PetscScalar*)data_fout,fout);CHKERRQ(ierr);
763b3a17365SAmlan Barua         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
764b3a17365SAmlan Barua       }
765b3a17365SAmlan Barua 
7663c3a9c75SAmlan Barua #else
767c92418ccSAmlan Barua       alloc_local = fftw_mpi_local_size(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start);
76811902ff2SHong Zhang //      printf("The value of alloc local is %d from process %d\n",alloc_local,rank);
76911902ff2SHong Zhang //      printf("The value of alloc local is %d",alloc_local);
77011902ff2SHong Zhang //      pdim=(ptrdiff_t *)calloc(ndim,sizeof(ptrdiff_t));
77111902ff2SHong Zhang //      for(i=0;i<ndim;i++)
77211902ff2SHong Zhang //         {
77311902ff2SHong Zhang //          pdim[i]=dim[i];printf("%d",pdim[i]);
77411902ff2SHong Zhang //         }
77511902ff2SHong Zhang //      alloc_local = fftw_mpi_local_size(ndim,pdim,comm,&local_n0,&local_0_start);
77611902ff2SHong Zhang //      printf("The quantity n is %d",n);
7776882372aSHong Zhang       if (fin) {
7786882372aSHong Zhang         data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
7796882372aSHong Zhang         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
7806882372aSHong Zhang         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
7816882372aSHong Zhang       }
7826882372aSHong Zhang       if (fout) {
7836882372aSHong Zhang         data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
7846882372aSHong Zhang         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
7856882372aSHong Zhang         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
7866882372aSHong Zhang       }
7873c3a9c75SAmlan Barua #endif
788b2b77a04SHong Zhang       break;
789b2b77a04SHong Zhang     }
790b2b77a04SHong Zhang   }
79154dd5118SAmlan Barua //  if (fin){
79254dd5118SAmlan Barua //    ierr = PetscLayoutReference(A->cmap,&(*fin)->map);CHKERRQ(ierr);
79354dd5118SAmlan Barua //  }
79454dd5118SAmlan Barua //  if (fout){
79554dd5118SAmlan Barua //    ierr = PetscLayoutReference(A->rmap,&(*fout)->map);CHKERRQ(ierr);
79654dd5118SAmlan Barua //  }
797b2b77a04SHong Zhang   PetscFunctionReturn(0);
798b2b77a04SHong Zhang }
799b2b77a04SHong Zhang 
800b2b77a04SHong Zhang #undef __FUNCT__
8013c3a9c75SAmlan Barua #define __FUNCT__ "InputTransformFFT"
8023c3a9c75SAmlan Barua PetscErrorCode InputTransformFFT(Mat A,Vec x,Vec y)
8033c3a9c75SAmlan Barua {
8043c3a9c75SAmlan Barua   PetscErrorCode ierr;
8053c3a9c75SAmlan Barua   PetscFunctionBegin;
8063c3a9c75SAmlan Barua   ierr = PetscTryMethod(A,"InputTransformFFT_C",(Mat,Vec,Vec),(A,x,y));CHKERRQ(ierr);
8073c3a9c75SAmlan Barua   PetscFunctionReturn(0);
8083c3a9c75SAmlan Barua }
80954efbe56SHong Zhang 
810b2b77a04SHong Zhang /*
8113c3a9c75SAmlan Barua       InputTransformFFT_FFTW - Copies the user data to the vector that goes into FFTW block
8123c3a9c75SAmlan Barua   Input A, x, y
8133c3a9c75SAmlan Barua   A - FFTW matrix
81454dd5118SAmlan Barua   x - user data
815b2b77a04SHong Zhang   Options Database Keys:
816b2b77a04SHong Zhang + -mat_fftw_plannerflags - set FFTW planner flags
817b2b77a04SHong Zhang 
818b2b77a04SHong Zhang    Level: intermediate
819b2b77a04SHong Zhang 
820b2b77a04SHong Zhang */
8213c3a9c75SAmlan Barua 
8223c3a9c75SAmlan Barua EXTERN_C_BEGIN
8233c3a9c75SAmlan Barua #undef __FUNCT__
8243c3a9c75SAmlan Barua #define __FUNCT__ "InputTransformFFT_FTTW"
8253c3a9c75SAmlan Barua PetscErrorCode InputTransformFFT_FFTW(Mat A,Vec x,Vec y)
8263c3a9c75SAmlan Barua {
8273c3a9c75SAmlan Barua   PetscErrorCode ierr;
8283c3a9c75SAmlan Barua   MPI_Comm       comm=((PetscObject)A)->comm;
8293c3a9c75SAmlan Barua   Mat_FFT        *fft = (Mat_FFT*)A->data;
8303c3a9c75SAmlan Barua   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
8313c3a9c75SAmlan Barua   PetscInt       N=fft->N, N1, n1 ,NM;
832b223cc97SAmlan Barua   PetscInt       ndim=fft->ndim,*dim=fft->dim;//n=fft->n;
833b223cc97SAmlan Barua   PetscInt       low, *indx1, *indx2, tempindx, tempindx1;
834e5d7f247SAmlan Barua   PetscInt       i,j,k,rank,size,partial_dim;
8353c3a9c75SAmlan Barua   ptrdiff_t      alloc_local,local_n0,local_0_start;
836e5d7f247SAmlan Barua   ptrdiff_t      local_n1,local_1_start,temp;
8373c3a9c75SAmlan Barua   VecScatter     vecscat;
8383c3a9c75SAmlan Barua   IS             list1,list2;
8393c3a9c75SAmlan Barua 
840f76f76beSAmlan Barua   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
841f76f76beSAmlan Barua   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
8423c3a9c75SAmlan Barua   ierr = VecGetOwnershipRange(y,&low,PETSC_NULL);
843f76f76beSAmlan Barua   printf("Local ownership starts at %d\n",low);
8443c3a9c75SAmlan Barua 
845e81bb599SAmlan Barua   if (size==1)
846e81bb599SAmlan Barua     {
8477536937bSAmlan Barua /*     switch (ndim){
848e81bb599SAmlan Barua      case 1:
849e81bb599SAmlan Barua           ierr = PetscMalloc(sizeof(PetscInt)*dim[0],&indx1);CHKERRQ(ierr);
850e81bb599SAmlan Barua           for (i=0;i<dim[0];i++)
851e81bb599SAmlan Barua              {
852e81bb599SAmlan Barua               indx1[i] = i;
853e81bb599SAmlan Barua              }
854e81bb599SAmlan Barua           ierr = ISCreateGeneral(comm,dim[0],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
855e81bb599SAmlan Barua           ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
856e81bb599SAmlan Barua           ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
857e81bb599SAmlan Barua           ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
858e81bb599SAmlan Barua           ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
859b223cc97SAmlan Barua           ierr = ISDestroy(&list1);CHKERRQ(ierr);
860b223cc97SAmlan Barua           ierr = PetscFree(indx1);CHKERRQ(ierr);
861e81bb599SAmlan Barua      break;
862e81bb599SAmlan Barua 
863e81bb599SAmlan Barua      case 2:
864e81bb599SAmlan Barua           ierr = PetscMalloc(sizeof(PetscInt)*dim[0]*dim[1],&indx1);CHKERRQ(ierr);
865e81bb599SAmlan Barua           for (i=0;i<dim[0];i++){
866e81bb599SAmlan Barua              for (j=0;j<dim[1];j++){
867e81bb599SAmlan Barua                 indx1[i*dim[1]+j] = i*dim[1] + j;
868e81bb599SAmlan Barua              }
869e81bb599SAmlan Barua           }
870e81bb599SAmlan Barua           ierr = ISCreateGeneral(comm,dim[0]*dim[1],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
871e81bb599SAmlan Barua           ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
872e81bb599SAmlan Barua           ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
873e81bb599SAmlan Barua           ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
874e81bb599SAmlan Barua           ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
875b223cc97SAmlan Barua           ierr = ISDestroy(&list1);CHKERRQ(ierr);
876b223cc97SAmlan Barua           ierr = PetscFree(indx1);CHKERRQ(ierr);
877e81bb599SAmlan Barua           break;
878e81bb599SAmlan Barua      case 3:
879e81bb599SAmlan Barua           ierr = PetscMalloc(sizeof(PetscInt)*dim[0]*dim[1]*dim[2],&indx1);CHKERRQ(ierr);
880e81bb599SAmlan Barua           for (i=0;i<dim[0];i++){
881e81bb599SAmlan Barua              for (j=0;j<dim[1];j++){
882e81bb599SAmlan Barua                 for (k=0;k<dim[2];k++){
883e81bb599SAmlan Barua                    indx1[i*dim[1]*dim[2]+j*dim[2]+k] = i*dim[1]*dim[2] + j*dim[2] + k;
884e81bb599SAmlan Barua                 }
885e81bb599SAmlan Barua              }
886e81bb599SAmlan Barua           }
887e81bb599SAmlan Barua           ierr = ISCreateGeneral(comm,dim[0]*dim[1]*dim[2],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
888e81bb599SAmlan Barua           ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
889e81bb599SAmlan Barua           ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
890e81bb599SAmlan Barua           ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
891e81bb599SAmlan Barua           ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
892b223cc97SAmlan Barua           ierr = ISDestroy(&list1);CHKERRQ(ierr);
893b223cc97SAmlan Barua           ierr = PetscFree(indx1);CHKERRQ(ierr);
894e81bb599SAmlan Barua           break;
895e81bb599SAmlan Barua      default:
8967536937bSAmlan Barua */
8976971673cSAmlan Barua           ierr = ISCreateStride(PETSC_COMM_SELF,N,0,1,&list1);
8986971673cSAmlan Barua           ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
8996971673cSAmlan Barua           ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9006971673cSAmlan Barua           ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9016971673cSAmlan Barua           ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
902b223cc97SAmlan Barua           ierr = ISDestroy(&list1);CHKERRQ(ierr);
9036971673cSAmlan Barua           //ierr = ISDestroy(list1);CHKERRQ(ierr);
9047536937bSAmlan Barua  //         break;
9057536937bSAmlan Barua   //    }
906e81bb599SAmlan Barua     }
907e81bb599SAmlan Barua 
908e81bb599SAmlan Barua  else{
909e81bb599SAmlan Barua 
9103c3a9c75SAmlan Barua  switch (ndim){
9113c3a9c75SAmlan Barua  case 1:
91264657d84SAmlan Barua #if defined(PETSC_USE_COMPLEX)
91364657d84SAmlan Barua       alloc_local = fftw_mpi_local_size_1d(dim[0],comm,FFTW_FORWARD,FFTW_ESTIMATE,&local_n0,&local_0_start,&local_n1,&local_1_start);
91464657d84SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0,local_0_start,1,&list1);
91564657d84SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0,low,1,&list2);
91664657d84SAmlan Barua       //ierr = ISView(list1,PETSC_VIEWER_STDOUT_WORLD);
91764657d84SAmlan Barua       //ierr = ISView(list2,PETSC_VIEWER_STDOUT_WORLD);
91864657d84SAmlan Barua       ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr);
91964657d84SAmlan Barua       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
92064657d84SAmlan Barua       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
92164657d84SAmlan Barua       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
92264657d84SAmlan Barua       ierr = ISDestroy(&list1);CHKERRQ(ierr);
92364657d84SAmlan Barua       ierr = ISDestroy(&list2);CHKERRQ(ierr);
92464657d84SAmlan Barua       break;
92564657d84SAmlan Barua #else
926e5d7f247SAmlan Barua   SETERRQ(comm,PETSC_ERR_SUP,"FFTW does not support parallel 1D real transform");
92764657d84SAmlan Barua #endif
9283c3a9c75SAmlan Barua   break;
9293c3a9c75SAmlan Barua  case 2:
930bd59e6c6SAmlan Barua #if defined(PETSC_USE_COMPLEX)
931bd59e6c6SAmlan Barua       //PetscInt my_value;
932bd59e6c6SAmlan Barua       alloc_local = fftw_mpi_local_size_2d(dim[0],dim[1],comm,&local_n0,&local_0_start);
933bd59e6c6SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1],local_0_start*dim[1],1,&list1);
934bd59e6c6SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1],low,1,&list2);
935bd59e6c6SAmlan Barua       //ierr = ISView(list1,PETSC_VIEWER_STDOUT_WORLD);
936bd59e6c6SAmlan Barua       //ierr = ISView(list2,PETSC_VIEWER_STDOUT_WORLD);
937bd59e6c6SAmlan Barua       ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr);
938bd59e6c6SAmlan Barua       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
939bd59e6c6SAmlan Barua       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
940bd59e6c6SAmlan Barua       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
941bd59e6c6SAmlan Barua       ierr = ISDestroy(&list1);CHKERRQ(ierr);
942bd59e6c6SAmlan Barua       ierr = ISDestroy(&list2);CHKERRQ(ierr);
943bd59e6c6SAmlan Barua       break;
944bd59e6c6SAmlan Barua #else
9453c3a9c75SAmlan Barua       alloc_local =  fftw_mpi_local_size_2d_transposed(dim[0],dim[1]/2+1,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
9463c3a9c75SAmlan Barua       N1 = 2*dim[0]*(dim[1]/2+1); n1 = 2*local_n0*(dim[1]/2+1);
9473c3a9c75SAmlan Barua 
9485b3e41ffSAmlan Barua       ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx1);CHKERRQ(ierr);
9495b3e41ffSAmlan Barua       ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx2);CHKERRQ(ierr);
9505b3e41ffSAmlan Barua       printf("Val local_0_start = %ld\n",local_0_start);
9513c3a9c75SAmlan Barua 
9523c3a9c75SAmlan Barua       if (dim[1]%2==0)
9533c3a9c75SAmlan Barua         NM = dim[1]+2;
9543c3a9c75SAmlan Barua       else
9553c3a9c75SAmlan Barua         NM = dim[1]+1;
9563c3a9c75SAmlan Barua 
9573c3a9c75SAmlan Barua       for (i=0;i<local_n0;i++){
9583c3a9c75SAmlan Barua          for (j=0;j<dim[1];j++){
9593c3a9c75SAmlan Barua             tempindx = i*dim[1] + j;
9603c3a9c75SAmlan Barua             tempindx1 = i*NM + j;
9615b3e41ffSAmlan Barua             indx1[tempindx]=local_0_start*dim[1]+tempindx;
9623c3a9c75SAmlan Barua             indx2[tempindx]=low+tempindx1;
9635b3e41ffSAmlan Barua            // printf("Val tempindx1 = %d\n",tempindx1);
9645b3e41ffSAmlan Barua   //          printf("index1 %d from proc %d is \n",indx1[tempindx],rank);
9655b3e41ffSAmlan Barua   //          printf("index2 %d from proc %d is \n",indx2[tempindx],rank);
9665b3e41ffSAmlan Barua   //          printf("-------------------------\n",indx2[tempindx],rank);
9673c3a9c75SAmlan Barua         }
9683c3a9c75SAmlan Barua      }
9693c3a9c75SAmlan Barua 
9703c3a9c75SAmlan Barua       ierr = ISCreateGeneral(comm,local_n0*dim[1],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
9713c3a9c75SAmlan Barua       ierr = ISCreateGeneral(comm,local_n0*dim[1],indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr);
9723c3a9c75SAmlan Barua 
973f76f76beSAmlan Barua       ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr);
974f76f76beSAmlan Barua       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
975f76f76beSAmlan Barua       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
976f76f76beSAmlan Barua       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
977b223cc97SAmlan Barua       ierr = ISDestroy(&list1);CHKERRQ(ierr);
978b223cc97SAmlan Barua       ierr = ISDestroy(&list2);CHKERRQ(ierr);
979b223cc97SAmlan Barua       ierr = PetscFree(indx1);CHKERRQ(ierr);
980b223cc97SAmlan Barua       ierr = PetscFree(indx2);CHKERRQ(ierr);
9813c3a9c75SAmlan Barua       break;
982bd59e6c6SAmlan Barua #endif
9833c3a9c75SAmlan Barua 
9843c3a9c75SAmlan Barua  case 3:
985bd59e6c6SAmlan Barua       #if defined(PETSC_USE_COMPLEX)
986bd59e6c6SAmlan Barua       alloc_local = fftw_mpi_local_size_3d(dim[0],dim[1],dim[2],comm,&local_n0,&local_0_start);
987bd59e6c6SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1]*dim[2],local_0_start*dim[1]*dim[2],1,&list1);
988bd59e6c6SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1]*dim[2],low,1,&list2);
989bd59e6c6SAmlan Barua       //ierr = ISView(list1,PETSC_VIEWER_STDOUT_WORLD);
990bd59e6c6SAmlan Barua       //ierr = ISView(list2,PETSC_VIEWER_STDOUT_WORLD);
991bd59e6c6SAmlan Barua       ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr);
992bd59e6c6SAmlan Barua       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
993bd59e6c6SAmlan Barua       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
994bd59e6c6SAmlan Barua       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
995bd59e6c6SAmlan Barua       ierr = ISDestroy(&list1);CHKERRQ(ierr);
996bd59e6c6SAmlan Barua       ierr = ISDestroy(&list2);CHKERRQ(ierr);
997bd59e6c6SAmlan Barua       break;
998bd59e6c6SAmlan Barua #else
99951d1eed7SAmlan Barua       alloc_local =  fftw_mpi_local_size_3d_transposed(dim[0],dim[1],dim[2]/2+1,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
100051d1eed7SAmlan Barua       N1 = 2*dim[0]*dim[1]*(dim[2]/2+1); n1 = 2*local_n0*dim[1]*(dim[2]/2+1);
100151d1eed7SAmlan Barua 
100251d1eed7SAmlan Barua       ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx1);CHKERRQ(ierr);
100351d1eed7SAmlan Barua       ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx2);CHKERRQ(ierr);
100451d1eed7SAmlan Barua       printf("Val local_0_start = %ld\n",local_0_start);
100551d1eed7SAmlan Barua 
100651d1eed7SAmlan Barua       if (dim[2]%2==0)
100751d1eed7SAmlan Barua         NM = dim[2]+2;
100851d1eed7SAmlan Barua       else
100951d1eed7SAmlan Barua         NM = dim[2]+1;
101051d1eed7SAmlan Barua 
101151d1eed7SAmlan Barua       for (i=0;i<local_n0;i++){
101251d1eed7SAmlan Barua          for (j=0;j<dim[1];j++){
101351d1eed7SAmlan Barua             for (k=0;k<dim[2];k++){
101451d1eed7SAmlan Barua                tempindx = i*dim[1]*dim[2] + j*dim[2] + k;
101551d1eed7SAmlan Barua                tempindx1 = i*dim[1]*NM + j*NM + k;
101651d1eed7SAmlan Barua                indx1[tempindx]=local_0_start*dim[1]*dim[2]+tempindx;
101751d1eed7SAmlan Barua                indx2[tempindx]=low+tempindx1;
101851d1eed7SAmlan Barua             }
101951d1eed7SAmlan Barua            // printf("Val tempindx1 = %d\n",tempindx1);
102051d1eed7SAmlan Barua            // printf("index1 %d from proc %d is \n",indx1[tempindx],rank);
102151d1eed7SAmlan Barua            // printf("index2 %d from proc %d is \n",indx2[tempindx],rank);
102251d1eed7SAmlan Barua            // printf("-------------------------\n",indx2[tempindx],rank);
102351d1eed7SAmlan Barua         }
102451d1eed7SAmlan Barua      }
102551d1eed7SAmlan Barua 
102651d1eed7SAmlan Barua       ierr = ISCreateGeneral(comm,local_n0*dim[1]*dim[2],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
102751d1eed7SAmlan Barua       ierr = ISCreateGeneral(comm,local_n0*dim[1]*dim[2],indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr);
102851d1eed7SAmlan Barua 
102951d1eed7SAmlan Barua       ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr);
103051d1eed7SAmlan Barua       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
103151d1eed7SAmlan Barua       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
103251d1eed7SAmlan Barua       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1033b223cc97SAmlan Barua       ierr = ISDestroy(&list1);CHKERRQ(ierr);
1034b223cc97SAmlan Barua       ierr = ISDestroy(&list2);CHKERRQ(ierr);
1035b223cc97SAmlan Barua       ierr = PetscFree(indx1);CHKERRQ(ierr);
1036b223cc97SAmlan Barua       ierr = PetscFree(indx2);CHKERRQ(ierr);
10373c3a9c75SAmlan Barua       break;
1038bd59e6c6SAmlan Barua #endif
10393c3a9c75SAmlan Barua 
10403c3a9c75SAmlan Barua  default:
1041bd59e6c6SAmlan Barua #if defined(PETSC_USE_COMPLEX)
1042bd59e6c6SAmlan Barua       alloc_local = fftw_mpi_local_size(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start);
1043bd59e6c6SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*(fftw->partial_dim),local_0_start*(fftw->partial_dim),1,&list1);
1044bd59e6c6SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*(fftw->partial_dim),low,1,&list2);
1045bd59e6c6SAmlan Barua       //ierr = ISView(list1,PETSC_VIEWER_STDOUT_WORLD);
1046bd59e6c6SAmlan Barua       //ierr = ISView(list2,PETSC_VIEWER_STDOUT_WORLD);
1047bd59e6c6SAmlan Barua       ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr);
1048bd59e6c6SAmlan Barua       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1049bd59e6c6SAmlan Barua       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1050bd59e6c6SAmlan Barua       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1051bd59e6c6SAmlan Barua       ierr = ISDestroy(&list1);CHKERRQ(ierr);
1052bd59e6c6SAmlan Barua       ierr = ISDestroy(&list2);CHKERRQ(ierr);
1053bd59e6c6SAmlan Barua 
1054bd59e6c6SAmlan Barua 
1055bd59e6c6SAmlan Barua #else
1056e5d7f247SAmlan Barua       temp = (fftw->dim_fftw)[fftw->ndim_fftw-1];
1057e5d7f247SAmlan Barua       printf("The value of temp is %ld\n",temp);
1058e5d7f247SAmlan Barua       (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp/2 + 1;
1059e5d7f247SAmlan Barua       alloc_local = fftw_mpi_local_size_transposed(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
1060e5d7f247SAmlan Barua       N1 = 2*N*(PetscInt)((fftw->dim_fftw)[fftw->ndim_fftw-1])/((PetscInt) temp);
1061e5d7f247SAmlan Barua       (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp;
1062e5d7f247SAmlan Barua 
1063e5d7f247SAmlan Barua       partial_dim = fftw->partial_dim;
1064e5d7f247SAmlan Barua       printf("The value of partial dim is %d\n",partial_dim);
1065e5d7f247SAmlan Barua 
1066e5d7f247SAmlan Barua       ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*partial_dim,&indx1);CHKERRQ(ierr);
1067e5d7f247SAmlan Barua       ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*partial_dim,&indx2);CHKERRQ(ierr);
1068e5d7f247SAmlan Barua       printf("Val local_0_start = %ld\n",local_0_start);
1069e5d7f247SAmlan Barua 
1070e5d7f247SAmlan Barua       if (dim[ndim-1]%2==0)
1071ba6e06d0SAmlan Barua         NM = 2;
1072e5d7f247SAmlan Barua       else
1073ba6e06d0SAmlan Barua         NM = 1;
1074e5d7f247SAmlan Barua 
10756971673cSAmlan Barua       j = low;
10766971673cSAmlan Barua       for (i=0,k=1; i<((PetscInt)local_n0)*partial_dim;i++,k++)
10776971673cSAmlan Barua          {
10786971673cSAmlan Barua           indx1[i] = local_0_start*partial_dim + i;
10796971673cSAmlan Barua           indx2[i] = j;
1080ba6e06d0SAmlan Barua           //printf("The values are %d and %d from %d\n",indx1[i],indx2[i],rank);
10816971673cSAmlan Barua           if (k%dim[ndim-1]==0)
10826971673cSAmlan Barua             { j+=NM;}
10836971673cSAmlan Barua           j++;
10846971673cSAmlan Barua          }
10856971673cSAmlan Barua       ierr = ISCreateGeneral(comm,local_n0*partial_dim,indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
10866971673cSAmlan Barua       ierr = ISCreateGeneral(comm,local_n0*partial_dim,indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr);
10876971673cSAmlan Barua       ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr);
10886971673cSAmlan Barua       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10896971673cSAmlan Barua       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10906971673cSAmlan Barua       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1091b223cc97SAmlan Barua       ierr = ISDestroy(&list1);CHKERRQ(ierr);
1092b223cc97SAmlan Barua       ierr = ISDestroy(&list2);CHKERRQ(ierr);
1093b223cc97SAmlan Barua       ierr = PetscFree(indx1);CHKERRQ(ierr);
1094b223cc97SAmlan Barua       ierr = PetscFree(indx2);CHKERRQ(ierr);
10953c3a9c75SAmlan Barua       break;
1096bd59e6c6SAmlan Barua #endif
10973c3a9c75SAmlan Barua   }
1098e81bb599SAmlan Barua  }
10993c3a9c75SAmlan Barua 
11003c3a9c75SAmlan Barua  return 0;
11013c3a9c75SAmlan Barua }
11023c3a9c75SAmlan Barua EXTERN_C_END
11033c3a9c75SAmlan Barua 
11043c3a9c75SAmlan Barua #undef __FUNCT__
11053c3a9c75SAmlan Barua #define __FUNCT__ "OutputTransformFFT"
11063c3a9c75SAmlan Barua PetscErrorCode OutputTransformFFT(Mat A,Vec x,Vec y)
11073c3a9c75SAmlan Barua {
11083c3a9c75SAmlan Barua   PetscErrorCode ierr;
11093c3a9c75SAmlan Barua   PetscFunctionBegin;
11103c3a9c75SAmlan Barua   ierr = PetscTryMethod(A,"OutputTransformFFT_C",(Mat,Vec,Vec),(A,x,y));CHKERRQ(ierr);
11113c3a9c75SAmlan Barua   PetscFunctionReturn(0);
11123c3a9c75SAmlan Barua }
111354efbe56SHong Zhang 
11145b3e41ffSAmlan Barua /*
11155b3e41ffSAmlan Barua       OutputTransformFFT_FFTW - Copies the FFTW output to the PETSc vector that user can use
11165b3e41ffSAmlan Barua   Input A, x, y
11175b3e41ffSAmlan Barua   A - FFTW matrix
11185b3e41ffSAmlan Barua   x - FFTW vector
11195b3e41ffSAmlan Barua   y - PETSc vector that user can read
11205b3e41ffSAmlan Barua   Options Database Keys:
11215b3e41ffSAmlan Barua + -mat_fftw_plannerflags - set FFTW planner flags
11225b3e41ffSAmlan Barua 
11235b3e41ffSAmlan Barua    Level: intermediate
11245b3e41ffSAmlan Barua 
11253c3a9c75SAmlan Barua */
11263c3a9c75SAmlan Barua 
11273c3a9c75SAmlan Barua EXTERN_C_BEGIN
11283c3a9c75SAmlan Barua #undef __FUNCT__
11295b3e41ffSAmlan Barua #define __FUNCT__ "OutputTransformFFT_FTTW"
11305b3e41ffSAmlan Barua PetscErrorCode OutputTransformFFT_FFTW(Mat A,Vec x,Vec y)
11315b3e41ffSAmlan Barua {
11325b3e41ffSAmlan Barua   PetscErrorCode ierr;
11335b3e41ffSAmlan Barua   MPI_Comm       comm=((PetscObject)A)->comm;
11345b3e41ffSAmlan Barua   Mat_FFT        *fft = (Mat_FFT*)A->data;
11355b3e41ffSAmlan Barua   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
11365b3e41ffSAmlan Barua   PetscInt       N=fft->N, N1, n1 ,NM;
1137b223cc97SAmlan Barua   PetscInt       ndim=fft->ndim,*dim=fft->dim;//n=fft->n;
1138b223cc97SAmlan Barua   PetscInt       low, *indx1, *indx2, tempindx, tempindx1;
1139ba6e06d0SAmlan Barua   PetscInt       i,j,k,rank,size,partial_dim;
11405b3e41ffSAmlan Barua   ptrdiff_t      alloc_local,local_n0,local_0_start;
1141ba6e06d0SAmlan Barua   ptrdiff_t      local_n1,local_1_start,temp;
11425b3e41ffSAmlan Barua   VecScatter     vecscat;
11435b3e41ffSAmlan Barua   IS             list1,list2;
11445b3e41ffSAmlan Barua 
11455b3e41ffSAmlan Barua   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
11465b3e41ffSAmlan Barua   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1147cfe1ae98SAmlan Barua   ierr = VecGetOwnershipRange(x,&low,PETSC_NULL);
11485b3e41ffSAmlan Barua 
1149e81bb599SAmlan Barua   if (size==1){
11507536937bSAmlan Barua /*
11515b3e41ffSAmlan Barua     switch (ndim){
11525b3e41ffSAmlan Barua     case 1:
1153e81bb599SAmlan Barua            ierr = PetscMalloc(sizeof(PetscInt)*dim[0],&indx1);CHKERRQ(ierr);
1154e81bb599SAmlan Barua           for (i=0;i<dim[0];i++)
1155e81bb599SAmlan Barua              {
1156e81bb599SAmlan Barua               indx1[i] = i;
1157e81bb599SAmlan Barua              }
1158e81bb599SAmlan Barua           ierr = ISCreateGeneral(comm,dim[0],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
1159e81bb599SAmlan Barua           ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
1160e81bb599SAmlan Barua           ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1161e81bb599SAmlan Barua           ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1162e81bb599SAmlan Barua           ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1163b223cc97SAmlan Barua           ierr = ISDestroy(&list1);CHKERRQ(ierr);
1164b223cc97SAmlan Barua           ierr = PetscFree(indx1);CHKERRQ(ierr);
1165e81bb599SAmlan Barua           break;
1166e81bb599SAmlan Barua 
1167e81bb599SAmlan Barua     case 2:
1168e81bb599SAmlan Barua          ierr = PetscMalloc(sizeof(PetscInt)*dim[0]*dim[1],&indx1);CHKERRQ(ierr);
1169e81bb599SAmlan Barua           for (i=0;i<dim[0];i++){
1170e81bb599SAmlan Barua              for (j=0;j<dim[1];j++){
1171e81bb599SAmlan Barua                 indx1[i*dim[1]+j] = i*dim[1] + j;
1172e81bb599SAmlan Barua              }
1173e81bb599SAmlan Barua           }
1174e81bb599SAmlan Barua          ierr = ISCreateGeneral(comm,dim[0]*dim[1],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
1175e81bb599SAmlan Barua          ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
1176e81bb599SAmlan Barua          ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1177e81bb599SAmlan Barua          ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1178e81bb599SAmlan Barua          ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1179b223cc97SAmlan Barua          ierr = ISDestroy(&list1);CHKERRQ(ierr);
1180b223cc97SAmlan Barua          ierr = PetscFree(indx1);CHKERRQ(ierr);
1181e81bb599SAmlan Barua          break;
1182e81bb599SAmlan Barua     case 3:
1183e81bb599SAmlan Barua          ierr = PetscMalloc(sizeof(PetscInt)*dim[0]*dim[1]*dim[2],&indx1);CHKERRQ(ierr);
1184e81bb599SAmlan Barua          for (i=0;i<dim[0];i++){
1185e81bb599SAmlan Barua             for (j=0;j<dim[1];j++){
1186e81bb599SAmlan Barua                for (k=0;k<dim[2];k++){
1187e81bb599SAmlan Barua                   indx1[i*dim[1]*dim[2]+j*dim[2]+k] = i*dim[1]*dim[2] + j*dim[2] + k;
1188e81bb599SAmlan Barua                }
1189e81bb599SAmlan Barua             }
1190e81bb599SAmlan Barua          }
1191e81bb599SAmlan Barua          ierr = ISCreateGeneral(comm,dim[0]*dim[1]*dim[2],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
1192e81bb599SAmlan Barua          ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
1193e81bb599SAmlan Barua          ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1194e81bb599SAmlan Barua          ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1195e81bb599SAmlan Barua          ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1196b223cc97SAmlan Barua          ierr = ISDestroy(&list1);CHKERRQ(ierr);
1197b223cc97SAmlan Barua          ierr = PetscFree(indx1);CHKERRQ(ierr);
1198e81bb599SAmlan Barua          break;
1199e81bb599SAmlan Barua     default:
12007536937bSAmlan Barua */
12016971673cSAmlan Barua          ierr = ISCreateStride(comm,N,0,1,&list1);
12026971673cSAmlan Barua          //ierr = ISView(list1,PETSC_VIEWER_STDOUT_SELF);
12036971673cSAmlan Barua          ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
12046971673cSAmlan Barua          ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12056971673cSAmlan Barua          ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12066971673cSAmlan Barua          ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1207b223cc97SAmlan Barua          ierr = ISDestroy(&list1);CHKERRQ(ierr);
12087536937bSAmlan Barua   //       break;
12097536937bSAmlan Barua    // }
1210e81bb599SAmlan Barua   }
1211e81bb599SAmlan Barua   else{
1212e81bb599SAmlan Barua 
1213e81bb599SAmlan Barua  switch (ndim){
1214e81bb599SAmlan Barua  case 1:
121564657d84SAmlan Barua #if defined(PETSC_USE_COMPLEX)
121664657d84SAmlan Barua       alloc_local = fftw_mpi_local_size_1d(dim[0],comm,FFTW_BACKWARD,FFTW_ESTIMATE,&local_n0,&local_0_start,&local_n1,&local_1_start);
121764657d84SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0,local_0_start,1,&list1);
121864657d84SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0,low,1,&list2);
121964657d84SAmlan Barua       //ierr = ISView(list1,PETSC_VIEWER_STDOUT_WORLD);
122064657d84SAmlan Barua       //ierr = ISView(list2,PETSC_VIEWER_STDOUT_WORLD);
122164657d84SAmlan Barua       ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr);
122264657d84SAmlan Barua       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
122364657d84SAmlan Barua       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
122464657d84SAmlan Barua       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
122564657d84SAmlan Barua       ierr = ISDestroy(&list1);CHKERRQ(ierr);
122664657d84SAmlan Barua       ierr = ISDestroy(&list2);CHKERRQ(ierr);
122764657d84SAmlan Barua       break;
122864657d84SAmlan Barua 
122964657d84SAmlan Barua #else
12306a506ed8SAmlan Barua       SETERRQ(comm,PETSC_ERR_SUP,"No support for real parallel 1D FFT");
123164657d84SAmlan Barua #endif
12325b3e41ffSAmlan Barua   break;
12335b3e41ffSAmlan Barua  case 2:
1234bd59e6c6SAmlan Barua #if defined(PETSC_USE_COMPLEX)
1235bd59e6c6SAmlan Barua       alloc_local = fftw_mpi_local_size_2d(dim[0],dim[1],comm,&local_n0,&local_0_start);
1236bd59e6c6SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1],local_0_start*dim[1],1,&list1);
1237bd59e6c6SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1],low,1,&list2);
1238bd59e6c6SAmlan Barua       ierr = VecScatterCreate(x,list2,y,list1,&vecscat);CHKERRQ(ierr);
1239bd59e6c6SAmlan Barua       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1240bd59e6c6SAmlan Barua       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1241bd59e6c6SAmlan Barua       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1242bd59e6c6SAmlan Barua       ierr = ISDestroy(&list1);CHKERRQ(ierr);
1243bd59e6c6SAmlan Barua       ierr = ISDestroy(&list2);CHKERRQ(ierr);
1244bd59e6c6SAmlan Barua       break;
1245bd59e6c6SAmlan Barua #else
12465b3e41ffSAmlan Barua       alloc_local =  fftw_mpi_local_size_2d_transposed(dim[0],dim[1]/2+1,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
12475b3e41ffSAmlan Barua       N1 = 2*dim[0]*(dim[1]/2+1); n1 = 2*local_n0*(dim[1]/2+1);
12485b3e41ffSAmlan Barua 
12495b3e41ffSAmlan Barua      ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx1);CHKERRQ(ierr);
12505b3e41ffSAmlan Barua      ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx2);CHKERRQ(ierr);
12515b3e41ffSAmlan Barua      printf("Val local_0_start = %ld\n",local_0_start);
12525b3e41ffSAmlan Barua 
12535b3e41ffSAmlan Barua      if (dim[1]%2==0)
12545b3e41ffSAmlan Barua       NM = dim[1]+2;
12555b3e41ffSAmlan Barua     else
12565b3e41ffSAmlan Barua       NM = dim[1]+1;
12575b3e41ffSAmlan Barua 
12585b3e41ffSAmlan Barua 
12595b3e41ffSAmlan Barua 
12605b3e41ffSAmlan Barua      for (i=0;i<local_n0;i++){
12615b3e41ffSAmlan Barua         for (j=0;j<dim[1];j++){
12625b3e41ffSAmlan Barua             tempindx = i*dim[1] + j;
12635b3e41ffSAmlan Barua             tempindx1 = i*NM + j;
12645b3e41ffSAmlan Barua             indx1[tempindx]=local_0_start*dim[1]+tempindx;
12655b3e41ffSAmlan Barua             indx2[tempindx]=low+tempindx1;
1266cfe1ae98SAmlan Barua        //     printf("Val tempindx1 = %d\n",tempindx1);
1267cfe1ae98SAmlan Barua        //     printf("index1 %d from proc %d is \n",indx1[tempindx],rank);
1268cfe1ae98SAmlan Barua        //     printf("index2 %d from proc %d is \n",indx2[tempindx],rank);
1269cfe1ae98SAmlan Barua        //     printf("-------------------------\n",indx2[tempindx],rank);
12705b3e41ffSAmlan Barua         }
12715b3e41ffSAmlan Barua      }
12725b3e41ffSAmlan Barua 
12735b3e41ffSAmlan Barua      ierr = ISCreateGeneral(comm,local_n0*dim[1],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
12745b3e41ffSAmlan Barua      ierr = ISCreateGeneral(comm,local_n0*dim[1],indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr);
12755b3e41ffSAmlan Barua 
12765b3e41ffSAmlan Barua      ierr = VecScatterCreate(x,list2,y,list1,&vecscat);CHKERRQ(ierr);
12775b3e41ffSAmlan Barua      ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12785b3e41ffSAmlan Barua      ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
12795b3e41ffSAmlan Barua      ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1280b223cc97SAmlan Barua      ierr = ISDestroy(&list1);CHKERRQ(ierr);
1281b223cc97SAmlan Barua      ierr = ISDestroy(&list2);CHKERRQ(ierr);
1282b223cc97SAmlan Barua      ierr = PetscFree(indx1);CHKERRQ(ierr);
1283b223cc97SAmlan Barua      ierr = PetscFree(indx2);CHKERRQ(ierr);
12845b3e41ffSAmlan Barua   break;
1285bd59e6c6SAmlan Barua #endif
12865b3e41ffSAmlan Barua  case 3:
1287bd59e6c6SAmlan Barua #if defined(PETSC_USE_COMPLEX)
1288bd59e6c6SAmlan Barua       alloc_local = fftw_mpi_local_size_3d(dim[0],dim[1],dim[2],comm,&local_n0,&local_0_start);
1289bd59e6c6SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1]*dim[2],local_0_start*dim[1]*dim[2],1,&list1);
1290bd59e6c6SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1]*dim[2],low,1,&list2);
1291bd59e6c6SAmlan Barua       //ierr = ISView(list1,PETSC_VIEWER_STDOUT_WORLD);
1292bd59e6c6SAmlan Barua       //ierr = ISView(list2,PETSC_VIEWER_STDOUT_WORLD);
1293bd59e6c6SAmlan Barua       ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr);
1294bd59e6c6SAmlan Barua       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1295bd59e6c6SAmlan Barua       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1296bd59e6c6SAmlan Barua       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1297bd59e6c6SAmlan Barua       ierr = ISDestroy(&list1);CHKERRQ(ierr);
1298bd59e6c6SAmlan Barua       ierr = ISDestroy(&list2);CHKERRQ(ierr);
1299bd59e6c6SAmlan Barua       break;
1300bd59e6c6SAmlan Barua #else
1301bd59e6c6SAmlan Barua 
130251d1eed7SAmlan Barua       alloc_local =  fftw_mpi_local_size_3d_transposed(dim[0],dim[1],dim[2]/2+1,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
130351d1eed7SAmlan Barua       N1 = 2*dim[0]*dim[1]*(dim[2]/2+1); n1 = 2*local_n0*dim[1]*(dim[2]/2+1);
130451d1eed7SAmlan Barua 
130551d1eed7SAmlan Barua      ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx1);CHKERRQ(ierr);
130651d1eed7SAmlan Barua      ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx2);CHKERRQ(ierr);
130751d1eed7SAmlan Barua      printf("Val local_0_start = %ld\n",local_0_start);
130851d1eed7SAmlan Barua 
130951d1eed7SAmlan Barua      if (dim[2]%2==0)
131051d1eed7SAmlan Barua       NM = dim[2]+2;
131151d1eed7SAmlan Barua     else
131251d1eed7SAmlan Barua       NM = dim[2]+1;
131351d1eed7SAmlan Barua 
131451d1eed7SAmlan Barua      for (i=0;i<local_n0;i++){
131551d1eed7SAmlan Barua         for (j=0;j<dim[1];j++){
131651d1eed7SAmlan Barua            for (k=0;k<dim[2];k++){
131751d1eed7SAmlan Barua               tempindx = i*dim[1]*dim[2] + j*dim[2] + k;
131851d1eed7SAmlan Barua               tempindx1 = i*dim[1]*NM + j*NM + k;
131951d1eed7SAmlan Barua               indx1[tempindx]=local_0_start*dim[1]*dim[2]+tempindx;
132051d1eed7SAmlan Barua               indx2[tempindx]=low+tempindx1;
132151d1eed7SAmlan Barua            }
132251d1eed7SAmlan Barua         //    printf("Val tempindx1 = %d\n",tempindx1);
132351d1eed7SAmlan Barua         //    printf("index1 %d from proc %d is \n",indx1[tempindx],rank);
132451d1eed7SAmlan Barua         //    printf("index2 %d from proc %d is \n",indx2[tempindx],rank);
132551d1eed7SAmlan Barua         //    printf("-------------------------\n",indx2[tempindx],rank);
132651d1eed7SAmlan Barua         }
132751d1eed7SAmlan Barua      }
132851d1eed7SAmlan Barua 
132951d1eed7SAmlan Barua      ierr = ISCreateGeneral(comm,local_n0*dim[1]*dim[2],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
133051d1eed7SAmlan Barua      ierr = ISCreateGeneral(comm,local_n0*dim[1]*dim[2],indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr);
133151d1eed7SAmlan Barua 
133251d1eed7SAmlan Barua      ierr = VecScatterCreate(x,list2,y,list1,&vecscat);CHKERRQ(ierr);
133351d1eed7SAmlan Barua      ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
133451d1eed7SAmlan Barua      ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
133551d1eed7SAmlan Barua      ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1336b223cc97SAmlan Barua      ierr = ISDestroy(&list1);CHKERRQ(ierr);
1337b223cc97SAmlan Barua      ierr = ISDestroy(&list2);CHKERRQ(ierr);
1338b223cc97SAmlan Barua      ierr = PetscFree(indx1);CHKERRQ(ierr);
1339b223cc97SAmlan Barua      ierr = PetscFree(indx2);CHKERRQ(ierr);
13405b3e41ffSAmlan Barua   break;
1341bd59e6c6SAmlan Barua #endif
13425b3e41ffSAmlan Barua   default:
1343bd59e6c6SAmlan Barua #if defined(PETSC_USE_COMPLEX)
1344bd59e6c6SAmlan Barua       alloc_local = fftw_mpi_local_size(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start);
1345bd59e6c6SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*(fftw->partial_dim),local_0_start*(fftw->partial_dim),1,&list1);
1346bd59e6c6SAmlan Barua       ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*(fftw->partial_dim),low,1,&list2);
1347bd59e6c6SAmlan Barua       //ierr = ISView(list1,PETSC_VIEWER_STDOUT_WORLD);
1348bd59e6c6SAmlan Barua       //ierr = ISView(list2,PETSC_VIEWER_STDOUT_WORLD);
1349bd59e6c6SAmlan Barua       ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr);
1350bd59e6c6SAmlan Barua       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1351bd59e6c6SAmlan Barua       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1352bd59e6c6SAmlan Barua       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1353bd59e6c6SAmlan Barua       ierr = ISDestroy(&list1);CHKERRQ(ierr);
1354bd59e6c6SAmlan Barua       ierr = ISDestroy(&list2);CHKERRQ(ierr);
1355bd59e6c6SAmlan Barua #else
1356ba6e06d0SAmlan Barua      temp = (fftw->dim_fftw)[fftw->ndim_fftw-1];
1357ba6e06d0SAmlan Barua      printf("The value of temp is %ld\n",temp);
1358ba6e06d0SAmlan Barua      (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp/2 + 1;
1359ba6e06d0SAmlan Barua      alloc_local = fftw_mpi_local_size_transposed(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
1360ba6e06d0SAmlan Barua      N1 = 2*N*(PetscInt)((fftw->dim_fftw)[fftw->ndim_fftw-1])/((PetscInt) temp);
1361ba6e06d0SAmlan Barua      (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp;
1362ba6e06d0SAmlan Barua 
1363ba6e06d0SAmlan Barua      partial_dim = fftw->partial_dim;
1364ba6e06d0SAmlan Barua      printf("The value of partial dim is %d\n",partial_dim);
1365ba6e06d0SAmlan Barua 
1366ba6e06d0SAmlan Barua      ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*partial_dim,&indx1);CHKERRQ(ierr);
1367ba6e06d0SAmlan Barua      ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*partial_dim,&indx2);CHKERRQ(ierr);
1368ba6e06d0SAmlan Barua      printf("Val local_0_start = %ld\n",local_0_start);
1369ba6e06d0SAmlan Barua 
1370ba6e06d0SAmlan Barua      if (dim[ndim-1]%2==0)
1371ba6e06d0SAmlan Barua        NM = 2;
1372ba6e06d0SAmlan Barua      else
1373ba6e06d0SAmlan Barua        NM = 1;
1374ba6e06d0SAmlan Barua 
1375ba6e06d0SAmlan Barua      j = low;
1376ba6e06d0SAmlan Barua      for (i=0,k=1; i<((PetscInt)local_n0)*partial_dim;i++,k++)
1377ba6e06d0SAmlan Barua         {
1378ba6e06d0SAmlan Barua          indx1[i] = local_0_start*partial_dim + i;
1379ba6e06d0SAmlan Barua          indx2[i] = j;
1380ba6e06d0SAmlan Barua          //printf("The values are %d and %d from %d\n",indx1[i],indx2[i],rank);
1381ba6e06d0SAmlan Barua          if (k%dim[ndim-1]==0)
1382ba6e06d0SAmlan Barua            { j+=NM;}
1383ba6e06d0SAmlan Barua          j++;
1384ba6e06d0SAmlan Barua         }
1385ba6e06d0SAmlan Barua      ierr = ISCreateGeneral(comm,local_n0*partial_dim,indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
1386ba6e06d0SAmlan Barua      ierr = ISCreateGeneral(comm,local_n0*partial_dim,indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr);
1387ba6e06d0SAmlan Barua 
1388ba6e06d0SAmlan Barua       //ISView(list1,PETSC_VIEWER_STDOUT_SELF);
1389ba6e06d0SAmlan Barua 
1390ba6e06d0SAmlan Barua 
1391ba6e06d0SAmlan Barua      ierr = VecScatterCreate(x,list2,y,list1,&vecscat);CHKERRQ(ierr);
1392ba6e06d0SAmlan Barua      ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1393ba6e06d0SAmlan Barua      ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1394ba6e06d0SAmlan Barua      ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1395b223cc97SAmlan Barua      ierr = ISDestroy(&list1);CHKERRQ(ierr);
1396b223cc97SAmlan Barua      ierr = ISDestroy(&list2);CHKERRQ(ierr);
1397b223cc97SAmlan Barua      ierr = PetscFree(indx1);CHKERRQ(ierr);
1398b223cc97SAmlan Barua      ierr = PetscFree(indx2);CHKERRQ(ierr);
1399ba6e06d0SAmlan Barua 
14005b3e41ffSAmlan Barua      break;
1401bd59e6c6SAmlan Barua #endif
14025b3e41ffSAmlan Barua  }
1403e81bb599SAmlan Barua  }
14045b3e41ffSAmlan Barua  return 0;
14055b3e41ffSAmlan Barua }
14065b3e41ffSAmlan Barua EXTERN_C_END
14075b3e41ffSAmlan Barua 
14085b3e41ffSAmlan Barua EXTERN_C_BEGIN
14095b3e41ffSAmlan Barua #undef __FUNCT__
14103c3a9c75SAmlan Barua #define __FUNCT__ "MatCreate_FFTW"
14113c3a9c75SAmlan Barua /*
14123c3a9c75SAmlan Barua       MatCreate_FFTW - Creates a matrix object that provides FFT
14133c3a9c75SAmlan Barua   via the external package FFTW
14143c3a9c75SAmlan Barua   Options Database Keys:
14153c3a9c75SAmlan Barua + -mat_fftw_plannerflags - set FFTW planner flags
14163c3a9c75SAmlan Barua 
14173c3a9c75SAmlan Barua    Level: intermediate
14183c3a9c75SAmlan Barua 
14193c3a9c75SAmlan Barua */
14203c3a9c75SAmlan Barua 
1421b2b77a04SHong Zhang PetscErrorCode MatCreate_FFTW(Mat A)
1422b2b77a04SHong Zhang {
1423b2b77a04SHong Zhang   PetscErrorCode ierr;
1424b2b77a04SHong Zhang   MPI_Comm       comm=((PetscObject)A)->comm;
1425b2b77a04SHong Zhang   Mat_FFT        *fft=(Mat_FFT*)A->data;
1426b2b77a04SHong Zhang   Mat_FFTW       *fftw;
1427b2b77a04SHong Zhang   PetscInt       n=fft->n,N=fft->N,ndim=fft->ndim,*dim = fft->dim;
1428b2b77a04SHong Zhang   const char     *p_flags[]={"FFTW_ESTIMATE","FFTW_MEASURE","FFTW_PATIENT","FFTW_EXHAUSTIVE"};
1429b2b77a04SHong Zhang   PetscBool      flg;
1430b3a17365SAmlan Barua   PetscInt       p_flag,partial_dim=1,ctr,N1;
143111902ff2SHong Zhang   PetscMPIInt    size,rank;
1432b3a17365SAmlan Barua   ptrdiff_t      *pdim, temp;
14334a16a297SSean Farley   ptrdiff_t      local_n1,local_1_start,local_1_end;
1434b2b77a04SHong Zhang 
1435b2b77a04SHong Zhang   PetscFunctionBegin;
1436b2b77a04SHong Zhang   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
143711902ff2SHong Zhang   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
143854dd5118SAmlan Barua   ierr = MPI_Barrier(PETSC_COMM_WORLD);
1439c92418ccSAmlan Barua 
144011902ff2SHong Zhang   pdim = (ptrdiff_t *)calloc(ndim,sizeof(ptrdiff_t));
144111902ff2SHong Zhang   pdim[0] = dim[0];
144211902ff2SHong Zhang   for(ctr=1;ctr<ndim;ctr++)
144311902ff2SHong Zhang       {
14446882372aSHong Zhang           partial_dim *= dim[ctr];
144511902ff2SHong Zhang           pdim[ctr] = dim[ctr];
14466882372aSHong Zhang       }
14473c3a9c75SAmlan Barua 
1448b2b77a04SHong Zhang   if (size == 1) {
1449e81bb599SAmlan Barua #if defined(PETSC_USE_COMPLEX)
1450b2b77a04SHong Zhang     ierr = MatSetSizes(A,N,N,N,N);CHKERRQ(ierr);
1451b2b77a04SHong Zhang     n = N;
1452e81bb599SAmlan Barua #else
1453e81bb599SAmlan Barua     int tot_dim = N*2*(dim[ndim-1]/2+1)/dim[ndim-1];
1454e81bb599SAmlan Barua     ierr = MatSetSizes(A,tot_dim,tot_dim,tot_dim,tot_dim);CHKERRQ(ierr);
1455e81bb599SAmlan Barua     n = tot_dim;
1456e81bb599SAmlan Barua #endif
1457e81bb599SAmlan Barua 
1458b2b77a04SHong Zhang   } else {
1459b223cc97SAmlan Barua     ptrdiff_t alloc_local,local_n0,local_0_start;//local_n1,local_1_end;
1460b2b77a04SHong Zhang     switch (ndim){
1461b2b77a04SHong Zhang     case 1:
14623c3a9c75SAmlan Barua #if !defined(PETSC_USE_COMPLEX)
14633c3a9c75SAmlan Barua   SETERRQ(comm,PETSC_ERR_SUP,"FFTW does not support parallel 1D real transform");
1464e5d7f247SAmlan Barua #else
14656882372aSHong Zhang       alloc_local = fftw_mpi_local_size_1d(dim[0],comm,FFTW_FORWARD,FFTW_ESTIMATE,&local_n0,&local_0_start,&local_n1,&local_1_end);
14666882372aSHong Zhang       n = (PetscInt)local_n0;
14676882372aSHong Zhang       ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr);
1468e5d7f247SAmlan Barua #endif
14693c3a9c75SAmlan Barua //      PetscObjectComposeFunctionDynamic((PetscObject)A,"MatGetVecs1DC_C","MatGetVecs1DC_FFTW",MatGetVecs1DC_FFTW);
1470b2b77a04SHong Zhang       break;
1471b2b77a04SHong Zhang     case 2:
14725b3e41ffSAmlan Barua #if defined(PETSC_USE_COMPLEX)
1473b2b77a04SHong Zhang       alloc_local = fftw_mpi_local_size_2d(dim[0],dim[1],comm,&local_n0,&local_0_start);
1474b2b77a04SHong Zhang       /*
1475b2b77a04SHong Zhang        PetscMPIInt    rank;
1476b2b77a04SHong Zhang        PetscSynchronizedPrintf(comm,"[%d] MatCreateSeqFFTW: local_n0, local_0_start %d %d, N %d,dim %d, %d\n",rank,(PetscInt)local_n0*dim[1],(PetscInt)local_0_start,m,dim[0],dim[1]);
1477b2b77a04SHong Zhang        PetscSynchronizedFlush(comm);
1478b2b77a04SHong Zhang        */
1479b2b77a04SHong Zhang       n = (PetscInt)local_n0*dim[1];
1480b2b77a04SHong Zhang       ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr);
14815b3e41ffSAmlan Barua #else
14825b3e41ffSAmlan Barua       alloc_local = fftw_mpi_local_size_2d_transposed(dim[0],dim[1]/2+1,PETSC_COMM_WORLD,&local_n0,&local_0_start,&local_n1,&local_1_start);
14835b3e41ffSAmlan Barua       n = 2*(PetscInt)local_n0*(dim[1]/2+1);
14845b3e41ffSAmlan Barua       ierr = MatSetSizes(A,n,n,2*dim[0]*(dim[1]/2+1),2*dim[0]*(dim[1]/2+1));
14855b3e41ffSAmlan Barua #endif
1486b2b77a04SHong Zhang       break;
1487b2b77a04SHong Zhang     case 3:
148811902ff2SHong Zhang //      printf("The value of alloc local is %d",alloc_local);
148951d1eed7SAmlan Barua #if defined(PETSC_USE_COMPLEX)
149051d1eed7SAmlan Barua       alloc_local = fftw_mpi_local_size_3d(dim[0],dim[1],dim[2],comm,&local_n0,&local_0_start);
14916882372aSHong Zhang       n = (PetscInt)local_n0*dim[1]*dim[2];
14926882372aSHong Zhang       ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr);
149351d1eed7SAmlan Barua #else
149451d1eed7SAmlan Barua       printf("The code comes here\n");
149551d1eed7SAmlan Barua       alloc_local = fftw_mpi_local_size_3d_transposed(dim[0],dim[1],dim[2]/2+1,PETSC_COMM_WORLD,&local_n0,&local_0_start,&local_n1,&local_1_start);
149651d1eed7SAmlan Barua       n = 2*(PetscInt)local_n0*dim[1]*(dim[2]/2+1);
149751d1eed7SAmlan Barua       ierr = MatSetSizes(A,n,n,2*dim[0]*dim[1]*(dim[2]/2+1),2*dim[0]*dim[1]*(dim[2]/2+1));
149851d1eed7SAmlan Barua #endif
1499b2b77a04SHong Zhang       break;
1500b2b77a04SHong Zhang     default:
1501b3a17365SAmlan Barua #if defined(PETSC_USE_COMPLEX)
150211902ff2SHong Zhang       alloc_local = fftw_mpi_local_size(ndim,pdim,comm,&local_n0,&local_0_start);
1503c92418ccSAmlan Barua //      printf("The value of alloc local is %ld from process %d\n",alloc_local,rank);
150411902ff2SHong Zhang //      alloc_local = fftw_mpi_local_size(ndim,dim,comm,&local_n0,&local_0_start);
15056882372aSHong Zhang       n = (PetscInt)local_n0*partial_dim;
150611902ff2SHong Zhang //      printf("New partial dimension is %d %d %d",n,N,ndim);
15076882372aSHong Zhang       ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr);
1508b3a17365SAmlan Barua #else
1509b3a17365SAmlan Barua       temp = pdim[ndim-1];
1510b3a17365SAmlan Barua       pdim[ndim-1]= temp/2 + 1;
1511b3a17365SAmlan Barua       printf("For Multi dim case temp = %ld, pdim[ndim-1] = %ld\n",temp,pdim[ndim-1]);
1512b3a17365SAmlan Barua       alloc_local = fftw_mpi_local_size_transposed(ndim,pdim,PETSC_COMM_WORLD,&local_n0,&local_0_start,&local_n1,&local_1_start);
1513b3a17365SAmlan Barua       n = 2*(PetscInt)local_n0*partial_dim*pdim[ndim-1]/temp;
1514b3a17365SAmlan Barua       N1 = 2*N*(PetscInt)pdim[ndim-1]/((PetscInt) temp);
1515b3a17365SAmlan Barua       pdim[ndim-1] = temp;
1516b3a17365SAmlan Barua       printf("For Multi dim case n = %d, N1  = %d\n",n,N1);
1517b3a17365SAmlan Barua       ierr = MatSetSizes(A,n,n,N1,N1);CHKERRQ(ierr);
1518b3a17365SAmlan Barua #endif
1519b2b77a04SHong Zhang       break;
1520b2b77a04SHong Zhang     }
1521b2b77a04SHong Zhang   }
1522b2b77a04SHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,MATFFTW);CHKERRQ(ierr);
1523b2b77a04SHong Zhang   ierr = PetscNewLog(A,Mat_FFTW,&fftw);CHKERRQ(ierr);
1524b2b77a04SHong Zhang   fft->data = (void*)fftw;
1525b2b77a04SHong Zhang 
1526b2b77a04SHong Zhang   fft->n           = n;
1527c92418ccSAmlan Barua   fftw->ndim_fftw  = (ptrdiff_t)ndim; // This is dimension of fft
1528e5d7f247SAmlan Barua   fftw->partial_dim = partial_dim;
1529c92418ccSAmlan Barua   ierr = PetscMalloc(ndim*sizeof(ptrdiff_t), (ptrdiff_t *)&(fftw->dim_fftw));CHKERRQ(ierr);
1530b1301b2eSHong Zhang   for(ctr=0;ctr<ndim;ctr++) (fftw->dim_fftw)[ctr]=dim[ctr];
1531c92418ccSAmlan Barua 
1532b2b77a04SHong Zhang   fftw->p_forward  = 0;
1533b2b77a04SHong Zhang   fftw->p_backward = 0;
1534b2b77a04SHong Zhang   fftw->p_flag     = FFTW_ESTIMATE;
1535b2b77a04SHong Zhang 
1536b2b77a04SHong Zhang   if (size == 1){
1537b2b77a04SHong Zhang     A->ops->mult          = MatMult_SeqFFTW;
1538b2b77a04SHong Zhang     A->ops->multtranspose = MatMultTranspose_SeqFFTW;
1539b2b77a04SHong Zhang   } else {
1540b2b77a04SHong Zhang     A->ops->mult          = MatMult_MPIFFTW;
1541b2b77a04SHong Zhang     A->ops->multtranspose = MatMultTranspose_MPIFFTW;
1542b2b77a04SHong Zhang   }
1543b2b77a04SHong Zhang   fft->matdestroy          = MatDestroy_FFTW;
15446a506ed8SAmlan Barua // if(ndim=1 && size>1) and also if it is complex then getvecs should be attached to MatGetVecs_FFTW1D
1545b2b77a04SHong Zhang   A->ops->getvecs       = MatGetVecs_FFTW;
1546b2b77a04SHong Zhang   A->assembled          = PETSC_TRUE;
15474f7415efSAmlan Barua   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatGetVecsFFT_C","MatGetVecsFFT_FFTW",MatGetVecsFFT_FFTW);
15483c3a9c75SAmlan Barua   PetscObjectComposeFunctionDynamic((PetscObject)A,"InputTransformFFT_C","InputTransformFFT_FFTW",InputTransformFFT_FFTW);
15495b3e41ffSAmlan Barua   PetscObjectComposeFunctionDynamic((PetscObject)A,"OutputTransformFFT_C","OutputTransformFFT_FFTW",OutputTransformFFT_FFTW);
1550b2b77a04SHong Zhang 
1551b2b77a04SHong Zhang   /* get runtime options */
1552b2b77a04SHong Zhang   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"FFTW Options","Mat");CHKERRQ(ierr);
1553b2b77a04SHong Zhang     ierr = PetscOptionsEList("-mat_fftw_plannerflags","Planner Flags","None",p_flags,4,p_flags[0],&p_flag,&flg);CHKERRQ(ierr);
1554b2b77a04SHong Zhang     if (flg) {fftw->p_flag = (unsigned)p_flag;}
1555b2b77a04SHong Zhang   PetscOptionsEnd();
1556b2b77a04SHong Zhang   PetscFunctionReturn(0);
1557b2b77a04SHong Zhang }
1558b2b77a04SHong Zhang EXTERN_C_END
15593c3a9c75SAmlan Barua 
15603c3a9c75SAmlan Barua 
15613c3a9c75SAmlan Barua 
15623c3a9c75SAmlan Barua 
1563