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