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> 13ad0e0ea3SStefano Zampini int main(int argc,char **args) 14c4762a1bSJed Brown { 15c4762a1bSJed Brown typedef enum {RANDOM, CONSTANT, TANH, NUM_FUNCS} FuncType; 16c4762a1bSJed Brown const char *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"}; 17c4762a1bSJed Brown Mat A, AA; 18c4762a1bSJed Brown PetscMPIInt size; 19c4762a1bSJed Brown PetscInt N,i, stencil=1,dof=1; 20c4762a1bSJed Brown PetscInt dim[3] = {10,10,10}, ndim = 3; 21c4762a1bSJed Brown Vec coords,x,y,z,xx,yy,zz; 22c4762a1bSJed Brown PetscReal h[3]; 23c4762a1bSJed Brown PetscScalar s; 24c4762a1bSJed Brown PetscRandom rdm; 25c4762a1bSJed Brown PetscReal norm, enorm; 26c4762a1bSJed Brown PetscInt func; 27c4762a1bSJed Brown FuncType function = TANH; 28c4762a1bSJed Brown DM da, coordsda; 29c4762a1bSJed Brown PetscBool view_x = PETSC_FALSE, view_y = PETSC_FALSE, view_z = PETSC_FALSE; 30c4762a1bSJed Brown 31*327415f7SBarry Smith PetscFunctionBeginUser; 329566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 3408401ef6SPierre Jolivet PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_SUP, "This is a uniprocessor example only!"); 35d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "USFFT Options", "ex27"); 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-function", "Function type", "ex27", funcNames, NUM_FUNCS, funcNames[function], &func, NULL)); 37c4762a1bSJed Brown function = (FuncType) func; 38d0609cedSBarry Smith PetscOptionsEnd(); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_x",&view_x,NULL)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_y",&view_y,NULL)); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_z",&view_z,NULL)); 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL,NULL,"-dim",dim,&ndim,NULL)); 43c4762a1bSJed Brown 44d0609cedSBarry Smith PetscCall(DMDACreate3d(PETSC_COMM_SELF,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,dim[0], dim[1], dim[2], 45d0609cedSBarry Smith PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE,dof, stencil,NULL, NULL, NULL,&da)); 469566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 479566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* Coordinates */ 509566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &coordsda)); 519566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(coordsda, &coords)); 529566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) coords, "Grid coordinates")); 53c4762a1bSJed Brown for (i = 0, N = 1; i < 3; i++) { 54c4762a1bSJed Brown h[i] = 1.0/dim[i]; 55c4762a1bSJed Brown PetscScalar *a; 569566063dSJacob Faibussowitsch PetscCall(VecGetArray(coords, &a)); 57c4762a1bSJed Brown PetscInt j,k,n = 0; 58c4762a1bSJed Brown for (i = 0; i < 3; ++i) { 59c4762a1bSJed Brown for (j = 0; j < dim[i]; ++j) { 60c4762a1bSJed Brown for (k = 0; k < 3; ++k) { 61c4762a1bSJed Brown a[n] = j*h[i]; /* coordinate along the j-th point in the i-th dimension */ 62c4762a1bSJed Brown ++n; 63c4762a1bSJed Brown } 64c4762a1bSJed Brown } 65c4762a1bSJed Brown } 669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coords, &a)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown } 699566063dSJacob Faibussowitsch PetscCall(DMSetCoordinates(da, coords)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* Work vectors */ 729566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &x)); 739566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) x, "Real space vector")); 749566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &xx)); 759566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) xx, "Real space vector")); 769566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &y)); 779566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) y, "USFFT frequency space vector")); 789566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &yy)); 799566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) yy, "FFTW frequency space vector")); 809566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &z)); 819566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) z, "USFFT reconstructed vector")); 829566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(da, &zz)); 839566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) zz, "FFTW reconstructed vector")); 84c4762a1bSJed Brown 8563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%3-" PetscInt_FMT ": USFFT on vector of ")); 86c4762a1bSJed Brown for (i = 0, N = 1; i < 3; i++) { 879566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "dim[%d] = %d ",i,dim[i])); 88c4762a1bSJed Brown N *= dim[i]; 89c4762a1bSJed Brown } 909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "; total size %d \n",N)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown if (function == RANDOM) { 939566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 949566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 959566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 969566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 97c4762a1bSJed Brown } else if (function == CONSTANT) { 989566063dSJacob Faibussowitsch PetscCall(VecSet(x, 1.0)); 99c4762a1bSJed Brown } else if (function == TANH) { 100c4762a1bSJed Brown PetscScalar *a; 1019566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &a)); 102c4762a1bSJed Brown PetscInt j,k = 0; 103c4762a1bSJed Brown for (i = 0; i < 3; ++i) { 104c4762a1bSJed Brown for (j = 0; j < dim[i]; ++j) { 105c4762a1bSJed Brown a[k] = tanh((j - dim[i]/2.0)*(10.0/dim[i])); 106c4762a1bSJed Brown ++k; 107c4762a1bSJed Brown } 108c4762a1bSJed Brown } 1099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &a)); 110c4762a1bSJed Brown } 1111baa6e33SBarry Smith if (view_x) PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 1129566063dSJacob Faibussowitsch PetscCall(VecCopy(x,xx)); 113c4762a1bSJed Brown 1149566063dSJacob Faibussowitsch PetscCall(VecNorm(x,NORM_2,&norm)); 1159566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|x|_2 = %g\n",norm)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* create USFFT object */ 1189566063dSJacob Faibussowitsch PetscCall(MatCreateSeqUSFFT(coords,da,&A)); 119c4762a1bSJed Brown /* create FFTW object */ 1209566063dSJacob Faibussowitsch PetscCall(MatCreateSeqFFTW(PETSC_COMM_SELF,3,dim,&AA)); 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* apply USFFT and FFTW FORWARD "preemptively", so the fftw_plans can be reused on different vectors */ 1239566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,z)); 1249566063dSJacob Faibussowitsch PetscCall(MatMult(AA,xx,zz)); 125c4762a1bSJed Brown /* Now apply USFFT and FFTW forward several (3) times */ 126c4762a1bSJed Brown for (i=0; i<3; ++i) { 1279566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,y)); 1289566063dSJacob Faibussowitsch PetscCall(MatMult(AA,xx,yy)); 1299566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,y,z)); 1309566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(AA,yy,zz)); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown if (view_y) { 1349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "y = \n")); 1359566063dSJacob Faibussowitsch PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD)); 1369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "yy = \n")); 1379566063dSJacob Faibussowitsch PetscCall(VecView(yy, PETSC_VIEWER_STDOUT_WORLD)); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown 140c4762a1bSJed Brown if (view_z) { 1419566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "z = \n")); 1429566063dSJacob Faibussowitsch PetscCall(VecView(z, PETSC_VIEWER_STDOUT_WORLD)); 1439566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "zz = \n")); 1449566063dSJacob Faibussowitsch PetscCall(VecView(zz, PETSC_VIEWER_STDOUT_WORLD)); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* compare x and z. USFFT computes an unnormalized DFT, thus z = N*x */ 148c4762a1bSJed Brown s = 1.0/(PetscReal)N; 1499566063dSJacob Faibussowitsch PetscCall(VecScale(z,s)); 1509566063dSJacob Faibussowitsch PetscCall(VecAXPY(x,-1.0,z)); 1519566063dSJacob Faibussowitsch PetscCall(VecNorm(x,NORM_1,&enorm)); 1529566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|x-z| = %g\n",enorm)); 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* compare xx and zz. FFTW computes an unnormalized DFT, thus zz = N*x */ 155c4762a1bSJed Brown s = 1.0/(PetscReal)N; 1569566063dSJacob Faibussowitsch PetscCall(VecScale(zz,s)); 1579566063dSJacob Faibussowitsch PetscCall(VecAXPY(xx,-1.0,zz)); 1589566063dSJacob Faibussowitsch PetscCall(VecNorm(xx,NORM_1,&enorm)); 1599566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|xx-zz| = %g\n",enorm)); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* compare y and yy: USFFT and FFTW results*/ 1629566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&norm)); 1639566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,-1.0,yy)); 1649566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_1,&enorm)); 1659566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|y|_2 = %g\n",norm)); 1669566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|y-yy| = %g\n",enorm)); 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* compare z and zz: USFFT and FFTW results*/ 1699566063dSJacob Faibussowitsch PetscCall(VecNorm(z,NORM_2,&norm)); 1709566063dSJacob Faibussowitsch PetscCall(VecAXPY(z,-1.0,zz)); 1719566063dSJacob Faibussowitsch PetscCall(VecNorm(z,NORM_1,&enorm)); 1729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|z|_2 = %g\n",norm)); 1739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "|z-zz| = %g\n",enorm)); 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* free spaces */ 1769566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da,&x)); 1779566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da,&xx)); 1789566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da,&y)); 1799566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da,&yy)); 1809566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da,&z)); 1819566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(da,&zz)); 1829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coords)); 1839566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1849566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 185b122ec5aSJacob Faibussowitsch return 0; 186c4762a1bSJed Brown } 187