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