Actual source code: fft.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1: /*
  2:     Provides an interface to the FFT packages.
  3: */

  5:  #include <../src/mat/impls/fft/fft.h>

  7: PetscErrorCode MatDestroy_FFT(Mat A)
  8: {
 10:   Mat_FFT        *fft = (Mat_FFT*)A->data;

 13:   if (fft->matdestroy) {
 14:     (fft->matdestroy)(A);
 15:   }
 16:   PetscFree(fft->dim);
 17:   PetscFree(A->data);
 18:   PetscObjectChangeTypeName((PetscObject)A,0);
 19:   return(0);
 20: }

 22: /*@C
 23:       MatCreateFFT - Creates a matrix object that provides FFT via an external package

 25:    Collective

 27:    Input Parameter:
 28: +   comm - MPI communicator
 29: .   ndim - the ndim-dimensional transform
 30: .   dim - array of size ndim, dim[i] contains the vector length in the i-dimension
 31: -   type - package type, e.g., FFTW or MATSEQCUFFT

 33:    Output Parameter:
 34: .   A  - the matrix

 36:    Options Database Keys:
 37: .   -mat_fft_type - set FFT type fft or seqcufft

 39:    Note: this serves as a base class for all FFT marix classes, currently MATFFTW or MATSEQCUFFT

 41:    Level: intermediate

 43: .seealso: MatCreateVecsFFTW()
 44: @*/
 45: PetscErrorCode MatCreateFFT(MPI_Comm comm,PetscInt ndim,const PetscInt dim[],MatType mattype,Mat *A)
 46: {
 48:   PetscMPIInt    size;
 49:   Mat            FFT;
 50:   PetscInt       N,i;
 51:   Mat_FFT        *fft;

 54:   if (ndim < 1) SETERRQ1(comm,PETSC_ERR_USER,"ndim %d must be > 0",ndim);
 55:   MPI_Comm_size(comm, &size);

 57:   MatCreate(comm,&FFT);
 58:   PetscNewLog(FFT,&fft);
 59:   FFT->data = (void*)fft;
 60:   N         = 1;
 61:   for (i=0; i<ndim; i++) {
 62:     if (dim[i] < 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"dim[%d]=%d must be > 0",i,dim[i]);
 63:     N *= dim[i];
 64:   }

 66:   PetscMalloc1(ndim,&fft->dim);
 67:   PetscArraycpy(fft->dim,dim,ndim);

 69:   fft->ndim = ndim;
 70:   fft->n    = PETSC_DECIDE;
 71:   fft->N    = N;
 72:   fft->data = NULL;

 74:   MatSetType(FFT,mattype);

 76:   FFT->ops->destroy = MatDestroy_FFT;

 78:   /* get runtime options */
 79:   PetscOptionsBegin(PetscObjectComm((PetscObject)FFT),((PetscObject)FFT)->prefix,"FFT Options","Mat");
 80:   PetscOptionsEnd();

 82:   *A = FFT;
 83:   return(0);
 84: }