xref: /petsc/src/mat/impls/fft/fftw/fftw.c (revision 54dd51181056e65e4e0ae3624b0a65d971e8e32b)
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   fftw_plan   p_forward,p_backward;
15   unsigned    p_flag; /* planner flags, FFTW_ESTIMATE,FFTW_MEASURE, FFTW_PATIENT, FFTW_EXHAUSTIVE */
16   PetscScalar *finarray,*foutarray,*binarray,*boutarray; /* keep track of arrays becaue fftw plan should be
17                                                             executed for the arrays with which the plan was created */
18 } Mat_FFTW;
19 
20 extern PetscErrorCode MatMult_SeqFFTW(Mat,Vec,Vec);
21 extern PetscErrorCode MatMultTranspose_SeqFFTW(Mat,Vec,Vec);
22 extern PetscErrorCode MatMult_MPIFFTW(Mat,Vec,Vec);
23 extern PetscErrorCode MatMultTranspose_MPIFFTW(Mat,Vec,Vec);
24 extern PetscErrorCode MatDestroy_FFTW(Mat);
25 extern PetscErrorCode VecDestroy_MPIFFTW(Vec);
26 extern PetscErrorCode MatGetVecs_FFTW(Mat,Vec*,Vec*);
27 
28 #undef __FUNCT__
29 #define __FUNCT__ "MatMult_SeqFFTW"
30 PetscErrorCode MatMult_SeqFFTW(Mat A,Vec x,Vec y)
31 {
32   PetscErrorCode ierr;
33   Mat_FFT        *fft  = (Mat_FFT*)A->data;
34   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
35   PetscScalar    *x_array,*y_array;
36   PetscInt       ndim=fft->ndim,*dim=fft->dim;
37 
38   PetscFunctionBegin;
39 #if !defined(PETSC_USE_COMPLEX)
40 
41   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers");
42 #endif
43   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
44   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
45   if (!fftw->p_forward){ /* create a plan, then excute it */
46     switch (ndim){
47     case 1:
48       fftw->p_forward = fftw_plan_dft_1d(dim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_FORWARD,fftw->p_flag);
49       break;
50     case 2:
51       fftw->p_forward = fftw_plan_dft_2d(dim[0],dim[1],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_FORWARD,fftw->p_flag);
52       break;
53     case 3:
54       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);
55       break;
56     default:
57       fftw->p_forward = fftw_plan_dft(ndim,dim,(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_FORWARD,fftw->p_flag);
58       break;
59     }
60     fftw->finarray  = x_array;
61     fftw->foutarray = y_array;
62     /* Warning: if (fftw->p_flag!==FFTW_ESTIMATE) The data in the in/out arrays is overwritten!
63                 planning should be done before x is initialized! See FFTW manual sec2.1 or sec4 */
64     fftw_execute(fftw->p_forward);
65   } else { /* use existing plan */
66     if (fftw->finarray != x_array || fftw->foutarray != y_array){ /* use existing plan on new arrays */
67       fftw_execute_dft(fftw->p_forward,(fftw_complex*)x_array,(fftw_complex*)y_array);
68     } else {
69       fftw_execute(fftw->p_forward);
70     }
71   }
72   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
73   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
74   PetscFunctionReturn(0);
75 }
76 
77 #undef __FUNCT__
78 #define __FUNCT__ "MatMultTranspose_SeqFFTW"
79 PetscErrorCode MatMultTranspose_SeqFFTW(Mat A,Vec x,Vec y)
80 {
81   PetscErrorCode ierr;
82   Mat_FFT        *fft = (Mat_FFT*)A->data;
83   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
84   PetscScalar    *x_array,*y_array;
85   PetscInt       ndim=fft->ndim,*dim=fft->dim;
86 
87   PetscFunctionBegin;
88 #if !defined(PETSC_USE_COMPLEX)
89   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers");
90 #endif
91   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
92   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
93   if (!fftw->p_backward){ /* create a plan, then excute it */
94     switch (ndim){
95     case 1:
96       fftw->p_backward = fftw_plan_dft_1d(dim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_BACKWARD,fftw->p_flag);
97       break;
98     case 2:
99       fftw->p_backward = fftw_plan_dft_2d(dim[0],dim[1],(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_BACKWARD,fftw->p_flag);
100       break;
101     case 3:
102       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);
103       break;
104     default:
105       fftw->p_backward = fftw_plan_dft(ndim,dim,(fftw_complex*)x_array,(fftw_complex*)y_array,FFTW_BACKWARD,fftw->p_flag);
106       break;
107     }
108     fftw->binarray  = x_array;
109     fftw->boutarray = y_array;
110     fftw_execute(fftw->p_backward);CHKERRQ(ierr);
111   } else { /* use existing plan */
112     if (fftw->binarray != x_array || fftw->boutarray != y_array){ /* use existing plan on new arrays */
113       fftw_execute_dft(fftw->p_backward,(fftw_complex*)x_array,(fftw_complex*)y_array);
114     } else {
115       fftw_execute(fftw->p_backward);CHKERRQ(ierr);
116     }
117   }
118   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
119   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
120   PetscFunctionReturn(0);
121 }
122 
123 #undef __FUNCT__
124 #define __FUNCT__ "MatMult_MPIFFTW"
125 PetscErrorCode MatMult_MPIFFTW(Mat A,Vec x,Vec y)
126 {
127   PetscErrorCode ierr;
128   Mat_FFT        *fft  = (Mat_FFT*)A->data;
129   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
130   PetscScalar    *x_array,*y_array;
131   PetscInt       ndim=fft->ndim,*dim=fft->dim;
132   MPI_Comm       comm=((PetscObject)A)->comm;
133 // PetscInt ctr;
134 //  ptrdiff_t      ndim1=(ptrdiff_t) ndim,*pdim;
135 //  ndim1=(ptrdiff_t) ndim;
136 //  pdim = (ptrdiff_t *)calloc(ndim,sizeof(ptrdiff_t));
137 
138 //  for(ctr=0;ctr<ndim;ctr++)
139 //     {
140 //      pdim[ctr] = dim[ctr];
141 //     }
142 
143   PetscFunctionBegin;
144 #if !defined(PETSC_USE_COMPLEX)
145   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers");
146 #endif
147 //  pdim = (ptrdiff_t *)calloc(ndim,sizeof(ptrdiff_t));
148 //  for (ctr=0; ctr<ndim; ctr++) pdim[ctr] = dim[ctr];
149 
150   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
151   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
152   if (!fftw->p_forward){ /* create a plan, then excute it */
153     switch (ndim){
154     case 1:
155 #if defined(PETSC_USE_COMPLEX)
156       fftw->p_forward = fftw_mpi_plan_dft_1d(dim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_FORWARD,fftw->p_flag);
157 #endif
158       break;
159     case 2:
160 #if defined(PETSC_USE_COMPLEX)
161       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);
162 #else
163       fftw->p_forward = fftw_mpi_plan_dft_r2c_2d(dim[0],dim[1],(double *)x_array,(fftw_complex*)y_array,comm,FFTW_ESTIMATE);
164 #endif
165       break;
166     case 3:
167 #if defined(PETSC_USE_COMPLEX)
168       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);
169 #else
170       fftw->p_forward = fftw_mpi_plan_dft_r2c_3d(dim[0],dim[1],dim[3],(double *)x_array,(fftw_complex*)y_array,comm,FFTW_ESTIMATE);
171 #endif
172       break;
173     default:
174 #if defined(PETSC_USE_COMPLEX)
175       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);
176 #else
177       fftw->p_forward = fftw_mpi_plan_dft_r2c(fftw->ndim_fftw,fftw->dim_fftw,(double *)x_array,(fftw_complex*)y_array,comm,FFTW_ESTIMATE);
178 #endif
179  //     fftw->p_forward = fftw_mpi_plan_dft(ndim,dim,(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_FORWARD,fftw->p_flag);
180       break;
181     }
182     fftw->finarray  = x_array;
183     fftw->foutarray = y_array;
184     /* Warning: if (fftw->p_flag!==FFTW_ESTIMATE) The data in the in/out arrays is overwritten!
185                 planning should be done before x is initialized! See FFTW manual sec2.1 or sec4 */
186     fftw_execute(fftw->p_forward);
187   } else { /* use existing plan */
188     if (fftw->finarray != x_array || fftw->foutarray != y_array){ /* use existing plan on new arrays */
189       fftw_execute_dft(fftw->p_forward,(fftw_complex*)x_array,(fftw_complex*)y_array);
190     } else {
191       fftw_execute(fftw->p_forward);
192     }
193   }
194   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
195   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
196   PetscFunctionReturn(0);
197 }
198 
199 #undef __FUNCT__
200 #define __FUNCT__ "MatMultTranspose_MPIFFTW"
201 PetscErrorCode MatMultTranspose_MPIFFTW(Mat A,Vec x,Vec y)
202 {
203   PetscErrorCode ierr;
204   Mat_FFT        *fft  = (Mat_FFT*)A->data;
205   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
206   PetscScalar    *x_array,*y_array;
207   PetscInt       ndim=fft->ndim,*dim=fft->dim;
208   MPI_Comm       comm=((PetscObject)A)->comm;
209 //  PetscInt       ctr;
210 //  ptrdiff_t      ndim1=(ptrdiff_t)ndim,*pdim;
211 //  ndim1=(ptrdiff_t) ndim;
212 //  pdim = (ptrdiff_t *)calloc(ndim,sizeof(ptrdiff_t));
213 
214 //  for(ctr=0;ctr<ndim;ctr++)
215 //     {
216 //      pdim[ctr] = dim[ctr];
217 //     }
218 
219   PetscFunctionBegin;
220 //#if !defined(PETSC_USE_COMPLEX)
221 //  SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers");
222 //#endif
223 //  ierr = PetscMalloc(ndim*sizeof(ptrdiff_t), (ptrdiff_t *)&pdim);CHKERRQ(ierr);
224 // should pdim be a member of Mat_FFTW?
225 //  for (ctr=0; ctr<ndim; ctr++) pdim[ctr] = dim[ctr];
226 
227   ierr = VecGetArray(x,&x_array);CHKERRQ(ierr);
228   ierr = VecGetArray(y,&y_array);CHKERRQ(ierr);
229   if (!fftw->p_backward){ /* create a plan, then excute it */
230     switch (ndim){
231     case 1:
232 #if defined(PETSC_USE_COMPLEX)
233       fftw->p_backward = fftw_mpi_plan_dft_1d(dim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_BACKWARD,fftw->p_flag);
234 #endif
235       break;
236     case 2:
237 #if defined(PETSC_USE_COMPLEX)
238       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);
239 #else
240       fftw->p_backward = fftw_mpi_plan_dft_c2r_2d(dim[0],dim[1],(fftw_complex*)x_array,(double *)y_array,comm,FFTW_ESTIMATE);
241 #endif
242       break;
243     case 3:
244 #if defined(PETSC_USE_COMPLEX)
245       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);
246 #else
247       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);
248 #endif
249       break;
250     default:
251 #if defined(PETSC_USE_COMPLEX)
252       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);
253 #else
254       fftw->p_backward = fftw_mpi_plan_dft_c2r(fftw->ndim_fftw,fftw->dim_fftw,(fftw_complex*)x_array,(double *)y_array,comm,FFTW_ESTIMATE);
255 #endif
256 //      fftw->p_backward = fftw_mpi_plan_dft(ndim,dim,(fftw_complex*)x_array,(fftw_complex*)y_array,comm,FFTW_BACKWARD,fftw->p_flag);
257       break;
258     }
259     fftw->binarray  = x_array;
260     fftw->boutarray = y_array;
261     fftw_execute(fftw->p_backward);CHKERRQ(ierr);
262   } else { /* use existing plan */
263     if (fftw->binarray != x_array || fftw->boutarray != y_array){ /* use existing plan on new arrays */
264       fftw_execute_dft(fftw->p_backward,(fftw_complex*)x_array,(fftw_complex*)y_array);
265     } else {
266       fftw_execute(fftw->p_backward);CHKERRQ(ierr);
267     }
268   }
269   ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr);
270   ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr);
271   PetscFunctionReturn(0);
272 }
273 
274 #undef __FUNCT__
275 #define __FUNCT__ "MatDestroy_FFTW"
276 PetscErrorCode MatDestroy_FFTW(Mat A)
277 {
278   Mat_FFT        *fft = (Mat_FFT*)A->data;
279   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
280   PetscErrorCode ierr;
281 
282   PetscFunctionBegin;
283 #if !defined(PETSC_USE_COMPLEX)
284   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers");
285 #endif
286   fftw_destroy_plan(fftw->p_forward);
287   fftw_destroy_plan(fftw->p_backward);
288   ierr = PetscFree(fftw->dim_fftw);CHKERRQ(ierr);
289   ierr = PetscFree(fft->data);CHKERRQ(ierr);
290   PetscFunctionReturn(0);
291 }
292 
293 #include <../src/vec/vec/impls/mpi/pvecimpl.h>   /*I  "petscvec.h"   I*/
294 #undef __FUNCT__
295 #define __FUNCT__ "VecDestroy_MPIFFTW"
296 PetscErrorCode VecDestroy_MPIFFTW(Vec v)
297 {
298   PetscErrorCode  ierr;
299   PetscScalar     *array;
300 
301   PetscFunctionBegin;
302 #if !defined(PETSC_USE_COMPLEX)
303   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers");
304 #endif
305   ierr = VecGetArray(v,&array);CHKERRQ(ierr);
306   fftw_free((fftw_complex*)array);CHKERRQ(ierr);
307   ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
308   ierr = VecDestroy_MPI(v);CHKERRQ(ierr);
309   PetscFunctionReturn(0);
310 }
311 
312 #undef __FUNCT__
313 #define __FUNCT__ "MatGetVecs1DC_FFTW"
314 /*
315    MatGetVecs_FFTW1D - Get Vectors(s) compatible with matrix, i.e. with the
316      parallel layout determined by FFTW-1D
317 
318 */
319 PetscErrorCode  MatGetVecs_FFTW1D(Mat A,Vec *fin,Vec *fout,Vec *bin,Vec *bout)
320 {
321   PetscErrorCode ierr;
322   PetscMPIInt    size,rank;
323   MPI_Comm       comm=((PetscObject)A)->comm;
324   Mat_FFT        *fft = (Mat_FFT*)A->data;
325 //  Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
326   PetscInt       N=fft->N;
327   PetscInt       ndim=fft->ndim,*dim=fft->dim;
328   ptrdiff_t      f_alloc_local,f_local_n0,f_local_0_start;
329   ptrdiff_t      f_local_n1,f_local_1_end;
330   ptrdiff_t      b_alloc_local,b_local_n0,b_local_0_start;
331   ptrdiff_t      b_local_n1,b_local_1_end;
332   fftw_complex   *data_fin,*data_fout,*data_bin,*data_bout;
333 
334   PetscFunctionBegin;
335 #if !defined(PETSC_USE_COMPLEX)
336   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers");
337 #endif
338   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
339   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
340   if (size == 1){
341     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Works only for parallel 1D");
342   }
343   else {
344       if (ndim>1){
345         SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Works only for parallel 1D");}
346       else {
347           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);
348           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);
349           if (fin) {
350             data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*f_alloc_local);
351             ierr = VecCreateMPIWithArray(comm,f_local_n0,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
352             (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
353           }
354           if (fout) {
355             data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*f_alloc_local);
356             ierr = VecCreateMPIWithArray(comm,f_local_n1,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
357             (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
358           }
359           if (bin) {
360             data_bin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*b_alloc_local);
361             ierr = VecCreateMPIWithArray(comm,b_local_n0,N,(const PetscScalar*)data_bin,bin);CHKERRQ(ierr);
362             (*bin)->ops->destroy   = VecDestroy_MPIFFTW;
363           }
364           if (bout) {
365             data_bout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*b_alloc_local);
366             ierr = VecCreateMPIWithArray(comm,b_local_n1,N,(const PetscScalar*)data_bout,bout);CHKERRQ(ierr);
367             (*bout)->ops->destroy   = VecDestroy_MPIFFTW;
368           }
369   }
370   if (fin){
371     ierr = PetscLayoutReference(A->cmap,&(*fin)->map);CHKERRQ(ierr);
372   }
373   if (fout){
374     ierr = PetscLayoutReference(A->rmap,&(*fout)->map);CHKERRQ(ierr);
375   }
376   if (bin){
377     ierr = PetscLayoutReference(A->rmap,&(*bin)->map);CHKERRQ(ierr);
378   }
379   if (bout){
380     ierr = PetscLayoutReference(A->rmap,&(*bout)->map);CHKERRQ(ierr);
381   }
382   PetscFunctionReturn(0);
383 }
384 
385 
386 }
387 
388 #undef __FUNCT__
389 #define __FUNCT__ "MatGetVecs_FFTW"
390 /*
391    MatGetVecs_FFTW - Get vector(s) compatible with the matrix, i.e. with the
392      parallel layout determined by FFTW
393 
394    Collective on Mat
395 
396    Input Parameter:
397 .  mat - the matrix
398 
399    Output Parameter:
400 +   fin - (optional) input vector of forward FFTW
401 -   fout - (optional) output vector of forward FFTW
402 
403   Level: advanced
404 
405 .seealso: MatCreateFFTW()
406 */
407 PetscErrorCode  MatGetVecs_FFTW(Mat A,Vec *fin,Vec *fout)
408 {
409   PetscErrorCode ierr;
410   PetscMPIInt    size,rank;
411   MPI_Comm       comm=((PetscObject)A)->comm;
412   Mat_FFT        *fft = (Mat_FFT*)A->data;
413   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
414   PetscInt       N=fft->N, N1, n1,vsize;
415 
416   PetscFunctionBegin;
417 //#if !defined(PETSC_USE_COMPLEX)
418 //  SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not support for real numbers");
419 //#endif
420   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
421   PetscValidType(A,1);
422 
423   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
424   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
425   if (size == 1){ /* sequential case */
426     if (fin) {ierr = VecCreateSeq(PETSC_COMM_SELF,N,fin);CHKERRQ(ierr);}
427     if (fout){ierr = VecCreateSeq(PETSC_COMM_SELF,N,fout);CHKERRQ(ierr);}
428     printf("The code successfully comes at the end of the routine with one processor\n");
429   } else {        /* mpi case */
430     ptrdiff_t      alloc_local,local_n0,local_0_start;
431     ptrdiff_t      local_n1,local_1_end;
432     PetscInt       ndim=fft->ndim,*dim=fft->dim,n=fft->n;
433     fftw_complex   *data_fin,*data_fout;
434     double         *data_finr, *data_foutr;
435     ptrdiff_t local_1_start;
436 //    PetscInt ctr;
437 //    ptrdiff_t      ndim1,*pdim;
438 //    ndim1=(ptrdiff_t) ndim;
439 //    pdim = (ptrdiff_t *)calloc(ndim,sizeof(ptrdiff_t));
440 
441 //    for(ctr=0;ctr<ndim;ctr++)
442 //        {
443 //           pdim[ctr] = dim[ctr];
444 //       }
445 
446     switch (ndim){
447     case 1:
448       /* Get local size */
449       /* We need to write an error message here saying that one cannot call this routine when doing paralllel 1D complex FFTW */
450 //      SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Works only for parallel Multi-dimensional FFTW, Dimension>1. Check Documentation for MatGetVecs_FFTW1D routine");
451       alloc_local = fftw_mpi_local_size_1d(dim[0],comm,FFTW_FORWARD,FFTW_ESTIMATE,&local_n0,&local_0_start,&local_n1,&local_1_end);
452       if (fin) {
453         data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
454         ierr = VecCreateMPIWithArray(comm,local_n0,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
455         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
456       }
457       if (fout) {
458         data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
459         ierr = VecCreateMPIWithArray(comm,local_n1,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
460         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
461       }
462       break;
463     case 2:
464 #if !defined(PETSC_USE_COMPLEX)
465       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);
466       N1 = 2*dim[0]*(dim[1]/2+1); n1 = 2*local_n0*(dim[1]/2+1);
467       if (fin) {
468         data_finr=(double *)fftw_malloc(sizeof(double)*alloc_local*2);
469         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,(PetscInt)n1,N1,(PetscScalar*)data_finr,fin);CHKERRQ(ierr);
470         ierr = VecGetSize(*fin,&vsize);CHKERRQ(ierr);
471         //printf("The code comes here with vector size %d\n",vsize);
472         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
473       }
474       if (fout) {
475         data_fout=(fftw_complex *)fftw_malloc(sizeof(fftw_complex)*alloc_local);
476         ierr = VecCreateMPIWithArray(PETSC_COMM_WORLD,n1,N1,(PetscScalar*)data_fout,fout);CHKERRQ(ierr);
477         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
478       }
479       printf("Vector size from fftw.c is  given by %d\n",N1);
480 
481 #else
482       /* Get local size */
483      printf("Hope this does not come here");
484       alloc_local = fftw_mpi_local_size_2d(dim[0],dim[1],comm,&local_n0,&local_0_start);
485       if (fin) {
486         data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
487         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
488         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
489       }
490       if (fout) {
491         data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
492         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
493         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
494       }
495      printf("Hope this does not come here");
496 #endif
497       break;
498     case 3:
499       /* Get local size */
500 #if !defined(PETSC_USE_COMPLEX)
501   SETERRQ(comm,PETSC_ERR_SUP,"Not done yet");
502 #else
503       alloc_local = fftw_mpi_local_size_3d(dim[0],dim[1],dim[2],comm,&local_n0,&local_0_start);
504 //      printf("The quantity n is %d",n);
505       if (fin) {
506         data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
507         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
508         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
509       }
510       if (fout) {
511         data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
512         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
513         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
514       }
515 #endif
516       break;
517     default:
518       /* Get local size */
519 #if !defined(PETSC_USE_COMPLEX)
520   SETERRQ(comm,PETSC_ERR_SUP,"Not done yet");
521 #else
522       alloc_local = fftw_mpi_local_size(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start);
523 //      printf("The value of alloc local is %d from process %d\n",alloc_local,rank);
524 //      printf("The value of alloc local is %d",alloc_local);
525 //      pdim=(ptrdiff_t *)calloc(ndim,sizeof(ptrdiff_t));
526 //      for(i=0;i<ndim;i++)
527 //         {
528 //          pdim[i]=dim[i];printf("%d",pdim[i]);
529 //         }
530 //      alloc_local = fftw_mpi_local_size(ndim,pdim,comm,&local_n0,&local_0_start);
531 //      printf("The quantity n is %d",n);
532       if (fin) {
533         data_fin  = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
534         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fin,fin);CHKERRQ(ierr);
535         (*fin)->ops->destroy   = VecDestroy_MPIFFTW;
536       }
537       if (fout) {
538         data_fout = (fftw_complex*)fftw_malloc(sizeof(fftw_complex)*alloc_local);
539         ierr = VecCreateMPIWithArray(comm,n,N,(const PetscScalar*)data_fout,fout);CHKERRQ(ierr);
540         (*fout)->ops->destroy   = VecDestroy_MPIFFTW;
541       }
542 #endif
543       break;
544     }
545   }
546 //  if (fin){
547 //    ierr = PetscLayoutReference(A->cmap,&(*fin)->map);CHKERRQ(ierr);
548 //  }
549 //  if (fout){
550 //    ierr = PetscLayoutReference(A->rmap,&(*fout)->map);CHKERRQ(ierr);
551 //  }
552   PetscFunctionReturn(0);
553 }
554 
555 //EXTERN_C_BEGIN - Do we need this?
556 #undef __FUNCT__
557 #define __FUNCT__ "InputTransformFFT"
558 PetscErrorCode InputTransformFFT(Mat A,Vec x,Vec y)
559 {
560   PetscErrorCode ierr;
561   PetscFunctionBegin;
562   ierr = PetscTryMethod(A,"InputTransformFFT_C",(Mat,Vec,Vec),(A,x,y));CHKERRQ(ierr);
563   PetscFunctionReturn(0);
564 }
565 //EXTERN_C_END - Do we need this?
566 /*
567       InputTransformFFT_FFTW - Copies the user data to the vector that goes into FFTW block
568   Input A, x, y
569   A - FFTW matrix
570   x - user data
571   Options Database Keys:
572 + -mat_fftw_plannerflags - set FFTW planner flags
573 
574    Level: intermediate
575 
576 */
577 
578 EXTERN_C_BEGIN
579 #undef __FUNCT__
580 #define __FUNCT__ "InputTransformFFT_FTTW"
581 PetscErrorCode InputTransformFFT_FFTW(Mat A,Vec x,Vec y)
582 {
583   PetscErrorCode ierr;
584   MPI_Comm       comm=((PetscObject)A)->comm;
585   Mat_FFT        *fft = (Mat_FFT*)A->data;
586   Mat_FFTW       *fftw = (Mat_FFTW*)fft->data;
587   PetscInt       N=fft->N, N1, n1 ,NM;
588   PetscInt       ndim=fft->ndim,*dim=fft->dim,n=fft->n;
589   PetscInt       low, *indx1, *indx2, tempindx, tempindx1;
590   PetscInt       i,j;
591   ptrdiff_t      alloc_local,local_n0,local_0_start;
592   ptrdiff_t      local_n1,local_1_start;
593   VecScatter     vecscat;
594   IS             list1,list2;
595 
596   ierr = VecGetOwnershipRange(y,&low,PETSC_NULL);
597 
598  switch (ndim){
599  case 1:
600   SETERRQ(comm,PETSC_ERR_SUP,"Not Supported by FFTW");
601   break;
602  case 2:
603       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);
604       N1 = 2*dim[0]*(dim[1]/2+1); n1 = 2*local_n0*(dim[1]/2+1);
605 
606      ierr = PetscMalloc(sizeof(PetscInt)*local_n0*N1,&indx1);CHKERRQ(ierr);
607      ierr = PetscMalloc(sizeof(PetscInt)*local_n0*N1,&indx2);CHKERRQ(ierr);
608 
609      if (dim[1]%2==0)
610       NM = dim[1]+2;
611     else
612       NM = dim[1]+1;
613 
614      for (i=0;i<local_n0;i++){
615         for (j=0;j<dim[1];j++){
616             tempindx = i*dim[1] + j;
617             tempindx1 = i*NM + j;
618             indx1[tempindx]=local_0_start*N1+tempindx;
619             indx2[tempindx]=low+tempindx1;
620   //          printf("index3 %d from proc %d is \n",indx3[tempindx],rank);
621   //          printf("index4 %d from proc %d is \n",indx4[tempindx],rank);
622         }
623      }
624 
625      ierr = ISCreateGeneral(comm,local_n0*dim[1],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr);
626      ierr = ISCreateGeneral(comm,local_n0*dim[1],indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr);
627 
628      ierr = VecScatterCreate(x,list1,y,list2,&vecscat);
629      ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);
630      ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);
631      ierr = VecScatterDestroy(&vecscat);
632      break;
633 
634  case 3:
635   SETERRQ(comm,PETSC_ERR_SUP,"Not Done Yet");
636   break;
637 
638  default:
639   SETERRQ(comm,PETSC_ERR_SUP,"Not Done Yet");
640   break;
641  }
642 
643  return 0;
644 }
645 EXTERN_C_END
646 
647 /*
648 //EXTERN_C_BEGIN - Do we need this?
649 #undef __FUNCT__
650 #define __FUNCT__ "OutputTransformFFT"
651 PetscErrorCode OutputTransformFFT(Mat A,Vec x,Vec y)
652 {
653   PetscErrorCode ierr;
654   PetscFunctionBegin;
655   ierr = PetscTryMethod(A,"OutputTransformFFT_C",(Mat,Vec,Vec),(A,x,y));CHKERRQ(ierr);
656   PetscFunctionReturn(0);
657 }
658 //EXTERN_C_END - Do we need this?
659 */
660 
661 EXTERN_C_BEGIN
662 #undef __FUNCT__
663 #define __FUNCT__ "MatCreate_FFTW"
664 /*
665       MatCreate_FFTW - Creates a matrix object that provides FFT
666   via the external package FFTW
667   Options Database Keys:
668 + -mat_fftw_plannerflags - set FFTW planner flags
669 
670    Level: intermediate
671 
672 */
673 
674 PetscErrorCode MatCreate_FFTW(Mat A)
675 {
676   PetscErrorCode ierr;
677   MPI_Comm       comm=((PetscObject)A)->comm;
678   Mat_FFT        *fft=(Mat_FFT*)A->data;
679   Mat_FFTW       *fftw;
680   PetscInt       n=fft->n,N=fft->N,ndim=fft->ndim,*dim = fft->dim;
681   const char     *p_flags[]={"FFTW_ESTIMATE","FFTW_MEASURE","FFTW_PATIENT","FFTW_EXHAUSTIVE"};
682   PetscBool      flg;
683   PetscInt       p_flag,partial_dim=1,ctr;
684   PetscMPIInt    size,rank;
685   ptrdiff_t      *pdim;
686 
687   PetscFunctionBegin;
688 //#if !defined(PETSC_USE_COMPLEX)
689 //  SETERRQ(comm,PETSC_ERR_SUP,"not support for real numbers");
690 //#endif
691 
692   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
693   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
694       printf("The code is coming here\n");
695   ierr = MPI_Barrier(PETSC_COMM_WORLD);
696 
697   pdim = (ptrdiff_t *)calloc(ndim,sizeof(ptrdiff_t));
698   pdim[0] = dim[0];
699   for(ctr=1;ctr<ndim;ctr++)
700       {
701           partial_dim *= dim[ctr];
702           pdim[ctr] = dim[ctr];
703       }
704 //#if !defined(PETSC_USE_COMPLEX)
705 //  SETERRQ(comm,PETSC_ERR_SUP,"not support for real numbers");
706 //#endif
707 
708 //  printf("partial dimension is %d",partial_dim);
709   if (size == 1) {
710     ierr = MatSetSizes(A,N,N,N,N);CHKERRQ(ierr);
711     n = N;
712   } else {
713     ptrdiff_t alloc_local,local_n0,local_0_start,local_n1,local_1_end;
714     switch (ndim){
715     case 1:
716 #if !defined(PETSC_USE_COMPLEX)
717       printf("The code is coming here\n");
718   SETERRQ(comm,PETSC_ERR_SUP,"FFTW does not support parallel 1D real transform");
719 #endif
720       alloc_local = fftw_mpi_local_size_1d(dim[0],comm,FFTW_FORWARD,FFTW_ESTIMATE,&local_n0,&local_0_start,&local_n1,&local_1_end);
721       n = (PetscInt)local_n0;
722       ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr);
723 //      PetscObjectComposeFunctionDynamic((PetscObject)A,"MatGetVecs1DC_C","MatGetVecs1DC_FFTW",MatGetVecs1DC_FFTW);
724       break;
725     case 2:
726       printf("The code is coming here\n");
727       alloc_local = fftw_mpi_local_size_2d(dim[0],dim[1],comm,&local_n0,&local_0_start);
728       /*
729        PetscMPIInt    rank;
730        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]);
731        PetscSynchronizedFlush(comm);
732        */
733       n = (PetscInt)local_n0*dim[1];
734       ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr);
735       break;
736     case 3:
737 //      printf("The value of alloc local is %d",alloc_local);
738       n = (PetscInt)local_n0*dim[1]*dim[2];
739       ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr);
740       break;
741     default:
742       alloc_local = fftw_mpi_local_size(ndim,pdim,comm,&local_n0,&local_0_start);
743 //      printf("The value of alloc local is %ld from process %d\n",alloc_local,rank);
744 //      alloc_local = fftw_mpi_local_size(ndim,dim,comm,&local_n0,&local_0_start);
745       n = (PetscInt)local_n0*partial_dim;
746 //      printf("New partial dimension is %d %d %d",n,N,ndim);
747       ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr);
748       break;
749     }
750   }
751   ierr = PetscObjectChangeTypeName((PetscObject)A,MATFFTW);CHKERRQ(ierr);
752 
753   ierr = PetscNewLog(A,Mat_FFTW,&fftw);CHKERRQ(ierr);
754   fft->data = (void*)fftw;
755 
756   fft->n           = n;
757   fftw->ndim_fftw  = (ptrdiff_t)ndim; // This is dimension of fft
758   ierr = PetscMalloc(ndim*sizeof(ptrdiff_t), (ptrdiff_t *)&(fftw->dim_fftw));CHKERRQ(ierr);
759   for(ctr=0;ctr<ndim;ctr++) (fftw->dim_fftw)[ctr]=dim[ctr];
760 
761   fftw->p_forward  = 0;
762   fftw->p_backward = 0;
763   fftw->p_flag     = FFTW_ESTIMATE;
764 
765   if (size == 1){
766     A->ops->mult          = MatMult_SeqFFTW;
767     A->ops->multtranspose = MatMultTranspose_SeqFFTW;
768   } else {
769     A->ops->mult          = MatMult_MPIFFTW;
770     A->ops->multtranspose = MatMultTranspose_MPIFFTW;
771   }
772   fft->matdestroy          = MatDestroy_FFTW;
773   A->ops->getvecs       = MatGetVecs_FFTW;
774   A->assembled          = PETSC_TRUE;
775   printf("The code is coming here\n");
776 #if !defined(PETSC_USE_COMPLEX)
777   printf("The code is coming here\n");
778   PetscObjectComposeFunctionDynamic((PetscObject)A,"InputTransformFFT_C","InputTransformFFT_FFTW",InputTransformFFT_FFTW);
779 //  PetscObjectComposeFunctionDynamic((PetscObject)A,"OutputTransformFFT_C","OutputTransformFFT_FFTW",OutputTransformFFT_FFTW);
780 #endif
781 
782   /* get runtime options */
783   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"FFTW Options","Mat");CHKERRQ(ierr);
784     ierr = PetscOptionsEList("-mat_fftw_plannerflags","Planner Flags","None",p_flags,4,p_flags[0],&p_flag,&flg);CHKERRQ(ierr);
785     if (flg) {fftw->p_flag = (unsigned)p_flag;}
786   PetscOptionsEnd();
787   PetscFunctionReturn(0);
788 }
789 EXTERN_C_END
790 
791 
792 
793 
794