1dedccee8SHong Zhang /*
2dedccee8SHong Zhang Provides an interface to the FFT packages.
3dedccee8SHong Zhang */
4dedccee8SHong Zhang
5c6db04a5SJed Brown #include <../src/mat/impls/fft/fft.h> /*I "petscmat.h" I*/
6dedccee8SHong Zhang
MatDestroy_FFT(Mat A)766976f2fSJacob Faibussowitsch static PetscErrorCode MatDestroy_FFT(Mat A)
8d71ae5a4SJacob Faibussowitsch {
9dedccee8SHong Zhang Mat_FFT *fft = (Mat_FFT *)A->data;
10dedccee8SHong Zhang
11dedccee8SHong Zhang PetscFunctionBegin;
12f4f49eeaSPierre Jolivet if (fft->matdestroy) PetscCall(fft->matdestroy(A));
139566063dSJacob Faibussowitsch PetscCall(PetscFree(fft->dim));
149566063dSJacob Faibussowitsch PetscCall(PetscFree(A->data));
159566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
17dedccee8SHong Zhang }
18dedccee8SHong Zhang
19*5d83a8b1SBarry Smith /*@
20dedccee8SHong Zhang MatCreateFFT - Creates a matrix object that provides FFT via an external package
21dedccee8SHong Zhang
22d083f849SBarry Smith Collective
23dedccee8SHong Zhang
24d8d19677SJose E. Roman Input Parameters:
25dedccee8SHong Zhang + comm - MPI communicator
26dedccee8SHong Zhang . ndim - the ndim-dimensional transform
27dedccee8SHong Zhang . dim - array of size ndim, dim[i] contains the vector length in the i-dimension
28fe59aa6dSJacob Faibussowitsch - mattype - package type, e.g., `MATFFTW` or `MATSEQCUFFT`
29dedccee8SHong Zhang
30dedccee8SHong Zhang Output Parameter:
31dedccee8SHong Zhang . A - the matrix
32dedccee8SHong Zhang
332ef1f0ffSBarry Smith Options Database Key:
3475f45d78SBarry Smith . -mat_fft_type - set FFT type fft or seqcufft
3575f45d78SBarry Smith
36dedccee8SHong Zhang Level: intermediate
37dedccee8SHong Zhang
3820f4b53cSBarry Smith Note:
39baca6076SPierre Jolivet This serves as a base class for all FFT matrix classes, currently `MATFFTW` or `MATSEQCUFFT`
4020f4b53cSBarry Smith
411cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATFFTW`, `MATSEQCUFFT`, `MatCreateVecsFFTW()`
42dedccee8SHong Zhang @*/
MatCreateFFT(MPI_Comm comm,PetscInt ndim,const PetscInt dim[],MatType mattype,Mat * A)43d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateFFT(MPI_Comm comm, PetscInt ndim, const PetscInt dim[], MatType mattype, Mat *A)
44d71ae5a4SJacob Faibussowitsch {
45dedccee8SHong Zhang PetscMPIInt size;
46dedccee8SHong Zhang Mat FFT;
47dedccee8SHong Zhang PetscInt N, i;
48dedccee8SHong Zhang Mat_FFT *fft;
49dedccee8SHong Zhang
50dedccee8SHong Zhang PetscFunctionBegin;
514f572ea9SToby Isaac PetscAssertPointer(dim, 3);
524f572ea9SToby Isaac PetscAssertPointer(A, 5);
535f80ce2aSJacob Faibussowitsch PetscCheck(ndim >= 1, comm, PETSC_ERR_USER, "ndim %" PetscInt_FMT " must be > 0", ndim);
549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
55dedccee8SHong Zhang
569566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &FFT));
574dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&fft));
58dedccee8SHong Zhang FFT->data = (void *)fft;
59dedccee8SHong Zhang N = 1;
60dedccee8SHong Zhang for (i = 0; i < ndim; i++) {
615f80ce2aSJacob Faibussowitsch PetscCheck(dim[i] >= 1, PETSC_COMM_SELF, PETSC_ERR_USER, "dim[%" PetscInt_FMT "]=%" PetscInt_FMT " must be > 0", i, dim[i]);
62dedccee8SHong Zhang N *= dim[i];
63dedccee8SHong Zhang }
64dedccee8SHong Zhang
659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ndim, &fft->dim));
669566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(fft->dim, dim, ndim));
67dedccee8SHong Zhang
68dedccee8SHong Zhang fft->ndim = ndim;
69dedccee8SHong Zhang fft->n = PETSC_DECIDE;
70dedccee8SHong Zhang fft->N = N;
710298fd71SBarry Smith fft->data = NULL;
72dedccee8SHong Zhang
739566063dSJacob Faibussowitsch PetscCall(MatSetType(FFT, mattype));
7426fbe8dcSKarl Rupp
75dedccee8SHong Zhang FFT->ops->destroy = MatDestroy_FFT;
76dedccee8SHong Zhang
775f80ce2aSJacob Faibussowitsch /* get runtime options... what options? */
78d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)FFT);
79d0609cedSBarry Smith PetscOptionsEnd();
80dedccee8SHong Zhang
81dedccee8SHong Zhang *A = FFT;
823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
83dedccee8SHong Zhang }
84