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