Actual source code: matusfft.c
petsc-3.6.1 2015-08-06
2: /*
3: Provides an implementation of the Unevenly Sampled FFT algorithm as a Mat.
4: Testing examples can be found in ~/src/mat/examples/tests FIX: should these be moved to dm/da/examples/tests?
5: */
7: #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
8: #include <petscdmda.h> /*I "petscdmda.h" I*/ /* Unlike equispaced FFT, USFFT requires geometric information encoded by a DMDA */
9: #include <fftw3.h>
11: typedef struct {
12: PetscInt dim;
13: Vec sampleCoords;
14: PetscInt dof;
15: DM freqDA; /* frequency DMDA */
16: PetscInt *freqSizes; /* sizes of the frequency DMDA, one per each dim */
17: DM resampleDa; /* the Battle-Lemarie interpolant DMDA */
18: Vec resample; /* Vec of samples, one per dof per sample point */
19: fftw_plan p_forward,p_backward;
20: unsigned p_flag; /* planner flags, FFTW_ESTIMATE,FFTW_MEASURE, FFTW_PATIENT, FFTW_EXHAUSTIVE */
21: } Mat_USFFT;
26: PetscErrorCode MatApply_USFFT_Private(Mat A, fftw_plan *plan, int direction, Vec x,Vec y)
27: {
28: #if 0
30: PetscScalar *r_array, *y_array;
31: Mat_USFFT* = (Mat_USFFT*)(A->data);
32: #endif
35: #if 0
36: /* resample x to usfft->resample */
37: MatResample_USFFT_Private(A, x);
39: /* NB: for now we use outdim for both x and y; this will change once a full USFFT is implemented */
40: VecGetArray(usfft->resample,&r_array);
41: VecGetArray(y,&y_array);
42: if (!*plan) { /* create a plan then execute it*/
43: if (usfft->dof == 1) {
44: #if defined(PETSC_DEBUG_USFFT)
45: PetscPrintf(PetscObjectComm((PetscObject)A), "direction = %d, usfft->ndim = %d\n", direction, usfft->ndim);
46: for (int ii = 0; ii < usfft->ndim; ++ii) {
47: PetscPrintf(PetscObjectComm((PetscObject)A), "usfft->outdim[%d] = %d\n", ii, usfft->outdim[ii]);
48: }
49: #endif
51: switch (usfft->dim) {
52: case 1:
53: *plan = fftw_plan_dft_1d(usfft->outdim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,direction,usfft->p_flag);
54: break;
55: case 2:
56: *plan = fftw_plan_dft_2d(usfft->outdim[0],usfft->outdim[1],(fftw_complex*)x_array,(fftw_complex*)y_array,direction,usfft->p_flag);
57: break;
58: case 3:
59: *plan = fftw_plan_dft_3d(usfft->outdim[0],usfft->outdim[1],usfft->outdim[2],(fftw_complex*)x_array,(fftw_complex*)y_array,direction,usfft->p_flag);
60: break;
61: default:
62: *plan = fftw_plan_dft(usfft->ndim,usfft->outdim,(fftw_complex*)x_array,(fftw_complex*)y_array,direction,usfft->p_flag);
63: break;
64: }
65: fftw_execute(*plan);
66: } /* if (dof == 1) */
67: else { /* if (dof > 1) */
68: *plan = fftw_plan_many_dft(/*rank*/usfft->ndim, /*n*/usfft->outdim, /*howmany*/usfft->dof,
69: (fftw_complex*)x_array, /*nembed*/usfft->outdim, /*stride*/usfft->dof, /*dist*/1,
70: (fftw_complex*)y_array, /*nembed*/usfft->outdim, /*stride*/usfft->dof, /*dist*/1,
71: /*sign*/direction, /*flags*/usfft->p_flag);
72: fftw_execute(*plan);
73: } /* if (dof > 1) */
74: } /* if (!*plan) */
75: else { /* if (*plan) */
76: /* use existing plan */
77: fftw_execute_dft(*plan,(fftw_complex*)x_array,(fftw_complex*)y_array);
78: }
79: VecRestoreArray(y,&y_array);
80: VecRestoreArray(x,&x_array);
81: #endif
82: return(0);
83: } /* MatApply_USFFT_Private() */
85: #if 0
88: PetscErrorCode Mat_USFFT_ProjectOnBattleLemarie_Private(Vec x,double *r)
89: /* Project onto the Battle-Lemarie function centered around r */
90: {
92: PetscScalar *x_array, *y_array;
95: return(0);
96: } /* Mat_USFFT_ProjectOnBattleLemarie_Private() */
100: PetscErrorCode MatInterpolate_USFFT_Private(Vec x,Vec y)
101: {
103: PetscScalar *x_array, *y_array;
106: return(0);
107: } /* MatInterpolate_USFFT_Private() */
112: PetscErrorCode MatMult_SeqUSFFT(Mat A,Vec x,Vec y)
113: {
115: Mat_USFFT *usfft = (Mat_USFFT*)A->data;
118: /* NB: for now we use outdim for both x and y; this will change once a full USFFT is implemented */
119: MatApply_USFFT_Private(A, &usfft->p_forward, FFTW_FORWARD, x,y);
120: return(0);
121: }
125: PetscErrorCode MatMultTranspose_SeqUSFFT(Mat A,Vec x,Vec y)
126: {
128: Mat_USFFT *usfft = (Mat_USFFT*)A->data;
131: /* NB: for now we use outdim for both x and y; this will change once a full USFFT is implemented */
132: MatApply_USFFT_Private(usfft, &usfft->p_backward, FFTW_BACKWARD, x,y);
133: return(0);
134: }
138: PetscErrorCode MatDestroy_SeqUSFFT(Mat A)
139: {
140: Mat_USFFT *usfft = (Mat_USFFT*)A->data;
144: fftw_destroy_plan(usfft->p_forward);
145: fftw_destroy_plan(usfft->p_backward);
146: PetscFree(usfft->indim);
147: PetscFree(usfft->outdim);
148: PetscFree(usfft);
149: PetscObjectChangeTypeName((PetscObject)A,0);
150: return(0);
151: } /* MatDestroy_SeqUSFFT() */
156: /*@C
157: MatCreateSeqUSFFT - Creates a matrix object that provides sequential USFFT
158: via the external package FFTW
160: Collective on MPI_Comm
162: Input Parameter:
163: + da - geometry of the domain encoded by a DMDA
165: Output Parameter:
166: . A - the matrix
168: Options Database Keys:
169: + -mat_usfft_plannerflags - set the FFTW planner flags
171: Level: intermediate
173: @*/
174: PetscErrorCode MatCreateSeqUSFFT(Vec sampleCoords, DMDA freqDA, Mat *A)
175: {
177: Mat_USFFT *usfft;
178: PetscInt m,n,M,N,i;
179: const char *p_flags[]={"FFTW_ESTIMATE","FFTW_MEASURE","FFTW_PATIENT","FFTW_EXHAUSTIVE"};
180: PetscBool flg;
181: PetscInt p_flag;
182: PetscInt dof, dim, freqSizes[3];
183: MPI_Comm comm;
184: PetscInt size;
187: PetscObjectGetComm((PetscObject)inda, &comm);
188: MPI_Comm_size(comm, &size);
189: if (size > 1) SETERRQ(comm,PETSC_ERR_USER, "Parallel DMDA (in) not yet supported by USFFT");
190: PetscObjectGetComm((PetscObject)outda, &comm);
191: MPI_Comm_size(comm, &size);
192: if (size > 1) SETERRQ(comm,PETSC_ERR_USER, "Parallel DMDA (out) not yet supported by USFFT");
193: MatCreate(comm,A);
194: PetscNewLog(*A,&usfft);
195: (*A)->data = (void*)usfft;
196: usfft->inda = inda;
197: usfft->outda = outda;
198: /* inda */
199: DMDAGetInfo(usfft->inda, &ndim, dim+0, dim+1, dim+2, NULL, NULL, NULL, &dof, NULL, NULL, NULL);
200: if (ndim <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"ndim %d must be > 0",ndim);
201: if (dof <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"dof %d must be > 0",dof);
202: usfft->ndim = ndim;
203: usfft->dof = dof;
204: usfft->freqDA = freqDA;
205: /* NB: we reverse the freq and resample DMDA sizes, since the DMDA ordering (natural on x-y-z, with x varying the fastest)
206: is the order opposite of that assumed by FFTW: z varying the fastest */
207: PetscMalloc1(usfft->ndim+1,&usfft->indim);
208: for (i = usfft->ndim; i > 0; --i) usfft->indim[usfft->ndim-i] = dim[i-1];
210: /* outda */
211: DMDAGetInfo(usfft->outda, &ndim, dim+0, dim+1, dim+2, NULL, NULL, NULL, &dof, NULL, NULL, NULL);
212: if (ndim != usfft->ndim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"in and out DMDA dimensions must match: %d != %d",usfft->ndim, ndim);
213: if (dof != usfft->dof) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"in and out DMDA dof must match: %d != %d",usfft->dof, dof);
214: /* Store output dimensions */
215: /* NB: we reverse the DMDA dimensions, since the DMDA ordering (natural on x-y-z, with x varying the fastest)
216: is the order opposite of that assumed by FFTW: z varying the fastest */
217: PetscMalloc1(usfft->ndim+1,&usfft->outdim);
218: for (i = usfft->ndim; i > 0; --i) usfft->outdim[usfft->ndim-i] = dim[i-1];
220: /* TODO: Use the new form of DMDACreate() */
221: #if 0
222: DMDACreate(comm,usfft->dim, DMDA_NONPERIODIC, DMDA_STENCIL_STAR, usfft->freqSizes[0], usfft->freqSizes[1], usfft->freqSizes[2],
223: PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, dof, 0, NULL, NULL, NULL, 0, &(usfft->resampleDA));
224: #endif
225: DMDAGetVec(usfft->resampleDA, usfft->resample);
228: /* CONTINUE: Need to build the connectivity "Sieve" attaching sample points to the resample points they are close to */
230: /* CONTINUE: recalculate matrix sizes based on the connectivity "Sieve" */
231: /* mat sizes */
232: m = 1; n = 1;
233: for (i=0; i<usfft->ndim; i++) {
234: if (usfft->indim[i] <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"indim[%d]=%d must be > 0",i,usfft->indim[i]);
235: if (usfft->outdim[i] <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"outdim[%d]=%d must be > 0",i,usfft->outdim[i]);
236: n *= usfft->indim[i];
237: m *= usfft->outdim[i];
238: }
239: N = n*usfft->dof;
240: M = m*usfft->dof;
241: MatSetSizes(*A,M,N,M,N); /* "in size" is the number of columns, "out size" is the number of rows" */
242: PetscObjectChangeTypeName((PetscObject)*A,MATSEQUSFFT);
243: usfft->m = m; usfft->n = n; usfft->M = M; usfft->N = N;
244: /* FFTW */
245: usfft->p_forward = 0;
246: usfft->p_backward = 0;
247: usfft->p_flag = FFTW_ESTIMATE;
248: /* set Mat ops */
249: (*A)->ops->mult = MatMult_SeqUSFFT;
250: (*A)->ops->multtranspose = MatMultTranspose_SeqUSFFT;
251: (*A)->assembled = PETSC_TRUE;
252: (*A)->ops->destroy = MatDestroy_SeqUSFFT;
253: /* get runtime options */
254: PetscOptionsBegin(((PetscObject)(*A))->comm,((PetscObject)(*A))->prefix,"USFFT Options","Mat");
255: PetscOptionsEList("-mat_usfft_fftw_plannerflags","Planner Flags","None",p_flags,4,p_flags[0],&p_flag,&flg);
256: if (flg) usfft->p_flag = (unsigned)p_flag;
257: PetscOptionsEnd();
258: return(0);
259: } /* MatCreateSeqUSFFT() */
261: #endif