xref: /petsc/src/mat/impls/fft/fftw/fftw.c (revision 0cfd297d8d96f598482e3bb17bf4e35d55570ac3)
1 
2 /*
3     Provides an interface to the FFTW package.
4     Testing examples can be found in ~src/mat/examples/tests
5 */
6 
7 #include <../src/mat/impls/fft/fft.h>   /*I "petscmat.h" I*/
8 EXTERN_C_BEGIN
9 #include <fftw3-mpi.h>
10 EXTERN_C_END
11 
12 typedef struct {
13   ptrdiff_t   ndim_fftw,*dim_fftw;
14   PetscInt   partial_dim;
15   fftw_plan   p_forward,p_backward;
16   unsigned    p_flag; /* planner flags, FFTW_ESTIMATE,FFTW_MEASURE, FFTW_PATIENT, FFTW_EXHAUSTIVE */
17   PetscScalar *finarray,*foutarray,*binarray,*boutarray; /* keep track of arrays becaue fftw plan should be
18                                                             executed for the arrays with which the plan was created */
19 } Mat_FFTW;
20 
21 extern PetscErrorCode MatMult_SeqFFTW(Mat,Vec,Vec);
22 extern PetscErrorCode MatMultTranspose_SeqFFTW(Mat,Vec,Vec);
23 extern PetscErrorCode MatMult_MPIFFTW(Mat,Vec,Vec);
24 extern PetscErrorCode MatMultTranspose_MPIFFTW(Mat,Vec,Vec);
25 extern PetscErrorCode MatDestroy_FFTW(Mat);
26 extern PetscErrorCode VecDestroy_MPIFFTW(Vec);
27 extern PetscErrorCode MatGetVecs_FFTW(Mat,Vec*,Vec*);
28 
29 #undef __FUNCT__
30 #define __FUNCT__ "MatMult_SeqFFTW"
31 PetscErrorCode MatMult_SeqFFTW(Mat A,Vec x,Vec y)
32 {
33   PetscErrorCode ierr;
34   Mat_FFT        *fft  = (Mat_FFT*)A->data;
35   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
36   PetscScalar    *x_array,*y_array;
37   PetscInt       ndim=fft->ndim,*dim=fft->dim;
38 
39   PetscFunctionBegin;
40   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
41   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
42   if (!fftw->p_forward){ /* create a plan, then excute it */
43     switch (ndim){
44     case 1:
45 #if defined(PETSC_USE_COMPLEX)
46       fftw->p_forward = fftw_plan_dft_1d(dim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_FORWARD,fftw->p_flag);
47 #else
48       fftw->p_forward = fftw_plan_dft_r2c_1d(dim[0],(double *)x_array,(fftw_complex*)y_array,fftw->p_flag);
49 #endif
50       break;
51     case 2:
52 #if defined(PETSC_USE_COMPLEX)
53       fftw->p_forward = fftw_plan_dft_2d(dim[0],dim[1],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_FORWARD,fftw->p_flag);
54 #else
55       fftw->p_forward = fftw_plan_dft_r2c_2d(dim[0],dim[1],(double *)x_array,(fftw_complex*)y_array,fftw->p_flag);
56 #endif
57       break;
58     case 3:
59 #if defined(PETSC_USE_COMPLEX)
60       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);
61 #else
62       fftw->p_forward = fftw_plan_dft_r2c_3d(dim[0],dim[1],dim[2],(double *)x_array,(fftw_complex*)y_array,fftw->p_flag);
63 #endif
64       break;
65     default:
66 #if defined(PETSC_USE_COMPLEX)
67       fftw->p_forward = fftw_plan_dft(ndim,dim,(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_FORWARD,fftw->p_flag);
68 #else
69       fftw->p_forward = fftw_plan_dft_r2c(ndim,dim,(double *)x_array,(fftw_complex*)y_array,fftw->p_flag);
70 #endif
71       break;
72     }
73     fftw->finarray  = x_array;
74     fftw->foutarray = y_array;
75     /* Warning: if (fftw->p_flag!==FFTW_ESTIMATE) The data in the in/out arrays is overwritten!
76                 planning should be done before x is initialized! See FFTW manual sec2.1 or sec4 */
77     fftw_execute(fftw->p_forward);
78   } else { /* use existing plan */
79     if (fftw->finarray != x_array || fftw->foutarray != y_array){ /* use existing plan on new arrays */
80       fftw_execute_dft(fftw->p_forward,(fftw_complex*)x_array,(fftw_complex*)y_array);
81     } else {
82       fftw_execute(fftw->p_forward);
83     }
84   }
85   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
86   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
87   PetscFunctionReturn(0);
88 }
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "MatMultTranspose_SeqFFTW"
92 PetscErrorCode MatMultTranspose_SeqFFTW(Mat A,Vec x,Vec y)
93 {
94   PetscErrorCode ierr;
95   Mat_FFT        *fft = (Mat_FFT*)A->data;
96   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
97   PetscScalar    *x_array,*y_array;
98   PetscInt       ndim=fft->ndim,*dim=fft->dim;
99 
100   PetscFunctionBegin;
101   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
102   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
103   if (!fftw->p_backward){ /* create a plan, then excute it */
104     switch (ndim){
105     case 1:
106 #if defined(PETSC_USE_COMPLEX)
107       fftw->p_backward = fftw_plan_dft_1d(dim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_BACKWARD,fftw->p_flag);
108 #else
109       fftw->p_backward= fftw_plan_dft_c2r_1d(dim[0],(fftw_complex*)x_array,(double *)y_array,fftw->p_flag);
110 #endif
111       break;
112     case 2:
113 #if defined(PETSC_USE_COMPLEX)
114       fftw->p_backward = fftw_plan_dft_2d(dim[0],dim[1],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_BACKWARD,fftw->p_flag);
115 #else
116       fftw->p_backward= fftw_plan_dft_c2r_2d(dim[0],dim[1],(fftw_complex*)x_array,(double *)y_array,fftw->p_flag);
117 #endif
118       break;
119     case 3:
120 #if defined(PETSC_USE_COMPLEX)
121       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);
122 #else
123       fftw->p_backward= fftw_plan_dft_c2r_3d(dim[0],dim[1],dim[2],(fftw_complex*)x_array,(double *)y_array,fftw->p_flag);
124 #endif
125       break;
126     default:
127 #if defined(PETSC_USE_COMPLEX)
128       fftw->p_backward = fftw_plan_dft(ndim,dim,(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_BACKWARD,fftw->p_flag);
129 #else
130       fftw->p_backward= fftw_plan_dft_c2r(ndim,dim,(fftw_complex*)x_array,(double *)y_array,fftw->p_flag);
131 #endif
132       break;
133     }
134     fftw->binarray  = x_array;
135     fftw->boutarray = y_array;
136     fftw_execute(fftw->p_backward);CHKERRQ(ierr);
137   } else { /* use existing plan */
138     if (fftw->binarray != x_array || fftw->boutarray != y_array){ /* use existing plan on new arrays */
139       fftw_execute_dft(fftw->p_backward,(fftw_complex*)x_array,(fftw_complex*)y_array);
140     } else {
141       fftw_execute(fftw->p_backward);CHKERRQ(ierr);
142     }
143   }
144   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
145   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
146   PetscFunctionReturn(0);
147 }
148 
149 #undef __FUNCT__
150 #define __FUNCT__ "MatMult_MPIFFTW"
151 PetscErrorCode MatMult_MPIFFTW(Mat A,Vec x,Vec y)
152 {
153   PetscErrorCode ierr;
154   Mat_FFT        *fft  = (Mat_FFT*)A->data;
155   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
156   PetscScalar    *x_array,*y_array;
157   PetscInt       ndim=fft->ndim,*dim=fft->dim;
158   MPI_Comm       comm=((PetscObject)A)->comm;
159 
160   PetscFunctionBegin;
161   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
162   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
163   if (!fftw->p_forward){ /* create a plan, then excute it */
164     switch (ndim){
165     case 1:
166 #if defined(PETSC_USE_COMPLEX)
167       fftw->p_forward = fftw_mpi_plan_dft_1d(dim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_FORWARD,fftw->p_flag);
168 #else
169       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers yet");
170 #endif
171       break;
172     case 2:
173 #if defined(PETSC_USE_COMPLEX)
174       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);
175 #else
176       fftw->p_forward = fftw_mpi_plan_dft_r2c_2d(dim[0],dim[1],(double *)x_array,(fftw_complex*)y_array,comm,FFTW_ESTIMATE);
177 #endif
178       break;
179     case 3:
180 #if defined(PETSC_USE_COMPLEX)
181       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);
182 #else
183       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);
184 #endif
185       break;
186     default:
187 #if defined(PETSC_USE_COMPLEX)
188       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);
189 #else
190       fftw->p_forward = fftw_mpi_plan_dft_r2c(fftw->ndim_fftw,fftw->dim_fftw,(double *)x_array,(fftw_complex*)y_array,comm,FFTW_ESTIMATE);
191 #endif
192       break;
193     }
194     fftw->finarray  = x_array;
195     fftw->foutarray = y_array;
196     /* Warning: if (fftw->p_flag!==FFTW_ESTIMATE) The data in the in/out arrays is overwritten!
197                 planning should be done before x is initialized! See FFTW manual sec2.1 or sec4 */
198     fftw_execute(fftw->p_forward);
199   } else { /* use existing plan */
200     if (fftw->finarray != x_array || fftw->foutarray != y_array){ /* use existing plan on new arrays */
201       fftw_execute_dft(fftw->p_forward,(fftw_complex*)x_array,(fftw_complex*)y_array);
202     } else {
203       fftw_execute(fftw->p_forward);
204     }
205   }
206   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
207   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "MatMultTranspose_MPIFFTW"
213 PetscErrorCode MatMultTranspose_MPIFFTW(Mat A,Vec x,Vec y)
214 {
215   PetscErrorCode ierr;
216   Mat_FFT        *fft  = (Mat_FFT*)A->data;
217   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
218   PetscScalar    *x_array,*y_array;
219   PetscInt       ndim=fft->ndim,*dim=fft->dim;
220   MPI_Comm       comm=((PetscObject)A)->comm;
221 
222   PetscFunctionBegin;
223   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
224   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
225   if (!fftw->p_backward){ /* create a plan, then excute it */
226     switch (ndim){
227     case 1:
228 #if defined(PETSC_USE_COMPLEX)
229       fftw->p_backward = fftw_mpi_plan_dft_1d(dim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_BACKWARD,fftw->p_flag);
230 #else
231       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers yet");
232 #endif
233       break;
234     case 2:
235 #if defined(PETSC_USE_COMPLEX)
236       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);
237 #else
238       fftw->p_backward = fftw_mpi_plan_dft_c2r_2d(dim[0],dim[1],(fftw_complex*)x_array,(double *)y_array,comm,FFTW_ESTIMATE);
239 #endif
240       break;
241     case 3:
242 #if defined(PETSC_USE_COMPLEX)
243       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);
244 #else
245       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);
246 #endif
247       break;
248     default:
249 #if defined(PETSC_USE_COMPLEX)
250       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);
251 #else
252       fftw->p_backward = fftw_mpi_plan_dft_c2r(fftw->ndim_fftw,fftw->dim_fftw,(fftw_complex*)x_array,(double *)y_array,comm,FFTW_ESTIMATE);
253 #endif
254       break;
255     }
256     fftw->binarray  = x_array;
257     fftw->boutarray = y_array;
258     fftw_execute(fftw->p_backward);CHKERRQ(ierr);
259   } else { /* use existing plan */
260     if (fftw->binarray != x_array || fftw->boutarray != y_array){ /* use existing plan on new arrays */
261       fftw_execute_dft(fftw->p_backward,(fftw_complex*)x_array,(fftw_complex*)y_array);
262     } else {
263       fftw_execute(fftw->p_backward);CHKERRQ(ierr);
264     }
265   }
266   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
267   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
268   PetscFunctionReturn(0);
269 }
270 
271 #undef __FUNCT__
272 #define __FUNCT__ "MatDestroy_FFTW"
273 PetscErrorCode MatDestroy_FFTW(Mat A)
274 {
275   Mat_FFT        *fft = (Mat_FFT*)A->data;
276   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
277   PetscErrorCode ierr;
278 
279   PetscFunctionBegin;
280   fftw_destroy_plan(fftw->p_forward);
281   fftw_destroy_plan(fftw->p_backward);
282   ierr = PetscFree(fftw->dim_fftw);CHKERRQ(ierr);
283   ierr = PetscFree(fft->data);CHKERRQ(ierr);
284   PetscFunctionReturn(0);
285 }
286 
287 #include <../src/vec/vec/impls/mpi/pvecimpl.h>   /*I  "petscvec.h"   I*/
288 #undef __FUNCT__
289 #define __FUNCT__ "VecDestroy_MPIFFTW"
290 PetscErrorCode VecDestroy_MPIFFTW(Vec v)
291 {
292   PetscErrorCode  ierr;
293   PetscScalar     *array;
294 
295   PetscFunctionBegin;
296   ierr = VecGetArray(v,&array);CHKERRQ(ierr);
297   fftw_free((fftw_complex*)array);CHKERRQ(ierr);
298   ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
299   ierr = VecDestroy_MPI(v);CHKERRQ(ierr);
300   PetscFunctionReturn(0);
301 }
302 
303 #undef __FUNCT__
304 #define __FUNCT__ "MatGetVecs1DC_FFTW"
305 /*
306    MatGetVecs_FFTW1D - Get Vectors(s) compatible with matrix, i.e. with the
307      parallel layout determined by FFTW-1D
308 
309 */
310 PetscErrorCode  MatGetVecs_FFTW1D(Mat A,Vec *fin,Vec *fout,Vec *bout)
311 {
312   PetscErrorCode ierr;
313   PetscMPIInt    size,rank;
314   MPI_Comm       comm=((PetscObject)A)->comm;
315   Mat_FFT        *fft = (Mat_FFT*)A->data;
316 //  Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
317   PetscInt       N=fft->N;
318   PetscInt       ndim=fft->ndim,*dim=fft->dim;
319   ptrdiff_t      f_alloc_local,f_local_n0,f_local_0_start;
320   ptrdiff_t      f_local_n1,f_local_1_end;
321   ptrdiff_t      b_alloc_local,b_local_n0,b_local_0_start;
322   ptrdiff_t      b_local_n1,b_local_1_end;
323   fftw_complex   *data_fin,*data_fout,*data_bout;
324 
325   PetscFunctionBegin;
326 #if !defined(PETSC_USE_COMPLEX)
327   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers");
328 #endif
329   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
330   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
331   if (size == 1){
332     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Works only for parallel 1D");
333   }
334   else {
335       if (ndim>1){
336         SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Works only for parallel 1D");}
337       else {
338           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);
339           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);
340           if (fin) {
341             data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*f_alloc_local);
342             ierr = VecCreateMPIWithArray(comm,f_local_n0,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
343             (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
344           }
345           if (fout) {
346             data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*f_alloc_local);
347             ierr = VecCreateMPIWithArray(comm,f_local_n1,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
348             (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
349           }
350           if (bout) {
351             data_bout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*b_alloc_local);
352             ierr = VecCreateMPIWithArray(comm,b_local_n1,N,(const PetscScalar*)data_bout,bout);CHKERRQ(ierr);
353             (*bout)->ops->destroy   = VecDestroy_MPIFFTW;
354           }
355   }
356   if (fin){
357     ierr = PetscLayoutReference(A->cmap,&(*fin)->map);CHKERRQ(ierr);
358   }
359   if (fout){
360     ierr = PetscLayoutReference(A->rmap,&(*fout)->map);CHKERRQ(ierr);
361   }
362   if (bout){
363     ierr = PetscLayoutReference(A->rmap,&(*bout)->map);CHKERRQ(ierr);
364   }
365   PetscFunctionReturn(0);
366 }
367 
368 
369 }
370 
371 #undef __FUNCT__
372 #define __FUNCT__ "MatGetVecs_FFTW"
373 /*
374    MatGetVecs_FFTW - Get vector(s) compatible with the matrix, i.e. with the
375      parallel layout determined by FFTW
376 
377    Collective on Mat
378 
379    Input Parameter:
380 .  mat - the matrix
381 
382    Output Parameter:
383 +   fin - (optional) input vector of forward FFTW
384 -   fout - (optional) output vector of forward FFTW
385 
386   Level: advanced
387 
388 .seealso: MatCreateFFTW()
389 */
390 PetscErrorCode  MatGetVecs_FFTW(Mat A,Vec *fin,Vec *fout)
391 {
392   PetscErrorCode ierr;
393   PetscMPIInt    size,rank;
394   MPI_Comm       comm=((PetscObject)A)->comm;
395   Mat_FFT        *fft = (Mat_FFT*)A->data;
396   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
397   PetscInt       N=fft->N, N1, n1,vsize;
398   PetscInt       ndim=fft->ndim,*dim=fft->dim,n=fft->n;
399 
400   PetscFunctionBegin;
401   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
402   PetscValidType(A,1);
403 
404   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
405   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
406   if (size == 1){ /* sequential case */
407 #if defined(PETSC_USE_COMPLEX)
408     if (fin) {ierr = VecCreateSeq(PETSC_COMM_SELF,N,fin);CHKERRQ(ierr);}
409     if (fout){ierr = VecCreateSeq(PETSC_COMM_SELF,N,fout);CHKERRQ(ierr);}
410 #else
411     if (fin) {ierr = VecCreateSeq(PETSC_COMM_SELF,N*2*(dim[ndim-1]/2+1)/dim[ndim-1],fin);CHKERRQ(ierr);}
412     if (fout){ierr = VecCreateSeq(PETSC_COMM_SELF,2*N*(dim[ndim-1]/2+1)/dim[ndim-1],fout);CHKERRQ(ierr);}
413 #endif
414   } else {        /* mpi case */
415     ptrdiff_t      alloc_local,local_n0,local_0_start;
416     ptrdiff_t      local_n1,local_1_end;
417     fftw_complex   *data_fin,*data_fout;
418     double         *data_finr ;
419     ptrdiff_t      local_1_start,temp;
420 //    PetscInt ctr;
421 //    ptrdiff_t      ndim1,*pdim;
422 //    ndim1=(ptrdiff_t) ndim;
423 //    pdim = (ptrdiff_t *)calloc(ndim,sizeof(ptrdiff_t));
424 
425 //    for(ctr=0;ctr<ndim;ctr++)
426 //        {k
427 //           pdim[ctr] = dim[ctr];
428 //       }
429 
430 
431 
432     switch (ndim){
433     case 1:
434       /* Get local size */
435       /* We need to write an error message here saying that one cannot call this routine when doing parallel 1D real FFTW */
436 //      SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Works only for parallel Multi-dimensional FFTW, Dimension>1. Check Documentation for MatGetVecs_FFTW1D routine");
437       alloc_local = fftw_mpi_local_size_1d(dim[0],comm,FFTW_FORWARD,FFTW_ESTIMATE,&local_n0,&local_0_start,&local_n1,&local_1_end);
438       if (fin) {
439         data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
440         ierr = VecCreateMPIWithArray(comm,local_n0,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
441         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
442       }
443       if (fout) {
444         data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
445         ierr = VecCreateMPIWithArray(comm,local_n1,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
446         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
447       }
448       break;
449     case 2:
450 #if !defined(PETSC_USE_COMPLEX)
451       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);
452       N1 = 2*dim[0]*(dim[1]/2+1); n1 = 2*local_n0*(dim[1]/2+1);
453       if (fin) {
454         data_finr=(double *)fftw_malloc(sizeof(double)*alloc_local*2);
455         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,(PetscInt)n1,N1,(PetscScalar*)data_finr,fin);CHKERRQ(ierr);
456         ierr = VecGetSize(*fin,&vsize);CHKERRQ(ierr);
457         //printf("The code comes here with vector size %d\n",vsize);
458         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
459       }
460       if (fout) {
461         data_fout=(fftw_complex *)fftw_malloc(sizeof(fftw_complex)*alloc_local);
462         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n1,N1,(PetscScalar*)data_fout,fout);CHKERRQ(ierr);
463         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
464       }
465       printf("Vector size from fftw.c is  given by %d, %d\n",n1,N1);
466 
467 #else
468       /* Get local size */
469      printf("Hope this does not come here");
470       alloc_local = fftw_mpi_local_size_2d(dim[0],dim[1],comm,&local_n0,&local_0_start);
471       if (fin) {
472         data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
473         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
474         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
475       }
476       if (fout) {
477         data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
478         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
479         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
480       }
481      printf("Hope this does not come here");
482 #endif
483       break;
484     case 3:
485       /* Get local size */
486 #if !defined(PETSC_USE_COMPLEX)
487       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);
488       N1 = 2*dim[0]*dim[1]*(dim[2]/2+1); n1 = 2*local_n0*dim[1]*(dim[2]/2+1);
489       if (fin) {
490         data_finr=(double *)fftw_malloc(sizeof(double)*alloc_local*2);
491         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,(PetscInt)n1,N1,(PetscScalar*)data_finr,fin);CHKERRQ(ierr);
492         ierr = VecGetSize(*fin,&vsize);CHKERRQ(ierr);
493         //printf("The code comes here with vector size %d\n",vsize);
494         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
495       }
496       if (fout) {
497         data_fout=(fftw_complex *)fftw_malloc(sizeof(fftw_complex)*alloc_local);
498         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n1,N1,(PetscScalar*)data_fout,fout);CHKERRQ(ierr);
499         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
500       }
501       printf("Vector size from fftw.c is  given by %d, %d\n",n1,N1);
502 
503 
504 #else
505       alloc_local = fftw_mpi_local_size_3d(dim[0],dim[1],dim[2],comm,&local_n0,&local_0_start);
506 //      printf("The quantity n is %d",n);
507       if (fin) {
508         data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
509         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
510         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
511       }
512       if (fout) {
513         data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
514         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
515         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
516       }
517 #endif
518       break;
519     default:
520       /* Get local size */
521 #if !defined(PETSC_USE_COMPLEX)
522       temp = (fftw->dim_fftw)[fftw->ndim_fftw-1];
523       printf("The value of temp is %ld\n",temp);
524       (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp/2 + 1;
525       alloc_local = fftw_mpi_local_size_transposed(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
526       N1 = 2*N*(PetscInt)((fftw->dim_fftw)[fftw->ndim_fftw-1])/((PetscInt) temp);
527       (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp;
528       if (fin) {
529         data_finr=(double *)fftw_malloc(sizeof(double)*alloc_local*2);
530         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,(PetscInt)n,N1,(PetscScalar*)data_finr,fin);CHKERRQ(ierr);
531         ierr = VecGetSize(*fin,&vsize);CHKERRQ(ierr);
532         //printf("The code comes here with vector size %d\n",vsize);
533         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
534       }
535       if (fout) {
536         data_fout=(fftw_complex *)fftw_malloc(sizeof(fftw_complex)*alloc_local);
537         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n,N1,(PetscScalar*)data_fout,fout);CHKERRQ(ierr);
538         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
539       }
540 
541 #else
542       alloc_local = fftw_mpi_local_size(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start);
543 //      printf("The value of alloc local is %d from process %d\n",alloc_local,rank);
544 //      printf("The value of alloc local is %d",alloc_local);
545 //      pdim=(ptrdiff_t *)calloc(ndim,sizeof(ptrdiff_t));
546 //      for(i=0;i<ndim;i++)
547 //         {
548 //          pdim[i]=dim[i];printf("%d",pdim[i]);
549 //         }
550 //      alloc_local = fftw_mpi_local_size(ndim,pdim,comm,&local_n0,&local_0_start);
551 //      printf("The quantity n is %d",n);
552       if (fin) {
553         data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
554         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
555         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
556       }
557       if (fout) {
558         data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
559         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
560         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
561       }
562 #endif
563       break;
564     }
565   }
566 //  if (fin){
567 //    ierr = PetscLayoutReference(A->cmap,&(*fin)->map);CHKERRQ(ierr);
568 //  }
569 //  if (fout){
570 //    ierr = PetscLayoutReference(A->rmap,&(*fout)->map);CHKERRQ(ierr);
571 //  }
572   PetscFunctionReturn(0);
573 }
574 
575 #undef __FUNCT__
576 #define __FUNCT__ "InputTransformFFT"
577 PetscErrorCode InputTransformFFT(Mat A,Vec x,Vec y)
578 {
579   PetscErrorCode ierr;
580   PetscFunctionBegin;
581   ierr = PetscTryMethod(A,"InputTransformFFT_C",(Mat,Vec,Vec),(A,x,y));CHKERRQ(ierr);
582   PetscFunctionReturn(0);
583 }
584 
585 /*
586       InputTransformFFT_FFTW - Copies the user data to the vector that goes into FFTW block
587   Input A, x, y
588   A - FFTW matrix
589   x - user data
590   Options Database Keys:
591 + -mat_fftw_plannerflags - set FFTW planner flags
592 
593    Level: intermediate
594 
595 */
596 
597 EXTERN_C_BEGIN
598 #undef __FUNCT__
599 #define __FUNCT__ "InputTransformFFT_FTTW"
600 PetscErrorCode InputTransformFFT_FFTW(Mat A,Vec x,Vec y)
601 {
602   PetscErrorCode ierr;
603   MPI_Comm       comm=((PetscObject)A)->comm;
604   Mat_FFT        *fft = (Mat_FFT*)A->data;
605   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
606   PetscInt       N=fft->N, N1, n1 ,NM;
607   PetscInt       ndim=fft->ndim,*dim=fft->dim;//n=fft->n;
608   PetscInt       low, *indx1, *indx2, tempindx, tempindx1;
609   PetscInt       i,j,k,rank,size,partial_dim;
610   ptrdiff_t      alloc_local,local_n0,local_0_start;
611   ptrdiff_t      local_n1,local_1_start,temp;
612   VecScatter     vecscat;
613   IS             list1,list2;
614 
615   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
616   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
617   ierr = VecGetOwnershipRange(y,&low,PETSC_NULL);
618   printf("Local ownership starts at %d\n",low);
619 
620   if (size==1)
621     {
622      switch (ndim){
623      case 1:
624           ierr = PetscMalloc(sizeof(PetscInt)*dim[0],&indx1);CHKERRQ(ierr);
625           for (i=0;i<dim[0];i++)
626              {
627               indx1[i] = i;
628              }
629           ierr = ISCreateGeneral(comm,dim[0],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
630           ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
631           ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
632           ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
633           ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
634           ierr = ISDestroy(&list1);CHKERRQ(ierr);
635           ierr = PetscFree(indx1);CHKERRQ(ierr);
636      break;
637 
638      case 2:
639           ierr = PetscMalloc(sizeof(PetscInt)*dim[0]*dim[1],&indx1);CHKERRQ(ierr);
640           for (i=0;i<dim[0];i++){
641              for (j=0;j<dim[1];j++){
642                 indx1[i*dim[1]+j] = i*dim[1] + j;
643              }
644           }
645           ierr = ISCreateGeneral(comm,dim[0]*dim[1],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
646           ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
647           ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
648           ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
649           ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
650           ierr = ISDestroy(&list1);CHKERRQ(ierr);
651           ierr = PetscFree(indx1);CHKERRQ(ierr);
652           break;
653      case 3:
654           ierr = PetscMalloc(sizeof(PetscInt)*dim[0]*dim[1]*dim[2],&indx1);CHKERRQ(ierr);
655           for (i=0;i<dim[0];i++){
656              for (j=0;j<dim[1];j++){
657                 for (k=0;k<dim[2];k++){
658                    indx1[i*dim[1]*dim[2]+j*dim[2]+k] = i*dim[1]*dim[2] + j*dim[2] + k;
659                 }
660              }
661           }
662           ierr = ISCreateGeneral(comm,dim[0]*dim[1]*dim[2],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
663           ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
664           ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
665           ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
666           ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
667           ierr = ISDestroy(&list1);CHKERRQ(ierr);
668           ierr = PetscFree(indx1);CHKERRQ(ierr);
669           break;
670      default:
671           ierr = ISCreateStride(PETSC_COMM_SELF,N,0,1,&list1);
672           ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
673           ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
674           ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
675           ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
676           ierr = ISDestroy(&list1);CHKERRQ(ierr);
677           //ierr = ISDestroy(list1);CHKERRQ(ierr);
678           break;
679       }
680     }
681 
682  else{
683 
684  switch (ndim){
685  case 1:
686   SETERRQ(comm,PETSC_ERR_SUP,"FFTW does not support parallel 1D real transform");
687   break;
688  case 2:
689       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);
690       N1 = 2*dim[0]*(dim[1]/2+1); n1 = 2*local_n0*(dim[1]/2+1);
691 
692       ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx1);CHKERRQ(ierr);
693       ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx2);CHKERRQ(ierr);
694       printf("Val local_0_start = %ld\n",local_0_start);
695 
696       if (dim[1]%2==0)
697         NM = dim[1]+2;
698       else
699         NM = dim[1]+1;
700 
701       for (i=0;i<local_n0;i++){
702          for (j=0;j<dim[1];j++){
703             tempindx = i*dim[1] + j;
704             tempindx1 = i*NM + j;
705             indx1[tempindx]=local_0_start*dim[1]+tempindx;
706             indx2[tempindx]=low+tempindx1;
707            // printf("Val tempindx1 = %d\n",tempindx1);
708   //          printf("index1 %d from proc %d is \n",indx1[tempindx],rank);
709   //          printf("index2 %d from proc %d is \n",indx2[tempindx],rank);
710   //          printf("-------------------------\n",indx2[tempindx],rank);
711         }
712      }
713 
714       ierr = ISCreateGeneral(comm,local_n0*dim[1],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
715       ierr = ISCreateGeneral(comm,local_n0*dim[1],indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr);
716 
717       ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr);
718       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
719       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
720       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
721       ierr = ISDestroy(&list1);CHKERRQ(ierr);
722       ierr = ISDestroy(&list2);CHKERRQ(ierr);
723       ierr = PetscFree(indx1);CHKERRQ(ierr);
724       ierr = PetscFree(indx2);CHKERRQ(ierr);
725       break;
726 
727  case 3:
728       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);
729       N1 = 2*dim[0]*dim[1]*(dim[2]/2+1); n1 = 2*local_n0*dim[1]*(dim[2]/2+1);
730 
731       ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx1);CHKERRQ(ierr);
732       ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx2);CHKERRQ(ierr);
733       printf("Val local_0_start = %ld\n",local_0_start);
734 
735       if (dim[2]%2==0)
736         NM = dim[2]+2;
737       else
738         NM = dim[2]+1;
739 
740       for (i=0;i<local_n0;i++){
741          for (j=0;j<dim[1];j++){
742             for (k=0;k<dim[2];k++){
743                tempindx = i*dim[1]*dim[2] + j*dim[2] + k;
744                tempindx1 = i*dim[1]*NM + j*NM + k;
745                indx1[tempindx]=local_0_start*dim[1]*dim[2]+tempindx;
746                indx2[tempindx]=low+tempindx1;
747             }
748            // printf("Val tempindx1 = %d\n",tempindx1);
749            // printf("index1 %d from proc %d is \n",indx1[tempindx],rank);
750            // printf("index2 %d from proc %d is \n",indx2[tempindx],rank);
751            // printf("-------------------------\n",indx2[tempindx],rank);
752         }
753      }
754 
755       ierr = ISCreateGeneral(comm,local_n0*dim[1]*dim[2],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
756       ierr = ISCreateGeneral(comm,local_n0*dim[1]*dim[2],indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr);
757 
758       ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr);
759       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
760       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
761       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
762       ierr = ISDestroy(&list1);CHKERRQ(ierr);
763       ierr = ISDestroy(&list2);CHKERRQ(ierr);
764       ierr = PetscFree(indx1);CHKERRQ(ierr);
765       ierr = PetscFree(indx2);CHKERRQ(ierr);
766       break;
767 
768  default:
769       temp = (fftw->dim_fftw)[fftw->ndim_fftw-1];
770       printf("The value of temp is %ld\n",temp);
771       (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp/2 + 1;
772       alloc_local = fftw_mpi_local_size_transposed(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
773       N1 = 2*N*(PetscInt)((fftw->dim_fftw)[fftw->ndim_fftw-1])/((PetscInt) temp);
774       (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp;
775 
776       partial_dim = fftw->partial_dim;
777       printf("The value of partial dim is %d\n",partial_dim);
778 
779       ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*partial_dim,&indx1);CHKERRQ(ierr);
780       ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*partial_dim,&indx2);CHKERRQ(ierr);
781       printf("Val local_0_start = %ld\n",local_0_start);
782 
783       if (dim[ndim-1]%2==0)
784         NM = 2;
785       else
786         NM = 1;
787 
788       j = low;
789       for (i=0,k=1; i<((PetscInt)local_n0)*partial_dim;i++,k++)
790          {
791           indx1[i] = local_0_start*partial_dim + i;
792           indx2[i] = j;
793           //printf("The values are %d and %d from %d\n",indx1[i],indx2[i],rank);
794           if (k%dim[ndim-1]==0)
795             { j+=NM;}
796           j++;
797          }
798       ierr = ISCreateGeneral(comm,local_n0*partial_dim,indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
799       ierr = ISCreateGeneral(comm,local_n0*partial_dim,indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr);
800       ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr);
801       ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
802       ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
803       ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
804       ierr = ISDestroy(&list1);CHKERRQ(ierr);
805       ierr = ISDestroy(&list2);CHKERRQ(ierr);
806       ierr = PetscFree(indx1);CHKERRQ(ierr);
807       ierr = PetscFree(indx2);CHKERRQ(ierr);
808       break;
809   }
810  }
811 
812  return 0;
813 }
814 EXTERN_C_END
815 
816 #undef __FUNCT__
817 #define __FUNCT__ "OutputTransformFFT"
818 PetscErrorCode OutputTransformFFT(Mat A,Vec x,Vec y)
819 {
820   PetscErrorCode ierr;
821   PetscFunctionBegin;
822   ierr = PetscTryMethod(A,"OutputTransformFFT_C",(Mat,Vec,Vec),(A,x,y));CHKERRQ(ierr);
823   PetscFunctionReturn(0);
824 }
825 
826 /*
827       OutputTransformFFT_FFTW - Copies the FFTW output to the PETSc vector that user can use
828   Input A, x, y
829   A - FFTW matrix
830   x - FFTW vector
831   y - PETSc vector that user can read
832   Options Database Keys:
833 + -mat_fftw_plannerflags - set FFTW planner flags
834 
835    Level: intermediate
836 
837 */
838 
839 EXTERN_C_BEGIN
840 #undef __FUNCT__
841 #define __FUNCT__ "OutputTransformFFT_FTTW"
842 PetscErrorCode OutputTransformFFT_FFTW(Mat A,Vec x,Vec y)
843 {
844   PetscErrorCode ierr;
845   MPI_Comm       comm=((PetscObject)A)->comm;
846   Mat_FFT        *fft = (Mat_FFT*)A->data;
847   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
848   PetscInt       N=fft->N, N1, n1 ,NM;
849   PetscInt       ndim=fft->ndim,*dim=fft->dim;//n=fft->n;
850   PetscInt       low, *indx1, *indx2, tempindx, tempindx1;
851   PetscInt       i,j,k,rank,size,partial_dim;
852   ptrdiff_t      alloc_local,local_n0,local_0_start;
853   ptrdiff_t      local_n1,local_1_start,temp;
854   VecScatter     vecscat;
855   IS             list1,list2;
856 
857   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
858   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
859   ierr = VecGetOwnershipRange(x,&low,PETSC_NULL);
860 
861   if (size==1){
862     switch (ndim){
863     case 1:
864            ierr = PetscMalloc(sizeof(PetscInt)*dim[0],&indx1);CHKERRQ(ierr);
865           for (i=0;i<dim[0];i++)
866              {
867               indx1[i] = i;
868              }
869           ierr = ISCreateGeneral(comm,dim[0],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
870           ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
871           ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
872           ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
873           ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
874           ierr = ISDestroy(&list1);CHKERRQ(ierr);
875           ierr = PetscFree(indx1);CHKERRQ(ierr);
876           break;
877 
878     case 2:
879          ierr = PetscMalloc(sizeof(PetscInt)*dim[0]*dim[1],&indx1);CHKERRQ(ierr);
880           for (i=0;i<dim[0];i++){
881              for (j=0;j<dim[1];j++){
882                 indx1[i*dim[1]+j] = i*dim[1] + j;
883              }
884           }
885          ierr = ISCreateGeneral(comm,dim[0]*dim[1],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
886          ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
887          ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
888          ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
889          ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
890          ierr = ISDestroy(&list1);CHKERRQ(ierr);
891          ierr = PetscFree(indx1);CHKERRQ(ierr);
892          break;
893     case 3:
894          ierr = PetscMalloc(sizeof(PetscInt)*dim[0]*dim[1]*dim[2],&indx1);CHKERRQ(ierr);
895          for (i=0;i<dim[0];i++){
896             for (j=0;j<dim[1];j++){
897                for (k=0;k<dim[2];k++){
898                   indx1[i*dim[1]*dim[2]+j*dim[2]+k] = i*dim[1]*dim[2] + j*dim[2] + k;
899                }
900             }
901          }
902          ierr = ISCreateGeneral(comm,dim[0]*dim[1]*dim[2],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
903          ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
904          ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
905          ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
906          ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
907          ierr = ISDestroy(&list1);CHKERRQ(ierr);
908          ierr = PetscFree(indx1);CHKERRQ(ierr);
909          break;
910     default:
911          ierr = ISCreateStride(comm,N,0,1,&list1);
912          //ierr = ISView(list1,PETSC_VIEWER_STDOUT_SELF);
913          ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr);
914          ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
915          ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
916          ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
917          ierr = ISDestroy(&list1);CHKERRQ(ierr);
918          break;
919     }
920   }
921   else{
922 
923  switch (ndim){
924  case 1:
925   SETERRQ(comm,PETSC_ERR_SUP,"No support for real parallel 1D FFT");
926   break;
927  case 2:
928       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);
929       N1 = 2*dim[0]*(dim[1]/2+1); n1 = 2*local_n0*(dim[1]/2+1);
930 
931      ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx1);CHKERRQ(ierr);
932      ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx2);CHKERRQ(ierr);
933      printf("Val local_0_start = %ld\n",local_0_start);
934 
935      if (dim[1]%2==0)
936       NM = dim[1]+2;
937     else
938       NM = dim[1]+1;
939 
940 
941 
942      for (i=0;i<local_n0;i++){
943         for (j=0;j<dim[1];j++){
944             tempindx = i*dim[1] + j;
945             tempindx1 = i*NM + j;
946             indx1[tempindx]=local_0_start*dim[1]+tempindx;
947             indx2[tempindx]=low+tempindx1;
948        //     printf("Val tempindx1 = %d\n",tempindx1);
949        //     printf("index1 %d from proc %d is \n",indx1[tempindx],rank);
950        //     printf("index2 %d from proc %d is \n",indx2[tempindx],rank);
951        //     printf("-------------------------\n",indx2[tempindx],rank);
952         }
953      }
954 
955      ierr = ISCreateGeneral(comm,local_n0*dim[1],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
956      ierr = ISCreateGeneral(comm,local_n0*dim[1],indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr);
957 
958      ierr = VecScatterCreate(x,list2,y,list1,&vecscat);CHKERRQ(ierr);
959      ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
960      ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
961      ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
962      ierr = ISDestroy(&list1);CHKERRQ(ierr);
963      ierr = ISDestroy(&list2);CHKERRQ(ierr);
964      ierr = PetscFree(indx1);CHKERRQ(ierr);
965      ierr = PetscFree(indx2);CHKERRQ(ierr);
966   break;
967 
968  case 3:
969       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);
970       N1 = 2*dim[0]*dim[1]*(dim[2]/2+1); n1 = 2*local_n0*dim[1]*(dim[2]/2+1);
971 
972      ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx1);CHKERRQ(ierr);
973      ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx2);CHKERRQ(ierr);
974      printf("Val local_0_start = %ld\n",local_0_start);
975 
976      if (dim[2]%2==0)
977       NM = dim[2]+2;
978     else
979       NM = dim[2]+1;
980 
981      for (i=0;i<local_n0;i++){
982         for (j=0;j<dim[1];j++){
983            for (k=0;k<dim[2];k++){
984               tempindx = i*dim[1]*dim[2] + j*dim[2] + k;
985               tempindx1 = i*dim[1]*NM + j*NM + k;
986               indx1[tempindx]=local_0_start*dim[1]*dim[2]+tempindx;
987               indx2[tempindx]=low+tempindx1;
988            }
989         //    printf("Val tempindx1 = %d\n",tempindx1);
990         //    printf("index1 %d from proc %d is \n",indx1[tempindx],rank);
991         //    printf("index2 %d from proc %d is \n",indx2[tempindx],rank);
992         //    printf("-------------------------\n",indx2[tempindx],rank);
993         }
994      }
995 
996      ierr = ISCreateGeneral(comm,local_n0*dim[1]*dim[2],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
997      ierr = ISCreateGeneral(comm,local_n0*dim[1]*dim[2],indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr);
998 
999      ierr = VecScatterCreate(x,list2,y,list1,&vecscat);CHKERRQ(ierr);
1000      ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1001      ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1002      ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1003      ierr = ISDestroy(&list1);CHKERRQ(ierr);
1004      ierr = ISDestroy(&list2);CHKERRQ(ierr);
1005      ierr = PetscFree(indx1);CHKERRQ(ierr);
1006      ierr = PetscFree(indx2);CHKERRQ(ierr);
1007   break;
1008 
1009   default:
1010      temp = (fftw->dim_fftw)[fftw->ndim_fftw-1];
1011      printf("The value of temp is %ld\n",temp);
1012      (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp/2 + 1;
1013      alloc_local = fftw_mpi_local_size_transposed(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start,&local_n1,&local_1_start);
1014      N1 = 2*N*(PetscInt)((fftw->dim_fftw)[fftw->ndim_fftw-1])/((PetscInt) temp);
1015      (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp;
1016 
1017      partial_dim = fftw->partial_dim;
1018      printf("The value of partial dim is %d\n",partial_dim);
1019 
1020      ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*partial_dim,&indx1);CHKERRQ(ierr);
1021      ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*partial_dim,&indx2);CHKERRQ(ierr);
1022      printf("Val local_0_start = %ld\n",local_0_start);
1023 
1024      if (dim[ndim-1]%2==0)
1025        NM = 2;
1026      else
1027        NM = 1;
1028 
1029      j = low;
1030      for (i=0,k=1; i<((PetscInt)local_n0)*partial_dim;i++,k++)
1031         {
1032          indx1[i] = local_0_start*partial_dim + i;
1033          indx2[i] = j;
1034          //printf("The values are %d and %d from %d\n",indx1[i],indx2[i],rank);
1035          if (k%dim[ndim-1]==0)
1036            { j+=NM;}
1037          j++;
1038         }
1039      ierr = ISCreateGeneral(comm,local_n0*partial_dim,indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
1040      ierr = ISCreateGeneral(comm,local_n0*partial_dim,indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr);
1041 
1042       //ISView(list1,PETSC_VIEWER_STDOUT_SELF);
1043 
1044 
1045      ierr = VecScatterCreate(x,list2,y,list1,&vecscat);CHKERRQ(ierr);
1046      ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1047      ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1048      ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr);
1049      ierr = ISDestroy(&list1);CHKERRQ(ierr);
1050      ierr = ISDestroy(&list2);CHKERRQ(ierr);
1051      ierr = PetscFree(indx1);CHKERRQ(ierr);
1052      ierr = PetscFree(indx2);CHKERRQ(ierr);
1053 
1054   break;
1055  }
1056  }
1057  return 0;
1058 }
1059 EXTERN_C_END
1060 
1061 EXTERN_C_BEGIN
1062 #undef __FUNCT__
1063 #define __FUNCT__ "MatCreate_FFTW"
1064 /*
1065       MatCreate_FFTW - Creates a matrix object that provides FFT
1066   via the external package FFTW
1067   Options Database Keys:
1068 + -mat_fftw_plannerflags - set FFTW planner flags
1069 
1070    Level: intermediate
1071 
1072 */
1073 
1074 PetscErrorCode MatCreate_FFTW(Mat A)
1075 {
1076   PetscErrorCode ierr;
1077   MPI_Comm       comm=((PetscObject)A)->comm;
1078   Mat_FFT        *fft=(Mat_FFT*)A->data;
1079   Mat_FFTW       *fftw;
1080   PetscInt       n=fft->n,N=fft->N,ndim=fft->ndim,*dim = fft->dim;
1081   const char     *p_flags[]={"FFTW_ESTIMATE","FFTW_MEASURE","FFTW_PATIENT","FFTW_EXHAUSTIVE"};
1082   PetscBool      flg;
1083   PetscInt       p_flag,partial_dim=1,ctr,N1;
1084   PetscMPIInt    size,rank;
1085   ptrdiff_t      *pdim, temp;
1086   ptrdiff_t      local_n1,local_1_start;
1087 
1088   PetscFunctionBegin;
1089   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1090   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1091   ierr = MPI_Barrier(PETSC_COMM_WORLD);
1092 
1093   pdim = (ptrdiff_t *)calloc(ndim,sizeof(ptrdiff_t));
1094   pdim[0] = dim[0];
1095   for(ctr=1;ctr<ndim;ctr++)
1096       {
1097           partial_dim *= dim[ctr];
1098           pdim[ctr] = dim[ctr];
1099       }
1100 
1101   if (size == 1) {
1102 #if defined(PETSC_USE_COMPLEX)
1103     ierr = MatSetSizes(A,N,N,N,N);CHKERRQ(ierr);
1104     n = N;
1105 #else
1106     int tot_dim = N*2*(dim[ndim-1]/2+1)/dim[ndim-1];
1107     ierr = MatSetSizes(A,tot_dim,tot_dim,tot_dim,tot_dim);CHKERRQ(ierr);
1108     n = tot_dim;
1109 #endif
1110 
1111   } else {
1112     ptrdiff_t alloc_local,local_n0,local_0_start;//local_n1,local_1_end;
1113     switch (ndim){
1114     case 1:
1115 #if !defined(PETSC_USE_COMPLEX)
1116   SETERRQ(comm,PETSC_ERR_SUP,"FFTW does not support parallel 1D real transform");
1117 #else
1118       alloc_local = fftw_mpi_local_size_1d(dim[0],comm,FFTW_FORWARD,FFTW_ESTIMATE,&local_n0,&local_0_start,&local_n1,&local_1_end);
1119       n = (PetscInt)local_n0;
1120       ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr);
1121 #endif
1122 //      PetscObjectComposeFunctionDynamic((PetscObject)A,"MatGetVecs1DC_C","MatGetVecs1DC_FFTW",MatGetVecs1DC_FFTW);
1123       break;
1124     case 2:
1125 #if defined(PETSC_USE_COMPLEX)
1126       alloc_local = fftw_mpi_local_size_2d(dim[0],dim[1],comm,&local_n0,&local_0_start);
1127       /*
1128        PetscMPIInt    rank;
1129        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]);
1130        PetscSynchronizedFlush(comm);
1131        */
1132       n = (PetscInt)local_n0*dim[1];
1133       ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr);
1134 #else
1135       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);
1136       n = 2*(PetscInt)local_n0*(dim[1]/2+1);
1137       ierr = MatSetSizes(A,n,n,2*dim[0]*(dim[1]/2+1),2*dim[0]*(dim[1]/2+1));
1138 #endif
1139       break;
1140     case 3:
1141 //      printf("The value of alloc local is %d",alloc_local);
1142 #if defined(PETSC_USE_COMPLEX)
1143       alloc_local = fftw_mpi_local_size_3d(dim[0],dim[1],dim[2],comm,&local_n0,&local_0_start);
1144       n = (PetscInt)local_n0*dim[1]*dim[2];
1145       ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr);
1146 #else
1147       printf("The code comes here\n");
1148       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);
1149       n = 2*(PetscInt)local_n0*dim[1]*(dim[2]/2+1);
1150       ierr = MatSetSizes(A,n,n,2*dim[0]*dim[1]*(dim[2]/2+1),2*dim[0]*dim[1]*(dim[2]/2+1));
1151 #endif
1152       break;
1153     default:
1154 #if defined(PETSC_USE_COMPLEX)
1155       alloc_local = fftw_mpi_local_size(ndim,pdim,comm,&local_n0,&local_0_start);
1156 //      printf("The value of alloc local is %ld from process %d\n",alloc_local,rank);
1157 //      alloc_local = fftw_mpi_local_size(ndim,dim,comm,&local_n0,&local_0_start);
1158       n = (PetscInt)local_n0*partial_dim;
1159 //      printf("New partial dimension is %d %d %d",n,N,ndim);
1160       ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr);
1161 #else
1162       temp = pdim[ndim-1];
1163       pdim[ndim-1]= temp/2 + 1;
1164       printf("For Multi dim case temp = %ld, pdim[ndim-1] = %ld\n",temp,pdim[ndim-1]);
1165       alloc_local = fftw_mpi_local_size_transposed(ndim,pdim,PETSC_COMM_WORLD,&local_n0,&local_0_start,&local_n1,&local_1_start);
1166       n = 2*(PetscInt)local_n0*partial_dim*pdim[ndim-1]/temp;
1167       N1 = 2*N*(PetscInt)pdim[ndim-1]/((PetscInt) temp);
1168       pdim[ndim-1] = temp;
1169       printf("For Multi dim case n = %d, N1  = %d\n",n,N1);
1170       ierr = MatSetSizes(A,n,n,N1,N1);CHKERRQ(ierr);
1171 #endif
1172       break;
1173     }
1174   }
1175   ierr = PetscObjectChangeTypeName((PetscObject)A,MATFFTW);CHKERRQ(ierr);
1176   ierr = PetscNewLog(A,Mat_FFTW,&fftw);CHKERRQ(ierr);
1177   fft->data = (void*)fftw;
1178 
1179   fft->n           = n;
1180   fftw->ndim_fftw  = (ptrdiff_t)ndim; // This is dimension of fft
1181   fftw->partial_dim = partial_dim;
1182   ierr = PetscMalloc(ndim*sizeof(ptrdiff_t), (ptrdiff_t *)&(fftw->dim_fftw));CHKERRQ(ierr);
1183   for(ctr=0;ctr<ndim;ctr++) (fftw->dim_fftw)[ctr]=dim[ctr];
1184 
1185   fftw->p_forward  = 0;
1186   fftw->p_backward = 0;
1187   fftw->p_flag     = FFTW_ESTIMATE;
1188 
1189   if (size == 1){
1190     A->ops->mult          = MatMult_SeqFFTW;
1191     A->ops->multtranspose = MatMultTranspose_SeqFFTW;
1192   } else {
1193     A->ops->mult          = MatMult_MPIFFTW;
1194     A->ops->multtranspose = MatMultTranspose_MPIFFTW;
1195   }
1196   fft->matdestroy          = MatDestroy_FFTW;
1197 // if(ndim=1 && size>1) and also if it is complex then getvecs should be attached to MatGetVecs_FFTW1D
1198   A->ops->getvecs       = MatGetVecs_FFTW;
1199   A->assembled          = PETSC_TRUE;
1200 #if !defined(PETSC_USE_COMPLEX)
1201   PetscObjectComposeFunctionDynamic((PetscObject)A,"InputTransformFFT_C","InputTransformFFT_FFTW",InputTransformFFT_FFTW);
1202   PetscObjectComposeFunctionDynamic((PetscObject)A,"OutputTransformFFT_C","OutputTransformFFT_FFTW",OutputTransformFFT_FFTW);
1203 #endif
1204 
1205   /* get runtime options */
1206   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"FFTW Options","Mat");CHKERRQ(ierr);
1207     ierr = PetscOptionsEList("-mat_fftw_plannerflags","Planner Flags","None",p_flags,4,p_flags[0],&p_flag,&flg);CHKERRQ(ierr);
1208     if (flg) {fftw->p_flag = (unsigned)p_flag;}
1209   PetscOptionsEnd();
1210   PetscFunctionReturn(0);
1211 }
1212 EXTERN_C_END
1213 
1214 
1215 
1216 
1217