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