1 static char help[] = "Test sequential FFTW convolution\n\n"; 2 3 /* 4 Compiling the code: 5 This code uses the complex numbers, so configure must be given --with-scalar-type=complex to enable this 6 */ 7 8 #include <petscmat.h> 9 10 int main(int argc,char **args) 11 { 12 typedef enum {RANDOM, CONSTANT, TANH, NUM_FUNCS} FuncType; 13 const char *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"}; 14 Mat A; 15 PetscMPIInt size; 16 PetscInt n = 10,N,ndim=4,dim[4],DIM,i,j; 17 Vec w,x,y1,y2,z1,z2; 18 PetscScalar *a, *a2, *a3; 19 PetscScalar s; 20 PetscRandom rdm; 21 PetscReal enorm; 22 PetscInt func = 0; 23 FuncType function = RANDOM; 24 PetscBool view = PETSC_FALSE; 25 PetscErrorCode ierr; 26 27 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 28 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 29 PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 30 ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex112");PetscCall(ierr); 31 PetscCall(PetscOptionsEList("-function", "Function type", "ex121", funcNames, NUM_FUNCS, funcNames[function], &func, NULL)); 32 PetscCall(PetscOptionsBool("-vec_view draw", "View the functions", "ex112", view, &view, NULL)); 33 function = (FuncType) func; 34 ierr = PetscOptionsEnd();PetscCall(ierr); 35 36 for (DIM = 0; DIM < ndim; DIM++) { 37 dim[DIM] = n; /* size of transformation in DIM-dimension */ 38 } 39 PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 40 PetscCall(PetscRandomSetFromOptions(rdm)); 41 42 for (DIM = 1; DIM < 5; DIM++) { 43 /* create vectors of length N=n^DIM */ 44 for (i = 0, N = 1; i < DIM; i++) N *= dim[i]; 45 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n %d-D: FFTW on vector of size %d \n",DIM,N)); 46 PetscCall(VecCreateSeq(PETSC_COMM_SELF,N,&x)); 47 PetscCall(PetscObjectSetName((PetscObject) x, "Real space vector")); 48 PetscCall(VecDuplicate(x,&w)); 49 PetscCall(PetscObjectSetName((PetscObject) w, "Window vector")); 50 PetscCall(VecDuplicate(x,&y1)); 51 PetscCall(PetscObjectSetName((PetscObject) y1, "Frequency space vector")); 52 PetscCall(VecDuplicate(x,&y2)); 53 PetscCall(PetscObjectSetName((PetscObject) y2, "Frequency space window vector")); 54 PetscCall(VecDuplicate(x,&z1)); 55 PetscCall(PetscObjectSetName((PetscObject) z1, "Reconstructed convolution")); 56 PetscCall(VecDuplicate(x,&z2)); 57 PetscCall(PetscObjectSetName((PetscObject) z2, "Real space convolution")); 58 59 if (function == RANDOM) { 60 PetscCall(VecSetRandom(x, rdm)); 61 } else if (function == CONSTANT) { 62 PetscCall(VecSet(x, 1.0)); 63 } else if (function == TANH) { 64 PetscCall(VecGetArray(x, &a)); 65 for (i = 0; i < N; ++i) { 66 a[i] = tanh((i - N/2.0)*(10.0/N)); 67 } 68 PetscCall(VecRestoreArray(x, &a)); 69 } 70 if (view) PetscCall(VecView(x, PETSC_VIEWER_DRAW_WORLD)); 71 72 /* Create window function */ 73 PetscCall(VecGetArray(w, &a)); 74 for (i = 0; i < N; ++i) { 75 /* Step Function */ 76 a[i] = (i > N/4 && i < 3*N/4) ? 1.0 : 0.0; 77 /* Delta Function */ 78 /*a[i] = (i == N/2)? 1.0: 0.0; */ 79 } 80 PetscCall(VecRestoreArray(w, &a)); 81 if (view) PetscCall(VecView(w, PETSC_VIEWER_DRAW_WORLD)); 82 83 /* create FFTW object */ 84 PetscCall(MatCreateFFT(PETSC_COMM_SELF,DIM,dim,MATFFTW,&A)); 85 86 /* Convolve x with w*/ 87 PetscCall(MatMult(A,x,y1)); 88 PetscCall(MatMult(A,w,y2)); 89 PetscCall(VecPointwiseMult(y1, y1, y2)); 90 if (view && i == 0) PetscCall(VecView(y1, PETSC_VIEWER_DRAW_WORLD)); 91 PetscCall(MatMultTranspose(A,y1,z1)); 92 93 /* Compute the real space convolution */ 94 PetscCall(VecGetArray(x, &a)); 95 PetscCall(VecGetArray(w, &a2)); 96 PetscCall(VecGetArray(z2, &a3)); 97 for (i = 0; i < N; ++i) { 98 /* PetscInt checkInd = (i > N/2-1)? i-N/2: i+N/2;*/ 99 100 a3[i] = 0.0; 101 for (j = -N/2+1; j < N/2; ++j) { 102 PetscInt xpInd = (j < 0) ? N+j : j; 103 PetscInt diffInd = (i-j < 0) ? N-(j-i) : (i-j > N-1) ? i-j-N : i-j; 104 105 a3[i] += a[xpInd]*a2[diffInd]; 106 } 107 } 108 PetscCall(VecRestoreArray(x, &a)); 109 PetscCall(VecRestoreArray(w, &a2)); 110 PetscCall(VecRestoreArray(z2, &a3)); 111 112 /* compare z1 and z2. FFTW computes an unnormalized DFT, thus z1 = N*z2 */ 113 s = 1.0/(PetscReal)N; 114 PetscCall(VecScale(z1,s)); 115 if (view) PetscCall(VecView(z1, PETSC_VIEWER_DRAW_WORLD)); 116 if (view) PetscCall(VecView(z2, PETSC_VIEWER_DRAW_WORLD)); 117 PetscCall(VecAXPY(z1,-1.0,z2)); 118 PetscCall(VecNorm(z1,NORM_1,&enorm)); 119 if (enorm > 1.e-11) { 120 PetscCall(PetscPrintf(PETSC_COMM_SELF," Error norm of |z1 - z2| %g\n",(double)enorm)); 121 } 122 123 /* free spaces */ 124 PetscCall(VecDestroy(&x)); 125 PetscCall(VecDestroy(&y1)); 126 PetscCall(VecDestroy(&y2)); 127 PetscCall(VecDestroy(&z1)); 128 PetscCall(VecDestroy(&z2)); 129 PetscCall(VecDestroy(&w)); 130 PetscCall(MatDestroy(&A)); 131 } 132 PetscCall(PetscRandomDestroy(&rdm)); 133 PetscCall(PetscFinalize()); 134 return 0; 135 } 136 137 /*TEST 138 139 build: 140 requires: fftw complex 141 142 test: 143 output_file: output/ex121.out 144 TODO: Example or FFTW interface is broken 145 146 TEST*/ 147