1 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"; 2 3 /* 4 Compiling the code: 5 This code uses the complex numbers version of PETSc and the FFTW package, so configure 6 must be run to enable these. 7 8 */ 9 10 #define DOF 3 11 12 #include <petscmat.h> 13 #include <petscdm.h> 14 #include <petscdmda.h> 15 int main(int argc,char **args) 16 { 17 typedef enum {RANDOM, CONSTANT, TANH, NUM_FUNCS} FuncType; 18 const char *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"}; 19 Mat A, AA; 20 PetscMPIInt size; 21 PetscInt N,i, stencil=1,dof=3; 22 PetscInt dim[3] = {10,10,10}, ndim = 3; 23 Vec coords,x,y,z,xx, yy, zz; 24 Vec xxsplit[DOF], yysplit[DOF], zzsplit[DOF]; 25 PetscReal h[3]; 26 PetscScalar s; 27 PetscRandom rdm; 28 PetscReal norm, enorm; 29 PetscInt func,ii; 30 FuncType function = TANH; 31 DM da, da1, coordsda; 32 PetscBool view_x = PETSC_FALSE, view_y = PETSC_FALSE, view_z = PETSC_FALSE; 33 PetscErrorCode ierr; 34 35 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 36 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 37 PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP, "This is a uniprocessor example only!"); 38 ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "USFFT Options", "ex27");CHKERRQ(ierr); 39 CHKERRQ(PetscOptionsEList("-function", "Function type", "ex27", funcNames, NUM_FUNCS, funcNames[function], &func, NULL)); 40 function = (FuncType) func; 41 ierr = PetscOptionsEnd();CHKERRQ(ierr); 42 CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_x",&view_x,NULL)); 43 CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_y",&view_y,NULL)); 44 CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_z",&view_z,NULL)); 45 CHKERRQ(PetscOptionsGetIntArray(NULL,NULL,"-dim",dim,&ndim,NULL)); 46 47 /* DMDA with the correct fiber dimension */ 48 ierr = 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, 49 dof,stencil,NULL,NULL,NULL,&da);CHKERRQ(ierr); 50 CHKERRQ(DMSetFromOptions(da)); 51 CHKERRQ(DMSetUp(da)); 52 /* DMDA with fiber dimension 1 for split fields */ 53 ierr = 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, 54 1,stencil,NULL,NULL,NULL,&da1);CHKERRQ(ierr); 55 CHKERRQ(DMSetFromOptions(da1)); 56 CHKERRQ(DMSetUp(da1)); 57 58 /* Coordinates */ 59 CHKERRQ(DMGetCoordinateDM(da,&coordsda)); 60 CHKERRQ(DMGetGlobalVector(coordsda,&coords)); 61 CHKERRQ(PetscObjectSetName((PetscObject) coords,"Grid coordinates")); 62 for (i = 0, N = 1; i < 3; i++) { 63 h[i] = 1.0/dim[i]; 64 PetscScalar *a; 65 CHKERRQ(VecGetArray(coords, &a)); 66 PetscInt j,k,n = 0; 67 for (i = 0; i < 3; ++i) { 68 for (j = 0; j < dim[i]; ++j) { 69 for (k = 0; k < 3; ++k) { 70 a[n] = j*h[i]; /* coordinate along the j-th point in the i-th dimension */ 71 ++n; 72 } 73 } 74 } 75 CHKERRQ(VecRestoreArray(coords, &a)); 76 77 } 78 CHKERRQ(DMSetCoordinates(da, coords)); 79 CHKERRQ(VecDestroy(&coords)); 80 81 /* Work vectors */ 82 CHKERRQ(DMGetGlobalVector(da, &x)); 83 CHKERRQ(PetscObjectSetName((PetscObject) x, "Real space vector")); 84 CHKERRQ(DMGetGlobalVector(da, &xx)); 85 CHKERRQ(PetscObjectSetName((PetscObject) xx, "Real space vector")); 86 CHKERRQ(DMGetGlobalVector(da, &y)); 87 CHKERRQ(PetscObjectSetName((PetscObject) y, "USFFT frequency space vector")); 88 CHKERRQ(DMGetGlobalVector(da, &yy)); 89 CHKERRQ(PetscObjectSetName((PetscObject) yy, "FFTW frequency space vector")); 90 CHKERRQ(DMGetGlobalVector(da, &z)); 91 CHKERRQ(PetscObjectSetName((PetscObject) z, "USFFT reconstructed vector")); 92 CHKERRQ(DMGetGlobalVector(da, &zz)); 93 CHKERRQ(PetscObjectSetName((PetscObject) zz, "FFTW reconstructed vector")); 94 /* Split vectors for FFTW */ 95 for (ii = 0; ii < 3; ++ii) { 96 CHKERRQ(DMGetGlobalVector(da1, &xxsplit[ii])); 97 CHKERRQ(PetscObjectSetName((PetscObject) xxsplit[ii], "Real space split vector")); 98 CHKERRQ(DMGetGlobalVector(da1, &yysplit[ii])); 99 CHKERRQ(PetscObjectSetName((PetscObject) yysplit[ii], "FFTW frequency space split vector")); 100 CHKERRQ(DMGetGlobalVector(da1, &zzsplit[ii])); 101 CHKERRQ(PetscObjectSetName((PetscObject) zzsplit[ii], "FFTW reconstructed split vector")); 102 } 103 104 CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "%3-D: USFFT on vector of ")); 105 for (i = 0, N = 1; i < 3; i++) { 106 CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "dim[%d] = %d ",i,dim[i])); 107 N *= dim[i]; 108 } 109 CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "; total size %d \n",N)); 110 111 if (function == RANDOM) { 112 CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 113 CHKERRQ(PetscRandomSetFromOptions(rdm)); 114 CHKERRQ(VecSetRandom(x, rdm)); 115 CHKERRQ(PetscRandomDestroy(&rdm)); 116 } else if (function == CONSTANT) { 117 CHKERRQ(VecSet(x, 1.0)); 118 } else if (function == TANH) { 119 PetscScalar *a; 120 CHKERRQ(VecGetArray(x, &a)); 121 PetscInt j,k = 0; 122 for (i = 0; i < 3; ++i) { 123 for (j = 0; j < dim[i]; ++j) { 124 a[k] = tanh((j - dim[i]/2.0)*(10.0/dim[i])); 125 ++k; 126 } 127 } 128 CHKERRQ(VecRestoreArray(x, &a)); 129 } 130 if (view_x) { 131 CHKERRQ(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 132 } 133 CHKERRQ(VecCopy(x,xx)); 134 /* Split xx */ 135 CHKERRQ(VecStrideGatherAll(xx,xxsplit, INSERT_VALUES)); /*YES! 'Gather' means 'split' (or maybe 'scatter'?)! */ 136 137 CHKERRQ(VecNorm(x,NORM_2,&norm)); 138 CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|x|_2 = %g\n",norm)); 139 140 /* create USFFT object */ 141 CHKERRQ(MatCreateSeqUSFFT(da,da,&A)); 142 /* create FFTW object */ 143 CHKERRQ(MatCreateSeqFFTW(PETSC_COMM_SELF,3,dim,&AA)); 144 145 /* apply USFFT and FFTW FORWARD "preemptively", so the fftw_plans can be reused on different vectors */ 146 CHKERRQ(MatMult(A,x,z)); 147 for (ii = 0; ii < 3; ++ii) { 148 CHKERRQ(MatMult(AA,xxsplit[ii],zzsplit[ii])); 149 } 150 /* Now apply USFFT and FFTW forward several (3) times */ 151 for (i=0; i<3; ++i) { 152 CHKERRQ(MatMult(A,x,y)); 153 for (ii = 0; ii < 3; ++ii) { 154 CHKERRQ(MatMult(AA,xxsplit[ii],yysplit[ii])); 155 } 156 CHKERRQ(MatMultTranspose(A,y,z)); 157 for (ii = 0; ii < 3; ++ii) { 158 CHKERRQ(MatMult(AA,yysplit[ii],zzsplit[ii])); 159 } 160 } 161 /* Unsplit yy */ 162 CHKERRQ(VecStrideScatterAll(yysplit, yy, INSERT_VALUES)); /*YES! 'Scatter' means 'collect' (or maybe 'gather'?)! */ 163 /* Unsplit zz */ 164 CHKERRQ(VecStrideScatterAll(zzsplit, zz, INSERT_VALUES)); /*YES! 'Scatter' means 'collect' (or maybe 'gather'?)! */ 165 166 if (view_y) { 167 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "y = \n")); 168 CHKERRQ(VecView(y, PETSC_VIEWER_STDOUT_WORLD)); 169 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "yy = \n")); 170 CHKERRQ(VecView(yy, PETSC_VIEWER_STDOUT_WORLD)); 171 } 172 173 if (view_z) { 174 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "z = \n")); 175 CHKERRQ(VecView(z, PETSC_VIEWER_STDOUT_WORLD)); 176 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "zz = \n")); 177 CHKERRQ(VecView(zz, PETSC_VIEWER_STDOUT_WORLD)); 178 } 179 180 /* compare x and z. USFFT computes an unnormalized DFT, thus z = N*x */ 181 s = 1.0/(PetscReal)N; 182 CHKERRQ(VecScale(z,s)); 183 CHKERRQ(VecAXPY(x,-1.0,z)); 184 CHKERRQ(VecNorm(x,NORM_1,&enorm)); 185 CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|x-z| = %g\n",enorm)); 186 187 /* compare xx and zz. FFTW computes an unnormalized DFT, thus zz = N*x */ 188 s = 1.0/(PetscReal)N; 189 CHKERRQ(VecScale(zz,s)); 190 CHKERRQ(VecAXPY(xx,-1.0,zz)); 191 CHKERRQ(VecNorm(xx,NORM_1,&enorm)); 192 CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|xx-zz| = %g\n",enorm)); 193 194 /* compare y and yy: USFFT and FFTW results*/ 195 CHKERRQ(VecNorm(y,NORM_2,&norm)); 196 CHKERRQ(VecAXPY(y,-1.0,yy)); 197 CHKERRQ(VecNorm(y,NORM_1,&enorm)); 198 CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|y|_2 = %g\n",norm)); 199 CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|y-yy| = %g\n",enorm)); 200 201 /* compare z and zz: USFFT and FFTW results*/ 202 CHKERRQ(VecNorm(z,NORM_2,&norm)); 203 CHKERRQ(VecAXPY(z,-1.0,zz)); 204 CHKERRQ(VecNorm(z,NORM_1,&enorm)); 205 CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|z|_2 = %g\n",norm)); 206 CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|z-zz| = %g\n",enorm)); 207 208 /* free spaces */ 209 CHKERRQ(DMRestoreGlobalVector(da,&x)); 210 CHKERRQ(DMRestoreGlobalVector(da,&xx)); 211 CHKERRQ(DMRestoreGlobalVector(da,&y)); 212 CHKERRQ(DMRestoreGlobalVector(da,&yy)); 213 CHKERRQ(DMRestoreGlobalVector(da,&z)); 214 CHKERRQ(DMRestoreGlobalVector(da,&zz)); 215 216 ierr = PetscFinalize(); 217 return ierr; 218 } 219