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