1c4762a1bSJed Brown static char help[] = "Test sequential USFFT interface on a 3-dof field over a uniform DMDA and compares to the result of FFTW acting on a split version of the field\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 #define DOF 3 11c4762a1bSJed Brown 12c4762a1bSJed Brown #include <petscmat.h> 13c4762a1bSJed Brown #include <petscdm.h> 14c4762a1bSJed Brown #include <petscdmda.h> 159371c9d4SSatish Balay int main(int argc, char **args) { 169371c9d4SSatish Balay typedef enum { 179371c9d4SSatish Balay RANDOM, 189371c9d4SSatish Balay CONSTANT, 199371c9d4SSatish Balay TANH, 209371c9d4SSatish Balay NUM_FUNCS 219371c9d4SSatish Balay } FuncType; 22c4762a1bSJed Brown const char *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"}; 23c4762a1bSJed Brown Mat A, AA; 24c4762a1bSJed Brown PetscMPIInt size; 25c4762a1bSJed Brown PetscInt N, i, stencil = 1, dof = 3; 26c4762a1bSJed Brown PetscInt dim[3] = {10, 10, 10}, ndim = 3; 27c4762a1bSJed Brown Vec coords, x, y, z, xx, yy, zz; 28c4762a1bSJed Brown Vec xxsplit[DOF], yysplit[DOF], zzsplit[DOF]; 29c4762a1bSJed Brown PetscReal h[3]; 30c4762a1bSJed Brown PetscScalar s; 31c4762a1bSJed Brown PetscRandom rdm; 32c4762a1bSJed Brown PetscReal norm, enorm; 33c4762a1bSJed Brown PetscInt func, ii; 34c4762a1bSJed Brown FuncType function = TANH; 35c4762a1bSJed Brown DM da, da1, coordsda; 36c4762a1bSJed Brown PetscBool view_x = PETSC_FALSE, view_y = PETSC_FALSE, view_z = PETSC_FALSE; 37c4762a1bSJed Brown 38327415f7SBarry Smith PetscFunctionBeginUser; 399566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 4108401ef6SPierre Jolivet PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only!"); 42d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "USFFT Options", "ex27"); 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-function", "Function type", "ex27", funcNames, NUM_FUNCS, funcNames[function], &func, NULL)); 44c4762a1bSJed Brown function = (FuncType)func; 45d0609cedSBarry Smith PetscOptionsEnd(); 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_x", &view_x, NULL)); 479566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_y", &view_y, NULL)); 489566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_z", &view_z, NULL)); 499566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dim", dim, &ndim, NULL)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* DMDA with the correct fiber dimension */ 529371c9d4SSatish 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)); 539566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 549566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 55c4762a1bSJed Brown /* DMDA with fiber dimension 1 for split fields */ 569371c9d4SSatish 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, 1, stencil, NULL, NULL, NULL, &da1)); 579566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da1)); 589566063dSJacob Faibussowitsch PetscCall(DMSetUp(da1)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Coordinates */ 619566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &coordsda)); 629566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(coordsda, &coords)); 639566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coords, "Grid coordinates")); 64c4762a1bSJed Brown for (i = 0, N = 1; i < 3; i++) { 65c4762a1bSJed Brown h[i] = 1.0 / dim[i]; 66c4762a1bSJed Brown PetscScalar *a; 679566063dSJacob Faibussowitsch PetscCall(VecGetArray(coords, &a)); 68c4762a1bSJed Brown PetscInt j, k, n = 0; 69c4762a1bSJed Brown for (i = 0; i < 3; ++i) { 70c4762a1bSJed Brown for (j = 0; j < dim[i]; ++j) { 71c4762a1bSJed Brown for (k = 0; k < 3; ++k) { 72c4762a1bSJed Brown a[n] = j * h[i]; /* coordinate along the j-th point in the i-th dimension */ 73c4762a1bSJed Brown ++n; 74c4762a1bSJed Brown } 75c4762a1bSJed Brown } 76c4762a1bSJed Brown } 779566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coords, &a)); 78c4762a1bSJed Brown } 799566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(da, coords)); 809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coords)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* Work vectors */ 839566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &x)); 849566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Real space vector")); 859566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &xx)); 869566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)xx, "Real space vector")); 879566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &y)); 889566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y, "USFFT frequency space vector")); 899566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &yy)); 909566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)yy, "FFTW frequency space vector")); 919566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &z)); 929566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)z, "USFFT reconstructed vector")); 939566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &zz)); 949566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)zz, "FFTW reconstructed vector")); 95c4762a1bSJed Brown /* Split vectors for FFTW */ 96c4762a1bSJed Brown for (ii = 0; ii < 3; ++ii) { 979566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da1, &xxsplit[ii])); 989566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)xxsplit[ii], "Real space split vector")); 999566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da1, &yysplit[ii])); 1009566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)yysplit[ii], "FFTW frequency space split vector")); 1019566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da1, &zzsplit[ii])); 1029566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)zzsplit[ii], "FFTW reconstructed split vector")); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 10563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%3-" PetscInt_FMT ": USFFT on vector of ")); 106c4762a1bSJed Brown for (i = 0, N = 1; i < 3; i++) { 1079566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "dim[%d] = %d ", i, dim[i])); 108c4762a1bSJed Brown N *= dim[i]; 109c4762a1bSJed Brown } 1109566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "; total size %d \n", N)); 111c4762a1bSJed Brown 112c4762a1bSJed Brown if (function == RANDOM) { 1139566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 1149566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 1159566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 1169566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 117c4762a1bSJed Brown } else if (function == CONSTANT) { 1189566063dSJacob Faibussowitsch PetscCall(VecSet(x, 1.0)); 119c4762a1bSJed Brown } else if (function == TANH) { 120c4762a1bSJed Brown PetscScalar *a; 1219566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &a)); 122c4762a1bSJed Brown PetscInt j, k = 0; 123c4762a1bSJed Brown for (i = 0; i < 3; ++i) { 124c4762a1bSJed Brown for (j = 0; j < dim[i]; ++j) { 125c4762a1bSJed Brown a[k] = tanh((j - dim[i] / 2.0) * (10.0 / dim[i])); 126c4762a1bSJed Brown ++k; 127c4762a1bSJed Brown } 128c4762a1bSJed Brown } 1299566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &a)); 130c4762a1bSJed Brown } 1311baa6e33SBarry Smith if (view_x) PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 1329566063dSJacob Faibussowitsch PetscCall(VecCopy(x, xx)); 133c4762a1bSJed Brown /* Split xx */ 1349566063dSJacob Faibussowitsch PetscCall(VecStrideGatherAll(xx, xxsplit, INSERT_VALUES)); /*YES! 'Gather' means 'split' (or maybe 'scatter'?)! */ 135c4762a1bSJed Brown 1369566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &norm)); 1379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|x|_2 = %g\n", norm)); 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* create USFFT object */ 1409566063dSJacob Faibussowitsch PetscCall(MatCreateSeqUSFFT(da, da, &A)); 141c4762a1bSJed Brown /* create FFTW object */ 1429566063dSJacob Faibussowitsch PetscCall(MatCreateSeqFFTW(PETSC_COMM_SELF, 3, dim, &AA)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* apply USFFT and FFTW FORWARD "preemptively", so the fftw_plans can be reused on different vectors */ 1459566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, z)); 146*48a46eb9SPierre Jolivet for (ii = 0; ii < 3; ++ii) PetscCall(MatMult(AA, xxsplit[ii], zzsplit[ii])); 147c4762a1bSJed Brown /* Now apply USFFT and FFTW forward several (3) times */ 148c4762a1bSJed Brown for (i = 0; i < 3; ++i) { 1499566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, y)); 150*48a46eb9SPierre Jolivet for (ii = 0; ii < 3; ++ii) PetscCall(MatMult(AA, xxsplit[ii], yysplit[ii])); 1519566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, y, z)); 152*48a46eb9SPierre Jolivet for (ii = 0; ii < 3; ++ii) PetscCall(MatMult(AA, yysplit[ii], zzsplit[ii])); 153c4762a1bSJed Brown } 154c4762a1bSJed Brown /* Unsplit yy */ 1559566063dSJacob Faibussowitsch PetscCall(VecStrideScatterAll(yysplit, yy, INSERT_VALUES)); /*YES! 'Scatter' means 'collect' (or maybe 'gather'?)! */ 156c4762a1bSJed Brown /* Unsplit zz */ 1579566063dSJacob Faibussowitsch PetscCall(VecStrideScatterAll(zzsplit, zz, INSERT_VALUES)); /*YES! 'Scatter' means 'collect' (or maybe 'gather'?)! */ 158c4762a1bSJed Brown 159c4762a1bSJed Brown if (view_y) { 1609566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "y = \n")); 1619566063dSJacob Faibussowitsch PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD)); 1629566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "yy = \n")); 1639566063dSJacob Faibussowitsch PetscCall(VecView(yy, PETSC_VIEWER_STDOUT_WORLD)); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown if (view_z) { 1679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "z = \n")); 1689566063dSJacob Faibussowitsch PetscCall(VecView(z, PETSC_VIEWER_STDOUT_WORLD)); 1699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "zz = \n")); 1709566063dSJacob Faibussowitsch PetscCall(VecView(zz, PETSC_VIEWER_STDOUT_WORLD)); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* compare x and z. USFFT computes an unnormalized DFT, thus z = N*x */ 174c4762a1bSJed Brown s = 1.0 / (PetscReal)N; 1759566063dSJacob Faibussowitsch PetscCall(VecScale(z, s)); 1769566063dSJacob Faibussowitsch PetscCall(VecAXPY(x, -1.0, z)); 1779566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_1, &enorm)); 1789566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|x-z| = %g\n", enorm)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* compare xx and zz. FFTW computes an unnormalized DFT, thus zz = N*x */ 181c4762a1bSJed Brown s = 1.0 / (PetscReal)N; 1829566063dSJacob Faibussowitsch PetscCall(VecScale(zz, s)); 1839566063dSJacob Faibussowitsch PetscCall(VecAXPY(xx, -1.0, zz)); 1849566063dSJacob Faibussowitsch PetscCall(VecNorm(xx, NORM_1, &enorm)); 1859566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|xx-zz| = %g\n", enorm)); 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* compare y and yy: USFFT and FFTW results*/ 1889566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_2, &norm)); 1899566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, yy)); 1909566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_1, &enorm)); 1919566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|y|_2 = %g\n", norm)); 1929566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|y-yy| = %g\n", enorm)); 193c4762a1bSJed Brown 194c4762a1bSJed Brown /* compare z and zz: USFFT and FFTW results*/ 1959566063dSJacob Faibussowitsch PetscCall(VecNorm(z, NORM_2, &norm)); 1969566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, -1.0, zz)); 1979566063dSJacob Faibussowitsch PetscCall(VecNorm(z, NORM_1, &enorm)); 1989566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|z|_2 = %g\n", norm)); 1999566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|z-zz| = %g\n", enorm)); 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* free spaces */ 2029566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &x)); 2039566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &xx)); 2049566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &y)); 2059566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &yy)); 2069566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &z)); 2079566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da, &zz)); 208c4762a1bSJed Brown 2099566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 210b122ec5aSJacob Faibussowitsch return 0; 211c4762a1bSJed Brown } 212