Actual source code: matusfft.c
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 <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_USFFT *usfft, fftw_plan *plan, int direction, Vec x,Vec y)
27: {
28: #if 0
30: PetscScalar *r_array, *y_array;
31: #endif
34: #if 0
35: /* resample x to usfft->resample */
36: MatResample_USFFT_Private(Mat_USFFT *usfft, x);
38: /* NB: for now we use outdim for both x and y; this will change once a full USFFT is implemented */
39: VecGetArray(usfft->resample,&r_array);
40: VecGetArray(y,&y_array);
41: if (!*plan){ /* create a plan then execute it*/
42: if(usfft->dof == 1) {
43: #ifdef PETSC_DEBUG_USFFT
44: PetscPrintf(PETSC_COMM_WORLD, "direction = %d, usfft->ndim = %d\n", direction, usfft->ndim);
45: for(int ii = 0; ii < usfft->ndim; ++ii) {
46: PetscPrintf(PETSC_COMM_WORLD, "usfft->outdim[%d] = %d\n", ii, usfft->outdim[ii]);
47: }
48: #endif
50: switch (usfft->dim){
51: case 1:
52: *plan = fftw_plan_dft_1d(usfft->outdim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,direction,usfft->p_flag);
53: break;
54: case 2:
55: *plan = fftw_plan_dft_2d(usfft->outdim[0],usfft->outdim[1],(fftw_complex*)x_array,(fftw_complex*)y_array,direction,usfft->p_flag);
56: break;
57: case 3:
58: *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);
59: break;
60: default:
61: *plan = fftw_plan_dft(usfft->ndim,usfft->outdim,(fftw_complex*)x_array,(fftw_complex*)y_array,direction,usfft->p_flag);
62: break;
63: }
64: fftw_execute(*plan);
65: }/* if(dof == 1) */
66: else { /* if(dof > 1) */
67: *plan = fftw_plan_many_dft(/*rank*/usfft->ndim, /*n*/usfft->outdim, /*howmany*/usfft->dof,
68: (fftw_complex*)x_array, /*nembed*/usfft->outdim, /*stride*/usfft->dof, /*dist*/1,
69: (fftw_complex*)y_array, /*nembed*/usfft->outdim, /*stride*/usfft->dof, /*dist*/1,
70: /*sign*/direction, /*flags*/usfft->p_flag);
71: fftw_execute(*plan);
72: }/* if(dof > 1) */
73: }/* if(!*plan) */
74: else { /* if(*plan) */
75: /* use existing plan */
76: fftw_execute_dft(*plan,(fftw_complex*)x_array,(fftw_complex*)y_array);
77: }
78: VecRestoreArray(y,&y_array);
79: VecRestoreArray(x,&x_array);
80: #endif
81: return(0);
82: }/* MatApply_USFFT_Private() */
84: #if 0
87: PetscErrorCode Mat_USFFT_ProjectOnBattleLemarie_Private(Vec x,double *r)
88: /* Project onto the Battle-Lemarie function centered around r */
89: {
91: PetscScalar *x_array, *y_array;
94: return(0);
95: }/* Mat_USFFT_ProjectOnBattleLemarie_Private() */
99: PetscErrorCode MatInterpolate_USFFT_Private(Vec x,Vec y)
100: {
102: PetscScalar *x_array, *y_array;
105: return(0);
106: }/* MatInterpolate_USFFT_Private() */
111: PetscErrorCode MatMult_SeqUSFFT(Mat A,Vec x,Vec y)
112: {
114: Mat_USFFT *usfft = (Mat_USFFT*)A->data;
117: /* NB: for now we use outdim for both x and y; this will change once a full USFFT is implemented */
118: MatApply_USFFT_Private(usfft, &usfft->p_forward, FFTW_FORWARD, x,y);
119: return(0);
120: }
124: PetscErrorCode MatMultTranspose_SeqUSFFT(Mat A,Vec x,Vec y)
125: {
127: Mat_USFFT *usfft = (Mat_USFFT*)A->data;
129: /* NB: for now we use outdim for both x and y; this will change once a full USFFT is implemented */
130: MatApply_USFFT_Private(usfft, &usfft->p_backward, FFTW_BACKWARD, x,y);
131: return(0);
132: }
136: PetscErrorCode MatDestroy_SeqUSFFT(Mat A)
137: {
138: Mat_USFFT *usfft = (Mat_USFFT*)A->data;
142: fftw_destroy_plan(usfft->p_forward);
143: fftw_destroy_plan(usfft->p_backward);
144: PetscFree(usfft->indim);
145: PetscFree(usfft->outdim);
146: PetscFree(usfft);
147: PetscObjectChangeTypeName((PetscObject)A,0);
148: return(0);
149: }/* MatDestroy_SeqUSFFT() */
154: /*@C
155: MatCreateSeqUSFFT - Creates a matrix object that provides sequential USFFT
156: via the external package FFTW
158: Collective on MPI_Comm
160: Input Parameter:
161: + da - geometry of the domain encoded by a DMDA
163: Output Parameter:
164: . A - the matrix
166: Options Database Keys:
167: + -mat_usfft_plannerflags - set the FFTW planner flags
169: Level: intermediate
170:
171: @*/
172: PetscErrorCode MatCreateSeqUSFFT(Vec sampleCoords, DMDA freqDA, Mat* A)
173: {
175: Mat_USFFT *usfft;
176: PetscInt m,n,M,N,i;
177: const char *p_flags[]={"FFTW_ESTIMATE","FFTW_MEASURE","FFTW_PATIENT","FFTW_EXHAUSTIVE"};
178: PetscBool flg;
179: PetscInt p_flag;
180: PetscInt dof, dim, freqSizes[3];
181: MPI_Comm comm;
182: PetscInt size;
185: PetscObjectGetComm((PetscObject)inda, &comm);
186: MPI_Comm_size(comm, &size);
187: if (size > 1) SETERRQ(comm,PETSC_ERR_USER, "Parallel DMDA (in) not yet supported by USFFT");
188: PetscObjectGetComm((PetscObject)outda, &comm);
189: MPI_Comm_size(comm, &size);
190: if (size > 1) SETERRQ(comm,PETSC_ERR_USER, "Parallel DMDA (out) not yet supported by USFFT");
191: MatCreate(comm,A);
192: PetscNewLog(*A,Mat_USFFT,&usfft);
193: (*A)->data = (void*)usfft;
194: usfft->inda = inda;
195: usfft->outda = outda;
196: /* inda */
197: DMDAGetInfo(usfft->inda, &ndim, dim+0, dim+1, dim+2, PETSC_NULL, PETSC_NULL, PETSC_NULL, &dof, PETSC_NULL, PETSC_NULL, PETSC_NULL);
198: if (ndim <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"ndim %d must be > 0",ndim);
199: if (dof <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"dof %d must be > 0",dof);
200: usfft->ndim = ndim;
201: usfft->dof = dof;
202: usfft->freqDA = freqDA;
203: /* NB: we reverse the freq and resample DMDA sizes, since the DMDA ordering (natural on x-y-z, with x varying the fastest)
204: is the order opposite of that assumed by FFTW: z varying the fastest */
205: PetscMalloc((usfft->ndim+1)*sizeof(PetscInt),&usfft->indim);
206: for(i = usfft->ndim; i > 0; --i) {
207: usfft->indim[usfft->ndim-i] = dim[i-1];
208: }
209: /* outda */
210: DMDAGetInfo(usfft->outda, &ndim, dim+0, dim+1, dim+2, PETSC_NULL, PETSC_NULL, PETSC_NULL, &dof, PETSC_NULL, PETSC_NULL, PETSC_NULL);
211: if (ndim != usfft->ndim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"in and out DMDA dimensions must match: %d != %d",usfft->ndim, ndim);
212: if (dof != usfft->dof) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"in and out DMDA dof must match: %d != %d",usfft->dof, dof);
213: /* Store output dimensions */
214: /* NB: we reverse the DMDA dimensions, since the DMDA ordering (natural on x-y-z, with x varying the fastest)
215: is the order opposite of that assumed by FFTW: z varying the fastest */
216: PetscMalloc((usfft->ndim+1)*sizeof(PetscInt),&usfft->outdim);
217: for(i = usfft->ndim; i > 0; --i) {
218: usfft->outdim[usfft->ndim-i] = dim[i-1];
219: }
221: /* TODO: Use the new form of DMDACreate() */
222: #if 0
223: DMDACreate(comm,usfft->dim, DMDA_NONPERIODIC, DMDA_STENCIL_STAR, usfft->freqSizes[0], usfft->freqSizes[1], usfft->freqSizes[2],
224: PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, dof, 0, PETSC_NULL, PETSC_NULL, PETSC_NULL, 0, &(usfft->resampleDA));
225: #endif
226: DMDAGetVec(usfft->resampleDA, usfft->resample);
229: /* CONTINUE: Need to build the connectivity "Sieve" attaching sample points to the resample points they are close to */
231: /* CONTINUE: recalculate matrix sizes based on the connectivity "Sieve" */
232: /* mat sizes */
233: m = 1; n = 1;
234: for (i=0; i<usfft->ndim; i++){
235: if (usfft->indim[i] <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"indim[%d]=%d must be > 0",i,usfft->indim[i]);
236: if (usfft->outdim[i] <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"outdim[%d]=%d must be > 0",i,usfft->outdim[i]);
237: n *= usfft->indim[i];
238: m *= usfft->outdim[i];
239: }
240: N = n*usfft->dof;
241: M = m*usfft->dof;
242: MatSetSizes(*A,M,N,M,N); /* "in size" is the number of columns, "out size" is the number of rows" */
243: PetscObjectChangeTypeName((PetscObject)*A,MATSEQUSFFT);
244: usfft->m = m; usfft->n = n; usfft->M = M; usfft->N = N;
245: /* FFTW */
246: usfft->p_forward = 0;
247: usfft->p_backward = 0;
248: usfft->p_flag = FFTW_ESTIMATE;
249: /* set Mat ops */
250: (*A)->ops->mult = MatMult_SeqUSFFT;
251: (*A)->ops->multtranspose = MatMultTranspose_SeqUSFFT;
252: (*A)->assembled = PETSC_TRUE;
253: (*A)->ops->destroy = MatDestroy_SeqUSFFT;
254: /* get runtime options */
255: PetscOptionsBegin(((PetscObject)(*A))->comm,((PetscObject)(*A))->prefix,"USFFT Options","Mat");
256: PetscOptionsEList("-mat_usfft_fftw_plannerflags","Planner Flags","None",p_flags,4,p_flags[0],&p_flag,&flg);
257: if (flg) {usfft->p_flag = (unsigned)p_flag;}
258: PetscOptionsEnd();
259: return(0);
260: }/* MatCreateSeqUSFFT() */
262: #endif