xref: /petsc/src/dm/tests/ex27.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
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 int main(int 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   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
33   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
34   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP, "This is a uniprocessor example only!");
35   ierr     = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "USFFT Options", "ex27");PetscCall(ierr);
36   PetscCall(PetscOptionsEList("-function", "Function type", "ex27", funcNames, NUM_FUNCS, funcNames[function], &func, NULL));
37   function = (FuncType) func;
38   ierr     = PetscOptionsEnd();PetscCall(ierr);
39   PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_x",&view_x,NULL));
40   PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_y",&view_y,NULL));
41   PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_z",&view_z,NULL));
42   PetscCall(PetscOptionsGetIntArray(NULL,NULL,"-dim",dim,&ndim,NULL));
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);PetscCall(ierr);
46   PetscCall(DMSetFromOptions(da));
47   PetscCall(DMSetUp(da));
48 
49   /* Coordinates */
50   PetscCall(DMGetCoordinateDM(da, &coordsda));
51   PetscCall(DMGetGlobalVector(coordsda, &coords));
52   PetscCall(PetscObjectSetName((PetscObject) coords, "Grid coordinates"));
53   for (i = 0, N = 1; i < 3; i++) {
54     h[i] = 1.0/dim[i];
55     PetscScalar *a;
56     PetscCall(VecGetArray(coords, &a));
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     PetscCall(VecRestoreArray(coords, &a));
67 
68   }
69   PetscCall(DMSetCoordinates(da, coords));
70 
71   /* Work vectors */
72   PetscCall(DMGetGlobalVector(da, &x));
73   PetscCall(PetscObjectSetName((PetscObject) x, "Real space vector"));
74   PetscCall(DMGetGlobalVector(da, &xx));
75   PetscCall(PetscObjectSetName((PetscObject) xx, "Real space vector"));
76   PetscCall(DMGetGlobalVector(da, &y));
77   PetscCall(PetscObjectSetName((PetscObject) y, "USFFT frequency space vector"));
78   PetscCall(DMGetGlobalVector(da, &yy));
79   PetscCall(PetscObjectSetName((PetscObject) yy, "FFTW frequency space vector"));
80   PetscCall(DMGetGlobalVector(da, &z));
81   PetscCall(PetscObjectSetName((PetscObject) z, "USFFT reconstructed vector"));
82   PetscCall(DMGetGlobalVector(da, &zz));
83   PetscCall(PetscObjectSetName((PetscObject) zz, "FFTW reconstructed vector"));
84 
85   PetscCall(PetscPrintf(PETSC_COMM_SELF, "%3-D: USFFT on vector of "));
86   for (i = 0, N = 1; i < 3; i++) {
87     PetscCall(PetscPrintf(PETSC_COMM_SELF, "dim[%d] = %d ",i,dim[i]));
88     N   *= dim[i];
89   }
90   PetscCall(PetscPrintf(PETSC_COMM_SELF, "; total size %d \n",N));
91 
92   if (function == RANDOM) {
93     PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
94     PetscCall(PetscRandomSetFromOptions(rdm));
95     PetscCall(VecSetRandom(x, rdm));
96     PetscCall(PetscRandomDestroy(&rdm));
97   } else if (function == CONSTANT) {
98     PetscCall(VecSet(x, 1.0));
99   } else if (function == TANH) {
100     PetscScalar *a;
101     PetscCall(VecGetArray(x, &a));
102     PetscInt j,k = 0;
103     for (i = 0; i < 3; ++i) {
104       for (j = 0; j < dim[i]; ++j) {
105         a[k] = tanh((j - dim[i]/2.0)*(10.0/dim[i]));
106         ++k;
107       }
108     }
109     PetscCall(VecRestoreArray(x, &a));
110   }
111   if (view_x) {
112     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
113   }
114   PetscCall(VecCopy(x,xx));
115 
116   PetscCall(VecNorm(x,NORM_2,&norm));
117   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|x|_2 = %g\n",norm));
118 
119   /* create USFFT object */
120   PetscCall(MatCreateSeqUSFFT(coords,da,&A));
121   /* create FFTW object */
122   PetscCall(MatCreateSeqFFTW(PETSC_COMM_SELF,3,dim,&AA));
123 
124   /* apply USFFT and FFTW FORWARD "preemptively", so the fftw_plans can be reused on different vectors */
125   PetscCall(MatMult(A,x,z));
126   PetscCall(MatMult(AA,xx,zz));
127   /* Now apply USFFT and FFTW forward several (3) times */
128   for (i=0; i<3; ++i) {
129     PetscCall(MatMult(A,x,y));
130     PetscCall(MatMult(AA,xx,yy));
131     PetscCall(MatMultTranspose(A,y,z));
132     PetscCall(MatMultTranspose(AA,yy,zz));
133   }
134 
135   if (view_y) {
136     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "y = \n"));
137     PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
138     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "yy = \n"));
139     PetscCall(VecView(yy, PETSC_VIEWER_STDOUT_WORLD));
140   }
141 
142   if (view_z) {
143     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "z = \n"));
144     PetscCall(VecView(z, PETSC_VIEWER_STDOUT_WORLD));
145     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "zz = \n"));
146     PetscCall(VecView(zz, PETSC_VIEWER_STDOUT_WORLD));
147   }
148 
149   /* compare x and z. USFFT computes an unnormalized DFT, thus z = N*x */
150   s    = 1.0/(PetscReal)N;
151   PetscCall(VecScale(z,s));
152   PetscCall(VecAXPY(x,-1.0,z));
153   PetscCall(VecNorm(x,NORM_1,&enorm));
154   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|x-z| = %g\n",enorm));
155 
156   /* compare xx and zz. FFTW computes an unnormalized DFT, thus zz = N*x */
157   s    = 1.0/(PetscReal)N;
158   PetscCall(VecScale(zz,s));
159   PetscCall(VecAXPY(xx,-1.0,zz));
160   PetscCall(VecNorm(xx,NORM_1,&enorm));
161   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|xx-zz| = %g\n",enorm));
162 
163   /* compare y and yy: USFFT and FFTW results*/
164   PetscCall(VecNorm(y,NORM_2,&norm));
165   PetscCall(VecAXPY(y,-1.0,yy));
166   PetscCall(VecNorm(y,NORM_1,&enorm));
167   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|y|_2 = %g\n",norm));
168   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|y-yy| = %g\n",enorm));
169 
170   /* compare z and zz: USFFT and FFTW results*/
171   PetscCall(VecNorm(z,NORM_2,&norm));
172   PetscCall(VecAXPY(z,-1.0,zz));
173   PetscCall(VecNorm(z,NORM_1,&enorm));
174   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|z|_2 = %g\n",norm));
175   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|z-zz| = %g\n",enorm));
176 
177   /* free spaces */
178   PetscCall(DMRestoreGlobalVector(da,&x));
179   PetscCall(DMRestoreGlobalVector(da,&xx));
180   PetscCall(DMRestoreGlobalVector(da,&y));
181   PetscCall(DMRestoreGlobalVector(da,&yy));
182   PetscCall(DMRestoreGlobalVector(da,&z));
183   PetscCall(DMRestoreGlobalVector(da,&zz));
184   PetscCall(VecDestroy(&coords));
185   PetscCall(DMDestroy(&da));
186   PetscCall(PetscFinalize());
187   return 0;
188 }
189