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 typedef enum { 12 RANDOM, 13 CONSTANT, 14 TANH, 15 NUM_FUNCS 16 } FuncType; 17 const char *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"}; 18 Mat A; 19 PetscMPIInt size; 20 PetscInt n = 10, N, ndim = 4, dim[4], DIM, i, j; 21 Vec w, x, y1, y2, z1, z2; 22 PetscScalar *a, *a2, *a3; 23 PetscScalar s; 24 PetscRandom rdm; 25 PetscReal enorm; 26 PetscInt func = 0; 27 FuncType function = RANDOM; 28 PetscBool view = PETSC_FALSE; 29 30 PetscFunctionBeginUser; 31 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 32 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 33 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 34 PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex112"); 35 PetscCall(PetscOptionsEList("-function", "Function type", "ex121", funcNames, NUM_FUNCS, funcNames[function], &func, NULL)); 36 PetscCall(PetscOptionsBool("-vec_view draw", "View the functions", "ex112", view, &view, NULL)); 37 function = (FuncType)func; 38 PetscOptionsEnd(); 39 40 for (DIM = 0; DIM < ndim; DIM++) { dim[DIM] = n; /* size of transformation in DIM-dimension */ } 41 PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 42 PetscCall(PetscRandomSetFromOptions(rdm)); 43 44 for (DIM = 1; DIM < 5; DIM++) { 45 /* create vectors of length N=n^DIM */ 46 for (i = 0, N = 1; i < DIM; i++) N *= dim[i]; 47 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n %d-D: FFTW on vector of size %d \n", DIM, N)); 48 PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &x)); 49 PetscCall(PetscObjectSetName((PetscObject)x, "Real space vector")); 50 PetscCall(VecDuplicate(x, &w)); 51 PetscCall(PetscObjectSetName((PetscObject)w, "Window vector")); 52 PetscCall(VecDuplicate(x, &y1)); 53 PetscCall(PetscObjectSetName((PetscObject)y1, "Frequency space vector")); 54 PetscCall(VecDuplicate(x, &y2)); 55 PetscCall(PetscObjectSetName((PetscObject)y2, "Frequency space window vector")); 56 PetscCall(VecDuplicate(x, &z1)); 57 PetscCall(PetscObjectSetName((PetscObject)z1, "Reconstructed convolution")); 58 PetscCall(VecDuplicate(x, &z2)); 59 PetscCall(PetscObjectSetName((PetscObject)z2, "Real space convolution")); 60 61 if (function == RANDOM) { 62 PetscCall(VecSetRandom(x, rdm)); 63 } else if (function == CONSTANT) { 64 PetscCall(VecSet(x, 1.0)); 65 } else if (function == TANH) { 66 PetscCall(VecGetArray(x, &a)); 67 for (i = 0; i < N; ++i) { a[i] = tanh((i - N / 2.0) * (10.0 / N)); } 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) { PetscCall(PetscPrintf(PETSC_COMM_SELF, " Error norm of |z1 - z2| %g\n", (double)enorm)); } 120 121 /* free spaces */ 122 PetscCall(VecDestroy(&x)); 123 PetscCall(VecDestroy(&y1)); 124 PetscCall(VecDestroy(&y2)); 125 PetscCall(VecDestroy(&z1)); 126 PetscCall(VecDestroy(&z2)); 127 PetscCall(VecDestroy(&w)); 128 PetscCall(MatDestroy(&A)); 129 } 130 PetscCall(PetscRandomDestroy(&rdm)); 131 PetscCall(PetscFinalize()); 132 return 0; 133 } 134 135 /*TEST 136 137 build: 138 requires: fftw complex 139 140 test: 141 output_file: output/ex121.out 142 TODO: Example or FFTW interface is broken 143 144 TEST*/ 145