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