Actual source code: matusfft.c

  1: #define PETSCMAT_DLL

  3: /*
  4:     Provides an implementation of the Unevenly Sampled FFT algorithm as a Mat.
  5:     Testing examples can be found in ~/src/mat/examples/tests FIX: should these be moved to dm/da/examples/tests?
  6: */

 8:  #include private/matimpl.h
  9: #include "petscda.h"                  /*I "petscda.h"  I*/ /* Unlike equispaced FFT, USFFT requires geometric information encoded by a DA */
 10: #include "fftw3.h"

 12: typedef struct {
 13:   PetscInt       dim;
 14:   Vec            sampleCoords;
 15:   PetscInt       dof;
 16:   DA             freqDA;       /* frequency DA */
 17:   PetscInt       *freqSizes;   // sizes of the frequency DA, one per each dim
 18:   DA             resampleDa;   /* the Battle-Lemarie interpolant DA */
 19:   Vec            resample;     // Vec of samples, one per dof per sample point
 20:   fftw_plan      p_forward,p_backward;
 21:   unsigned       p_flag; /* planner flags, FFTW_ESTIMATE,FFTW_MEASURE, FFTW_PATIENT, FFTW_EXHAUSTIVE */
 22: } Mat_USFFT;


 27: PetscErrorCode MatApply_USFFT_Private(Mat_USFFT *usfft, fftw_plan *plan, int direction, Vec x,Vec y)
 28: {
 30:   PetscScalar    *r_array, *y_array;

 33: #if 0
 34:   // resample x to usfft->resample
 35:   MatResample_USFFT_Private(Mat_USFFT *usfft, x);

 37:   /* NB: for now we use outdim for both x and y; this will change once a full USFFT is implemented */
 38:   VecGetArray(usfft->resample,&r_array);
 39:   VecGetArray(y,&y_array);
 40:   if (!*plan){ /* create a plan then execute it*/
 41:     if(usfft->dof == 1) {
 42: #ifdef PETSC_DEBUG_USFFT
 43:       PetscPrintf(PETSC_COMM_WORLD, "direction = %d, usfft->ndim = %d\n", direction, usfft->ndim);
 44:       for(int ii = 0; ii < usfft->ndim; ++ii) {
 45:         PetscPrintf(PETSC_COMM_WORLD, "usfft->outdim[%d] = %d\n", ii, usfft->outdim[ii]);
 46:       }
 47: #endif 

 49:       switch (usfft->dim){
 50:       case 1:
 51:         *plan = fftw_plan_dft_1d(usfft->outdim[0],(fftw_complex*)x_array,(fftw_complex*)y_array,direction,usfft->p_flag);
 52:         break;
 53:       case 2:
 54:         *plan = fftw_plan_dft_2d(usfft->outdim[0],usfft->outdim[1],(fftw_complex*)x_array,(fftw_complex*)y_array,direction,usfft->p_flag);
 55:         break;
 56:       case 3:
 57:         *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);
 58:         break;
 59:       default:
 60:         *plan = fftw_plan_dft(usfft->ndim,usfft->outdim,(fftw_complex*)x_array,(fftw_complex*)y_array,direction,usfft->p_flag);
 61:         break;
 62:       }
 63:       fftw_execute(*plan);
 64:     }/* if(dof == 1) */
 65:     else { /* if(dof > 1) */
 66:       *plan = fftw_plan_many_dft(/*rank*/usfft->ndim, /*n*/usfft->outdim, /*howmany*/usfft->dof,
 67:                                  (fftw_complex*)x_array, /*nembed*/usfft->outdim, /*stride*/usfft->dof, /*dist*/1,
 68:                                  (fftw_complex*)y_array, /*nembed*/usfft->outdim, /*stride*/usfft->dof, /*dist*/1,
 69:                                  /*sign*/direction, /*flags*/usfft->p_flag);
 70:       fftw_execute(*plan);
 71:     }/* if(dof > 1) */
 72:   }/* if(!*plan) */
 73:   else { /* if(*plan) */
 74:     /* use existing plan */
 75:     fftw_execute_dft(*plan,(fftw_complex*)x_array,(fftw_complex*)y_array);
 76:   }
 77:   VecRestoreArray(y,&y_array);
 78:   VecRestoreArray(x,&x_array);
 79: #endif
 80:   return(0);
 81: }/* MatApply_USFFT_Private() */

 83: #if 0
 86: PetscErrorCode Mat_USFFT_ProjectOnBattleLemarie_Private(Vec x,double *r)
 87: /* Project onto the Battle-Lemarie function centered around r */
 88: {
 90:   PetscScalar    *x_array, *y_array;

 93:   return(0);
 94: }/* Mat_USFFT_ProjectOnBattleLemarie_Private() */

 98: PetscErrorCode MatInterpolate_USFFT_Private(Vec x,Vec y)
 99: {
101:   PetscScalar    *x_array, *y_array;

104:   return(0);
105: }/* MatInterpolate_USFFT_Private() */


110: PetscErrorCode MatMult_SeqUSFFT(Mat A,Vec x,Vec y)
111: {
113:   Mat_USFFT       *usfft = (Mat_USFFT*)A->data;

116:   /* NB: for now we use outdim for both x and y; this will change once a full USFFT is implemented */
117:   MatApply_USFFT_Private(usfft, &usfft->p_forward, FFTW_FORWARD, x,y);
118:   return(0);
119: }

123: PetscErrorCode MatMultTranspose_SeqUSFFT(Mat A,Vec x,Vec y)
124: {
126:   Mat_USFFT       *usfft = (Mat_USFFT*)A->data;
128:   /* NB: for now we use outdim for both x and y; this will change once a full USFFT is implemented */
129:   MatApply_USFFT_Private(usfft, &usfft->p_backward, FFTW_BACKWARD, x,y);
130:   return(0);
131: }

135: PetscErrorCode MatDestroy_SeqUSFFT(Mat A)
136: {
137:   Mat_USFFT       *usfft = (Mat_USFFT*)A->data;

141:   fftw_destroy_plan(usfft->p_forward);
142:   fftw_destroy_plan(usfft->p_backward);
143:   PetscFree(usfft->indim);
144:   PetscFree(usfft->outdim);
145:   PetscFree(usfft);
146:   PetscObjectChangeTypeName((PetscObject)A,0);
147:   return(0);
148: }/* MatDestroy_SeqUSFFT() */


153: /*@C
154:       MatCreateSeqUSFFT - Creates a matrix object that provides sequential USFFT
155:   via the external package FFTW

157:    Collective on MPI_Comm

159:    Input Parameter:
160: +   da - geometry of the domain encoded by a DA

162:    Output Parameter:
163: .   A  - the matrix

165:   Options Database Keys:
166: + -mat_usfft_plannerflags - set the FFTW planner flags

168:    Level: intermediate
169:    
170: @*/
171: PetscErrorCode  MatCreateSeqUSFFT(Vec sampleCoords, DA freqDA, Mat* A)
172: {
174:   Mat_USFFT      *usfft;
175:   PetscInt       m,n,M,N,i;
176:   const char     *p_flags[]={"FFTW_ESTIMATE","FFTW_MEASURE","FFTW_PATIENT","FFTW_EXHAUSTIVE"};
177:   PetscTruth     flg;
178:   PetscInt       p_flag;
179:   PetscInt       dof, dim, freqSizes[3];
180:   MPI_Comm       comm;
181:   PetscInt       size;

184:   PetscObjectGetComm((PetscObject)inda, &comm);
185:   MPI_Comm_size(comm, &size);
186:   if (size > 1) SETERRQ(PETSC_ERR_USER, "Parallel DA (in) not yet supported by USFFT");
187:   PetscObjectGetComm((PetscObject)outda, &comm);
188:   MPI_Comm_size(comm, &size);
189:   if (size > 1) SETERRQ(PETSC_ERR_USER, "Parallel DA (out) not yet supported by USFFT");
190:   MatCreate(comm,A);
191:   PetscNewLog(*A,Mat_USFFT,&usfft);
192:   (*A)->data = (void*)usfft;
193:   usfft->inda = inda;
194:   usfft->outda = outda;
195:   /* inda */
196:   DAGetInfo(usfft->inda, &ndim, dim+0, dim+1, dim+2, PETSC_NULL, PETSC_NULL, PETSC_NULL, &dof, PETSC_NULL, PETSC_NULL, PETSC_NULL);
197:   if (ndim <= 0) SETERRQ1(PETSC_ERR_USER,"ndim %d must be > 0",ndim);
198:   if (dof <= 0) SETERRQ1(PETSC_ERR_USER,"dof %d must be > 0",dof);
199:   usfft->ndim = ndim;
200:   usfft->dof = dof;
201:   usfft->freqDA     = freqDA;
202:   /* NB: we reverse the freq and resample DA sizes, since the DA ordering (natural on x-y-z, with x varying the fastest) 
203:      is the order opposite of that assumed by FFTW: z varying the fastest */
204:   PetscMalloc((usfft->ndim+1)*sizeof(PetscInt),&usfft->indim);
205:   for(i = usfft->ndim; i > 0; --i) {
206:     usfft->indim[usfft->ndim-i] = dim[i-1];
207:   }
208:   /* outda */
209:   DAGetInfo(usfft->outda, &ndim, dim+0, dim+1, dim+2, PETSC_NULL, PETSC_NULL, PETSC_NULL, &dof, PETSC_NULL, PETSC_NULL, PETSC_NULL);
210:   if (ndim != usfft->ndim) SETERRQ2(PETSC_ERR_USER,"in and out DA dimensions must match: %d != %d",usfft->ndim, ndim);
211:   if (dof != usfft->dof) SETERRQ2(PETSC_ERR_USER,"in and out DA dof must match: %d != %d",usfft->dof, dof);
212:   /* Store output dimensions */
213:   /* NB: we reverse the DA dimensions, since the DA ordering (natural on x-y-z, with x varying the fastest) 
214:      is the order opposite of that assumed by FFTW: z varying the fastest */
215:   PetscMalloc((usfft->ndim+1)*sizeof(PetscInt),&usfft->outdim);
216:   for(i = usfft->ndim; i > 0; --i) {
217:     usfft->outdim[usfft->ndim-i] = dim[i-1];
218:   }

220:   // TODO: Use the new form of DACreate()
221:   //DACreate(comm,usfft->dim, DA_NONPERIODIC, DA_STENCIL_STAR, usfft->freqSizes[0], usfft->freqSizes[1], usfft->freqSizes[2],
222:   //                PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, dof, 0, PETSC_NULL, PETSC_NULL, PETSC_NULL,  0, &(usfft->resampleDA));
223:   DAGetVec(usfft->resampleDA, usfft->resample);


226:   // CONTINUE: Need to build the connectivity "Sieve" attaching sample points to the resample points they are close to

228:   // CONTINUE: recalculate matrix sizes based on the connectivity "Sieve"
229:   /* mat sizes */
230:   m = 1; n = 1;
231:   for (i=0; i<usfft->ndim; i++){
232:     if (usfft->indim[i] <= 0) SETERRQ2(PETSC_ERR_USER,"indim[%d]=%d must be > 0",i,usfft->indim[i]);
233:     if (usfft->outdim[i] <= 0) SETERRQ2(PETSC_ERR_USER,"outdim[%d]=%d must be > 0",i,usfft->outdim[i]);
234:     n *= usfft->indim[i];
235:     m *= usfft->outdim[i];
236:   }
237:   N = n*usfft->dof;
238:   M = m*usfft->dof;
239:   MatSetSizes(*A,M,N,M,N);  /* "in size" is the number of columns, "out size" is the number of rows" */
240:   PetscObjectChangeTypeName((PetscObject)*A,MATSEQUSFFT);
241:   usfft->m = m; usfft->n = n; usfft->M = M; usfft->N = N;
242:   /* FFTW */
243:   usfft->p_forward  = 0;
244:   usfft->p_backward = 0;
245:   usfft->p_flag     = FFTW_ESTIMATE;
246:   /* set Mat ops */
247:   (*A)->ops->mult          = MatMult_SeqUSFFT;
248:   (*A)->ops->multtranspose = MatMultTranspose_SeqUSFFT;
249:   (*A)->assembled          = PETSC_TRUE;
250:   (*A)->ops->destroy       = MatDestroy_SeqUSFFT;
251:   /* get runtime options */
252:   PetscOptionsBegin(((PetscObject)(*A))->comm,((PetscObject)(*A))->prefix,"USFFT Options","Mat");
253:   PetscOptionsEList("-mat_usfft_fftw_plannerflags","Planner Flags","None",p_flags,4,p_flags[0],&p_flag,&flg);
254:   if (flg) {usfft->p_flag = (unsigned)p_flag;}
255:   PetscOptionsEnd();
256:   return(0);
257: }/* MatCreateSeqUSFFT() */

259: #endif