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 */ 672 ierr = ISCreateStride(PETSC_COMM_SELF,N,0,1,&list1); 673 ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr); 674 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 675 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 676 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 677 ierr = ISDestroy(&list1);CHKERRQ(ierr); 678 //ierr = ISDestroy(list1);CHKERRQ(ierr); 679 // break; 680 // } 681 } 682 683 else{ 684 685 switch (ndim){ 686 case 1: 687 SETERRQ(comm,PETSC_ERR_SUP,"FFTW does not support parallel 1D real transform"); 688 break; 689 case 2: 690 #if defined(PETSC_USE_COMPLEX) 691 //PetscInt my_value; 692 alloc_local = fftw_mpi_local_size_2d(dim[0],dim[1],comm,&local_n0,&local_0_start); 693 ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1],local_0_start*dim[1],1,&list1); 694 ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1],low,1,&list2); 695 //ierr = ISView(list1,PETSC_VIEWER_STDOUT_WORLD); 696 //ierr = ISView(list2,PETSC_VIEWER_STDOUT_WORLD); 697 ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr); 698 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 699 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 700 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 701 ierr = ISDestroy(&list1);CHKERRQ(ierr); 702 ierr = ISDestroy(&list2);CHKERRQ(ierr); 703 break; 704 #else 705 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); 706 N1 = 2*dim[0]*(dim[1]/2+1); n1 = 2*local_n0*(dim[1]/2+1); 707 708 ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx1);CHKERRQ(ierr); 709 ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx2);CHKERRQ(ierr); 710 printf("Val local_0_start = %ld\n",local_0_start); 711 712 if (dim[1]%2==0) 713 NM = dim[1]+2; 714 else 715 NM = dim[1]+1; 716 717 for (i=0;i<local_n0;i++){ 718 for (j=0;j<dim[1];j++){ 719 tempindx = i*dim[1] + j; 720 tempindx1 = i*NM + j; 721 indx1[tempindx]=local_0_start*dim[1]+tempindx; 722 indx2[tempindx]=low+tempindx1; 723 // printf("Val tempindx1 = %d\n",tempindx1); 724 // printf("index1 %d from proc %d is \n",indx1[tempindx],rank); 725 // printf("index2 %d from proc %d is \n",indx2[tempindx],rank); 726 // printf("-------------------------\n",indx2[tempindx],rank); 727 } 728 } 729 730 ierr = ISCreateGeneral(comm,local_n0*dim[1],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr); 731 ierr = ISCreateGeneral(comm,local_n0*dim[1],indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr); 732 733 ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr); 734 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 735 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 736 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 737 ierr = ISDestroy(&list1);CHKERRQ(ierr); 738 ierr = ISDestroy(&list2);CHKERRQ(ierr); 739 ierr = PetscFree(indx1);CHKERRQ(ierr); 740 ierr = PetscFree(indx2);CHKERRQ(ierr); 741 break; 742 #endif 743 744 case 3: 745 #if defined(PETSC_USE_COMPLEX) 746 alloc_local = fftw_mpi_local_size_3d(dim[0],dim[1],dim[2],comm,&local_n0,&local_0_start); 747 ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1]*dim[2],local_0_start*dim[1]*dim[2],1,&list1); 748 ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1]*dim[2],low,1,&list2); 749 //ierr = ISView(list1,PETSC_VIEWER_STDOUT_WORLD); 750 //ierr = ISView(list2,PETSC_VIEWER_STDOUT_WORLD); 751 ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr); 752 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 753 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 754 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 755 ierr = ISDestroy(&list1);CHKERRQ(ierr); 756 ierr = ISDestroy(&list2);CHKERRQ(ierr); 757 break; 758 #else 759 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); 760 N1 = 2*dim[0]*dim[1]*(dim[2]/2+1); n1 = 2*local_n0*dim[1]*(dim[2]/2+1); 761 762 ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx1);CHKERRQ(ierr); 763 ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx2);CHKERRQ(ierr); 764 printf("Val local_0_start = %ld\n",local_0_start); 765 766 if (dim[2]%2==0) 767 NM = dim[2]+2; 768 else 769 NM = dim[2]+1; 770 771 for (i=0;i<local_n0;i++){ 772 for (j=0;j<dim[1];j++){ 773 for (k=0;k<dim[2];k++){ 774 tempindx = i*dim[1]*dim[2] + j*dim[2] + k; 775 tempindx1 = i*dim[1]*NM + j*NM + k; 776 indx1[tempindx]=local_0_start*dim[1]*dim[2]+tempindx; 777 indx2[tempindx]=low+tempindx1; 778 } 779 // printf("Val tempindx1 = %d\n",tempindx1); 780 // printf("index1 %d from proc %d is \n",indx1[tempindx],rank); 781 // printf("index2 %d from proc %d is \n",indx2[tempindx],rank); 782 // printf("-------------------------\n",indx2[tempindx],rank); 783 } 784 } 785 786 ierr = ISCreateGeneral(comm,local_n0*dim[1]*dim[2],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr); 787 ierr = ISCreateGeneral(comm,local_n0*dim[1]*dim[2],indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr); 788 789 ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr); 790 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 791 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 792 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 793 ierr = ISDestroy(&list1);CHKERRQ(ierr); 794 ierr = ISDestroy(&list2);CHKERRQ(ierr); 795 ierr = PetscFree(indx1);CHKERRQ(ierr); 796 ierr = PetscFree(indx2);CHKERRQ(ierr); 797 break; 798 #endif 799 800 default: 801 #if defined(PETSC_USE_COMPLEX) 802 alloc_local = fftw_mpi_local_size(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start); 803 ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*(fftw->partial_dim),local_0_start*(fftw->partial_dim),1,&list1); 804 ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*(fftw->partial_dim),low,1,&list2); 805 //ierr = ISView(list1,PETSC_VIEWER_STDOUT_WORLD); 806 //ierr = ISView(list2,PETSC_VIEWER_STDOUT_WORLD); 807 ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr); 808 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 809 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 810 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 811 ierr = ISDestroy(&list1);CHKERRQ(ierr); 812 ierr = ISDestroy(&list2);CHKERRQ(ierr); 813 814 815 #else 816 temp = (fftw->dim_fftw)[fftw->ndim_fftw-1]; 817 printf("The value of temp is %ld\n",temp); 818 (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp/2 + 1; 819 alloc_local = fftw_mpi_local_size_transposed(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start,&local_n1,&local_1_start); 820 N1 = 2*N*(PetscInt)((fftw->dim_fftw)[fftw->ndim_fftw-1])/((PetscInt) temp); 821 (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp; 822 823 partial_dim = fftw->partial_dim; 824 printf("The value of partial dim is %d\n",partial_dim); 825 826 ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*partial_dim,&indx1);CHKERRQ(ierr); 827 ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*partial_dim,&indx2);CHKERRQ(ierr); 828 printf("Val local_0_start = %ld\n",local_0_start); 829 830 if (dim[ndim-1]%2==0) 831 NM = 2; 832 else 833 NM = 1; 834 835 j = low; 836 for (i=0,k=1; i<((PetscInt)local_n0)*partial_dim;i++,k++) 837 { 838 indx1[i] = local_0_start*partial_dim + i; 839 indx2[i] = j; 840 //printf("The values are %d and %d from %d\n",indx1[i],indx2[i],rank); 841 if (k%dim[ndim-1]==0) 842 { j+=NM;} 843 j++; 844 } 845 ierr = ISCreateGeneral(comm,local_n0*partial_dim,indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr); 846 ierr = ISCreateGeneral(comm,local_n0*partial_dim,indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr); 847 ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr); 848 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 849 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 850 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 851 ierr = ISDestroy(&list1);CHKERRQ(ierr); 852 ierr = ISDestroy(&list2);CHKERRQ(ierr); 853 ierr = PetscFree(indx1);CHKERRQ(ierr); 854 ierr = PetscFree(indx2);CHKERRQ(ierr); 855 break; 856 #endif 857 } 858 } 859 860 return 0; 861 } 862 EXTERN_C_END 863 864 #undef __FUNCT__ 865 #define __FUNCT__ "OutputTransformFFT" 866 PetscErrorCode OutputTransformFFT(Mat A,Vec x,Vec y) 867 { 868 PetscErrorCode ierr; 869 PetscFunctionBegin; 870 ierr = PetscTryMethod(A,"OutputTransformFFT_C",(Mat,Vec,Vec),(A,x,y));CHKERRQ(ierr); 871 PetscFunctionReturn(0); 872 } 873 874 /* 875 OutputTransformFFT_FFTW - Copies the FFTW output to the PETSc vector that user can use 876 Input A, x, y 877 A - FFTW matrix 878 x - FFTW vector 879 y - PETSc vector that user can read 880 Options Database Keys: 881 + -mat_fftw_plannerflags - set FFTW planner flags 882 883 Level: intermediate 884 885 */ 886 887 EXTERN_C_BEGIN 888 #undef __FUNCT__ 889 #define __FUNCT__ "OutputTransformFFT_FTTW" 890 PetscErrorCode OutputTransformFFT_FFTW(Mat A,Vec x,Vec y) 891 { 892 PetscErrorCode ierr; 893 MPI_Comm comm=((PetscObject)A)->comm; 894 Mat_FFT *fft = (Mat_FFT*)A->data; 895 Mat_FFTW *fftw = (Mat_FFTW*)fft->data; 896 PetscInt N=fft->N, N1, n1 ,NM; 897 PetscInt ndim=fft->ndim,*dim=fft->dim;//n=fft->n; 898 PetscInt low, *indx1, *indx2, tempindx, tempindx1; 899 PetscInt i,j,k,rank,size,partial_dim; 900 ptrdiff_t alloc_local,local_n0,local_0_start; 901 ptrdiff_t local_n1,local_1_start,temp; 902 VecScatter vecscat; 903 IS list1,list2; 904 905 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 906 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 907 ierr = VecGetOwnershipRange(x,&low,PETSC_NULL); 908 909 if (size==1){ 910 /* 911 switch (ndim){ 912 case 1: 913 ierr = PetscMalloc(sizeof(PetscInt)*dim[0],&indx1);CHKERRQ(ierr); 914 for (i=0;i<dim[0];i++) 915 { 916 indx1[i] = i; 917 } 918 ierr = ISCreateGeneral(comm,dim[0],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr); 919 ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr); 920 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 921 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 922 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 923 ierr = ISDestroy(&list1);CHKERRQ(ierr); 924 ierr = PetscFree(indx1);CHKERRQ(ierr); 925 break; 926 927 case 2: 928 ierr = PetscMalloc(sizeof(PetscInt)*dim[0]*dim[1],&indx1);CHKERRQ(ierr); 929 for (i=0;i<dim[0];i++){ 930 for (j=0;j<dim[1];j++){ 931 indx1[i*dim[1]+j] = i*dim[1] + j; 932 } 933 } 934 ierr = ISCreateGeneral(comm,dim[0]*dim[1],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr); 935 ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr); 936 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 937 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 938 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 939 ierr = ISDestroy(&list1);CHKERRQ(ierr); 940 ierr = PetscFree(indx1);CHKERRQ(ierr); 941 break; 942 case 3: 943 ierr = PetscMalloc(sizeof(PetscInt)*dim[0]*dim[1]*dim[2],&indx1);CHKERRQ(ierr); 944 for (i=0;i<dim[0];i++){ 945 for (j=0;j<dim[1];j++){ 946 for (k=0;k<dim[2];k++){ 947 indx1[i*dim[1]*dim[2]+j*dim[2]+k] = i*dim[1]*dim[2] + j*dim[2] + k; 948 } 949 } 950 } 951 ierr = ISCreateGeneral(comm,dim[0]*dim[1]*dim[2],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr); 952 ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr); 953 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 954 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 955 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 956 ierr = ISDestroy(&list1);CHKERRQ(ierr); 957 ierr = PetscFree(indx1);CHKERRQ(ierr); 958 break; 959 default: 960 */ 961 ierr = ISCreateStride(comm,N,0,1,&list1); 962 //ierr = ISView(list1,PETSC_VIEWER_STDOUT_SELF); 963 ierr = VecScatterCreate(x,list1,y,list1,&vecscat);CHKERRQ(ierr); 964 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 965 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 966 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 967 ierr = ISDestroy(&list1);CHKERRQ(ierr); 968 // break; 969 // } 970 } 971 else{ 972 973 switch (ndim){ 974 case 1: 975 SETERRQ(comm,PETSC_ERR_SUP,"No support for real parallel 1D FFT"); 976 break; 977 case 2: 978 #if defined(PETSC_USE_COMPLEX) 979 alloc_local = fftw_mpi_local_size_2d(dim[0],dim[1],comm,&local_n0,&local_0_start); 980 ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1],local_0_start*dim[1],1,&list1); 981 ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1],low,1,&list2); 982 ierr = VecScatterCreate(x,list2,y,list1,&vecscat);CHKERRQ(ierr); 983 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 984 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 985 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 986 ierr = ISDestroy(&list1);CHKERRQ(ierr); 987 ierr = ISDestroy(&list2);CHKERRQ(ierr); 988 break; 989 #else 990 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); 991 N1 = 2*dim[0]*(dim[1]/2+1); n1 = 2*local_n0*(dim[1]/2+1); 992 993 ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx1);CHKERRQ(ierr); 994 ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx2);CHKERRQ(ierr); 995 printf("Val local_0_start = %ld\n",local_0_start); 996 997 if (dim[1]%2==0) 998 NM = dim[1]+2; 999 else 1000 NM = dim[1]+1; 1001 1002 1003 1004 for (i=0;i<local_n0;i++){ 1005 for (j=0;j<dim[1];j++){ 1006 tempindx = i*dim[1] + j; 1007 tempindx1 = i*NM + j; 1008 indx1[tempindx]=local_0_start*dim[1]+tempindx; 1009 indx2[tempindx]=low+tempindx1; 1010 // printf("Val tempindx1 = %d\n",tempindx1); 1011 // printf("index1 %d from proc %d is \n",indx1[tempindx],rank); 1012 // printf("index2 %d from proc %d is \n",indx2[tempindx],rank); 1013 // printf("-------------------------\n",indx2[tempindx],rank); 1014 } 1015 } 1016 1017 ierr = ISCreateGeneral(comm,local_n0*dim[1],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr); 1018 ierr = ISCreateGeneral(comm,local_n0*dim[1],indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr); 1019 1020 ierr = VecScatterCreate(x,list2,y,list1,&vecscat);CHKERRQ(ierr); 1021 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1022 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1023 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 1024 ierr = ISDestroy(&list1);CHKERRQ(ierr); 1025 ierr = ISDestroy(&list2);CHKERRQ(ierr); 1026 ierr = PetscFree(indx1);CHKERRQ(ierr); 1027 ierr = PetscFree(indx2);CHKERRQ(ierr); 1028 break; 1029 #endif 1030 case 3: 1031 #if defined(PETSC_USE_COMPLEX) 1032 alloc_local = fftw_mpi_local_size_3d(dim[0],dim[1],dim[2],comm,&local_n0,&local_0_start); 1033 ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1]*dim[2],local_0_start*dim[1]*dim[2],1,&list1); 1034 ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*dim[1]*dim[2],low,1,&list2); 1035 //ierr = ISView(list1,PETSC_VIEWER_STDOUT_WORLD); 1036 //ierr = ISView(list2,PETSC_VIEWER_STDOUT_WORLD); 1037 ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr); 1038 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1039 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1040 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 1041 ierr = ISDestroy(&list1);CHKERRQ(ierr); 1042 ierr = ISDestroy(&list2);CHKERRQ(ierr); 1043 break; 1044 #else 1045 1046 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); 1047 N1 = 2*dim[0]*dim[1]*(dim[2]/2+1); n1 = 2*local_n0*dim[1]*(dim[2]/2+1); 1048 1049 ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx1);CHKERRQ(ierr); 1050 ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*N1,&indx2);CHKERRQ(ierr); 1051 printf("Val local_0_start = %ld\n",local_0_start); 1052 1053 if (dim[2]%2==0) 1054 NM = dim[2]+2; 1055 else 1056 NM = dim[2]+1; 1057 1058 for (i=0;i<local_n0;i++){ 1059 for (j=0;j<dim[1];j++){ 1060 for (k=0;k<dim[2];k++){ 1061 tempindx = i*dim[1]*dim[2] + j*dim[2] + k; 1062 tempindx1 = i*dim[1]*NM + j*NM + k; 1063 indx1[tempindx]=local_0_start*dim[1]*dim[2]+tempindx; 1064 indx2[tempindx]=low+tempindx1; 1065 } 1066 // printf("Val tempindx1 = %d\n",tempindx1); 1067 // printf("index1 %d from proc %d is \n",indx1[tempindx],rank); 1068 // printf("index2 %d from proc %d is \n",indx2[tempindx],rank); 1069 // printf("-------------------------\n",indx2[tempindx],rank); 1070 } 1071 } 1072 1073 ierr = ISCreateGeneral(comm,local_n0*dim[1]*dim[2],indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr); 1074 ierr = ISCreateGeneral(comm,local_n0*dim[1]*dim[2],indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr); 1075 1076 ierr = VecScatterCreate(x,list2,y,list1,&vecscat);CHKERRQ(ierr); 1077 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1078 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1079 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 1080 ierr = ISDestroy(&list1);CHKERRQ(ierr); 1081 ierr = ISDestroy(&list2);CHKERRQ(ierr); 1082 ierr = PetscFree(indx1);CHKERRQ(ierr); 1083 ierr = PetscFree(indx2);CHKERRQ(ierr); 1084 break; 1085 #endif 1086 default: 1087 #if defined(PETSC_USE_COMPLEX) 1088 alloc_local = fftw_mpi_local_size(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start); 1089 ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*(fftw->partial_dim),local_0_start*(fftw->partial_dim),1,&list1); 1090 ierr = ISCreateStride(PETSC_COMM_WORLD,local_n0*(fftw->partial_dim),low,1,&list2); 1091 //ierr = ISView(list1,PETSC_VIEWER_STDOUT_WORLD); 1092 //ierr = ISView(list2,PETSC_VIEWER_STDOUT_WORLD); 1093 ierr = VecScatterCreate(x,list1,y,list2,&vecscat);CHKERRQ(ierr); 1094 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1095 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1096 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 1097 ierr = ISDestroy(&list1);CHKERRQ(ierr); 1098 ierr = ISDestroy(&list2);CHKERRQ(ierr); 1099 #else 1100 temp = (fftw->dim_fftw)[fftw->ndim_fftw-1]; 1101 printf("The value of temp is %ld\n",temp); 1102 (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp/2 + 1; 1103 alloc_local = fftw_mpi_local_size_transposed(fftw->ndim_fftw,fftw->dim_fftw,comm,&local_n0,&local_0_start,&local_n1,&local_1_start); 1104 N1 = 2*N*(PetscInt)((fftw->dim_fftw)[fftw->ndim_fftw-1])/((PetscInt) temp); 1105 (fftw->dim_fftw)[fftw->ndim_fftw-1] = temp; 1106 1107 partial_dim = fftw->partial_dim; 1108 printf("The value of partial dim is %d\n",partial_dim); 1109 1110 ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*partial_dim,&indx1);CHKERRQ(ierr); 1111 ierr = PetscMalloc(sizeof(PetscInt)*((PetscInt)local_n0)*partial_dim,&indx2);CHKERRQ(ierr); 1112 printf("Val local_0_start = %ld\n",local_0_start); 1113 1114 if (dim[ndim-1]%2==0) 1115 NM = 2; 1116 else 1117 NM = 1; 1118 1119 j = low; 1120 for (i=0,k=1; i<((PetscInt)local_n0)*partial_dim;i++,k++) 1121 { 1122 indx1[i] = local_0_start*partial_dim + i; 1123 indx2[i] = j; 1124 //printf("The values are %d and %d from %d\n",indx1[i],indx2[i],rank); 1125 if (k%dim[ndim-1]==0) 1126 { j+=NM;} 1127 j++; 1128 } 1129 ierr = ISCreateGeneral(comm,local_n0*partial_dim,indx1,PETSC_COPY_VALUES,&list1);CHKERRQ(ierr); 1130 ierr = ISCreateGeneral(comm,local_n0*partial_dim,indx2,PETSC_COPY_VALUES,&list2);CHKERRQ(ierr); 1131 1132 //ISView(list1,PETSC_VIEWER_STDOUT_SELF); 1133 1134 1135 ierr = VecScatterCreate(x,list2,y,list1,&vecscat);CHKERRQ(ierr); 1136 ierr = VecScatterBegin(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1137 ierr = VecScatterEnd(vecscat,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1138 ierr = VecScatterDestroy(&vecscat);CHKERRQ(ierr); 1139 ierr = ISDestroy(&list1);CHKERRQ(ierr); 1140 ierr = ISDestroy(&list2);CHKERRQ(ierr); 1141 ierr = PetscFree(indx1);CHKERRQ(ierr); 1142 ierr = PetscFree(indx2);CHKERRQ(ierr); 1143 1144 break; 1145 #endif 1146 } 1147 } 1148 return 0; 1149 } 1150 EXTERN_C_END 1151 1152 EXTERN_C_BEGIN 1153 #undef __FUNCT__ 1154 #define __FUNCT__ "MatCreate_FFTW" 1155 /* 1156 MatCreate_FFTW - Creates a matrix object that provides FFT 1157 via the external package FFTW 1158 Options Database Keys: 1159 + -mat_fftw_plannerflags - set FFTW planner flags 1160 1161 Level: intermediate 1162 1163 */ 1164 1165 PetscErrorCode MatCreate_FFTW(Mat A) 1166 { 1167 PetscErrorCode ierr; 1168 MPI_Comm comm=((PetscObject)A)->comm; 1169 Mat_FFT *fft=(Mat_FFT*)A->data; 1170 Mat_FFTW *fftw; 1171 PetscInt n=fft->n,N=fft->N,ndim=fft->ndim,*dim = fft->dim; 1172 const char *p_flags[]={"FFTW_ESTIMATE","FFTW_MEASURE","FFTW_PATIENT","FFTW_EXHAUSTIVE"}; 1173 PetscBool flg; 1174 PetscInt p_flag,partial_dim=1,ctr,N1; 1175 PetscMPIInt size,rank; 1176 ptrdiff_t *pdim, temp; 1177 ptrdiff_t local_n1,local_1_start,local_1_end; 1178 1179 PetscFunctionBegin; 1180 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 1181 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1182 ierr = MPI_Barrier(PETSC_COMM_WORLD); 1183 1184 pdim = (ptrdiff_t *)calloc(ndim,sizeof(ptrdiff_t)); 1185 pdim[0] = dim[0]; 1186 for(ctr=1;ctr<ndim;ctr++) 1187 { 1188 partial_dim *= dim[ctr]; 1189 pdim[ctr] = dim[ctr]; 1190 } 1191 1192 if (size == 1) { 1193 #if defined(PETSC_USE_COMPLEX) 1194 ierr = MatSetSizes(A,N,N,N,N);CHKERRQ(ierr); 1195 n = N; 1196 #else 1197 int tot_dim = N*2*(dim[ndim-1]/2+1)/dim[ndim-1]; 1198 ierr = MatSetSizes(A,tot_dim,tot_dim,tot_dim,tot_dim);CHKERRQ(ierr); 1199 n = tot_dim; 1200 #endif 1201 1202 } else { 1203 ptrdiff_t alloc_local,local_n0,local_0_start;//local_n1,local_1_end; 1204 switch (ndim){ 1205 case 1: 1206 #if !defined(PETSC_USE_COMPLEX) 1207 SETERRQ(comm,PETSC_ERR_SUP,"FFTW does not support parallel 1D real transform"); 1208 #else 1209 alloc_local = fftw_mpi_local_size_1d(dim[0],comm,FFTW_FORWARD,FFTW_ESTIMATE,&local_n0,&local_0_start,&local_n1,&local_1_end); 1210 n = (PetscInt)local_n0; 1211 ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr); 1212 #endif 1213 // PetscObjectComposeFunctionDynamic((PetscObject)A,"MatGetVecs1DC_C","MatGetVecs1DC_FFTW",MatGetVecs1DC_FFTW); 1214 break; 1215 case 2: 1216 #if defined(PETSC_USE_COMPLEX) 1217 alloc_local = fftw_mpi_local_size_2d(dim[0],dim[1],comm,&local_n0,&local_0_start); 1218 /* 1219 PetscMPIInt rank; 1220 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]); 1221 PetscSynchronizedFlush(comm); 1222 */ 1223 n = (PetscInt)local_n0*dim[1]; 1224 ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr); 1225 #else 1226 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); 1227 n = 2*(PetscInt)local_n0*(dim[1]/2+1); 1228 ierr = MatSetSizes(A,n,n,2*dim[0]*(dim[1]/2+1),2*dim[0]*(dim[1]/2+1)); 1229 #endif 1230 break; 1231 case 3: 1232 // printf("The value of alloc local is %d",alloc_local); 1233 #if defined(PETSC_USE_COMPLEX) 1234 alloc_local = fftw_mpi_local_size_3d(dim[0],dim[1],dim[2],comm,&local_n0,&local_0_start); 1235 n = (PetscInt)local_n0*dim[1]*dim[2]; 1236 ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr); 1237 #else 1238 printf("The code comes here\n"); 1239 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); 1240 n = 2*(PetscInt)local_n0*dim[1]*(dim[2]/2+1); 1241 ierr = MatSetSizes(A,n,n,2*dim[0]*dim[1]*(dim[2]/2+1),2*dim[0]*dim[1]*(dim[2]/2+1)); 1242 #endif 1243 break; 1244 default: 1245 #if defined(PETSC_USE_COMPLEX) 1246 alloc_local = fftw_mpi_local_size(ndim,pdim,comm,&local_n0,&local_0_start); 1247 // printf("The value of alloc local is %ld from process %d\n",alloc_local,rank); 1248 // alloc_local = fftw_mpi_local_size(ndim,dim,comm,&local_n0,&local_0_start); 1249 n = (PetscInt)local_n0*partial_dim; 1250 // printf("New partial dimension is %d %d %d",n,N,ndim); 1251 ierr = MatSetSizes(A,n,n,N,N);CHKERRQ(ierr); 1252 #else 1253 temp = pdim[ndim-1]; 1254 pdim[ndim-1]= temp/2 + 1; 1255 printf("For Multi dim case temp = %ld, pdim[ndim-1] = %ld\n",temp,pdim[ndim-1]); 1256 alloc_local = fftw_mpi_local_size_transposed(ndim,pdim,PETSC_COMM_WORLD,&local_n0,&local_0_start,&local_n1,&local_1_start); 1257 n = 2*(PetscInt)local_n0*partial_dim*pdim[ndim-1]/temp; 1258 N1 = 2*N*(PetscInt)pdim[ndim-1]/((PetscInt) temp); 1259 pdim[ndim-1] = temp; 1260 printf("For Multi dim case n = %d, N1 = %d\n",n,N1); 1261 ierr = MatSetSizes(A,n,n,N1,N1);CHKERRQ(ierr); 1262 #endif 1263 break; 1264 } 1265 } 1266 ierr = PetscObjectChangeTypeName((PetscObject)A,MATFFTW);CHKERRQ(ierr); 1267 ierr = PetscNewLog(A,Mat_FFTW,&fftw);CHKERRQ(ierr); 1268 fft->data = (void*)fftw; 1269 1270 fft->n = n; 1271 fftw->ndim_fftw = (ptrdiff_t)ndim; // This is dimension of fft 1272 fftw->partial_dim = partial_dim; 1273 ierr = PetscMalloc(ndim*sizeof(ptrdiff_t), (ptrdiff_t *)&(fftw->dim_fftw));CHKERRQ(ierr); 1274 for(ctr=0;ctr<ndim;ctr++) (fftw->dim_fftw)[ctr]=dim[ctr]; 1275 1276 fftw->p_forward = 0; 1277 fftw->p_backward = 0; 1278 fftw->p_flag = FFTW_ESTIMATE; 1279 1280 if (size == 1){ 1281 A->ops->mult = MatMult_SeqFFTW; 1282 A->ops->multtranspose = MatMultTranspose_SeqFFTW; 1283 } else { 1284 A->ops->mult = MatMult_MPIFFTW; 1285 A->ops->multtranspose = MatMultTranspose_MPIFFTW; 1286 } 1287 fft->matdestroy = MatDestroy_FFTW; 1288 // if(ndim=1 && size>1) and also if it is complex then getvecs should be attached to MatGetVecs_FFTW1D 1289 A->ops->getvecs = MatGetVecs_FFTW; 1290 A->assembled = PETSC_TRUE; 1291 PetscObjectComposeFunctionDynamic((PetscObject)A,"InputTransformFFT_C","InputTransformFFT_FFTW",InputTransformFFT_FFTW); 1292 PetscObjectComposeFunctionDynamic((PetscObject)A,"OutputTransformFFT_C","OutputTransformFFT_FFTW",OutputTransformFFT_FFTW); 1293 1294 /* get runtime options */ 1295 ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"FFTW Options","Mat");CHKERRQ(ierr); 1296 ierr = PetscOptionsEList("-mat_fftw_plannerflags","Planner Flags","None",p_flags,4,p_flags[0],&p_flag,&flg);CHKERRQ(ierr); 1297 if (flg) {fftw->p_flag = (unsigned)p_flag;} 1298 PetscOptionsEnd(); 1299 PetscFunctionReturn(0); 1300 } 1301 EXTERN_C_END 1302 1303 1304 1305 1306