1c4762a1bSJed Brown static char help[] = "Test sequential USFFT interface on a uniform DMDA and compares the result to FFTW\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Compiling the code: 5c4762a1bSJed Brown This code uses the complex numbers version of PETSc and the FFTW package, so configure 6c4762a1bSJed Brown must be run to enable these. 7c4762a1bSJed Brown 8c4762a1bSJed Brown */ 9c4762a1bSJed Brown 10c4762a1bSJed Brown #include <petscmat.h> 11c4762a1bSJed Brown #include <petscdm.h> 12c4762a1bSJed Brown #include <petscdmda.h> 13d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 14d71ae5a4SJacob Faibussowitsch { 159371c9d4SSatish Balay typedef enum { 169371c9d4SSatish Balay RANDOM, 179371c9d4SSatish Balay CONSTANT, 189371c9d4SSatish Balay TANH, 199371c9d4SSatish Balay NUM_FUNCS 209371c9d4SSatish Balay } FuncType; 21c4762a1bSJed Brown const char *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"}; 22c4762a1bSJed Brown Mat A, AA; 23c4762a1bSJed Brown PetscMPIInt size; 24c4762a1bSJed Brown PetscInt N, i, stencil = 1, dof = 1; 25c4762a1bSJed Brown PetscInt dim[3] = {10, 10, 10}, ndim = 3; 26c4762a1bSJed Brown Vec coords, x, y, z, xx, yy, zz; 27c4762a1bSJed Brown PetscReal h[3]; 28c4762a1bSJed Brown PetscScalar s; 29c4762a1bSJed Brown PetscRandom rdm; 30c4762a1bSJed Brown PetscReal norm, enorm; 31c4762a1bSJed Brown PetscInt func; 32c4762a1bSJed Brown FuncType function = TANH; 33c4762a1bSJed Brown DM da, coordsda; 34c4762a1bSJed Brown PetscBool view_x = PETSC_FALSE, view_y = PETSC_FALSE, view_z = PETSC_FALSE; 35c4762a1bSJed Brown 36327415f7SBarry Smith PetscFunctionBeginUser; 37c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 39*bd158744SPierre Jolivet PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 40d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "USFFT Options", "ex27"); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-function", "Function type", "ex27", funcNames, NUM_FUNCS, funcNames[function], &func, NULL)); 42c4762a1bSJed Brown function = (FuncType)func; 43d0609cedSBarry Smith PetscOptionsEnd(); 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_x", &view_x, NULL)); 459566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_y", &view_y, NULL)); 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_z", &view_z, NULL)); 479566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dim", dim, &ndim, NULL)); 48c4762a1bSJed Brown 499371c9d4SSatish Balay PetscCall(DMDACreate3d(PETSC_COMM_SELF, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, dim[0], dim[1], dim[2], PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, dof, stencil, NULL, NULL, NULL, &da)); 509566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 519566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Coordinates */ 549566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &coordsda)); 559566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(coordsda, &coords)); 569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coords, "Grid coordinates")); 57c4762a1bSJed Brown for (i = 0, N = 1; i < 3; i++) { 58c4762a1bSJed Brown h[i] = 1.0 / dim[i]; 59c4762a1bSJed Brown PetscScalar *a; 609566063dSJacob Faibussowitsch PetscCall(VecGetArray(coords, &a)); 61c4762a1bSJed Brown PetscInt j, k, n = 0; 62c4762a1bSJed Brown for (i = 0; i < 3; ++i) { 63c4762a1bSJed Brown for (j = 0; j < dim[i]; ++j) { 64c4762a1bSJed Brown for (k = 0; k < 3; ++k) { 65c4762a1bSJed Brown a[n] = j * h[i]; /* coordinate along the j-th point in the i-th dimension */ 66c4762a1bSJed Brown ++n; 67c4762a1bSJed Brown } 68c4762a1bSJed Brown } 69c4762a1bSJed Brown } 709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coords, &a)); 71c4762a1bSJed Brown } 729566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(da, coords)); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* Work vectors */ 759566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &x)); 769566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Real space vector")); 779566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &xx)); 789566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)xx, "Real space vector")); 799566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &y)); 809566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y, "USFFT frequency space vector")); 819566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &yy)); 829566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)yy, "FFTW frequency space vector")); 839566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &z)); 849566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)z, "USFFT reconstructed vector")); 859566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &zz)); 869566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)zz, "FFTW reconstructed vector")); 87c4762a1bSJed Brown 8863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%3-" PetscInt_FMT ": USFFT on vector of ")); 89c4762a1bSJed Brown for (i = 0, N = 1; i < 3; i++) { 909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "dim[%d] = %d ", i, dim[i])); 91c4762a1bSJed Brown N *= dim[i]; 92c4762a1bSJed Brown } 939566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "; total size %d \n", N)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown if (function == RANDOM) { 969566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 979566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 989566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 999566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 100c4762a1bSJed Brown } else if (function == CONSTANT) { 1019566063dSJacob Faibussowitsch PetscCall(VecSet(x, 1.0)); 102c4762a1bSJed Brown } else if (function == TANH) { 103c4762a1bSJed Brown PetscScalar *a; 1049566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &a)); 105c4762a1bSJed Brown PetscInt j, k = 0; 106c4762a1bSJed Brown for (i = 0; i < 3; ++i) { 107c4762a1bSJed Brown for (j = 0; j < dim[i]; ++j) { 108c4762a1bSJed Brown a[k] = tanh((j - dim[i] / 2.0) * (10.0 / dim[i])); 109c4762a1bSJed Brown ++k; 110c4762a1bSJed Brown } 111c4762a1bSJed Brown } 1129566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &a)); 113c4762a1bSJed Brown } 1141baa6e33SBarry Smith if (view_x) PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 1159566063dSJacob Faibussowitsch PetscCall(VecCopy(x, xx)); 116c4762a1bSJed Brown 1179566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &norm)); 1189566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|x|_2 = %g\n", norm)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* create USFFT object */ 1219566063dSJacob Faibussowitsch PetscCall(MatCreateSeqUSFFT(coords, da, &A)); 122c4762a1bSJed Brown /* create FFTW object */ 1239566063dSJacob Faibussowitsch PetscCall(MatCreateSeqFFTW(PETSC_COMM_SELF, 3, dim, &AA)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* apply USFFT and FFTW FORWARD "preemptively", so the fftw_plans can be reused on different vectors */ 1269566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, z)); 1279566063dSJacob Faibussowitsch PetscCall(MatMult(AA, xx, zz)); 128c4762a1bSJed Brown /* Now apply USFFT and FFTW forward several (3) times */ 129c4762a1bSJed Brown for (i = 0; i < 3; ++i) { 1309566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, y)); 1319566063dSJacob Faibussowitsch PetscCall(MatMult(AA, xx, yy)); 1329566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, y, z)); 1339566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(AA, yy, zz)); 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 136c4762a1bSJed Brown if (view_y) { 1379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "y = \n")); 1389566063dSJacob Faibussowitsch PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD)); 1399566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "yy = \n")); 1409566063dSJacob Faibussowitsch PetscCall(VecView(yy, PETSC_VIEWER_STDOUT_WORLD)); 141c4762a1bSJed Brown } 142c4762a1bSJed Brown 143c4762a1bSJed Brown if (view_z) { 1449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "z = \n")); 1459566063dSJacob Faibussowitsch PetscCall(VecView(z, PETSC_VIEWER_STDOUT_WORLD)); 1469566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "zz = \n")); 1479566063dSJacob Faibussowitsch PetscCall(VecView(zz, PETSC_VIEWER_STDOUT_WORLD)); 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* compare x and z. USFFT computes an unnormalized DFT, thus z = N*x */ 151c4762a1bSJed Brown s = 1.0 / (PetscReal)N; 1529566063dSJacob Faibussowitsch PetscCall(VecScale(z, s)); 1539566063dSJacob Faibussowitsch PetscCall(VecAXPY(x, -1.0, z)); 1549566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_1, &enorm)); 1559566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|x-z| = %g\n", enorm)); 156c4762a1bSJed Brown 157c4762a1bSJed Brown /* compare xx and zz. FFTW computes an unnormalized DFT, thus zz = N*x */ 158c4762a1bSJed Brown s = 1.0 / (PetscReal)N; 1599566063dSJacob Faibussowitsch PetscCall(VecScale(zz, s)); 1609566063dSJacob Faibussowitsch PetscCall(VecAXPY(xx, -1.0, zz)); 1619566063dSJacob Faibussowitsch PetscCall(VecNorm(xx, NORM_1, &enorm)); 1629566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|xx-zz| = %g\n", enorm)); 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* compare y and yy: USFFT and FFTW results*/ 1659566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 1669566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, yy)); 1679566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_1, &enorm)); 1689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|y|_2 = %g\n", norm)); 1699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|y-yy| = %g\n", enorm)); 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* compare z and zz: USFFT and FFTW results*/ 1729566063dSJacob Faibussowitsch PetscCall(VecNorm(z, NORM_2, &norm)); 1739566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, -1.0, zz)); 1749566063dSJacob Faibussowitsch PetscCall(VecNorm(z, NORM_1, &enorm)); 1759566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|z|_2 = %g\n", norm)); 1769566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|z-zz| = %g\n", enorm)); 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* free spaces */ 1799566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &x)); 1809566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &xx)); 1819566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &y)); 1829566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &yy)); 1839566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &z)); 1849566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &zz)); 1859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coords)); 1869566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1879566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 188b122ec5aSJacob Faibussowitsch return 0; 189c4762a1bSJed Brown } 190