xref: /petsc/src/dm/tests/ex28.c (revision 2e65eb737b5b07432530db55f6b2a145ebc548b2)
1 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";
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 #define DOF 3
11 
12 #include <petscmat.h>
13 #include <petscdm.h>
14 #include <petscdmda.h>
15 int main(int argc,char **args)
16 {
17   typedef enum {RANDOM, CONSTANT, TANH, NUM_FUNCS} FuncType;
18   const char     *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"};
19   Mat            A, AA;
20   PetscMPIInt    size;
21   PetscInt       N,i, stencil=1,dof=3;
22   PetscInt       dim[3] = {10,10,10}, ndim = 3;
23   Vec            coords,x,y,z,xx, yy, zz;
24   Vec            xxsplit[DOF], yysplit[DOF], zzsplit[DOF];
25   PetscReal      h[3];
26   PetscScalar    s;
27   PetscRandom    rdm;
28   PetscReal      norm, enorm;
29   PetscInt       func,ii;
30   FuncType       function = TANH;
31   DM             da, da1, coordsda;
32   PetscBool      view_x = PETSC_FALSE, view_y = PETSC_FALSE, view_z = PETSC_FALSE;
33 
34   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
35   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
36   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_SUP, "This is a uniprocessor example only!");
37   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "USFFT Options", "ex27");
38   PetscCall(PetscOptionsEList("-function", "Function type", "ex27", funcNames, NUM_FUNCS, funcNames[function], &func, NULL));
39   function = (FuncType) func;
40   PetscOptionsEnd();
41   PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_x",&view_x,NULL));
42   PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_y",&view_y,NULL));
43   PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_z",&view_z,NULL));
44   PetscCall(PetscOptionsGetIntArray(NULL,NULL,"-dim",dim,&ndim,NULL));
45 
46   /* DMDA with the correct fiber dimension */
47   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,
48                          dof,stencil,NULL,NULL,NULL,&da));
49   PetscCall(DMSetFromOptions(da));
50   PetscCall(DMSetUp(da));
51   /* DMDA with fiber dimension 1 for split fields */
52   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,
53                          1,stencil,NULL,NULL,NULL,&da1));
54   PetscCall(DMSetFromOptions(da1));
55   PetscCall(DMSetUp(da1));
56 
57   /* Coordinates */
58   PetscCall(DMGetCoordinateDM(da,&coordsda));
59   PetscCall(DMGetGlobalVector(coordsda,&coords));
60   PetscCall(PetscObjectSetName((PetscObject) coords,"Grid coordinates"));
61   for (i = 0, N = 1; i < 3; i++) {
62     h[i] = 1.0/dim[i];
63     PetscScalar *a;
64     PetscCall(VecGetArray(coords, &a));
65     PetscInt j,k,n = 0;
66     for (i = 0; i < 3; ++i) {
67       for (j = 0; j < dim[i]; ++j) {
68         for (k = 0; k < 3; ++k) {
69           a[n] = j*h[i]; /* coordinate along the j-th point in the i-th dimension */
70           ++n;
71         }
72       }
73     }
74     PetscCall(VecRestoreArray(coords, &a));
75 
76   }
77   PetscCall(DMSetCoordinates(da, coords));
78   PetscCall(VecDestroy(&coords));
79 
80   /* Work vectors */
81   PetscCall(DMGetGlobalVector(da, &x));
82   PetscCall(PetscObjectSetName((PetscObject) x, "Real space vector"));
83   PetscCall(DMGetGlobalVector(da, &xx));
84   PetscCall(PetscObjectSetName((PetscObject) xx, "Real space vector"));
85   PetscCall(DMGetGlobalVector(da, &y));
86   PetscCall(PetscObjectSetName((PetscObject) y, "USFFT frequency space vector"));
87   PetscCall(DMGetGlobalVector(da, &yy));
88   PetscCall(PetscObjectSetName((PetscObject) yy, "FFTW frequency space vector"));
89   PetscCall(DMGetGlobalVector(da, &z));
90   PetscCall(PetscObjectSetName((PetscObject) z, "USFFT reconstructed vector"));
91   PetscCall(DMGetGlobalVector(da, &zz));
92   PetscCall(PetscObjectSetName((PetscObject) zz, "FFTW reconstructed vector"));
93   /* Split vectors for FFTW */
94   for (ii = 0; ii < 3; ++ii) {
95     PetscCall(DMGetGlobalVector(da1, &xxsplit[ii]));
96     PetscCall(PetscObjectSetName((PetscObject) xxsplit[ii], "Real space split vector"));
97     PetscCall(DMGetGlobalVector(da1, &yysplit[ii]));
98     PetscCall(PetscObjectSetName((PetscObject) yysplit[ii], "FFTW frequency space split vector"));
99     PetscCall(DMGetGlobalVector(da1, &zzsplit[ii]));
100     PetscCall(PetscObjectSetName((PetscObject) zzsplit[ii], "FFTW reconstructed split vector"));
101   }
102 
103   PetscCall(PetscPrintf(PETSC_COMM_SELF, "%3-" PetscInt_FMT ": USFFT on vector of "));
104   for (i = 0, N = 1; i < 3; i++) {
105     PetscCall(PetscPrintf(PETSC_COMM_SELF, "dim[%d] = %d ",i,dim[i]));
106     N   *= dim[i];
107   }
108   PetscCall(PetscPrintf(PETSC_COMM_SELF, "; total size %d \n",N));
109 
110   if (function == RANDOM) {
111     PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
112     PetscCall(PetscRandomSetFromOptions(rdm));
113     PetscCall(VecSetRandom(x, rdm));
114     PetscCall(PetscRandomDestroy(&rdm));
115   } else if (function == CONSTANT) {
116     PetscCall(VecSet(x, 1.0));
117   } else if (function == TANH) {
118     PetscScalar *a;
119     PetscCall(VecGetArray(x, &a));
120     PetscInt j,k = 0;
121     for (i = 0; i < 3; ++i) {
122       for (j = 0; j < dim[i]; ++j) {
123         a[k] = tanh((j - dim[i]/2.0)*(10.0/dim[i]));
124         ++k;
125       }
126     }
127     PetscCall(VecRestoreArray(x, &a));
128   }
129   if (view_x) {
130     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
131   }
132   PetscCall(VecCopy(x,xx));
133   /* Split xx */
134   PetscCall(VecStrideGatherAll(xx,xxsplit, INSERT_VALUES)); /*YES! 'Gather' means 'split' (or maybe 'scatter'?)! */
135 
136   PetscCall(VecNorm(x,NORM_2,&norm));
137   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|x|_2 = %g\n",norm));
138 
139   /* create USFFT object */
140   PetscCall(MatCreateSeqUSFFT(da,da,&A));
141   /* create FFTW object */
142   PetscCall(MatCreateSeqFFTW(PETSC_COMM_SELF,3,dim,&AA));
143 
144   /* apply USFFT and FFTW FORWARD "preemptively", so the fftw_plans can be reused on different vectors */
145   PetscCall(MatMult(A,x,z));
146   for (ii = 0; ii < 3; ++ii) {
147     PetscCall(MatMult(AA,xxsplit[ii],zzsplit[ii]));
148   }
149   /* Now apply USFFT and FFTW forward several (3) times */
150   for (i=0; i<3; ++i) {
151     PetscCall(MatMult(A,x,y));
152     for (ii = 0; ii < 3; ++ii) {
153       PetscCall(MatMult(AA,xxsplit[ii],yysplit[ii]));
154     }
155     PetscCall(MatMultTranspose(A,y,z));
156     for (ii = 0; ii < 3; ++ii) {
157       PetscCall(MatMult(AA,yysplit[ii],zzsplit[ii]));
158     }
159   }
160   /* Unsplit yy */
161   PetscCall(VecStrideScatterAll(yysplit, yy, INSERT_VALUES)); /*YES! 'Scatter' means 'collect' (or maybe 'gather'?)! */
162   /* Unsplit zz */
163   PetscCall(VecStrideScatterAll(zzsplit, zz, INSERT_VALUES)); /*YES! 'Scatter' means 'collect' (or maybe 'gather'?)! */
164 
165   if (view_y) {
166     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "y = \n"));
167     PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
168     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "yy = \n"));
169     PetscCall(VecView(yy, PETSC_VIEWER_STDOUT_WORLD));
170   }
171 
172   if (view_z) {
173     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "z = \n"));
174     PetscCall(VecView(z, PETSC_VIEWER_STDOUT_WORLD));
175     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "zz = \n"));
176     PetscCall(VecView(zz, PETSC_VIEWER_STDOUT_WORLD));
177   }
178 
179   /* compare x and z. USFFT computes an unnormalized DFT, thus z = N*x */
180   s    = 1.0/(PetscReal)N;
181   PetscCall(VecScale(z,s));
182   PetscCall(VecAXPY(x,-1.0,z));
183   PetscCall(VecNorm(x,NORM_1,&enorm));
184   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|x-z| = %g\n",enorm));
185 
186   /* compare xx and zz. FFTW computes an unnormalized DFT, thus zz = N*x */
187   s    = 1.0/(PetscReal)N;
188   PetscCall(VecScale(zz,s));
189   PetscCall(VecAXPY(xx,-1.0,zz));
190   PetscCall(VecNorm(xx,NORM_1,&enorm));
191   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|xx-zz| = %g\n",enorm));
192 
193   /* compare y and yy: USFFT and FFTW results*/
194   PetscCall(VecNorm(y,NORM_2,&norm));
195   PetscCall(VecAXPY(y,-1.0,yy));
196   PetscCall(VecNorm(y,NORM_1,&enorm));
197   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|y|_2 = %g\n",norm));
198   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|y-yy| = %g\n",enorm));
199 
200   /* compare z and zz: USFFT and FFTW results*/
201   PetscCall(VecNorm(z,NORM_2,&norm));
202   PetscCall(VecAXPY(z,-1.0,zz));
203   PetscCall(VecNorm(z,NORM_1,&enorm));
204   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|z|_2 = %g\n",norm));
205   PetscCall(PetscPrintf(PETSC_COMM_SELF, "|z-zz| = %g\n",enorm));
206 
207   /* free spaces */
208   PetscCall(DMRestoreGlobalVector(da,&x));
209   PetscCall(DMRestoreGlobalVector(da,&xx));
210   PetscCall(DMRestoreGlobalVector(da,&y));
211   PetscCall(DMRestoreGlobalVector(da,&yy));
212   PetscCall(DMRestoreGlobalVector(da,&z));
213   PetscCall(DMRestoreGlobalVector(da,&zz));
214 
215   PetscCall(PetscFinalize());
216   return 0;
217 }
218