1c4762a1bSJed Brown static char help[] = "Test sequential FFTW convolution\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Compiling the code: 5c4762a1bSJed Brown This code uses the complex numbers, so configure must be given --with-scalar-type=complex to enable this 6c4762a1bSJed Brown */ 7c4762a1bSJed Brown 8c4762a1bSJed Brown #include <petscmat.h> 9c4762a1bSJed Brown 109371c9d4SSatish Balay int main(int argc, char **args) { 119371c9d4SSatish Balay typedef enum { 129371c9d4SSatish Balay RANDOM, 139371c9d4SSatish Balay CONSTANT, 149371c9d4SSatish Balay TANH, 159371c9d4SSatish Balay NUM_FUNCS 169371c9d4SSatish Balay } FuncType; 17c4762a1bSJed Brown const char *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"}; 18c4762a1bSJed Brown Mat A; 19c4762a1bSJed Brown PetscMPIInt size; 20c4762a1bSJed Brown PetscInt n = 10, N, ndim = 4, dim[4], DIM, i, j; 21c4762a1bSJed Brown Vec w, x, y1, y2, z1, z2; 22c4762a1bSJed Brown PetscScalar *a, *a2, *a3; 23c4762a1bSJed Brown PetscScalar s; 24c4762a1bSJed Brown PetscRandom rdm; 25c4762a1bSJed Brown PetscReal enorm; 26c4762a1bSJed Brown PetscInt func = 0; 27c4762a1bSJed Brown FuncType function = RANDOM; 28c4762a1bSJed Brown PetscBool view = PETSC_FALSE; 29c4762a1bSJed Brown 30327415f7SBarry Smith PetscFunctionBeginUser; 319566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 33be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 34d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex112"); 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-function", "Function type", "ex121", funcNames, NUM_FUNCS, funcNames[function], &func, NULL)); 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-vec_view draw", "View the functions", "ex112", view, &view, NULL)); 37c4762a1bSJed Brown function = (FuncType)func; 38d0609cedSBarry Smith PetscOptionsEnd(); 39c4762a1bSJed Brown 409371c9d4SSatish Balay for (DIM = 0; DIM < ndim; DIM++) { dim[DIM] = n; /* size of transformation in DIM-dimension */ } 419566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 429566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 43c4762a1bSJed Brown 44c4762a1bSJed Brown for (DIM = 1; DIM < 5; DIM++) { 45c4762a1bSJed Brown /* create vectors of length N=n^DIM */ 46c4762a1bSJed Brown for (i = 0, N = 1; i < DIM; i++) N *= dim[i]; 479566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n %d-D: FFTW on vector of size %d \n", DIM, N)); 489566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &x)); 499566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Real space vector")); 509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &w)); 519566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)w, "Window vector")); 529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y1)); 539566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y1, "Frequency space vector")); 549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y2)); 559566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y2, "Frequency space window vector")); 569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &z1)); 579566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)z1, "Reconstructed convolution")); 589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &z2)); 599566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)z2, "Real space convolution")); 60c4762a1bSJed Brown 61c4762a1bSJed Brown if (function == RANDOM) { 629566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 63c4762a1bSJed Brown } else if (function == CONSTANT) { 649566063dSJacob Faibussowitsch PetscCall(VecSet(x, 1.0)); 65c4762a1bSJed Brown } else if (function == TANH) { 669566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &a)); 679371c9d4SSatish Balay for (i = 0; i < N; ++i) { a[i] = tanh((i - N / 2.0) * (10.0 / N)); } 689566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &a)); 69c4762a1bSJed Brown } 709566063dSJacob Faibussowitsch if (view) PetscCall(VecView(x, PETSC_VIEWER_DRAW_WORLD)); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* Create window function */ 739566063dSJacob Faibussowitsch PetscCall(VecGetArray(w, &a)); 74c4762a1bSJed Brown for (i = 0; i < N; ++i) { 75c4762a1bSJed Brown /* Step Function */ 76c4762a1bSJed Brown a[i] = (i > N / 4 && i < 3 * N / 4) ? 1.0 : 0.0; 77c4762a1bSJed Brown /* Delta Function */ 78c4762a1bSJed Brown /*a[i] = (i == N/2)? 1.0: 0.0; */ 79c4762a1bSJed Brown } 809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w, &a)); 819566063dSJacob Faibussowitsch if (view) PetscCall(VecView(w, PETSC_VIEWER_DRAW_WORLD)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* create FFTW object */ 849566063dSJacob Faibussowitsch PetscCall(MatCreateFFT(PETSC_COMM_SELF, DIM, dim, MATFFTW, &A)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Convolve x with w*/ 879566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, y1)); 889566063dSJacob Faibussowitsch PetscCall(MatMult(A, w, y2)); 899566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y1, y1, y2)); 909566063dSJacob Faibussowitsch if (view && i == 0) PetscCall(VecView(y1, PETSC_VIEWER_DRAW_WORLD)); 919566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, y1, z1)); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* Compute the real space convolution */ 949566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &a)); 959566063dSJacob Faibussowitsch PetscCall(VecGetArray(w, &a2)); 969566063dSJacob Faibussowitsch PetscCall(VecGetArray(z2, &a3)); 97c4762a1bSJed Brown for (i = 0; i < N; ++i) { 98c4762a1bSJed Brown /* PetscInt checkInd = (i > N/2-1)? i-N/2: i+N/2;*/ 99c4762a1bSJed Brown 100c4762a1bSJed Brown a3[i] = 0.0; 101c4762a1bSJed Brown for (j = -N / 2 + 1; j < N / 2; ++j) { 102c4762a1bSJed Brown PetscInt xpInd = (j < 0) ? N + j : j; 103c4762a1bSJed Brown PetscInt diffInd = (i - j < 0) ? N - (j - i) : (i - j > N - 1) ? i - j - N : i - j; 104c4762a1bSJed Brown 105c4762a1bSJed Brown a3[i] += a[xpInd] * a2[diffInd]; 106c4762a1bSJed Brown } 107c4762a1bSJed Brown } 1089566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &a)); 1099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w, &a2)); 1109566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(z2, &a3)); 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* compare z1 and z2. FFTW computes an unnormalized DFT, thus z1 = N*z2 */ 113c4762a1bSJed Brown s = 1.0 / (PetscReal)N; 1149566063dSJacob Faibussowitsch PetscCall(VecScale(z1, s)); 1159566063dSJacob Faibussowitsch if (view) PetscCall(VecView(z1, PETSC_VIEWER_DRAW_WORLD)); 1169566063dSJacob Faibussowitsch if (view) PetscCall(VecView(z2, PETSC_VIEWER_DRAW_WORLD)); 1179566063dSJacob Faibussowitsch PetscCall(VecAXPY(z1, -1.0, z2)); 1189566063dSJacob Faibussowitsch PetscCall(VecNorm(z1, NORM_1, &enorm)); 119*48a46eb9SPierre Jolivet if (enorm > 1.e-11) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Error norm of |z1 - z2| %g\n", (double)enorm)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* free spaces */ 1229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y1)); 1249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y2)); 1259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z1)); 1269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z2)); 1279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&w)); 1289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 129c4762a1bSJed Brown } 1309566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 1319566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 132b122ec5aSJacob Faibussowitsch return 0; 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown /*TEST 136c4762a1bSJed Brown 137c4762a1bSJed Brown build: 138c4762a1bSJed Brown requires: fftw complex 139c4762a1bSJed Brown 140c4762a1bSJed Brown test: 141c4762a1bSJed Brown output_file: output/ex121.out 142c4762a1bSJed Brown TODO: Example or FFTW interface is broken 143c4762a1bSJed Brown 144c4762a1bSJed Brown TEST*/ 145