xref: /petsc/src/dm/tests/ex27.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
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   typedef enum {
15     RANDOM,
16     CONSTANT,
17     TANH,
18     NUM_FUNCS
19   } FuncType;
20   const char *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"};
21   Mat         A, AA;
22   PetscMPIInt size;
23   PetscInt    N, i, stencil = 1, dof = 1;
24   PetscInt    dim[3] = {10, 10, 10}, ndim = 3;
25   Vec         coords, x, y, z, xx, yy, zz;
26   PetscReal   h[3];
27   PetscScalar s;
28   PetscRandom rdm;
29   PetscReal   norm, enorm;
30   PetscInt    func;
31   FuncType    function = TANH;
32   DM          da, coordsda;
33   PetscBool   view_x = PETSC_FALSE, view_y = PETSC_FALSE, view_z = PETSC_FALSE;
34 
35   PetscFunctionBeginUser;
36   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
37   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
38   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only!");
39   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "USFFT Options", "ex27");
40   PetscCall(PetscOptionsEList("-function", "Function type", "ex27", funcNames, NUM_FUNCS, funcNames[function], &func, NULL));
41   function = (FuncType)func;
42   PetscOptionsEnd();
43   PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_x", &view_x, NULL));
44   PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_y", &view_y, NULL));
45   PetscCall(PetscOptionsGetBool(NULL, NULL, "-view_z", &view_z, NULL));
46   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dim", dim, &ndim, NULL));
47 
48   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, dof, stencil, NULL, NULL, NULL, &da));
49   PetscCall(DMSetFromOptions(da));
50   PetscCall(DMSetUp(da));
51 
52   /* Coordinates */
53   PetscCall(DMGetCoordinateDM(da, &coordsda));
54   PetscCall(DMGetGlobalVector(coordsda, &coords));
55   PetscCall(PetscObjectSetName((PetscObject)coords, "Grid coordinates"));
56   for (i = 0, N = 1; i < 3; i++) {
57     h[i] = 1.0 / dim[i];
58     PetscScalar *a;
59     PetscCall(VecGetArray(coords, &a));
60     PetscInt j, k, n = 0;
61     for (i = 0; i < 3; ++i) {
62       for (j = 0; j < dim[i]; ++j) {
63         for (k = 0; k < 3; ++k) {
64           a[n] = j * h[i]; /* coordinate along the j-th point in the i-th dimension */
65           ++n;
66         }
67       }
68     }
69     PetscCall(VecRestoreArray(coords, &a));
70   }
71   PetscCall(DMSetCoordinates(da, coords));
72 
73   /* Work vectors */
74   PetscCall(DMGetGlobalVector(da, &x));
75   PetscCall(PetscObjectSetName((PetscObject)x, "Real space vector"));
76   PetscCall(DMGetGlobalVector(da, &xx));
77   PetscCall(PetscObjectSetName((PetscObject)xx, "Real space vector"));
78   PetscCall(DMGetGlobalVector(da, &y));
79   PetscCall(PetscObjectSetName((PetscObject)y, "USFFT frequency space vector"));
80   PetscCall(DMGetGlobalVector(da, &yy));
81   PetscCall(PetscObjectSetName((PetscObject)yy, "FFTW frequency space vector"));
82   PetscCall(DMGetGlobalVector(da, &z));
83   PetscCall(PetscObjectSetName((PetscObject)z, "USFFT reconstructed vector"));
84   PetscCall(DMGetGlobalVector(da, &zz));
85   PetscCall(PetscObjectSetName((PetscObject)zz, "FFTW reconstructed vector"));
86 
87   PetscCall(PetscPrintf(PETSC_COMM_SELF, "%3-" PetscInt_FMT ": USFFT on vector of "));
88   for (i = 0, N = 1; i < 3; i++) {
89     PetscCall(PetscPrintf(PETSC_COMM_SELF, "dim[%d] = %d ", i, dim[i]));
90     N *= dim[i];
91   }
92   PetscCall(PetscPrintf(PETSC_COMM_SELF, "; total size %d \n", N));
93 
94   if (function == RANDOM) {
95     PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
96     PetscCall(PetscRandomSetFromOptions(rdm));
97     PetscCall(VecSetRandom(x, rdm));
98     PetscCall(PetscRandomDestroy(&rdm));
99   } else if (function == CONSTANT) {
100     PetscCall(VecSet(x, 1.0));
101   } else if (function == TANH) {
102     PetscScalar *a;
103     PetscCall(VecGetArray(x, &a));
104     PetscInt j, k = 0;
105     for (i = 0; i < 3; ++i) {
106       for (j = 0; j < dim[i]; ++j) {
107         a[k] = tanh((j - dim[i] / 2.0) * (10.0 / dim[i]));
108         ++k;
109       }
110     }
111     PetscCall(VecRestoreArray(x, &a));
112   }
113   if (view_x) PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
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