1*c4762a1bSJed Brown static char help[] = "Test sequential FFTW convolution\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown /* 4*c4762a1bSJed Brown Compiling the code: 5*c4762a1bSJed Brown This code uses the complex numbers, so configure must be given --with-scalar-type=complex to enable this 6*c4762a1bSJed Brown */ 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown #include <petscmat.h> 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown PetscInt main(PetscInt argc,char **args) 11*c4762a1bSJed Brown { 12*c4762a1bSJed Brown typedef enum {RANDOM, CONSTANT, TANH, NUM_FUNCS} FuncType; 13*c4762a1bSJed Brown const char *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"}; 14*c4762a1bSJed Brown Mat A; 15*c4762a1bSJed Brown PetscMPIInt size; 16*c4762a1bSJed Brown PetscInt n = 10,N,ndim=4,dim[4],DIM,i,j; 17*c4762a1bSJed Brown Vec w,x,y1,y2,z1,z2; 18*c4762a1bSJed Brown PetscScalar *a, *a2, *a3; 19*c4762a1bSJed Brown PetscScalar s; 20*c4762a1bSJed Brown PetscRandom rdm; 21*c4762a1bSJed Brown PetscReal enorm; 22*c4762a1bSJed Brown PetscInt func = 0; 23*c4762a1bSJed Brown FuncType function = RANDOM; 24*c4762a1bSJed Brown PetscBool view = PETSC_FALSE; 25*c4762a1bSJed Brown PetscErrorCode ierr; 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 28*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr); 29*c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP, "This is a uniprocessor example only!"); 30*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex112");CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = PetscOptionsEList("-function", "Function type", "ex121", funcNames, NUM_FUNCS, funcNames[function], &func, NULL);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = PetscOptionsBool("-vec_view draw", "View the functions", "ex112", view, &view, NULL);CHKERRQ(ierr); 33*c4762a1bSJed Brown function = (FuncType) func; 34*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown for (DIM = 0; DIM < ndim; DIM++) { 37*c4762a1bSJed Brown dim[DIM] = n; /* size of transformation in DIM-dimension */ 38*c4762a1bSJed Brown } 39*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_SELF, &rdm);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 41*c4762a1bSJed Brown 42*c4762a1bSJed Brown for (DIM = 1; DIM < 5; DIM++) { 43*c4762a1bSJed Brown /* create vectors of length N=n^DIM */ 44*c4762a1bSJed Brown for (i = 0, N = 1; i < DIM; i++) N *= dim[i]; 45*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "\n %d-D: FFTW on vector of size %d \n",DIM,N);CHKERRQ(ierr); 46*c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,N,&x);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) x, "Real space vector");CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = VecDuplicate(x,&w);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) w, "Window vector");CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = VecDuplicate(x,&y1);CHKERRQ(ierr); 51*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) y1, "Frequency space vector");CHKERRQ(ierr); 52*c4762a1bSJed Brown ierr = VecDuplicate(x,&y2);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) y2, "Frequency space window vector");CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = VecDuplicate(x,&z1);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) z1, "Reconstructed convolution");CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = VecDuplicate(x,&z2);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) z2, "Real space convolution");CHKERRQ(ierr); 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown if (function == RANDOM) { 60*c4762a1bSJed Brown ierr = VecSetRandom(x, rdm);CHKERRQ(ierr); 61*c4762a1bSJed Brown } else if (function == CONSTANT) { 62*c4762a1bSJed Brown ierr = VecSet(x, 1.0);CHKERRQ(ierr); 63*c4762a1bSJed Brown } else if (function == TANH) { 64*c4762a1bSJed Brown ierr = VecGetArray(x, &a);CHKERRQ(ierr); 65*c4762a1bSJed Brown for (i = 0; i < N; ++i) { 66*c4762a1bSJed Brown a[i] = tanh((i - N/2.0)*(10.0/N)); 67*c4762a1bSJed Brown } 68*c4762a1bSJed Brown ierr = VecRestoreArray(x, &a);CHKERRQ(ierr); 69*c4762a1bSJed Brown } 70*c4762a1bSJed Brown if (view) {ierr = VecView(x, PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);} 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown /* Create window function */ 73*c4762a1bSJed Brown ierr = VecGetArray(w, &a);CHKERRQ(ierr); 74*c4762a1bSJed Brown for (i = 0; i < N; ++i) { 75*c4762a1bSJed Brown /* Step Function */ 76*c4762a1bSJed Brown a[i] = (i > N/4 && i < 3*N/4) ? 1.0 : 0.0; 77*c4762a1bSJed Brown /* Delta Function */ 78*c4762a1bSJed Brown /*a[i] = (i == N/2)? 1.0: 0.0; */ 79*c4762a1bSJed Brown } 80*c4762a1bSJed Brown ierr = VecRestoreArray(w, &a);CHKERRQ(ierr); 81*c4762a1bSJed Brown if (view) {ierr = VecView(w, PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);} 82*c4762a1bSJed Brown 83*c4762a1bSJed Brown /* create FFTW object */ 84*c4762a1bSJed Brown ierr = MatCreateFFT(PETSC_COMM_SELF,DIM,dim,MATFFTW,&A);CHKERRQ(ierr); 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown /* Convolve x with w*/ 87*c4762a1bSJed Brown ierr = MatMult(A,x,y1);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = MatMult(A,w,y2);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = VecPointwiseMult(y1, y1, y2);CHKERRQ(ierr); 90*c4762a1bSJed Brown if (view && i == 0) {ierr = VecView(y1, PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);} 91*c4762a1bSJed Brown ierr = MatMultTranspose(A,y1,z1);CHKERRQ(ierr); 92*c4762a1bSJed Brown 93*c4762a1bSJed Brown /* Compute the real space convolution */ 94*c4762a1bSJed Brown ierr = VecGetArray(x, &a);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = VecGetArray(w, &a2);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = VecGetArray(z2, &a3);CHKERRQ(ierr); 97*c4762a1bSJed Brown for (i = 0; i < N; ++i) { 98*c4762a1bSJed Brown /* PetscInt checkInd = (i > N/2-1)? i-N/2: i+N/2;*/ 99*c4762a1bSJed Brown 100*c4762a1bSJed Brown /*if (!(i%100)) PetscPrintf(PETSC_COMM_WORLD, "Finished convolution row %d\n", i);*/ 101*c4762a1bSJed Brown a3[i] = 0.0; 102*c4762a1bSJed Brown for (j = -N/2+1; j < N/2; ++j) { 103*c4762a1bSJed Brown PetscInt xpInd = (j < 0) ? N+j : j; 104*c4762a1bSJed Brown PetscInt diffInd = (i-j < 0) ? N-(j-i) : (i-j > N-1) ? i-j-N : i-j; 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown a3[i] += a[xpInd]*a2[diffInd]; 107*c4762a1bSJed Brown } 108*c4762a1bSJed Brown /*if (PetscAbsScalar(a3[i]) > PetscAbsScalar(a[checkInd])+0.1) PetscPrintf(PETSC_COMM_WORLD, "Invalid convolution at row %d\n", i);*/ 109*c4762a1bSJed Brown } 110*c4762a1bSJed Brown ierr = VecRestoreArray(x, &a);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = VecRestoreArray(w, &a2);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = VecRestoreArray(z2, &a3);CHKERRQ(ierr); 113*c4762a1bSJed Brown 114*c4762a1bSJed Brown /* compare z1 and z2. FFTW computes an unnormalized DFT, thus z1 = N*z2 */ 115*c4762a1bSJed Brown s = 1.0/(PetscReal)N; 116*c4762a1bSJed Brown ierr = VecScale(z1,s);CHKERRQ(ierr); 117*c4762a1bSJed Brown if (view) {ierr = VecView(z1, PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);} 118*c4762a1bSJed Brown if (view) {ierr = VecView(z2, PETSC_VIEWER_DRAW_WORLD);CHKERRQ(ierr);} 119*c4762a1bSJed Brown ierr = VecAXPY(z1,-1.0,z2);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = VecNorm(z1,NORM_1,&enorm);CHKERRQ(ierr); 121*c4762a1bSJed Brown if (enorm > 1.e-11) { 122*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," Error norm of |z1 - z2| %g\n",(double)enorm);CHKERRQ(ierr); 123*c4762a1bSJed Brown } 124*c4762a1bSJed Brown 125*c4762a1bSJed Brown /* free spaces */ 126*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 127*c4762a1bSJed Brown ierr = VecDestroy(&y1);CHKERRQ(ierr); 128*c4762a1bSJed Brown ierr = VecDestroy(&y2);CHKERRQ(ierr); 129*c4762a1bSJed Brown ierr = VecDestroy(&z1);CHKERRQ(ierr); 130*c4762a1bSJed Brown ierr = VecDestroy(&z2);CHKERRQ(ierr); 131*c4762a1bSJed Brown ierr = VecDestroy(&w);CHKERRQ(ierr); 132*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 133*c4762a1bSJed Brown } 134*c4762a1bSJed Brown ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 135*c4762a1bSJed Brown ierr = PetscFinalize(); 136*c4762a1bSJed Brown return ierr; 137*c4762a1bSJed Brown } 138*c4762a1bSJed Brown 139*c4762a1bSJed Brown 140*c4762a1bSJed Brown /*TEST 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown build: 143*c4762a1bSJed Brown requires: fftw complex 144*c4762a1bSJed Brown 145*c4762a1bSJed Brown test: 146*c4762a1bSJed Brown output_file: output/ex121.out 147*c4762a1bSJed Brown TODO: Example or FFTW interface is broken 148*c4762a1bSJed Brown 149*c4762a1bSJed Brown TEST*/ 150