xref: /petsc/src/mat/tests/ex228.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Test duplication/destruction of FFTW vecs \n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown  Compiling the code:
5c4762a1bSJed Brown    This code uses the FFTW interface.
6c4762a1bSJed Brown    Use one of the options below to configure:
7c4762a1bSJed Brown    --with-fftw-dir=/.... or --download-fftw
8c4762a1bSJed Brown  Usage:
9c4762a1bSJed Brown    mpiexec -np <np> ./ex228
10c4762a1bSJed Brown */
11c4762a1bSJed Brown 
12c4762a1bSJed Brown #include <petscmat.h>
13c4762a1bSJed Brown int main(int argc,char **args)
14c4762a1bSJed Brown {
15c4762a1bSJed Brown   Mat            A;         /* FFT Matrix */
16c4762a1bSJed Brown   Vec            x,y,z;     /* Work vectors */
17c4762a1bSJed Brown   Vec            x1,y1,z1;  /* Duplicate vectors */
18c4762a1bSJed Brown   PetscInt       i,k;       /* for iterating over dimensions */
19c4762a1bSJed Brown   PetscRandom    rdm;       /* for creating random input */
20c4762a1bSJed Brown   PetscScalar    a;         /* used to scale output */
21c4762a1bSJed Brown   PetscReal      enorm;     /* norm for sanity check */
22c4762a1bSJed Brown   PetscInt       n=10,N=1;  /* FFT dimension params */
23c4762a1bSJed Brown   PetscInt       DIM,dim[5];/* FFT params */
24c4762a1bSJed Brown 
25*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* To create random input vector */
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF, &rdm));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rdm));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   /* Iterate over dimensions, use PETSc-FFTW interface */
33c4762a1bSJed Brown   for (i=1; i<5; i++) {
34c4762a1bSJed Brown     DIM = i;
35c4762a1bSJed Brown     N = 1;
36c4762a1bSJed Brown     for (k=0; k<i; k++){dim[k] = n; N*=n;}
37c4762a1bSJed Brown 
385f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\n %" PetscInt_FMT " dimensions: FFTW on vector of size %" PetscInt_FMT " \n",DIM,N));
39c4762a1bSJed Brown 
40c4762a1bSJed Brown     /* create FFTW object */
415f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateFFT(PETSC_COMM_SELF,DIM,dim,MATFFTW,&A));
42c4762a1bSJed Brown     /* create vectors of length N */
435f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecsFFTW(A,&x,&y,&z));
44c4762a1bSJed Brown 
455f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) x, "Real space vector"));
465f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) y, "Frequency space vector"));
475f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) z, "Reconstructed vector"));
48c4762a1bSJed Brown 
49c4762a1bSJed Brown     /* Test vector duplication*/
505f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(x,&x1));
515f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(y,&y1));
525f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(z,&z1));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown     /* Set values of space vector x, copy to duplicate */
555f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetRandom(x,rdm));
565f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(x,x1));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown     /* Apply FFTW_FORWARD and FFTW_BACKWARD */
595f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(A,x,y));
605f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTranspose(A,y,z));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown     /* Apply FFTW_FORWARD and FFTW_BACKWARD for duplicate vecs */
635f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(A,x1,y1));
645f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTranspose(A,y1,z1));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown     /* Compare x and z1. FFTW computes an unnormalized DFT, thus z1 = N*x */
67c4762a1bSJed Brown     a    = 1.0/(PetscReal)N;
685f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(z1,a));
695f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(z1,-1.0,x));
705f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(z1,NORM_1,&enorm));
715f80ce2aSJacob Faibussowitsch     if (enorm > 1.e-9) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"  Error norm of |x - z1| %g\n",enorm));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown     /* free spaces */
745f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&x1));
755f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&y1));
765f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&z1));
775f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&x));
785f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&y));
795f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&z));
805f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&A));
81c4762a1bSJed Brown   }
82c4762a1bSJed Brown 
835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rdm));
84*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
85*b122ec5aSJacob Faibussowitsch   return 0;
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
88c4762a1bSJed Brown /*TEST
89c4762a1bSJed Brown 
90c4762a1bSJed Brown     build:
91c4762a1bSJed Brown       requires: fftw complex
92c4762a1bSJed Brown 
93c4762a1bSJed Brown     test:
94c4762a1bSJed Brown       suffix: 2
95c4762a1bSJed Brown       nsize : 4
96c4762a1bSJed Brown       args: -mat_fftw_plannerflags FFTW_ESTIMATE -n 16
97c4762a1bSJed Brown 
98c4762a1bSJed Brown     test:
99c4762a1bSJed Brown       suffix: 3
100c4762a1bSJed Brown       nsize : 2
101c4762a1bSJed Brown       args: -mat_fftw_plannerflags FFTW_MEASURE -n 12
102c4762a1bSJed Brown 
103c4762a1bSJed Brown     test:
104c4762a1bSJed Brown       suffix: 4
105c4762a1bSJed Brown       nsize : 2
106c4762a1bSJed Brown       args: -mat_fftw_plannerflags FFTW_PATIENT -n 10
107c4762a1bSJed Brown 
108c4762a1bSJed Brown     test:
109c4762a1bSJed Brown       suffix: 5
110c4762a1bSJed Brown       nsize : 1
111c4762a1bSJed Brown       args: -mat_fftw_plannerflags FFTW_EXHAUSTIVE -n 5
112c4762a1bSJed Brown 
113c4762a1bSJed Brown TEST*/
114