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