1c4762a1bSJed Brown static char help[] = "Test sequential r2c/c2r FFTW without PETSc interface \n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Compiling the code: 5c4762a1bSJed Brown This code uses the real numbers version of PETSc 6c4762a1bSJed Brown */ 7c4762a1bSJed Brown 8c4762a1bSJed Brown #include <petscmat.h> 9c4762a1bSJed Brown #include <fftw3.h> 10c4762a1bSJed Brown 11*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 12*d71ae5a4SJacob Faibussowitsch { 139371c9d4SSatish Balay typedef enum { 149371c9d4SSatish Balay RANDOM, 159371c9d4SSatish Balay CONSTANT, 169371c9d4SSatish Balay TANH, 179371c9d4SSatish Balay NUM_FUNCS 189371c9d4SSatish Balay } FuncType; 19c4762a1bSJed Brown const char *funcNames[NUM_FUNCS] = {"random", "constant", "tanh"}; 20c4762a1bSJed Brown PetscMPIInt size; 21c4762a1bSJed Brown int n = 10, N, Ny, ndim = 4, i, dim[4], DIM; 22c4762a1bSJed Brown Vec x, y, z; 23c4762a1bSJed Brown PetscScalar s; 24c4762a1bSJed Brown PetscRandom rdm; 25c4762a1bSJed Brown PetscReal enorm; 26c4762a1bSJed Brown PetscInt func = RANDOM; 27c4762a1bSJed Brown FuncType function = RANDOM; 28c4762a1bSJed Brown PetscBool view = PETSC_FALSE; 29c4762a1bSJed Brown PetscScalar *x_array, *y_array, *z_array; 30c4762a1bSJed Brown fftw_plan fplan, bplan; 31c4762a1bSJed Brown 32327415f7SBarry Smith PetscFunctionBeginUser; 339566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 34c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 35c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This example requires real numbers"); 36c4762a1bSJed Brown #endif 37c4762a1bSJed Brown 389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 39be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 40d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FFTW Options", "ex142"); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-function", "Function type", "ex142", funcNames, NUM_FUNCS, funcNames[function], &func, NULL)); 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-vec_view draw", "View the functions", "ex142", view, &view, NULL)); 43c4762a1bSJed Brown function = (FuncType)func; 44d0609cedSBarry Smith PetscOptionsEnd(); 45c4762a1bSJed Brown 469371c9d4SSatish Balay for (DIM = 0; DIM < ndim; DIM++) { dim[DIM] = n; /* size of real space vector in DIM-dimension */ } 479566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rdm)); 489566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 49c4762a1bSJed Brown 50c4762a1bSJed Brown for (DIM = 1; DIM < 5; DIM++) { 51c4762a1bSJed Brown /* create vectors of length N=dim[0]*dim[1]* ...*dim[DIM-1] */ 52c4762a1bSJed Brown /*----------------------------------------------------------*/ 53c4762a1bSJed Brown N = Ny = 1; 54ad540459SPierre Jolivet for (i = 0; i < DIM - 1; i++) N *= dim[i]; 559371c9d4SSatish Balay Ny = N; 569371c9d4SSatish Balay Ny *= 2 * (dim[DIM - 1] / 2 + 1); /* add padding elements to output vector y */ 57c4762a1bSJed Brown N *= dim[DIM - 1]; 58c4762a1bSJed Brown 599566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n %d-D: FFTW on vector of size %d \n", DIM, N)); 609566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &x)); 619566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x, "Real space vector")); 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, Ny, &y)); 649566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)y, "Frequency space vector")); 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &z)); 679566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)z, "Reconstructed vector")); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Set fftw plan */ 70c4762a1bSJed Brown /*----------------------------------*/ 719566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_array)); 729566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &y_array)); 739566063dSJacob Faibussowitsch PetscCall(VecGetArray(z, &z_array)); 74c4762a1bSJed Brown 75c4762a1bSJed Brown unsigned int flags = FFTW_ESTIMATE; /*or FFTW_MEASURE */ 76c4762a1bSJed Brown /* The data in the in/out arrays is overwritten during FFTW_MEASURE planning, so such planning 77c4762a1bSJed Brown should be done before the input is initialized by the user. */ 789566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "DIM: %d, N %d, Ny %d\n", DIM, N, Ny)); 79c4762a1bSJed Brown 80c4762a1bSJed Brown switch (DIM) { 81c4762a1bSJed Brown case 1: 82c4762a1bSJed Brown fplan = fftw_plan_dft_r2c_1d(dim[0], (double *)x_array, (fftw_complex *)y_array, flags); 83c4762a1bSJed Brown bplan = fftw_plan_dft_c2r_1d(dim[0], (fftw_complex *)y_array, (double *)z_array, flags); 84c4762a1bSJed Brown break; 85c4762a1bSJed Brown case 2: 86c4762a1bSJed Brown fplan = fftw_plan_dft_r2c_2d(dim[0], dim[1], (double *)x_array, (fftw_complex *)y_array, flags); 87c4762a1bSJed Brown bplan = fftw_plan_dft_c2r_2d(dim[0], dim[1], (fftw_complex *)y_array, (double *)z_array, flags); 88c4762a1bSJed Brown break; 89c4762a1bSJed Brown case 3: 90c4762a1bSJed Brown fplan = fftw_plan_dft_r2c_3d(dim[0], dim[1], dim[2], (double *)x_array, (fftw_complex *)y_array, flags); 91c4762a1bSJed Brown bplan = fftw_plan_dft_c2r_3d(dim[0], dim[1], dim[2], (fftw_complex *)y_array, (double *)z_array, flags); 92c4762a1bSJed Brown break; 93c4762a1bSJed Brown default: 94c4762a1bSJed Brown fplan = fftw_plan_dft_r2c(DIM, (int *)dim, (double *)x_array, (fftw_complex *)y_array, flags); 95c4762a1bSJed Brown bplan = fftw_plan_dft_c2r(DIM, (int *)dim, (fftw_complex *)y_array, (double *)z_array, flags); 96c4762a1bSJed Brown break; 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 999566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_array)); 1009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &y_array)); 1019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(z, &z_array)); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Initialize Real space vector x: 104c4762a1bSJed Brown The data in the in/out arrays is overwritten during FFTW_MEASURE planning, so planning 105c4762a1bSJed Brown should be done before the input is initialized by the user. 106c4762a1bSJed Brown --------------------------------------------------------*/ 107c4762a1bSJed Brown if (function == RANDOM) { 1089566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x, rdm)); 109c4762a1bSJed Brown } else if (function == CONSTANT) { 1109566063dSJacob Faibussowitsch PetscCall(VecSet(x, 1.0)); 111c4762a1bSJed Brown } else if (function == TANH) { 1129566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_array)); 113ad540459SPierre Jolivet for (i = 0; i < N; ++i) x_array[i] = tanh((i - N / 2.0) * (10.0 / N)); 1149566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_array)); 115c4762a1bSJed Brown } 1161baa6e33SBarry Smith if (view) PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* FFT - also test repeated transformation */ 119c4762a1bSJed Brown /*-------------------------------------------*/ 1209566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_array)); 1219566063dSJacob Faibussowitsch PetscCall(VecGetArray(y, &y_array)); 1229566063dSJacob Faibussowitsch PetscCall(VecGetArray(z, &z_array)); 123c4762a1bSJed Brown for (i = 0; i < 4; i++) { 124c4762a1bSJed Brown /* FFTW_FORWARD */ 125c4762a1bSJed Brown fftw_execute(fplan); 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* FFTW_BACKWARD: destroys its input array 'y_array' even for out-of-place transforms! */ 128c4762a1bSJed Brown fftw_execute(bplan); 129c4762a1bSJed Brown } 1309566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_array)); 1319566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y, &y_array)); 1329566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(z, &z_array)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* Compare x and z. FFTW computes an unnormalized DFT, thus z = N*x */ 135c4762a1bSJed Brown /*------------------------------------------------------------------*/ 136c4762a1bSJed Brown s = 1.0 / (PetscReal)N; 1379566063dSJacob Faibussowitsch PetscCall(VecScale(z, s)); 1389566063dSJacob Faibussowitsch if (view) PetscCall(VecView(x, PETSC_VIEWER_DRAW_WORLD)); 1399566063dSJacob Faibussowitsch if (view) PetscCall(VecView(z, PETSC_VIEWER_DRAW_WORLD)); 1409566063dSJacob Faibussowitsch PetscCall(VecAXPY(z, -1.0, x)); 1419566063dSJacob Faibussowitsch PetscCall(VecNorm(z, NORM_1, &enorm)); 14248a46eb9SPierre Jolivet if (enorm > 1.e-11) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Error norm of |x - z| %g\n", (double)enorm)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* free spaces */ 145c4762a1bSJed Brown fftw_destroy_plan(fplan); 146c4762a1bSJed Brown fftw_destroy_plan(bplan); 1479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&z)); 150c4762a1bSJed Brown } 1519566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 1529566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 153b122ec5aSJacob Faibussowitsch return 0; 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156c4762a1bSJed Brown /*TEST 157c4762a1bSJed Brown 158c4762a1bSJed Brown build: 159c4762a1bSJed Brown requires: fftw !complex 160c4762a1bSJed Brown 161c4762a1bSJed Brown test: 162c4762a1bSJed Brown output_file: output/ex142.out 163c4762a1bSJed Brown 164c4762a1bSJed Brown TEST*/ 165