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