Actual source code: fft.c
1: /*
2: Provides an interface to the FFT packages.
3: */
5: #include <../src/mat/impls/fft/fft.h>
7: static PetscErrorCode MatDestroy_FFT(Mat A)
8: {
9: Mat_FFT *fft = (Mat_FFT *)A->data;
11: PetscFunctionBegin;
12: if (fft->matdestroy) PetscCall(fft->matdestroy(A));
13: PetscCall(PetscFree(fft->dim));
14: PetscCall(PetscFree(A->data));
15: PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
16: PetscFunctionReturn(PETSC_SUCCESS);
17: }
19: /*@
20: MatCreateFFT - Creates a matrix object that provides FFT via an external package
22: Collective
24: Input Parameters:
25: + comm - MPI communicator
26: . ndim - the ndim-dimensional transform
27: . dim - array of size ndim, dim[i] contains the vector length in the i-dimension
28: - mattype - package type, e.g., `MATFFTW` or `MATSEQCUFFT`
30: Output Parameter:
31: . A - the matrix
33: Options Database Key:
34: . -mat_fft_type - set FFT type fft or seqcufft
36: Level: intermediate
38: Note:
39: This serves as a base class for all FFT matrix classes, currently `MATFFTW` or `MATSEQCUFFT`
41: .seealso: [](ch_matrices), `Mat`, `MATFFTW`, `MATSEQCUFFT`, `MatCreateVecsFFTW()`
42: @*/
43: PetscErrorCode MatCreateFFT(MPI_Comm comm, PetscInt ndim, const PetscInt dim[], MatType mattype, Mat *A)
44: {
45: PetscMPIInt size;
46: Mat FFT;
47: PetscInt N, i;
48: Mat_FFT *fft;
50: PetscFunctionBegin;
51: PetscAssertPointer(dim, 3);
52: PetscAssertPointer(A, 5);
53: PetscCheck(ndim >= 1, comm, PETSC_ERR_USER, "ndim %" PetscInt_FMT " must be > 0", ndim);
54: PetscCallMPI(MPI_Comm_size(comm, &size));
56: PetscCall(MatCreate(comm, &FFT));
57: PetscCall(PetscNew(&fft));
58: FFT->data = (void *)fft;
59: N = 1;
60: for (i = 0; i < ndim; i++) {
61: PetscCheck(dim[i] >= 1, PETSC_COMM_SELF, PETSC_ERR_USER, "dim[%" PetscInt_FMT "]=%" PetscInt_FMT " must be > 0", i, dim[i]);
62: N *= dim[i];
63: }
65: PetscCall(PetscMalloc1(ndim, &fft->dim));
66: PetscCall(PetscArraycpy(fft->dim, dim, ndim));
68: fft->ndim = ndim;
69: fft->n = PETSC_DECIDE;
70: fft->N = N;
71: fft->data = NULL;
73: PetscCall(MatSetType(FFT, mattype));
75: FFT->ops->destroy = MatDestroy_FFT;
77: /* get runtime options... what options? */
78: PetscObjectOptionsBegin((PetscObject)FFT);
79: PetscOptionsEnd();
81: *A = FFT;
82: PetscFunctionReturn(PETSC_SUCCESS);
83: }