xref: /petsc/src/dm/tests/ex28.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 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";
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 #define DOF 3
11c4762a1bSJed Brown 
12c4762a1bSJed Brown #include <petscmat.h>
13c4762a1bSJed Brown #include <petscdm.h>
14c4762a1bSJed Brown #include <petscdmda.h>
15ad0e0ea3SStefano Zampini int main(int argc,char **args)
16c4762a1bSJed Brown {
17c4762a1bSJed Brown   typedef enum {RANDOM, CONSTANT, TANH, NUM_FUNCS} FuncType;
18c4762a1bSJed Brown   const char     *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"};
19c4762a1bSJed Brown   Mat            A, AA;
20c4762a1bSJed Brown   PetscMPIInt    size;
21c4762a1bSJed Brown   PetscInt       N,i, stencil=1,dof=3;
22c4762a1bSJed Brown   PetscInt       dim[3] = {10,10,10}, ndim = 3;
23c4762a1bSJed Brown   Vec            coords,x,y,z,xx, yy, zz;
24c4762a1bSJed Brown   Vec            xxsplit[DOF], yysplit[DOF], zzsplit[DOF];
25c4762a1bSJed Brown   PetscReal      h[3];
26c4762a1bSJed Brown   PetscScalar    s;
27c4762a1bSJed Brown   PetscRandom    rdm;
28c4762a1bSJed Brown   PetscReal      norm, enorm;
29c4762a1bSJed Brown   PetscInt       func,ii;
30c4762a1bSJed Brown   FuncType       function = TANH;
31c4762a1bSJed Brown   DM             da, da1, coordsda;
32c4762a1bSJed Brown   PetscBool      view_x = PETSC_FALSE, view_y = PETSC_FALSE, view_z = PETSC_FALSE;
33c4762a1bSJed Brown   PetscErrorCode ierr;
34c4762a1bSJed Brown 
35*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
365f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
372c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP, "This is a uniprocessor example only!");
38c4762a1bSJed Brown   ierr     = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "USFFT Options", "ex27");CHKERRQ(ierr);
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEList("-function", "Function type", "ex27", funcNames, NUM_FUNCS, funcNames[function], &func, NULL));
40c4762a1bSJed Brown   function = (FuncType) func;
41c4762a1bSJed Brown   ierr     = PetscOptionsEnd();CHKERRQ(ierr);
425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_x",&view_x,NULL));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_y",&view_y,NULL));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_z",&view_z,NULL));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetIntArray(NULL,NULL,"-dim",dim,&ndim,NULL));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   /* DMDA with the correct fiber dimension */
48c4762a1bSJed Brown   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,
49c4762a1bSJed Brown                       dof,stencil,NULL,NULL,NULL,&da);CHKERRQ(ierr);
505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
52c4762a1bSJed Brown   /* DMDA with fiber dimension 1 for split fields */
53c4762a1bSJed Brown   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,
54c4762a1bSJed Brown                       1,stencil,NULL,NULL,NULL,&da1);CHKERRQ(ierr);
555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da1));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da1));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* Coordinates */
595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(da,&coordsda));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(coordsda,&coords));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) coords,"Grid coordinates"));
62c4762a1bSJed Brown   for (i = 0, N = 1; i < 3; i++) {
63c4762a1bSJed Brown     h[i] = 1.0/dim[i];
64c4762a1bSJed Brown     PetscScalar *a;
655f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(coords, &a));
66c4762a1bSJed Brown     PetscInt j,k,n = 0;
67c4762a1bSJed Brown     for (i = 0; i < 3; ++i) {
68c4762a1bSJed Brown       for (j = 0; j < dim[i]; ++j) {
69c4762a1bSJed Brown         for (k = 0; k < 3; ++k) {
70c4762a1bSJed Brown           a[n] = j*h[i]; /* coordinate along the j-th point in the i-th dimension */
71c4762a1bSJed Brown           ++n;
72c4762a1bSJed Brown         }
73c4762a1bSJed Brown       }
74c4762a1bSJed Brown     }
755f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(coords, &a));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   }
785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetCoordinates(da, coords));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&coords));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /* Work vectors */
825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(da, &x));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) x, "Real space vector"));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(da, &xx));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) xx, "Real space vector"));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(da, &y));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) y, "USFFT frequency space vector"));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(da, &yy));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) yy, "FFTW frequency space vector"));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(da, &z));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) z, "USFFT reconstructed vector"));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(da, &zz));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) zz, "FFTW reconstructed vector"));
94c4762a1bSJed Brown   /* Split vectors for FFTW */
95c4762a1bSJed Brown   for (ii = 0; ii < 3; ++ii) {
965f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetGlobalVector(da1, &xxsplit[ii]));
975f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) xxsplit[ii], "Real space split vector"));
985f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetGlobalVector(da1, &yysplit[ii]));
995f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) yysplit[ii], "FFTW frequency space split vector"));
1005f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetGlobalVector(da1, &zzsplit[ii]));
1015f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) zzsplit[ii], "FFTW reconstructed split vector"));
102c4762a1bSJed Brown   }
103c4762a1bSJed Brown 
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "%3-D: USFFT on vector of "));
105c4762a1bSJed Brown   for (i = 0, N = 1; i < 3; i++) {
1065f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "dim[%d] = %d ",i,dim[i]));
107c4762a1bSJed Brown     N   *= dim[i];
108c4762a1bSJed Brown   }
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "; total size %d \n",N));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   if (function == RANDOM) {
1125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
1135f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomSetFromOptions(rdm));
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(x, rdm));
1155f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomDestroy(&rdm));
116c4762a1bSJed Brown   } else if (function == CONSTANT) {
1175f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(x, 1.0));
118c4762a1bSJed Brown   } else if (function == TANH) {
119c4762a1bSJed Brown     PetscScalar *a;
1205f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(x, &a));
121c4762a1bSJed Brown     PetscInt j,k = 0;
122c4762a1bSJed Brown     for (i = 0; i < 3; ++i) {
123c4762a1bSJed Brown       for (j = 0; j < dim[i]; ++j) {
124c4762a1bSJed Brown         a[k] = tanh((j - dim[i]/2.0)*(10.0/dim[i]));
125c4762a1bSJed Brown         ++k;
126c4762a1bSJed Brown       }
127c4762a1bSJed Brown     }
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(x, &a));
129c4762a1bSJed Brown   }
130c4762a1bSJed Brown   if (view_x) {
1315f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
132c4762a1bSJed Brown   }
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(x,xx));
134c4762a1bSJed Brown   /* Split xx */
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecStrideGatherAll(xx,xxsplit, INSERT_VALUES)); /*YES! 'Gather' means 'split' (or maybe 'scatter'?)! */
136c4762a1bSJed Brown 
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(x,NORM_2,&norm));
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|x|_2 = %g\n",norm));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* create USFFT object */
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqUSFFT(da,da,&A));
142c4762a1bSJed Brown   /* create FFTW object */
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqFFTW(PETSC_COMM_SELF,3,dim,&AA));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* apply USFFT and FFTW FORWARD "preemptively", so the fftw_plans can be reused on different vectors */
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,x,z));
147c4762a1bSJed Brown   for (ii = 0; ii < 3; ++ii) {
1485f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(AA,xxsplit[ii],zzsplit[ii]));
149c4762a1bSJed Brown   }
150c4762a1bSJed Brown   /* Now apply USFFT and FFTW forward several (3) times */
151c4762a1bSJed Brown   for (i=0; i<3; ++i) {
1525f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(A,x,y));
153c4762a1bSJed Brown     for (ii = 0; ii < 3; ++ii) {
1545f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(AA,xxsplit[ii],yysplit[ii]));
155c4762a1bSJed Brown     }
1565f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTranspose(A,y,z));
157c4762a1bSJed Brown     for (ii = 0; ii < 3; ++ii) {
1585f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMult(AA,yysplit[ii],zzsplit[ii]));
159c4762a1bSJed Brown     }
160c4762a1bSJed Brown   }
161c4762a1bSJed Brown   /* Unsplit yy */
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecStrideScatterAll(yysplit, yy, INSERT_VALUES)); /*YES! 'Scatter' means 'collect' (or maybe 'gather'?)! */
163c4762a1bSJed Brown   /* Unsplit zz */
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecStrideScatterAll(zzsplit, zz, INSERT_VALUES)); /*YES! 'Scatter' means 'collect' (or maybe 'gather'?)! */
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   if (view_y) {
1675f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "y = \n"));
1685f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
1695f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "yy = \n"));
1705f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(yy, PETSC_VIEWER_STDOUT_WORLD));
171c4762a1bSJed Brown   }
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   if (view_z) {
1745f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "z = \n"));
1755f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(z, PETSC_VIEWER_STDOUT_WORLD));
1765f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "zz = \n"));
1775f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(zz, PETSC_VIEWER_STDOUT_WORLD));
178c4762a1bSJed Brown   }
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   /* compare x and z. USFFT computes an unnormalized DFT, thus z = N*x */
181c4762a1bSJed Brown   s    = 1.0/(PetscReal)N;
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(z,s));
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(x,-1.0,z));
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(x,NORM_1,&enorm));
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|x-z| = %g\n",enorm));
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /* compare xx and zz. FFTW computes an unnormalized DFT, thus zz = N*x */
188c4762a1bSJed Brown   s    = 1.0/(PetscReal)N;
1895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(zz,s));
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(xx,-1.0,zz));
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(xx,NORM_1,&enorm));
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|xx-zz| = %g\n",enorm));
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   /* compare y and yy: USFFT and FFTW results*/
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(y,NORM_2,&norm));
1965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(y,-1.0,yy));
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(y,NORM_1,&enorm));
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|y|_2 = %g\n",norm));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|y-yy| = %g\n",enorm));
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /* compare z and zz: USFFT and FFTW results*/
2025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(z,NORM_2,&norm));
2035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(z,-1.0,zz));
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(z,NORM_1,&enorm));
2055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|z|_2 = %g\n",norm));
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "|z-zz| = %g\n",enorm));
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   /* free spaces */
2095f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(da,&x));
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(da,&xx));
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(da,&y));
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(da,&yy));
2135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(da,&z));
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(da,&zz));
215c4762a1bSJed Brown 
216*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
217*b122ec5aSJacob Faibussowitsch   return 0;
218c4762a1bSJed Brown }
219