1c4762a1bSJed Brown static char help[] = "Poiseuille flow problem. Viscous, laminar flow in a 2D channel with parabolic velocity\n\ 2c4762a1bSJed Brown profile and linear pressure drop, exact solution of the 2D Stokes\n"; 3c4762a1bSJed Brown 454703344SBarry Smith /* 554703344SBarry Smith M A R I T I M E R E S E A R C H I N S T I T U T E N E T H E R L A N D S 654703344SBarry Smith author : Christiaan M. Klaij 754703344SBarry Smith 854703344SBarry Smith Poiseuille flow problem. 954703344SBarry Smith 1054703344SBarry Smith Viscous, laminar flow in a 2D channel with parabolic velocity 1154703344SBarry Smith profile and linear pressure drop, exact solution of the 2D Stokes 1254703344SBarry Smith equations. 1354703344SBarry Smith 1454703344SBarry Smith Discretized with the cell-centered finite-volume method on a 1554703344SBarry Smith Cartesian grid with co-located variables. Variables ordered as 1654703344SBarry Smith [u1...uN v1...vN p1...pN]^T. Matrix [A00 A01; A10, A11] solved with 1754703344SBarry Smith PCFIELDSPLIT. Lower factorization is used to mimic the Semi-Implicit 1854703344SBarry Smith Method for Pressure Linked Equations (SIMPLE) used as preconditioner 1954703344SBarry Smith instead of solver. 2054703344SBarry Smith 2154703344SBarry Smith Disclaimer: does not contain the pressure-weighed interpolation 2254703344SBarry Smith method needed to suppress spurious pressure modes in real-life 2354703344SBarry Smith problems. 2454703344SBarry Smith 2554703344SBarry Smith Usage: 2654703344SBarry Smith mpiexec -n 2 ./ex70 -nx 32 -ny 48 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_1_pc_type none 2754703344SBarry Smith 2854703344SBarry Smith Runs with PCFIELDSPLIT on 32x48 grid, no PC for the Schur 2954703344SBarry Smith complement because A11 is zero. FGMRES is needed because 3054703344SBarry Smith PCFIELDSPLIT is a variable preconditioner. 3154703344SBarry Smith 3254703344SBarry Smith mpiexec -n 2 ./ex70 -nx 32 -ny 48 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc 3354703344SBarry Smith 3454703344SBarry Smith Same as above but with user defined PC for the true Schur 3554703344SBarry Smith complement. PC based on the SIMPLE-type approximation (inverse of 3654703344SBarry Smith A00 approximated by inverse of its diagonal). 3754703344SBarry Smith 3854703344SBarry Smith mpiexec -n 2 ./ex70 -nx 32 -ny 48 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_ksp 3954703344SBarry Smith 4054703344SBarry Smith Replace the true Schur complement with a user defined Schur 4154703344SBarry Smith complement based on the SIMPLE-type approximation. Same matrix is 4254703344SBarry Smith used as PC. 4354703344SBarry Smith 4454703344SBarry Smith mpiexec -n 2 ./ex70 -nx 32 -ny 48 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_0_ksp_type gmres -fieldsplit_0_pc_type bjacobi -fieldsplit_1_pc_type jacobi -fieldsplit_1_inner_ksp_type preonly -fieldsplit_1_inner_pc_type jacobi -fieldsplit_1_upper_ksp_type preonly -fieldsplit_1_upper_pc_type jacobi 4554703344SBarry Smith 4654703344SBarry Smith Out-of-the-box SIMPLE-type preconditioning. The major advantage 4754703344SBarry Smith is that the user neither needs to provide the approximation of 4854703344SBarry Smith the Schur complement, nor the corresponding preconditioner. 4954703344SBarry Smith */ 50c4762a1bSJed Brown 51c4762a1bSJed Brown #include <petscksp.h> 52c4762a1bSJed Brown 53c4762a1bSJed Brown typedef struct { 54c4762a1bSJed Brown PetscBool userPC, userKSP, matsymmetric; /* user defined preconditioner and matrix for the Schur complement */ 55c4762a1bSJed Brown PetscInt nx, ny; /* nb of cells in x- and y-direction */ 56c4762a1bSJed Brown PetscReal hx, hy; /* mesh size in x- and y-direction */ 57c4762a1bSJed Brown Mat A; /* block matrix */ 58c4762a1bSJed Brown Mat subA[4]; /* the four blocks */ 59c4762a1bSJed Brown Mat myS; /* the approximation of the Schur complement */ 60c4762a1bSJed Brown Vec x, b, y; /* solution, rhs and temporary vector */ 61c4762a1bSJed Brown IS isg[2]; /* index sets of split "0" and "1" */ 62c4762a1bSJed Brown } Stokes; 63c4762a1bSJed Brown 64c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock00(Stokes *); /* setup the block Q */ 65c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock01(Stokes *); /* setup the block G */ 66c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock10(Stokes *); /* setup the block D (equal to the transpose of G) */ 67c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock11(Stokes *); /* setup the block C (equal to zero) */ 68c4762a1bSJed Brown 69c4762a1bSJed Brown PetscErrorCode StokesGetPosition(Stokes *, PetscInt, PetscInt *, PetscInt *); /* row number j*nx+i corresponds to position (i,j) in grid */ 70c4762a1bSJed Brown 71c4762a1bSJed Brown PetscErrorCode StokesStencilLaplacian(Stokes *, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscScalar *); /* stencil of the Laplacian operator */ 72c4762a1bSJed Brown PetscErrorCode StokesStencilGradientX(Stokes *, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscScalar *); /* stencil of the Gradient operator (x-component) */ 73c4762a1bSJed Brown PetscErrorCode StokesStencilGradientY(Stokes *, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscScalar *); /* stencil of the Gradient operator (y-component) */ 74c4762a1bSJed Brown 75c4762a1bSJed Brown PetscErrorCode StokesRhs(Stokes *); /* rhs vector */ 76c4762a1bSJed Brown PetscErrorCode StokesRhsMomX(Stokes *, PetscInt, PetscInt, PetscScalar *); /* right hand side of velocity (x-component) */ 77c4762a1bSJed Brown PetscErrorCode StokesRhsMomY(Stokes *, PetscInt, PetscInt, PetscScalar *); /* right hand side of velocity (y-component) */ 78c4762a1bSJed Brown PetscErrorCode StokesRhsMass(Stokes *, PetscInt, PetscInt, PetscScalar *); /* right hand side of pressure */ 79c4762a1bSJed Brown 80c4762a1bSJed Brown PetscErrorCode StokesSetupApproxSchur(Stokes *); /* approximation of the Schur complement */ 81c4762a1bSJed Brown 82c4762a1bSJed Brown PetscErrorCode StokesExactSolution(Stokes *); /* exact solution vector */ 83c4762a1bSJed Brown PetscErrorCode StokesWriteSolution(Stokes *); /* write solution to file */ 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* exact solution for the velocity (x-component, y-component is zero) */ 869371c9d4SSatish Balay PetscScalar StokesExactVelocityX(const PetscScalar y) { 87c4762a1bSJed Brown return 4.0 * y * (1.0 - y); 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* exact solution for the pressure */ 919371c9d4SSatish Balay PetscScalar StokesExactPressure(const PetscScalar x) { 92c4762a1bSJed Brown return 8.0 * (2.0 - x); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown 959371c9d4SSatish Balay PetscErrorCode StokesSetupPC(Stokes *s, KSP ksp) { 96c4762a1bSJed Brown KSP *subksp; 97c4762a1bSJed Brown PC pc; 98c4762a1bSJed Brown PetscInt n = 1; 99c4762a1bSJed Brown 100c4762a1bSJed Brown PetscFunctionBeginUser; 1019566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 1029566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, "0", s->isg[0])); 1039566063dSJacob Faibussowitsch PetscCall(PCFieldSplitSetIS(pc, "1", s->isg[1])); 1041baa6e33SBarry Smith if (s->userPC) PetscCall(PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_USER, s->myS)); 105c4762a1bSJed Brown if (s->userKSP) { 1069566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 1079566063dSJacob Faibussowitsch PetscCall(PCFieldSplitGetSubKSP(pc, &n, &subksp)); 1089566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(subksp[1], s->myS, s->myS)); 1099566063dSJacob Faibussowitsch PetscCall(PetscFree(subksp)); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown PetscFunctionReturn(0); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 1149371c9d4SSatish Balay PetscErrorCode StokesWriteSolution(Stokes *s) { 115c4762a1bSJed Brown PetscMPIInt size; 116c4762a1bSJed Brown PetscInt n, i, j; 117c4762a1bSJed Brown const PetscScalar *array; 118c4762a1bSJed Brown 119c4762a1bSJed Brown PetscFunctionBeginUser; 120c4762a1bSJed Brown /* write data (*warning* only works sequential) */ 1219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size)); 122c4762a1bSJed Brown if (size == 1) { 123c4762a1bSJed Brown PetscViewer viewer; 1249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(s->x, &array)); 1259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "solution.dat", &viewer)); 1269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "# x, y, u, v, p\n")); 127c4762a1bSJed Brown for (j = 0; j < s->ny; j++) { 128c4762a1bSJed Brown for (i = 0; i < s->nx; i++) { 129c4762a1bSJed Brown n = j * s->nx + i; 1309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%.12g %.12g %.12g %.12g %.12g\n", (double)(i * s->hx + s->hx / 2), (double)(j * s->hy + s->hy / 2), (double)PetscRealPart(array[n]), (double)PetscRealPart(array[n + s->nx * s->ny]), (double)PetscRealPart(array[n + 2 * s->nx * s->ny]))); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown } 1339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(s->x, &array)); 1349566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown PetscFunctionReturn(0); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 1399371c9d4SSatish Balay PetscErrorCode StokesSetupIndexSets(Stokes *s) { 140c4762a1bSJed Brown PetscFunctionBeginUser; 141c4762a1bSJed Brown /* the two index sets */ 1429566063dSJacob Faibussowitsch PetscCall(MatNestGetISs(s->A, s->isg, NULL)); 143c4762a1bSJed Brown PetscFunctionReturn(0); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 1469371c9d4SSatish Balay PetscErrorCode StokesSetupVectors(Stokes *s) { 147c4762a1bSJed Brown PetscFunctionBeginUser; 148c4762a1bSJed Brown /* solution vector x */ 1499566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &s->x)); 1509566063dSJacob Faibussowitsch PetscCall(VecSetSizes(s->x, PETSC_DECIDE, 3 * s->nx * s->ny)); 1519566063dSJacob Faibussowitsch PetscCall(VecSetType(s->x, VECMPI)); 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* exact solution y */ 1549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(s->x, &s->y)); 1559566063dSJacob Faibussowitsch PetscCall(StokesExactSolution(s)); 156c4762a1bSJed Brown 157c4762a1bSJed Brown /* rhs vector b */ 1589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(s->x, &s->b)); 1599566063dSJacob Faibussowitsch PetscCall(StokesRhs(s)); 160c4762a1bSJed Brown PetscFunctionReturn(0); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 1639371c9d4SSatish Balay PetscErrorCode StokesGetPosition(Stokes *s, PetscInt row, PetscInt *i, PetscInt *j) { 164c4762a1bSJed Brown PetscInt n; 165c4762a1bSJed Brown 166c4762a1bSJed Brown PetscFunctionBeginUser; 167c4762a1bSJed Brown /* cell number n=j*nx+i has position (i,j) in grid */ 168c4762a1bSJed Brown n = row % (s->nx * s->ny); 169c4762a1bSJed Brown *i = n % s->nx; 170c4762a1bSJed Brown *j = (n - (*i)) / s->nx; 171c4762a1bSJed Brown PetscFunctionReturn(0); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown 1749371c9d4SSatish Balay PetscErrorCode StokesExactSolution(Stokes *s) { 175c4762a1bSJed Brown PetscInt row, start, end, i, j; 176c4762a1bSJed Brown PetscScalar val; 177c4762a1bSJed Brown Vec y0, y1; 178c4762a1bSJed Brown 179c4762a1bSJed Brown PetscFunctionBeginUser; 180c4762a1bSJed Brown /* velocity part */ 1819566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(s->y, s->isg[0], &y0)); 1829566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(y0, &start, &end)); 183c4762a1bSJed Brown for (row = start; row < end; row++) { 1849566063dSJacob Faibussowitsch PetscCall(StokesGetPosition(s, row, &i, &j)); 185c4762a1bSJed Brown if (row < s->nx * s->ny) { 186c4762a1bSJed Brown val = StokesExactVelocityX(j * s->hy + s->hy / 2); 187c4762a1bSJed Brown } else { 188c4762a1bSJed Brown val = 0; 189c4762a1bSJed Brown } 1909566063dSJacob Faibussowitsch PetscCall(VecSetValue(y0, row, val, INSERT_VALUES)); 191c4762a1bSJed Brown } 1929566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(s->y, s->isg[0], &y0)); 193c4762a1bSJed Brown 194c4762a1bSJed Brown /* pressure part */ 1959566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(s->y, s->isg[1], &y1)); 1969566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(y1, &start, &end)); 197c4762a1bSJed Brown for (row = start; row < end; row++) { 1989566063dSJacob Faibussowitsch PetscCall(StokesGetPosition(s, row, &i, &j)); 199c4762a1bSJed Brown val = StokesExactPressure(i * s->hx + s->hx / 2); 2009566063dSJacob Faibussowitsch PetscCall(VecSetValue(y1, row, val, INSERT_VALUES)); 201c4762a1bSJed Brown } 2029566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(s->y, s->isg[1], &y1)); 203c4762a1bSJed Brown PetscFunctionReturn(0); 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 2069371c9d4SSatish Balay PetscErrorCode StokesRhs(Stokes *s) { 207c4762a1bSJed Brown PetscInt row, start, end, i, j; 208c4762a1bSJed Brown PetscScalar val; 209c4762a1bSJed Brown Vec b0, b1; 210c4762a1bSJed Brown 211c4762a1bSJed Brown PetscFunctionBeginUser; 212c4762a1bSJed Brown /* velocity part */ 2139566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(s->b, s->isg[0], &b0)); 2149566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(b0, &start, &end)); 215c4762a1bSJed Brown for (row = start; row < end; row++) { 2169566063dSJacob Faibussowitsch PetscCall(StokesGetPosition(s, row, &i, &j)); 217c4762a1bSJed Brown if (row < s->nx * s->ny) { 2189566063dSJacob Faibussowitsch PetscCall(StokesRhsMomX(s, i, j, &val)); 219c4762a1bSJed Brown } else { 2209566063dSJacob Faibussowitsch PetscCall(StokesRhsMomY(s, i, j, &val)); 221c4762a1bSJed Brown } 2229566063dSJacob Faibussowitsch PetscCall(VecSetValue(b0, row, val, INSERT_VALUES)); 223c4762a1bSJed Brown } 2249566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(s->b, s->isg[0], &b0)); 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* pressure part */ 2279566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(s->b, s->isg[1], &b1)); 2289566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(b1, &start, &end)); 229c4762a1bSJed Brown for (row = start; row < end; row++) { 2309566063dSJacob Faibussowitsch PetscCall(StokesGetPosition(s, row, &i, &j)); 2319566063dSJacob Faibussowitsch PetscCall(StokesRhsMass(s, i, j, &val)); 2329371c9d4SSatish Balay if (s->matsymmetric) { val = -1.0 * val; } 2339566063dSJacob Faibussowitsch PetscCall(VecSetValue(b1, row, val, INSERT_VALUES)); 234c4762a1bSJed Brown } 2359566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(s->b, s->isg[1], &b1)); 236c4762a1bSJed Brown PetscFunctionReturn(0); 237c4762a1bSJed Brown } 238c4762a1bSJed Brown 2399371c9d4SSatish Balay PetscErrorCode StokesSetupMatBlock00(Stokes *s) { 240c4762a1bSJed Brown PetscInt row, start, end, sz, i, j; 241c4762a1bSJed Brown PetscInt cols[5]; 242c4762a1bSJed Brown PetscScalar vals[5]; 243c4762a1bSJed Brown 244c4762a1bSJed Brown PetscFunctionBeginUser; 245c4762a1bSJed Brown /* A[0] is 2N-by-2N */ 2469566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &s->subA[0])); 2479566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(s->subA[0], "a00_")); 2489566063dSJacob Faibussowitsch PetscCall(MatSetSizes(s->subA[0], PETSC_DECIDE, PETSC_DECIDE, 2 * s->nx * s->ny, 2 * s->nx * s->ny)); 2499566063dSJacob Faibussowitsch PetscCall(MatSetType(s->subA[0], MATMPIAIJ)); 2509566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(s->subA[0], 5, NULL, 5, NULL)); 2519566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(s->subA[0], &start, &end)); 252c4762a1bSJed Brown 253c4762a1bSJed Brown for (row = start; row < end; row++) { 2549566063dSJacob Faibussowitsch PetscCall(StokesGetPosition(s, row, &i, &j)); 255c4762a1bSJed Brown /* first part: rows 0 to (nx*ny-1) */ 2569566063dSJacob Faibussowitsch PetscCall(StokesStencilLaplacian(s, i, j, &sz, cols, vals)); 257c4762a1bSJed Brown /* second part: rows (nx*ny) to (2*nx*ny-1) */ 258c4762a1bSJed Brown if (row >= s->nx * s->ny) { 259c4762a1bSJed Brown for (i = 0; i < sz; i++) cols[i] += s->nx * s->ny; 260c4762a1bSJed Brown } 261c4762a1bSJed Brown for (i = 0; i < sz; i++) vals[i] = -1.0 * vals[i]; /* dynamic viscosity coef mu=-1 */ 2629566063dSJacob Faibussowitsch PetscCall(MatSetValues(s->subA[0], 1, &row, sz, cols, vals, INSERT_VALUES)); 263c4762a1bSJed Brown } 2649566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(s->subA[0], MAT_FINAL_ASSEMBLY)); 2659566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(s->subA[0], MAT_FINAL_ASSEMBLY)); 266c4762a1bSJed Brown PetscFunctionReturn(0); 267c4762a1bSJed Brown } 268c4762a1bSJed Brown 2699371c9d4SSatish Balay PetscErrorCode StokesSetupMatBlock01(Stokes *s) { 270c4762a1bSJed Brown PetscInt row, start, end, sz, i, j; 271c4762a1bSJed Brown PetscInt cols[5]; 272c4762a1bSJed Brown PetscScalar vals[5]; 273c4762a1bSJed Brown 274c4762a1bSJed Brown PetscFunctionBeginUser; 275c4762a1bSJed Brown /* A[1] is 2N-by-N */ 2769566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &s->subA[1])); 2779566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(s->subA[1], "a01_")); 2789566063dSJacob Faibussowitsch PetscCall(MatSetSizes(s->subA[1], PETSC_DECIDE, PETSC_DECIDE, 2 * s->nx * s->ny, s->nx * s->ny)); 2799566063dSJacob Faibussowitsch PetscCall(MatSetType(s->subA[1], MATMPIAIJ)); 2809566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(s->subA[1], 5, NULL, 5, NULL)); 2819566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(s->subA[1], &start, &end)); 282c4762a1bSJed Brown 2839566063dSJacob Faibussowitsch PetscCall(MatSetOption(s->subA[1], MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); 284c4762a1bSJed Brown 285c4762a1bSJed Brown for (row = start; row < end; row++) { 2869566063dSJacob Faibussowitsch PetscCall(StokesGetPosition(s, row, &i, &j)); 287c4762a1bSJed Brown /* first part: rows 0 to (nx*ny-1) */ 288c4762a1bSJed Brown if (row < s->nx * s->ny) { 2899566063dSJacob Faibussowitsch PetscCall(StokesStencilGradientX(s, i, j, &sz, cols, vals)); 290c4762a1bSJed Brown } else { /* second part: rows (nx*ny) to (2*nx*ny-1) */ 2919566063dSJacob Faibussowitsch PetscCall(StokesStencilGradientY(s, i, j, &sz, cols, vals)); 292c4762a1bSJed Brown } 2939566063dSJacob Faibussowitsch PetscCall(MatSetValues(s->subA[1], 1, &row, sz, cols, vals, INSERT_VALUES)); 294c4762a1bSJed Brown } 2959566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(s->subA[1], MAT_FINAL_ASSEMBLY)); 2969566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(s->subA[1], MAT_FINAL_ASSEMBLY)); 297c4762a1bSJed Brown PetscFunctionReturn(0); 298c4762a1bSJed Brown } 299c4762a1bSJed Brown 3009371c9d4SSatish Balay PetscErrorCode StokesSetupMatBlock10(Stokes *s) { 301c4762a1bSJed Brown PetscFunctionBeginUser; 302c4762a1bSJed Brown /* A[2] is minus transpose of A[1] */ 3039566063dSJacob Faibussowitsch PetscCall(MatTranspose(s->subA[1], MAT_INITIAL_MATRIX, &s->subA[2])); 304*48a46eb9SPierre Jolivet if (!s->matsymmetric) PetscCall(MatScale(s->subA[2], -1.0)); 3059566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(s->subA[2], "a10_")); 306c4762a1bSJed Brown PetscFunctionReturn(0); 307c4762a1bSJed Brown } 308c4762a1bSJed Brown 3099371c9d4SSatish Balay PetscErrorCode StokesSetupMatBlock11(Stokes *s) { 310c4762a1bSJed Brown PetscFunctionBeginUser; 311c4762a1bSJed Brown /* A[3] is N-by-N null matrix */ 3129566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &s->subA[3])); 3139566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(s->subA[3], "a11_")); 3149566063dSJacob Faibussowitsch PetscCall(MatSetSizes(s->subA[3], PETSC_DECIDE, PETSC_DECIDE, s->nx * s->ny, s->nx * s->ny)); 3159566063dSJacob Faibussowitsch PetscCall(MatSetType(s->subA[3], MATMPIAIJ)); 3169566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(s->subA[3], 0, NULL, 0, NULL)); 3179566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(s->subA[3], MAT_FINAL_ASSEMBLY)); 3189566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(s->subA[3], MAT_FINAL_ASSEMBLY)); 319c4762a1bSJed Brown PetscFunctionReturn(0); 320c4762a1bSJed Brown } 321c4762a1bSJed Brown 3229371c9d4SSatish Balay PetscErrorCode StokesSetupApproxSchur(Stokes *s) { 323c4762a1bSJed Brown Vec diag; 324c4762a1bSJed Brown 325c4762a1bSJed Brown PetscFunctionBeginUser; 32654703344SBarry Smith /* Schur complement approximation: myS = A11 - A10 inv(DIAGFORM(A00)) A01 */ 327c4762a1bSJed Brown /* note: A11 is zero */ 328c4762a1bSJed Brown /* note: in real life this matrix would be build directly, */ 329c4762a1bSJed Brown /* i.e. without MatMatMult */ 330c4762a1bSJed Brown 331c4762a1bSJed Brown /* inverse of diagonal of A00 */ 3329566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &diag)); 3339566063dSJacob Faibussowitsch PetscCall(VecSetSizes(diag, PETSC_DECIDE, 2 * s->nx * s->ny)); 3349566063dSJacob Faibussowitsch PetscCall(VecSetType(diag, VECMPI)); 3359566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(s->subA[0], diag)); 3369566063dSJacob Faibussowitsch PetscCall(VecReciprocal(diag)); 337c4762a1bSJed Brown 33854703344SBarry Smith /* compute: - A10 inv(DIAGFORM(A00)) A01 */ 3399566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(s->subA[1], diag, NULL)); /* (*warning* overwrites subA[1]) */ 3409566063dSJacob Faibussowitsch PetscCall(MatMatMult(s->subA[2], s->subA[1], MAT_INITIAL_MATRIX, PETSC_DEFAULT, &s->myS)); 3419566063dSJacob Faibussowitsch PetscCall(MatScale(s->myS, -1.0)); 342c4762a1bSJed Brown 343c4762a1bSJed Brown /* restore A10 */ 3449566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(s->subA[0], diag)); 3459566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(s->subA[1], diag, NULL)); 3469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&diag)); 347c4762a1bSJed Brown PetscFunctionReturn(0); 348c4762a1bSJed Brown } 349c4762a1bSJed Brown 3509371c9d4SSatish Balay PetscErrorCode StokesSetupMatrix(Stokes *s) { 351c4762a1bSJed Brown PetscFunctionBeginUser; 3529566063dSJacob Faibussowitsch PetscCall(StokesSetupMatBlock00(s)); 3539566063dSJacob Faibussowitsch PetscCall(StokesSetupMatBlock01(s)); 3549566063dSJacob Faibussowitsch PetscCall(StokesSetupMatBlock10(s)); 3559566063dSJacob Faibussowitsch PetscCall(StokesSetupMatBlock11(s)); 3569566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, s->subA, &s->A)); 3579566063dSJacob Faibussowitsch PetscCall(StokesSetupApproxSchur(s)); 358c4762a1bSJed Brown PetscFunctionReturn(0); 359c4762a1bSJed Brown } 360c4762a1bSJed Brown 3619371c9d4SSatish Balay PetscErrorCode StokesStencilLaplacian(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals) { 362c4762a1bSJed Brown PetscInt p = j * s->nx + i, w = p - 1, e = p + 1, s2 = p - s->nx, n = p + s->nx; 363c4762a1bSJed Brown PetscScalar ae = s->hy / s->hx, aeb = 0; 364c4762a1bSJed Brown PetscScalar aw = s->hy / s->hx, awb = s->hy / (s->hx / 2); 365c4762a1bSJed Brown PetscScalar as = s->hx / s->hy, asb = s->hx / (s->hy / 2); 366c4762a1bSJed Brown PetscScalar an = s->hx / s->hy, anb = s->hx / (s->hy / 2); 367c4762a1bSJed Brown 368c4762a1bSJed Brown PetscFunctionBeginUser; 369c4762a1bSJed Brown if (i == 0 && j == 0) { /* south-west corner */ 370c4762a1bSJed Brown *sz = 3; 3719371c9d4SSatish Balay cols[0] = p; 3729371c9d4SSatish Balay vals[0] = -(ae + awb + asb + an); 3739371c9d4SSatish Balay cols[1] = e; 3749371c9d4SSatish Balay vals[1] = ae; 3759371c9d4SSatish Balay cols[2] = n; 3769371c9d4SSatish Balay vals[2] = an; 377c4762a1bSJed Brown } else if (i == 0 && j == s->ny - 1) { /* north-west corner */ 378c4762a1bSJed Brown *sz = 3; 3799371c9d4SSatish Balay cols[0] = s2; 3809371c9d4SSatish Balay vals[0] = as; 3819371c9d4SSatish Balay cols[1] = p; 3829371c9d4SSatish Balay vals[1] = -(ae + awb + as + anb); 3839371c9d4SSatish Balay cols[2] = e; 3849371c9d4SSatish Balay vals[2] = ae; 385c4762a1bSJed Brown } else if (i == s->nx - 1 && j == 0) { /* south-east corner */ 386c4762a1bSJed Brown *sz = 3; 3879371c9d4SSatish Balay cols[0] = w; 3889371c9d4SSatish Balay vals[0] = aw; 3899371c9d4SSatish Balay cols[1] = p; 3909371c9d4SSatish Balay vals[1] = -(aeb + aw + asb + an); 3919371c9d4SSatish Balay cols[2] = n; 3929371c9d4SSatish Balay vals[2] = an; 393c4762a1bSJed Brown } else if (i == s->nx - 1 && j == s->ny - 1) { /* north-east corner */ 394c4762a1bSJed Brown *sz = 3; 3959371c9d4SSatish Balay cols[0] = s2; 3969371c9d4SSatish Balay vals[0] = as; 3979371c9d4SSatish Balay cols[1] = w; 3989371c9d4SSatish Balay vals[1] = aw; 3999371c9d4SSatish Balay cols[2] = p; 4009371c9d4SSatish Balay vals[2] = -(aeb + aw + as + anb); 401c4762a1bSJed Brown } else if (i == 0) { /* west boundary */ 402c4762a1bSJed Brown *sz = 4; 4039371c9d4SSatish Balay cols[0] = s2; 4049371c9d4SSatish Balay vals[0] = as; 4059371c9d4SSatish Balay cols[1] = p; 4069371c9d4SSatish Balay vals[1] = -(ae + awb + as + an); 4079371c9d4SSatish Balay cols[2] = e; 4089371c9d4SSatish Balay vals[2] = ae; 4099371c9d4SSatish Balay cols[3] = n; 4109371c9d4SSatish Balay vals[3] = an; 411c4762a1bSJed Brown } else if (i == s->nx - 1) { /* east boundary */ 412c4762a1bSJed Brown *sz = 4; 4139371c9d4SSatish Balay cols[0] = s2; 4149371c9d4SSatish Balay vals[0] = as; 4159371c9d4SSatish Balay cols[1] = w; 4169371c9d4SSatish Balay vals[1] = aw; 4179371c9d4SSatish Balay cols[2] = p; 4189371c9d4SSatish Balay vals[2] = -(aeb + aw + as + an); 4199371c9d4SSatish Balay cols[3] = n; 4209371c9d4SSatish Balay vals[3] = an; 421c4762a1bSJed Brown } else if (j == 0) { /* south boundary */ 422c4762a1bSJed Brown *sz = 4; 4239371c9d4SSatish Balay cols[0] = w; 4249371c9d4SSatish Balay vals[0] = aw; 4259371c9d4SSatish Balay cols[1] = p; 4269371c9d4SSatish Balay vals[1] = -(ae + aw + asb + an); 4279371c9d4SSatish Balay cols[2] = e; 4289371c9d4SSatish Balay vals[2] = ae; 4299371c9d4SSatish Balay cols[3] = n; 4309371c9d4SSatish Balay vals[3] = an; 431c4762a1bSJed Brown } else if (j == s->ny - 1) { /* north boundary */ 432c4762a1bSJed Brown *sz = 4; 4339371c9d4SSatish Balay cols[0] = s2; 4349371c9d4SSatish Balay vals[0] = as; 4359371c9d4SSatish Balay cols[1] = w; 4369371c9d4SSatish Balay vals[1] = aw; 4379371c9d4SSatish Balay cols[2] = p; 4389371c9d4SSatish Balay vals[2] = -(ae + aw + as + anb); 4399371c9d4SSatish Balay cols[3] = e; 4409371c9d4SSatish Balay vals[3] = ae; 441c4762a1bSJed Brown } else { /* interior */ 442c4762a1bSJed Brown *sz = 5; 4439371c9d4SSatish Balay cols[0] = s2; 4449371c9d4SSatish Balay vals[0] = as; 4459371c9d4SSatish Balay cols[1] = w; 4469371c9d4SSatish Balay vals[1] = aw; 4479371c9d4SSatish Balay cols[2] = p; 4489371c9d4SSatish Balay vals[2] = -(ae + aw + as + an); 4499371c9d4SSatish Balay cols[3] = e; 4509371c9d4SSatish Balay vals[3] = ae; 4519371c9d4SSatish Balay cols[4] = n; 4529371c9d4SSatish Balay vals[4] = an; 453c4762a1bSJed Brown } 454c4762a1bSJed Brown PetscFunctionReturn(0); 455c4762a1bSJed Brown } 456c4762a1bSJed Brown 4579371c9d4SSatish Balay PetscErrorCode StokesStencilGradientX(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals) { 458c4762a1bSJed Brown PetscInt p = j * s->nx + i, w = p - 1, e = p + 1; 459c4762a1bSJed Brown PetscScalar ae = s->hy / 2, aeb = s->hy; 460c4762a1bSJed Brown PetscScalar aw = -s->hy / 2, awb = 0; 461c4762a1bSJed Brown 462c4762a1bSJed Brown PetscFunctionBeginUser; 463c4762a1bSJed Brown if (i == 0 && j == 0) { /* south-west corner */ 464c4762a1bSJed Brown *sz = 2; 4659371c9d4SSatish Balay cols[0] = p; 4669371c9d4SSatish Balay vals[0] = -(ae + awb); 4679371c9d4SSatish Balay cols[1] = e; 4689371c9d4SSatish Balay vals[1] = ae; 469c4762a1bSJed Brown } else if (i == 0 && j == s->ny - 1) { /* north-west corner */ 470c4762a1bSJed Brown *sz = 2; 4719371c9d4SSatish Balay cols[0] = p; 4729371c9d4SSatish Balay vals[0] = -(ae + awb); 4739371c9d4SSatish Balay cols[1] = e; 4749371c9d4SSatish Balay vals[1] = ae; 475c4762a1bSJed Brown } else if (i == s->nx - 1 && j == 0) { /* south-east corner */ 476c4762a1bSJed Brown *sz = 2; 4779371c9d4SSatish Balay cols[0] = w; 4789371c9d4SSatish Balay vals[0] = aw; 4799371c9d4SSatish Balay cols[1] = p; 4809371c9d4SSatish Balay vals[1] = -(aeb + aw); 481c4762a1bSJed Brown } else if (i == s->nx - 1 && j == s->ny - 1) { /* north-east corner */ 482c4762a1bSJed Brown *sz = 2; 4839371c9d4SSatish Balay cols[0] = w; 4849371c9d4SSatish Balay vals[0] = aw; 4859371c9d4SSatish Balay cols[1] = p; 4869371c9d4SSatish Balay vals[1] = -(aeb + aw); 487c4762a1bSJed Brown } else if (i == 0) { /* west boundary */ 488c4762a1bSJed Brown *sz = 2; 4899371c9d4SSatish Balay cols[0] = p; 4909371c9d4SSatish Balay vals[0] = -(ae + awb); 4919371c9d4SSatish Balay cols[1] = e; 4929371c9d4SSatish Balay vals[1] = ae; 493c4762a1bSJed Brown } else if (i == s->nx - 1) { /* east boundary */ 494c4762a1bSJed Brown *sz = 2; 4959371c9d4SSatish Balay cols[0] = w; 4969371c9d4SSatish Balay vals[0] = aw; 4979371c9d4SSatish Balay cols[1] = p; 4989371c9d4SSatish Balay vals[1] = -(aeb + aw); 499c4762a1bSJed Brown } else if (j == 0) { /* south boundary */ 500c4762a1bSJed Brown *sz = 3; 5019371c9d4SSatish Balay cols[0] = w; 5029371c9d4SSatish Balay vals[0] = aw; 5039371c9d4SSatish Balay cols[1] = p; 5049371c9d4SSatish Balay vals[1] = -(ae + aw); 5059371c9d4SSatish Balay cols[2] = e; 5069371c9d4SSatish Balay vals[2] = ae; 507c4762a1bSJed Brown } else if (j == s->ny - 1) { /* north boundary */ 508c4762a1bSJed Brown *sz = 3; 5099371c9d4SSatish Balay cols[0] = w; 5109371c9d4SSatish Balay vals[0] = aw; 5119371c9d4SSatish Balay cols[1] = p; 5129371c9d4SSatish Balay vals[1] = -(ae + aw); 5139371c9d4SSatish Balay cols[2] = e; 5149371c9d4SSatish Balay vals[2] = ae; 515c4762a1bSJed Brown } else { /* interior */ 516c4762a1bSJed Brown *sz = 3; 5179371c9d4SSatish Balay cols[0] = w; 5189371c9d4SSatish Balay vals[0] = aw; 5199371c9d4SSatish Balay cols[1] = p; 5209371c9d4SSatish Balay vals[1] = -(ae + aw); 5219371c9d4SSatish Balay cols[2] = e; 5229371c9d4SSatish Balay vals[2] = ae; 523c4762a1bSJed Brown } 524c4762a1bSJed Brown PetscFunctionReturn(0); 525c4762a1bSJed Brown } 526c4762a1bSJed Brown 5279371c9d4SSatish Balay PetscErrorCode StokesStencilGradientY(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals) { 528c4762a1bSJed Brown PetscInt p = j * s->nx + i, s2 = p - s->nx, n = p + s->nx; 529c4762a1bSJed Brown PetscScalar as = -s->hx / 2, asb = 0; 530c4762a1bSJed Brown PetscScalar an = s->hx / 2, anb = 0; 531c4762a1bSJed Brown 532c4762a1bSJed Brown PetscFunctionBeginUser; 533c4762a1bSJed Brown if (i == 0 && j == 0) { /* south-west corner */ 534c4762a1bSJed Brown *sz = 2; 5359371c9d4SSatish Balay cols[0] = p; 5369371c9d4SSatish Balay vals[0] = -(asb + an); 5379371c9d4SSatish Balay cols[1] = n; 5389371c9d4SSatish Balay vals[1] = an; 539c4762a1bSJed Brown } else if (i == 0 && j == s->ny - 1) { /* north-west corner */ 540c4762a1bSJed Brown *sz = 2; 5419371c9d4SSatish Balay cols[0] = s2; 5429371c9d4SSatish Balay vals[0] = as; 5439371c9d4SSatish Balay cols[1] = p; 5449371c9d4SSatish Balay vals[1] = -(as + anb); 545c4762a1bSJed Brown } else if (i == s->nx - 1 && j == 0) { /* south-east corner */ 546c4762a1bSJed Brown *sz = 2; 5479371c9d4SSatish Balay cols[0] = p; 5489371c9d4SSatish Balay vals[0] = -(asb + an); 5499371c9d4SSatish Balay cols[1] = n; 5509371c9d4SSatish Balay vals[1] = an; 551c4762a1bSJed Brown } else if (i == s->nx - 1 && j == s->ny - 1) { /* north-east corner */ 552c4762a1bSJed Brown *sz = 2; 5539371c9d4SSatish Balay cols[0] = s2; 5549371c9d4SSatish Balay vals[0] = as; 5559371c9d4SSatish Balay cols[1] = p; 5569371c9d4SSatish Balay vals[1] = -(as + anb); 557c4762a1bSJed Brown } else if (i == 0) { /* west boundary */ 558c4762a1bSJed Brown *sz = 3; 5599371c9d4SSatish Balay cols[0] = s2; 5609371c9d4SSatish Balay vals[0] = as; 5619371c9d4SSatish Balay cols[1] = p; 5629371c9d4SSatish Balay vals[1] = -(as + an); 5639371c9d4SSatish Balay cols[2] = n; 5649371c9d4SSatish Balay vals[2] = an; 565c4762a1bSJed Brown } else if (i == s->nx - 1) { /* east boundary */ 566c4762a1bSJed Brown *sz = 3; 5679371c9d4SSatish Balay cols[0] = s2; 5689371c9d4SSatish Balay vals[0] = as; 5699371c9d4SSatish Balay cols[1] = p; 5709371c9d4SSatish Balay vals[1] = -(as + an); 5719371c9d4SSatish Balay cols[2] = n; 5729371c9d4SSatish Balay vals[2] = an; 573c4762a1bSJed Brown } else if (j == 0) { /* south boundary */ 574c4762a1bSJed Brown *sz = 2; 5759371c9d4SSatish Balay cols[0] = p; 5769371c9d4SSatish Balay vals[0] = -(asb + an); 5779371c9d4SSatish Balay cols[1] = n; 5789371c9d4SSatish Balay vals[1] = an; 579c4762a1bSJed Brown } else if (j == s->ny - 1) { /* north boundary */ 580c4762a1bSJed Brown *sz = 2; 5819371c9d4SSatish Balay cols[0] = s2; 5829371c9d4SSatish Balay vals[0] = as; 5839371c9d4SSatish Balay cols[1] = p; 5849371c9d4SSatish Balay vals[1] = -(as + anb); 585c4762a1bSJed Brown } else { /* interior */ 586c4762a1bSJed Brown *sz = 3; 5879371c9d4SSatish Balay cols[0] = s2; 5889371c9d4SSatish Balay vals[0] = as; 5899371c9d4SSatish Balay cols[1] = p; 5909371c9d4SSatish Balay vals[1] = -(as + an); 5919371c9d4SSatish Balay cols[2] = n; 5929371c9d4SSatish Balay vals[2] = an; 593c4762a1bSJed Brown } 594c4762a1bSJed Brown PetscFunctionReturn(0); 595c4762a1bSJed Brown } 596c4762a1bSJed Brown 5979371c9d4SSatish Balay PetscErrorCode StokesRhsMomX(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val) { 598c4762a1bSJed Brown PetscScalar y = j * s->hy + s->hy / 2; 599c4762a1bSJed Brown PetscScalar awb = s->hy / (s->hx / 2); 600c4762a1bSJed Brown 601c4762a1bSJed Brown PetscFunctionBeginUser; 602c4762a1bSJed Brown if (i == 0) { /* west boundary */ 603c4762a1bSJed Brown *val = awb * StokesExactVelocityX(y); 604c4762a1bSJed Brown } else { 605c4762a1bSJed Brown *val = 0.0; 606c4762a1bSJed Brown } 607c4762a1bSJed Brown PetscFunctionReturn(0); 608c4762a1bSJed Brown } 609c4762a1bSJed Brown 6109371c9d4SSatish Balay PetscErrorCode StokesRhsMomY(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val) { 611c4762a1bSJed Brown PetscFunctionBeginUser; 612c4762a1bSJed Brown *val = 0.0; 613c4762a1bSJed Brown PetscFunctionReturn(0); 614c4762a1bSJed Brown } 615c4762a1bSJed Brown 6169371c9d4SSatish Balay PetscErrorCode StokesRhsMass(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val) { 617c4762a1bSJed Brown PetscScalar y = j * s->hy + s->hy / 2; 618c4762a1bSJed Brown PetscScalar aeb = s->hy; 619c4762a1bSJed Brown 620c4762a1bSJed Brown PetscFunctionBeginUser; 621c4762a1bSJed Brown if (i == 0) { /* west boundary */ 622c4762a1bSJed Brown *val = aeb * StokesExactVelocityX(y); 623c4762a1bSJed Brown } else { 624c4762a1bSJed Brown *val = 0.0; 625c4762a1bSJed Brown } 626c4762a1bSJed Brown PetscFunctionReturn(0); 627c4762a1bSJed Brown } 628c4762a1bSJed Brown 6299371c9d4SSatish Balay PetscErrorCode StokesCalcResidual(Stokes *s) { 630c4762a1bSJed Brown PetscReal val; 631c4762a1bSJed Brown Vec b0, b1; 632c4762a1bSJed Brown 633c4762a1bSJed Brown PetscFunctionBeginUser; 634c4762a1bSJed Brown /* residual Ax-b (*warning* overwrites b) */ 6359566063dSJacob Faibussowitsch PetscCall(VecScale(s->b, -1.0)); 6369566063dSJacob Faibussowitsch PetscCall(MatMultAdd(s->A, s->x, s->b, s->b)); 637c4762a1bSJed Brown 638c4762a1bSJed Brown /* residual velocity */ 6399566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(s->b, s->isg[0], &b0)); 6409566063dSJacob Faibussowitsch PetscCall(VecNorm(b0, NORM_2, &val)); 6419566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " residual u = %g\n", (double)val)); 6429566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(s->b, s->isg[0], &b0)); 643c4762a1bSJed Brown 644c4762a1bSJed Brown /* residual pressure */ 6459566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(s->b, s->isg[1], &b1)); 6469566063dSJacob Faibussowitsch PetscCall(VecNorm(b1, NORM_2, &val)); 6479566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " residual p = %g\n", (double)val)); 6489566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(s->b, s->isg[1], &b1)); 649c4762a1bSJed Brown 650c4762a1bSJed Brown /* total residual */ 6519566063dSJacob Faibussowitsch PetscCall(VecNorm(s->b, NORM_2, &val)); 6529566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " residual [u,p] = %g\n", (double)val)); 653c4762a1bSJed Brown PetscFunctionReturn(0); 654c4762a1bSJed Brown } 655c4762a1bSJed Brown 6569371c9d4SSatish Balay PetscErrorCode StokesCalcError(Stokes *s) { 657c4762a1bSJed Brown PetscScalar scale = PetscSqrtReal((double)s->nx * s->ny); 658c4762a1bSJed Brown PetscReal val; 659c4762a1bSJed Brown Vec y0, y1; 660c4762a1bSJed Brown 661c4762a1bSJed Brown PetscFunctionBeginUser; 662c4762a1bSJed Brown /* error y-x */ 6639566063dSJacob Faibussowitsch PetscCall(VecAXPY(s->y, -1.0, s->x)); 664c4762a1bSJed Brown 665c4762a1bSJed Brown /* error in velocity */ 6669566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(s->y, s->isg[0], &y0)); 6679566063dSJacob Faibussowitsch PetscCall(VecNorm(y0, NORM_2, &val)); 6689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " discretization error u = %g\n", (double)(PetscRealPart(val / scale)))); 6699566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(s->y, s->isg[0], &y0)); 670c4762a1bSJed Brown 671c4762a1bSJed Brown /* error in pressure */ 6729566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(s->y, s->isg[1], &y1)); 6739566063dSJacob Faibussowitsch PetscCall(VecNorm(y1, NORM_2, &val)); 6749566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " discretization error p = %g\n", (double)(PetscRealPart(val / scale)))); 6759566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(s->y, s->isg[1], &y1)); 676c4762a1bSJed Brown 677c4762a1bSJed Brown /* total error */ 6789566063dSJacob Faibussowitsch PetscCall(VecNorm(s->y, NORM_2, &val)); 6799566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " discretization error [u,p] = %g\n", (double)PetscRealPart((val / scale)))); 680c4762a1bSJed Brown PetscFunctionReturn(0); 681c4762a1bSJed Brown } 682c4762a1bSJed Brown 6839371c9d4SSatish Balay int main(int argc, char **argv) { 684c4762a1bSJed Brown Stokes s; 685c4762a1bSJed Brown KSP ksp; 686c4762a1bSJed Brown 687327415f7SBarry Smith PetscFunctionBeginUser; 6889566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 689c4762a1bSJed Brown s.nx = 4; 690c4762a1bSJed Brown s.ny = 6; 6919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nx", &s.nx, NULL)); 6929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-ny", &s.ny, NULL)); 693c4762a1bSJed Brown s.hx = 2.0 / s.nx; 694c4762a1bSJed Brown s.hy = 1.0 / s.ny; 695c4762a1bSJed Brown s.matsymmetric = PETSC_FALSE; 6969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-mat_set_symmetric", &s.matsymmetric, NULL)); 697c4762a1bSJed Brown s.userPC = s.userKSP = PETSC_FALSE; 6989566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-user_pc", &s.userPC)); 6999566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-user_ksp", &s.userKSP)); 700c4762a1bSJed Brown 7019566063dSJacob Faibussowitsch PetscCall(StokesSetupMatrix(&s)); 7029566063dSJacob Faibussowitsch PetscCall(StokesSetupIndexSets(&s)); 7039566063dSJacob Faibussowitsch PetscCall(StokesSetupVectors(&s)); 704c4762a1bSJed Brown 7059566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp)); 7069566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, s.A, s.A)); 7079566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 7089566063dSJacob Faibussowitsch PetscCall(StokesSetupPC(&s, ksp)); 7099566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, s.b, s.x)); 710c4762a1bSJed Brown 711c4762a1bSJed Brown /* don't trust, verify! */ 7129566063dSJacob Faibussowitsch PetscCall(StokesCalcResidual(&s)); 7139566063dSJacob Faibussowitsch PetscCall(StokesCalcError(&s)); 7149566063dSJacob Faibussowitsch PetscCall(StokesWriteSolution(&s)); 715c4762a1bSJed Brown 7169566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 7179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&s.subA[0])); 7189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&s.subA[1])); 7199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&s.subA[2])); 7209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&s.subA[3])); 7219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&s.A)); 7229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s.x)); 7239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s.b)); 7249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s.y)); 7259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&s.myS)); 7269566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 727b122ec5aSJacob Faibussowitsch return 0; 728c4762a1bSJed Brown } 729c4762a1bSJed Brown 730c4762a1bSJed Brown /*TEST 731c4762a1bSJed Brown 732c4762a1bSJed Brown test: 733c4762a1bSJed Brown nsize: 2 734c4762a1bSJed Brown args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_1_pc_type none 735c4762a1bSJed Brown 736c4762a1bSJed Brown test: 737c4762a1bSJed Brown suffix: 2 738c4762a1bSJed Brown nsize: 2 739c4762a1bSJed Brown args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc 740c4762a1bSJed Brown 741c4762a1bSJed Brown test: 742c4762a1bSJed Brown suffix: 3 743c4762a1bSJed Brown nsize: 2 744c4762a1bSJed Brown args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc 745c4762a1bSJed Brown 746c4762a1bSJed Brown test: 747c4762a1bSJed Brown suffix: 4 748c4762a1bSJed Brown nsize: 2 749c4762a1bSJed Brown args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_0_ksp_type gmres -fieldsplit_0_pc_type bjacobi -fieldsplit_1_pc_type jacobi -fieldsplit_1_inner_ksp_type preonly -fieldsplit_1_inner_pc_type jacobi -fieldsplit_1_upper_ksp_type preonly -fieldsplit_1_upper_pc_type jacobi 750c4762a1bSJed Brown 751c4762a1bSJed Brown test: 752c4762a1bSJed Brown suffix: 4_pcksp 753c4762a1bSJed Brown nsize: 2 754c4762a1bSJed Brown args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_0_ksp_type gmres -fieldsplit_0_pc_type bjacobi -fieldsplit_1_pc_type jacobi -fieldsplit_1_inner_ksp_type preonly -fieldsplit_1_upper_ksp_type preonly -fieldsplit_1_upper_pc_type jacobi 755c4762a1bSJed Brown 756c4762a1bSJed Brown test: 757c4762a1bSJed Brown suffix: 5 758c4762a1bSJed Brown nsize: 2 759c4762a1bSJed Brown args: -nx 4 -ny 8 -mat_set_symmetric -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type gkb -fieldsplit_0_ksp_type cg -fieldsplit_0_pc_type jacobi -fieldsplit_0_ksp_rtol 1e-10 760c4762a1bSJed Brown 761c4762a1bSJed Brown test: 762c4762a1bSJed Brown suffix: 6 763c4762a1bSJed Brown nsize: 2 764c4762a1bSJed Brown args: -nx 4 -ny 8 -mat_set_symmetric -ksp_type preonly -pc_type fieldsplit -pc_fieldsplit_type gkb -fieldsplit_0_ksp_type cg -fieldsplit_0_pc_type jacobi -fieldsplit_0_ksp_rtol 1e-10 765c4762a1bSJed Brown 766c4762a1bSJed Brown test: 767c4762a1bSJed Brown suffix: 7 768c4762a1bSJed Brown nsize: 2 769c4762a1bSJed Brown args: -nx 4 -ny 8 -mat_set_symmetric -ksp_type preonly -pc_type fieldsplit -pc_fieldsplit_type gkb -pc_fieldsplit_gkb_tol 1e-4 -pc_fieldsplit_gkb_nu 5 -fieldsplit_0_ksp_type cg -fieldsplit_0_pc_type jacobi -fieldsplit_0_ksp_rtol 1e-6 770c4762a1bSJed Brown 771c4762a1bSJed Brown TEST*/ 772