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 4c4762a1bSJed Brown /*---------------------------------------------------------------------------- */ 5c4762a1bSJed Brown /* 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 */ 6c4762a1bSJed Brown /*---------------------------------------------------------------------------- */ 7c4762a1bSJed Brown /* author : Christiaan M. Klaij */ 8c4762a1bSJed Brown /*---------------------------------------------------------------------------- */ 9c4762a1bSJed Brown /* */ 10c4762a1bSJed Brown /* Poiseuille flow problem. */ 11c4762a1bSJed Brown /* */ 12c4762a1bSJed Brown /* Viscous, laminar flow in a 2D channel with parabolic velocity */ 13c4762a1bSJed Brown /* profile and linear pressure drop, exact solution of the 2D Stokes */ 14c4762a1bSJed Brown /* equations. */ 15c4762a1bSJed Brown /* */ 16c4762a1bSJed Brown /* Discretized with the cell-centered finite-volume method on a */ 17c4762a1bSJed Brown /* Cartesian grid with co-located variables. Variables ordered as */ 18c4762a1bSJed Brown /* [u1...uN v1...vN p1...pN]^T. Matrix [A00 A01; A10, A11] solved with */ 19c4762a1bSJed Brown /* PCFIELDSPLIT. Lower factorization is used to mimick the Semi-Implicit */ 20c4762a1bSJed Brown /* Method for Pressure Linked Equations (SIMPLE) used as preconditioner */ 21c4762a1bSJed Brown /* instead of solver. */ 22c4762a1bSJed Brown /* */ 23c4762a1bSJed Brown /* Disclaimer: does not contain the pressure-weighed interpolation */ 24c4762a1bSJed Brown /* method needed to suppress spurious pressure modes in real-life */ 25c4762a1bSJed Brown /* problems. */ 26c4762a1bSJed Brown /* */ 27c4762a1bSJed Brown /* Usage: */ 28c4762a1bSJed Brown /* */ 29c4762a1bSJed Brown /* 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 */ 30c4762a1bSJed Brown /* */ 31c4762a1bSJed Brown /* Runs with PCFIELDSPLIT on 32x48 grid, no PC for the Schur */ 32c4762a1bSJed Brown /* complement because A11 is zero. FGMRES is needed because */ 33c4762a1bSJed Brown /* PCFIELDSPLIT is a variable preconditioner. */ 34c4762a1bSJed Brown /* */ 35c4762a1bSJed Brown /* 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 */ 36c4762a1bSJed Brown /* */ 37c4762a1bSJed Brown /* Same as above but with user defined PC for the true Schur */ 38c4762a1bSJed Brown /* complement. PC based on the SIMPLE-type approximation (inverse of */ 39c4762a1bSJed Brown /* A00 approximated by inverse of its diagonal). */ 40c4762a1bSJed Brown /* */ 41c4762a1bSJed Brown /* 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 */ 42c4762a1bSJed Brown /* */ 43c4762a1bSJed Brown /* Replace the true Schur complement with a user defined Schur */ 44c4762a1bSJed Brown /* complement based on the SIMPLE-type approximation. Same matrix is */ 45c4762a1bSJed Brown /* used as PC. */ 46c4762a1bSJed Brown /* */ 47c4762a1bSJed Brown /* 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 */ 48c4762a1bSJed Brown /* */ 49c4762a1bSJed Brown /* Out-of-the-box SIMPLE-type preconditioning. The major advantage */ 50c4762a1bSJed Brown /* is that the user neither needs to provide the approximation of */ 51c4762a1bSJed Brown /* the Schur complement, nor the corresponding preconditioner. */ 52c4762a1bSJed Brown /* */ 53c4762a1bSJed Brown /*---------------------------------------------------------------------------- */ 54c4762a1bSJed Brown 55c4762a1bSJed Brown 56c4762a1bSJed Brown #include <petscksp.h> 57c4762a1bSJed Brown 58c4762a1bSJed Brown typedef struct { 59c4762a1bSJed Brown PetscBool userPC, userKSP, matsymmetric; /* user defined preconditioner and matrix for the Schur complement */ 60c4762a1bSJed Brown PetscInt nx, ny; /* nb of cells in x- and y-direction */ 61c4762a1bSJed Brown PetscReal hx, hy; /* mesh size in x- and y-direction */ 62c4762a1bSJed Brown Mat A; /* block matrix */ 63c4762a1bSJed Brown Mat subA[4]; /* the four blocks */ 64c4762a1bSJed Brown Mat myS; /* the approximation of the Schur complement */ 65c4762a1bSJed Brown Vec x, b, y; /* solution, rhs and temporary vector */ 66c4762a1bSJed Brown IS isg[2]; /* index sets of split "0" and "1" */ 67c4762a1bSJed Brown } Stokes; 68c4762a1bSJed Brown 69c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock00(Stokes*); /* setup the block Q */ 70c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock01(Stokes*); /* setup the block G */ 71c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock10(Stokes*); /* setup the block D (equal to the transpose of G) */ 72c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock11(Stokes*); /* setup the block C (equal to zero) */ 73c4762a1bSJed Brown 74c4762a1bSJed Brown PetscErrorCode StokesGetPosition(Stokes*, PetscInt, PetscInt*, PetscInt*); /* row number j*nx+i corresponds to position (i,j) in grid */ 75c4762a1bSJed Brown 76c4762a1bSJed Brown PetscErrorCode StokesStencilLaplacian(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*); /* stencil of the Laplacian operator */ 77c4762a1bSJed Brown PetscErrorCode StokesStencilGradientX(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*); /* stencil of the Gradient operator (x-component) */ 78c4762a1bSJed Brown PetscErrorCode StokesStencilGradientY(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*); /* stencil of the Gradient operator (y-component) */ 79c4762a1bSJed Brown 80c4762a1bSJed Brown PetscErrorCode StokesRhs(Stokes*); /* rhs vector */ 81c4762a1bSJed Brown PetscErrorCode StokesRhsMomX(Stokes*, PetscInt, PetscInt, PetscScalar*); /* right hand side of velocity (x-component) */ 82c4762a1bSJed Brown PetscErrorCode StokesRhsMomY(Stokes*, PetscInt, PetscInt, PetscScalar*); /* right hand side of velocity (y-component) */ 83c4762a1bSJed Brown PetscErrorCode StokesRhsMass(Stokes*, PetscInt, PetscInt, PetscScalar*); /* right hand side of pressure */ 84c4762a1bSJed Brown 85c4762a1bSJed Brown PetscErrorCode StokesSetupApproxSchur(Stokes*); /* approximation of the Schur complement */ 86c4762a1bSJed Brown 87c4762a1bSJed Brown PetscErrorCode StokesExactSolution(Stokes*); /* exact solution vector */ 88c4762a1bSJed Brown PetscErrorCode StokesWriteSolution(Stokes*); /* write solution to file */ 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* exact solution for the velocity (x-component, y-component is zero) */ 91c4762a1bSJed Brown PetscScalar StokesExactVelocityX(const PetscScalar y) 92c4762a1bSJed Brown { 93c4762a1bSJed Brown return 4.0*y*(1.0-y); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* exact solution for the pressure */ 97c4762a1bSJed Brown PetscScalar StokesExactPressure(const PetscScalar x) 98c4762a1bSJed Brown { 99c4762a1bSJed Brown return 8.0*(2.0-x); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown PetscErrorCode StokesSetupPC(Stokes *s, KSP ksp) 103c4762a1bSJed Brown { 104c4762a1bSJed Brown KSP *subksp; 105c4762a1bSJed Brown PC pc; 106c4762a1bSJed Brown PetscInt n = 1; 107c4762a1bSJed Brown PetscErrorCode ierr; 108c4762a1bSJed Brown 109c4762a1bSJed Brown PetscFunctionBeginUser; 110c4762a1bSJed Brown ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = PCFieldSplitSetIS(pc, "0", s->isg[0]);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = PCFieldSplitSetIS(pc, "1", s->isg[1]);CHKERRQ(ierr); 113c4762a1bSJed Brown if (s->userPC) { 114c4762a1bSJed Brown ierr = PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_USER, s->myS);CHKERRQ(ierr); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown if (s->userKSP) { 117c4762a1bSJed Brown ierr = PCSetUp(pc);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = PCFieldSplitGetSubKSP(pc, &n, &subksp);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = KSPSetOperators(subksp[1], s->myS, s->myS);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = PetscFree(subksp);CHKERRQ(ierr); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown PetscFunctionReturn(0); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125c4762a1bSJed Brown PetscErrorCode StokesWriteSolution(Stokes *s) 126c4762a1bSJed Brown { 127c4762a1bSJed Brown PetscMPIInt size; 128c4762a1bSJed Brown PetscInt n,i,j; 129c4762a1bSJed Brown const PetscScalar *array; 130c4762a1bSJed Brown PetscErrorCode ierr; 131c4762a1bSJed Brown 132c4762a1bSJed Brown PetscFunctionBeginUser; 133c4762a1bSJed Brown /* write data (*warning* only works sequential) */ 134*ffc4695bSBarry Smith ierr = MPI_Comm_size(MPI_COMM_WORLD,&size);CHKERRMPI(ierr); 135c4762a1bSJed Brown /*ierr = PetscPrintf(PETSC_COMM_WORLD," number of processors = %D\n",size);CHKERRQ(ierr);*/ 136c4762a1bSJed Brown if (size == 1) { 137c4762a1bSJed Brown PetscViewer viewer; 138c4762a1bSJed Brown ierr = VecGetArrayRead(s->x, &array);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD, "solution.dat", &viewer);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = PetscViewerASCIIPrintf(viewer, "# x, y, u, v, p\n");CHKERRQ(ierr); 141c4762a1bSJed Brown for (j = 0; j < s->ny; j++) { 142c4762a1bSJed Brown for (i = 0; i < s->nx; i++) { 143c4762a1bSJed Brown n = j*s->nx+i; 144c4762a1bSJed Brown ierr = 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]));CHKERRQ(ierr); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown } 147c4762a1bSJed Brown ierr = VecRestoreArrayRead(s->x, &array);CHKERRQ(ierr); 148c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 149c4762a1bSJed Brown } 150c4762a1bSJed Brown PetscFunctionReturn(0); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153c4762a1bSJed Brown PetscErrorCode StokesSetupIndexSets(Stokes *s) 154c4762a1bSJed Brown { 155c4762a1bSJed Brown PetscErrorCode ierr; 156c4762a1bSJed Brown 157c4762a1bSJed Brown PetscFunctionBeginUser; 158c4762a1bSJed Brown /* the two index sets */ 159c4762a1bSJed Brown ierr = MatNestGetISs(s->A, s->isg, NULL);CHKERRQ(ierr); 160c4762a1bSJed Brown /* ISView(isg[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 161c4762a1bSJed Brown /* ISView(isg[1],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 162c4762a1bSJed Brown PetscFunctionReturn(0); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 165c4762a1bSJed Brown PetscErrorCode StokesSetupVectors(Stokes *s) 166c4762a1bSJed Brown { 167c4762a1bSJed Brown PetscErrorCode ierr; 168c4762a1bSJed Brown 169c4762a1bSJed Brown PetscFunctionBeginUser; 170c4762a1bSJed Brown /* solution vector x */ 171c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD, &s->x);CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = VecSetSizes(s->x, PETSC_DECIDE, 3*s->nx*s->ny);CHKERRQ(ierr); 173c4762a1bSJed Brown ierr = VecSetType(s->x, VECMPI);CHKERRQ(ierr); 174c4762a1bSJed Brown /* ierr = VecSetRandom(s->x, NULL);CHKERRQ(ierr); */ 175c4762a1bSJed Brown /* ierr = VecView(s->x, (PetscViewer) PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); */ 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* exact solution y */ 178c4762a1bSJed Brown ierr = VecDuplicate(s->x, &s->y);CHKERRQ(ierr); 179c4762a1bSJed Brown ierr = StokesExactSolution(s);CHKERRQ(ierr); 180c4762a1bSJed Brown /* ierr = VecView(s->y, (PetscViewer) PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); */ 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* rhs vector b */ 183c4762a1bSJed Brown ierr = VecDuplicate(s->x, &s->b);CHKERRQ(ierr); 184c4762a1bSJed Brown ierr = StokesRhs(s);CHKERRQ(ierr); 185c4762a1bSJed Brown /*ierr = VecView(s->b, (PetscViewer) PETSC_VIEWER_DEFAULT);CHKERRQ(ierr);*/ 186c4762a1bSJed Brown PetscFunctionReturn(0); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown PetscErrorCode StokesGetPosition(Stokes *s, PetscInt row, PetscInt *i, PetscInt *j) 190c4762a1bSJed Brown { 191c4762a1bSJed Brown PetscInt n; 192c4762a1bSJed Brown 193c4762a1bSJed Brown PetscFunctionBeginUser; 194c4762a1bSJed Brown /* cell number n=j*nx+i has position (i,j) in grid */ 195c4762a1bSJed Brown n = row%(s->nx*s->ny); 196c4762a1bSJed Brown *i = n%s->nx; 197c4762a1bSJed Brown *j = (n-(*i))/s->nx; 198c4762a1bSJed Brown PetscFunctionReturn(0); 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201c4762a1bSJed Brown PetscErrorCode StokesExactSolution(Stokes *s) 202c4762a1bSJed Brown { 203c4762a1bSJed Brown PetscInt row, start, end, i, j; 204c4762a1bSJed Brown PetscScalar val; 205c4762a1bSJed Brown Vec y0,y1; 206c4762a1bSJed Brown PetscErrorCode ierr; 207c4762a1bSJed Brown 208c4762a1bSJed Brown PetscFunctionBeginUser; 209c4762a1bSJed Brown /* velocity part */ 210c4762a1bSJed Brown ierr = VecGetSubVector(s->y, s->isg[0], &y0);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = VecGetOwnershipRange(y0, &start, &end);CHKERRQ(ierr); 212c4762a1bSJed Brown for (row = start; row < end; row++) { 213c4762a1bSJed Brown ierr = StokesGetPosition(s, row,&i,&j);CHKERRQ(ierr); 214c4762a1bSJed Brown if (row < s->nx*s->ny) { 215c4762a1bSJed Brown val = StokesExactVelocityX(j*s->hy+s->hy/2); 216c4762a1bSJed Brown } else { 217c4762a1bSJed Brown val = 0; 218c4762a1bSJed Brown } 219c4762a1bSJed Brown ierr = VecSetValue(y0, row, val, INSERT_VALUES);CHKERRQ(ierr); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown ierr = VecRestoreSubVector(s->y, s->isg[0], &y0);CHKERRQ(ierr); 222c4762a1bSJed Brown 223c4762a1bSJed Brown /* pressure part */ 224c4762a1bSJed Brown ierr = VecGetSubVector(s->y, s->isg[1], &y1);CHKERRQ(ierr); 225c4762a1bSJed Brown ierr = VecGetOwnershipRange(y1, &start, &end);CHKERRQ(ierr); 226c4762a1bSJed Brown for (row = start; row < end; row++) { 227c4762a1bSJed Brown ierr = StokesGetPosition(s, row, &i, &j);CHKERRQ(ierr); 228c4762a1bSJed Brown val = StokesExactPressure(i*s->hx+s->hx/2); 229c4762a1bSJed Brown ierr = VecSetValue(y1, row, val, INSERT_VALUES);CHKERRQ(ierr); 230c4762a1bSJed Brown } 231c4762a1bSJed Brown ierr = VecRestoreSubVector(s->y, s->isg[1], &y1);CHKERRQ(ierr); 232c4762a1bSJed Brown PetscFunctionReturn(0); 233c4762a1bSJed Brown } 234c4762a1bSJed Brown 235c4762a1bSJed Brown PetscErrorCode StokesRhs(Stokes *s) 236c4762a1bSJed Brown { 237c4762a1bSJed Brown PetscInt row, start, end, i, j; 238c4762a1bSJed Brown PetscScalar val; 239c4762a1bSJed Brown Vec b0,b1; 240c4762a1bSJed Brown PetscErrorCode ierr; 241c4762a1bSJed Brown 242c4762a1bSJed Brown PetscFunctionBeginUser; 243c4762a1bSJed Brown /* velocity part */ 244c4762a1bSJed Brown ierr = VecGetSubVector(s->b, s->isg[0], &b0);CHKERRQ(ierr); 245c4762a1bSJed Brown ierr = VecGetOwnershipRange(b0, &start, &end);CHKERRQ(ierr); 246c4762a1bSJed Brown for (row = start; row < end; row++) { 247c4762a1bSJed Brown ierr = StokesGetPosition(s, row, &i, &j);CHKERRQ(ierr); 248c4762a1bSJed Brown if (row < s->nx*s->ny) { 249c4762a1bSJed Brown ierr = StokesRhsMomX(s, i, j, &val);CHKERRQ(ierr); 250c4762a1bSJed Brown } else { 251c4762a1bSJed Brown ierr = StokesRhsMomY(s, i, j, &val);CHKERRQ(ierr); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown ierr = VecSetValue(b0, row, val, INSERT_VALUES);CHKERRQ(ierr); 254c4762a1bSJed Brown } 255c4762a1bSJed Brown ierr = VecRestoreSubVector(s->b, s->isg[0], &b0);CHKERRQ(ierr); 256c4762a1bSJed Brown 257c4762a1bSJed Brown /* pressure part */ 258c4762a1bSJed Brown ierr = VecGetSubVector(s->b, s->isg[1], &b1);CHKERRQ(ierr); 259c4762a1bSJed Brown ierr = VecGetOwnershipRange(b1, &start, &end);CHKERRQ(ierr); 260c4762a1bSJed Brown for (row = start; row < end; row++) { 261c4762a1bSJed Brown ierr = StokesGetPosition(s, row, &i, &j);CHKERRQ(ierr); 262c4762a1bSJed Brown ierr = StokesRhsMass(s, i, j, &val);CHKERRQ(ierr); 263c4762a1bSJed Brown if (s->matsymmetric) { 264c4762a1bSJed Brown val = -1.0*val; 265c4762a1bSJed Brown } 266c4762a1bSJed Brown ierr = VecSetValue(b1, row, val, INSERT_VALUES);CHKERRQ(ierr); 267c4762a1bSJed Brown } 268c4762a1bSJed Brown ierr = VecRestoreSubVector(s->b, s->isg[1], &b1);CHKERRQ(ierr); 269c4762a1bSJed Brown PetscFunctionReturn(0); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 272c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock00(Stokes *s) 273c4762a1bSJed Brown { 274c4762a1bSJed Brown PetscInt row, start, end, sz, i, j; 275c4762a1bSJed Brown PetscInt cols[5]; 276c4762a1bSJed Brown PetscScalar vals[5]; 277c4762a1bSJed Brown PetscErrorCode ierr; 278c4762a1bSJed Brown 279c4762a1bSJed Brown PetscFunctionBeginUser; 280c4762a1bSJed Brown /* A[0] is 2N-by-2N */ 281c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&s->subA[0]);CHKERRQ(ierr); 282c4762a1bSJed Brown ierr = MatSetOptionsPrefix(s->subA[0],"a00_");CHKERRQ(ierr); 283c4762a1bSJed Brown ierr = MatSetSizes(s->subA[0],PETSC_DECIDE,PETSC_DECIDE,2*s->nx*s->ny,2*s->nx*s->ny);CHKERRQ(ierr); 284c4762a1bSJed Brown ierr = MatSetType(s->subA[0],MATMPIAIJ);CHKERRQ(ierr); 285c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(s->subA[0],5,NULL,5,NULL);CHKERRQ(ierr); 286c4762a1bSJed Brown ierr = MatGetOwnershipRange(s->subA[0], &start, &end);CHKERRQ(ierr); 287c4762a1bSJed Brown 288c4762a1bSJed Brown for (row = start; row < end; row++) { 289c4762a1bSJed Brown ierr = StokesGetPosition(s, row, &i, &j);CHKERRQ(ierr); 290c4762a1bSJed Brown /* first part: rows 0 to (nx*ny-1) */ 291c4762a1bSJed Brown ierr = StokesStencilLaplacian(s, i, j, &sz, cols, vals);CHKERRQ(ierr); 292c4762a1bSJed Brown /* second part: rows (nx*ny) to (2*nx*ny-1) */ 293c4762a1bSJed Brown if (row >= s->nx*s->ny) { 294c4762a1bSJed Brown for (i = 0; i < sz; i++) cols[i] += s->nx*s->ny; 295c4762a1bSJed Brown } 296c4762a1bSJed Brown for (i = 0; i < sz; i++) vals[i] = -1.0*vals[i]; /* dynamic viscosity coef mu=-1 */ 297c4762a1bSJed Brown ierr = MatSetValues(s->subA[0], 1, &row, sz, cols, vals, INSERT_VALUES);CHKERRQ(ierr); 298c4762a1bSJed Brown } 299c4762a1bSJed Brown ierr = MatAssemblyBegin(s->subA[0], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 300c4762a1bSJed Brown ierr = MatAssemblyEnd(s->subA[0], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 301c4762a1bSJed Brown PetscFunctionReturn(0); 302c4762a1bSJed Brown } 303c4762a1bSJed Brown 304c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock01(Stokes *s) 305c4762a1bSJed Brown { 306c4762a1bSJed Brown PetscInt row, start, end, sz, i, j; 307c4762a1bSJed Brown PetscInt cols[5]; 308c4762a1bSJed Brown PetscScalar vals[5]; 309c4762a1bSJed Brown PetscErrorCode ierr; 310c4762a1bSJed Brown 311c4762a1bSJed Brown PetscFunctionBeginUser; 312c4762a1bSJed Brown /* A[1] is 2N-by-N */ 313c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD, &s->subA[1]);CHKERRQ(ierr); 314c4762a1bSJed Brown ierr = MatSetOptionsPrefix(s->subA[1],"a01_");CHKERRQ(ierr); 315c4762a1bSJed Brown ierr = MatSetSizes(s->subA[1],PETSC_DECIDE,PETSC_DECIDE,2*s->nx*s->ny,s->nx*s->ny);CHKERRQ(ierr); 316c4762a1bSJed Brown ierr = MatSetType(s->subA[1],MATMPIAIJ);CHKERRQ(ierr); 317c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(s->subA[1],5,NULL,5,NULL);CHKERRQ(ierr); 318c4762a1bSJed Brown ierr = MatGetOwnershipRange(s->subA[1],&start,&end);CHKERRQ(ierr); 319c4762a1bSJed Brown 320c4762a1bSJed Brown ierr = MatSetOption(s->subA[1],MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 321c4762a1bSJed Brown 322c4762a1bSJed Brown for (row = start; row < end; row++) { 323c4762a1bSJed Brown ierr = StokesGetPosition(s, row, &i, &j);CHKERRQ(ierr); 324c4762a1bSJed Brown /* first part: rows 0 to (nx*ny-1) */ 325c4762a1bSJed Brown if (row < s->nx*s->ny) { 326c4762a1bSJed Brown ierr = StokesStencilGradientX(s, i, j, &sz, cols, vals);CHKERRQ(ierr); 327c4762a1bSJed Brown } else { /* second part: rows (nx*ny) to (2*nx*ny-1) */ 328c4762a1bSJed Brown ierr = StokesStencilGradientY(s, i, j, &sz, cols, vals);CHKERRQ(ierr); 329c4762a1bSJed Brown } 330c4762a1bSJed Brown ierr = MatSetValues(s->subA[1], 1, &row, sz, cols, vals, INSERT_VALUES);CHKERRQ(ierr); 331c4762a1bSJed Brown } 332c4762a1bSJed Brown ierr = MatAssemblyBegin(s->subA[1], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 333c4762a1bSJed Brown ierr = MatAssemblyEnd(s->subA[1], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 334c4762a1bSJed Brown PetscFunctionReturn(0); 335c4762a1bSJed Brown } 336c4762a1bSJed Brown 337c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock10(Stokes *s) 338c4762a1bSJed Brown { 339c4762a1bSJed Brown PetscErrorCode ierr; 340c4762a1bSJed Brown 341c4762a1bSJed Brown PetscFunctionBeginUser; 342c4762a1bSJed Brown /* A[2] is minus transpose of A[1] */ 343c4762a1bSJed Brown ierr = MatTranspose(s->subA[1], MAT_INITIAL_MATRIX, &s->subA[2]);CHKERRQ(ierr); 344c4762a1bSJed Brown if (!s->matsymmetric) { 345c4762a1bSJed Brown ierr = MatScale(s->subA[2], -1.0);CHKERRQ(ierr); 346c4762a1bSJed Brown } 347c4762a1bSJed Brown ierr = MatSetOptionsPrefix(s->subA[2], "a10_");CHKERRQ(ierr); 348c4762a1bSJed Brown PetscFunctionReturn(0); 349c4762a1bSJed Brown } 350c4762a1bSJed Brown 351c4762a1bSJed Brown PetscErrorCode StokesSetupMatBlock11(Stokes *s) 352c4762a1bSJed Brown { 353c4762a1bSJed Brown PetscErrorCode ierr; 354c4762a1bSJed Brown 355c4762a1bSJed Brown PetscFunctionBeginUser; 356c4762a1bSJed Brown /* A[3] is N-by-N null matrix */ 357c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD, &s->subA[3]);CHKERRQ(ierr); 358c4762a1bSJed Brown ierr = MatSetOptionsPrefix(s->subA[3], "a11_");CHKERRQ(ierr); 359c4762a1bSJed Brown ierr = MatSetSizes(s->subA[3], PETSC_DECIDE, PETSC_DECIDE, s->nx*s->ny, s->nx*s->ny);CHKERRQ(ierr); 360c4762a1bSJed Brown ierr = MatSetType(s->subA[3], MATMPIAIJ);CHKERRQ(ierr); 361c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(s->subA[3], 0, NULL, 0, NULL);CHKERRQ(ierr); 362c4762a1bSJed Brown ierr = MatAssemblyBegin(s->subA[3], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 363c4762a1bSJed Brown ierr = MatAssemblyEnd(s->subA[3], MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 364c4762a1bSJed Brown PetscFunctionReturn(0); 365c4762a1bSJed Brown } 366c4762a1bSJed Brown 367c4762a1bSJed Brown PetscErrorCode StokesSetupApproxSchur(Stokes *s) 368c4762a1bSJed Brown { 369c4762a1bSJed Brown Vec diag; 370c4762a1bSJed Brown PetscErrorCode ierr; 371c4762a1bSJed Brown 372c4762a1bSJed Brown PetscFunctionBeginUser; 373c4762a1bSJed Brown /* Schur complement approximation: myS = A11 - A10 diag(A00)^(-1) A01 */ 374c4762a1bSJed Brown /* note: A11 is zero */ 375c4762a1bSJed Brown /* note: in real life this matrix would be build directly, */ 376c4762a1bSJed Brown /* i.e. without MatMatMult */ 377c4762a1bSJed Brown 378c4762a1bSJed Brown /* inverse of diagonal of A00 */ 379c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&diag);CHKERRQ(ierr); 380c4762a1bSJed Brown ierr = VecSetSizes(diag,PETSC_DECIDE,2*s->nx*s->ny);CHKERRQ(ierr); 381c4762a1bSJed Brown ierr = VecSetType(diag,VECMPI);CHKERRQ(ierr); 382c4762a1bSJed Brown ierr = MatGetDiagonal(s->subA[0],diag);CHKERRQ(ierr); 383c4762a1bSJed Brown ierr = VecReciprocal(diag);CHKERRQ(ierr); 384c4762a1bSJed Brown 385c4762a1bSJed Brown /* compute: - A10 diag(A00)^(-1) A01 */ 386c4762a1bSJed Brown ierr = MatDiagonalScale(s->subA[1],diag,NULL);CHKERRQ(ierr); /* (*warning* overwrites subA[1]) */ 387c4762a1bSJed Brown ierr = MatMatMult(s->subA[2],s->subA[1],MAT_INITIAL_MATRIX,PETSC_DEFAULT,&s->myS);CHKERRQ(ierr); 388c4762a1bSJed Brown ierr = MatScale(s->myS,-1.0);CHKERRQ(ierr); 389c4762a1bSJed Brown 390c4762a1bSJed Brown /* restore A10 */ 391c4762a1bSJed Brown ierr = MatGetDiagonal(s->subA[0],diag);CHKERRQ(ierr); 392c4762a1bSJed Brown ierr = MatDiagonalScale(s->subA[1],diag,NULL);CHKERRQ(ierr); 393c4762a1bSJed Brown ierr = VecDestroy(&diag);CHKERRQ(ierr); 394c4762a1bSJed Brown PetscFunctionReturn(0); 395c4762a1bSJed Brown } 396c4762a1bSJed Brown 397c4762a1bSJed Brown PetscErrorCode StokesSetupMatrix(Stokes *s) 398c4762a1bSJed Brown { 399c4762a1bSJed Brown PetscErrorCode ierr; 400c4762a1bSJed Brown 401c4762a1bSJed Brown PetscFunctionBeginUser; 402c4762a1bSJed Brown ierr = StokesSetupMatBlock00(s);CHKERRQ(ierr); 403c4762a1bSJed Brown ierr = StokesSetupMatBlock01(s);CHKERRQ(ierr); 404c4762a1bSJed Brown ierr = StokesSetupMatBlock10(s);CHKERRQ(ierr); 405c4762a1bSJed Brown ierr = StokesSetupMatBlock11(s);CHKERRQ(ierr); 406c4762a1bSJed Brown ierr = MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, s->subA, &s->A);CHKERRQ(ierr); 407c4762a1bSJed Brown ierr = StokesSetupApproxSchur(s);CHKERRQ(ierr); 408c4762a1bSJed Brown PetscFunctionReturn(0); 409c4762a1bSJed Brown } 410c4762a1bSJed Brown 411c4762a1bSJed Brown PetscErrorCode StokesStencilLaplacian(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals) 412c4762a1bSJed Brown { 413c4762a1bSJed Brown PetscInt p =j*s->nx+i, w=p-1, e=p+1, s2=p-s->nx, n=p+s->nx; 414c4762a1bSJed Brown PetscScalar ae=s->hy/s->hx, aeb=0; 415c4762a1bSJed Brown PetscScalar aw=s->hy/s->hx, awb=s->hy/(s->hx/2); 416c4762a1bSJed Brown PetscScalar as=s->hx/s->hy, asb=s->hx/(s->hy/2); 417c4762a1bSJed Brown PetscScalar an=s->hx/s->hy, anb=s->hx/(s->hy/2); 418c4762a1bSJed Brown 419c4762a1bSJed Brown PetscFunctionBeginUser; 420c4762a1bSJed Brown if (i==0 && j==0) { /* south-west corner */ 421c4762a1bSJed Brown *sz =3; 422c4762a1bSJed Brown cols[0]=p; vals[0]=-(ae+awb+asb+an); 423c4762a1bSJed Brown cols[1]=e; vals[1]=ae; 424c4762a1bSJed Brown cols[2]=n; vals[2]=an; 425c4762a1bSJed Brown } else if (i==0 && j==s->ny-1) { /* north-west corner */ 426c4762a1bSJed Brown *sz =3; 427c4762a1bSJed Brown cols[0]=s2; vals[0]=as; 428c4762a1bSJed Brown cols[1]=p; vals[1]=-(ae+awb+as+anb); 429c4762a1bSJed Brown cols[2]=e; vals[2]=ae; 430c4762a1bSJed Brown } else if (i==s->nx-1 && j==0) { /* south-east corner */ 431c4762a1bSJed Brown *sz =3; 432c4762a1bSJed Brown cols[0]=w; vals[0]=aw; 433c4762a1bSJed Brown cols[1]=p; vals[1]=-(aeb+aw+asb+an); 434c4762a1bSJed Brown cols[2]=n; vals[2]=an; 435c4762a1bSJed Brown } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */ 436c4762a1bSJed Brown *sz =3; 437c4762a1bSJed Brown cols[0]=s2; vals[0]=as; 438c4762a1bSJed Brown cols[1]=w; vals[1]=aw; 439c4762a1bSJed Brown cols[2]=p; vals[2]=-(aeb+aw+as+anb); 440c4762a1bSJed Brown } else if (i==0) { /* west boundary */ 441c4762a1bSJed Brown *sz =4; 442c4762a1bSJed Brown cols[0]=s2; vals[0]=as; 443c4762a1bSJed Brown cols[1]=p; vals[1]=-(ae+awb+as+an); 444c4762a1bSJed Brown cols[2]=e; vals[2]=ae; 445c4762a1bSJed Brown cols[3]=n; vals[3]=an; 446c4762a1bSJed Brown } else if (i==s->nx-1) { /* east boundary */ 447c4762a1bSJed Brown *sz =4; 448c4762a1bSJed Brown cols[0]=s2; vals[0]=as; 449c4762a1bSJed Brown cols[1]=w; vals[1]=aw; 450c4762a1bSJed Brown cols[2]=p; vals[2]=-(aeb+aw+as+an); 451c4762a1bSJed Brown cols[3]=n; vals[3]=an; 452c4762a1bSJed Brown } else if (j==0) { /* south boundary */ 453c4762a1bSJed Brown *sz =4; 454c4762a1bSJed Brown cols[0]=w; vals[0]=aw; 455c4762a1bSJed Brown cols[1]=p; vals[1]=-(ae+aw+asb+an); 456c4762a1bSJed Brown cols[2]=e; vals[2]=ae; 457c4762a1bSJed Brown cols[3]=n; vals[3]=an; 458c4762a1bSJed Brown } else if (j==s->ny-1) { /* north boundary */ 459c4762a1bSJed Brown *sz =4; 460c4762a1bSJed Brown cols[0]=s2; vals[0]=as; 461c4762a1bSJed Brown cols[1]=w; vals[1]=aw; 462c4762a1bSJed Brown cols[2]=p; vals[2]=-(ae+aw+as+anb); 463c4762a1bSJed Brown cols[3]=e; vals[3]=ae; 464c4762a1bSJed Brown } else { /* interior */ 465c4762a1bSJed Brown *sz =5; 466c4762a1bSJed Brown cols[0]=s2; vals[0]=as; 467c4762a1bSJed Brown cols[1]=w; vals[1]=aw; 468c4762a1bSJed Brown cols[2]=p; vals[2]=-(ae+aw+as+an); 469c4762a1bSJed Brown cols[3]=e; vals[3]=ae; 470c4762a1bSJed Brown cols[4]=n; vals[4]=an; 471c4762a1bSJed Brown } 472c4762a1bSJed Brown PetscFunctionReturn(0); 473c4762a1bSJed Brown } 474c4762a1bSJed Brown 475c4762a1bSJed Brown PetscErrorCode StokesStencilGradientX(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals) 476c4762a1bSJed Brown { 477c4762a1bSJed Brown PetscInt p =j*s->nx+i, w=p-1, e=p+1; 478c4762a1bSJed Brown PetscScalar ae= s->hy/2, aeb=s->hy; 479c4762a1bSJed Brown PetscScalar aw=-s->hy/2, awb=0; 480c4762a1bSJed Brown 481c4762a1bSJed Brown PetscFunctionBeginUser; 482c4762a1bSJed Brown if (i==0 && j==0) { /* south-west corner */ 483c4762a1bSJed Brown *sz =2; 484c4762a1bSJed Brown cols[0]=p; vals[0]=-(ae+awb); 485c4762a1bSJed Brown cols[1]=e; vals[1]=ae; 486c4762a1bSJed Brown } else if (i==0 && j==s->ny-1) { /* north-west corner */ 487c4762a1bSJed Brown *sz =2; 488c4762a1bSJed Brown cols[0]=p; vals[0]=-(ae+awb); 489c4762a1bSJed Brown cols[1]=e; vals[1]=ae; 490c4762a1bSJed Brown } else if (i==s->nx-1 && j==0) { /* south-east corner */ 491c4762a1bSJed Brown *sz =2; 492c4762a1bSJed Brown cols[0]=w; vals[0]=aw; 493c4762a1bSJed Brown cols[1]=p; vals[1]=-(aeb+aw); 494c4762a1bSJed Brown } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */ 495c4762a1bSJed Brown *sz =2; 496c4762a1bSJed Brown cols[0]=w; vals[0]=aw; 497c4762a1bSJed Brown cols[1]=p; vals[1]=-(aeb+aw); 498c4762a1bSJed Brown } else if (i==0) { /* west boundary */ 499c4762a1bSJed Brown *sz =2; 500c4762a1bSJed Brown cols[0]=p; vals[0]=-(ae+awb); 501c4762a1bSJed Brown cols[1]=e; vals[1]=ae; 502c4762a1bSJed Brown } else if (i==s->nx-1) { /* east boundary */ 503c4762a1bSJed Brown *sz =2; 504c4762a1bSJed Brown cols[0]=w; vals[0]=aw; 505c4762a1bSJed Brown cols[1]=p; vals[1]=-(aeb+aw); 506c4762a1bSJed Brown } else if (j==0) { /* south boundary */ 507c4762a1bSJed Brown *sz =3; 508c4762a1bSJed Brown cols[0]=w; vals[0]=aw; 509c4762a1bSJed Brown cols[1]=p; vals[1]=-(ae+aw); 510c4762a1bSJed Brown cols[2]=e; vals[2]=ae; 511c4762a1bSJed Brown } else if (j==s->ny-1) { /* north boundary */ 512c4762a1bSJed Brown *sz =3; 513c4762a1bSJed Brown cols[0]=w; vals[0]=aw; 514c4762a1bSJed Brown cols[1]=p; vals[1]=-(ae+aw); 515c4762a1bSJed Brown cols[2]=e; vals[2]=ae; 516c4762a1bSJed Brown } else { /* interior */ 517c4762a1bSJed Brown *sz =3; 518c4762a1bSJed Brown cols[0]=w; vals[0]=aw; 519c4762a1bSJed Brown cols[1]=p; vals[1]=-(ae+aw); 520c4762a1bSJed Brown cols[2]=e; vals[2]=ae; 521c4762a1bSJed Brown } 522c4762a1bSJed Brown PetscFunctionReturn(0); 523c4762a1bSJed Brown } 524c4762a1bSJed Brown 525c4762a1bSJed Brown PetscErrorCode StokesStencilGradientY(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals) 526c4762a1bSJed Brown { 527c4762a1bSJed Brown PetscInt p =j*s->nx+i, s2=p-s->nx, n=p+s->nx; 528c4762a1bSJed Brown PetscScalar as=-s->hx/2, asb=0; 529c4762a1bSJed Brown PetscScalar an= s->hx/2, anb=0; 530c4762a1bSJed Brown 531c4762a1bSJed Brown PetscFunctionBeginUser; 532c4762a1bSJed Brown if (i==0 && j==0) { /* south-west corner */ 533c4762a1bSJed Brown *sz =2; 534c4762a1bSJed Brown cols[0]=p; vals[0]=-(asb+an); 535c4762a1bSJed Brown cols[1]=n; vals[1]=an; 536c4762a1bSJed Brown } else if (i==0 && j==s->ny-1) { /* north-west corner */ 537c4762a1bSJed Brown *sz =2; 538c4762a1bSJed Brown cols[0]=s2; vals[0]=as; 539c4762a1bSJed Brown cols[1]=p; vals[1]=-(as+anb); 540c4762a1bSJed Brown } else if (i==s->nx-1 && j==0) { /* south-east corner */ 541c4762a1bSJed Brown *sz =2; 542c4762a1bSJed Brown cols[0]=p; vals[0]=-(asb+an); 543c4762a1bSJed Brown cols[1]=n; vals[1]=an; 544c4762a1bSJed Brown } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */ 545c4762a1bSJed Brown *sz =2; 546c4762a1bSJed Brown cols[0]=s2; vals[0]=as; 547c4762a1bSJed Brown cols[1]=p; vals[1]=-(as+anb); 548c4762a1bSJed Brown } else if (i==0) { /* west boundary */ 549c4762a1bSJed Brown *sz =3; 550c4762a1bSJed Brown cols[0]=s2; vals[0]=as; 551c4762a1bSJed Brown cols[1]=p; vals[1]=-(as+an); 552c4762a1bSJed Brown cols[2]=n; vals[2]=an; 553c4762a1bSJed Brown } else if (i==s->nx-1) { /* east boundary */ 554c4762a1bSJed Brown *sz =3; 555c4762a1bSJed Brown cols[0]=s2; vals[0]=as; 556c4762a1bSJed Brown cols[1]=p; vals[1]=-(as+an); 557c4762a1bSJed Brown cols[2]=n; vals[2]=an; 558c4762a1bSJed Brown } else if (j==0) { /* south boundary */ 559c4762a1bSJed Brown *sz =2; 560c4762a1bSJed Brown cols[0]=p; vals[0]=-(asb+an); 561c4762a1bSJed Brown cols[1]=n; vals[1]=an; 562c4762a1bSJed Brown } else if (j==s->ny-1) { /* north boundary */ 563c4762a1bSJed Brown *sz =2; 564c4762a1bSJed Brown cols[0]=s2; vals[0]=as; 565c4762a1bSJed Brown cols[1]=p; vals[1]=-(as+anb); 566c4762a1bSJed Brown } else { /* interior */ 567c4762a1bSJed Brown *sz =3; 568c4762a1bSJed Brown cols[0]=s2; vals[0]=as; 569c4762a1bSJed Brown cols[1]=p; vals[1]=-(as+an); 570c4762a1bSJed Brown cols[2]=n; vals[2]=an; 571c4762a1bSJed Brown } 572c4762a1bSJed Brown PetscFunctionReturn(0); 573c4762a1bSJed Brown } 574c4762a1bSJed Brown 575c4762a1bSJed Brown PetscErrorCode StokesRhsMomX(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val) 576c4762a1bSJed Brown { 577c4762a1bSJed Brown PetscScalar y = j*s->hy+s->hy/2; 578c4762a1bSJed Brown PetscScalar awb = s->hy/(s->hx/2); 579c4762a1bSJed Brown 580c4762a1bSJed Brown PetscFunctionBeginUser; 581c4762a1bSJed Brown if (i == 0) { /* west boundary */ 582c4762a1bSJed Brown *val = awb*StokesExactVelocityX(y); 583c4762a1bSJed Brown } else { 584c4762a1bSJed Brown *val = 0.0; 585c4762a1bSJed Brown } 586c4762a1bSJed Brown PetscFunctionReturn(0); 587c4762a1bSJed Brown } 588c4762a1bSJed Brown 589c4762a1bSJed Brown PetscErrorCode StokesRhsMomY(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val) 590c4762a1bSJed Brown { 591c4762a1bSJed Brown PetscFunctionBeginUser; 592c4762a1bSJed Brown *val = 0.0; 593c4762a1bSJed Brown PetscFunctionReturn(0); 594c4762a1bSJed Brown } 595c4762a1bSJed Brown 596c4762a1bSJed Brown PetscErrorCode StokesRhsMass(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val) 597c4762a1bSJed Brown { 598c4762a1bSJed Brown PetscScalar y = j*s->hy+s->hy/2; 599c4762a1bSJed Brown PetscScalar aeb = s->hy; 600c4762a1bSJed Brown 601c4762a1bSJed Brown PetscFunctionBeginUser; 602c4762a1bSJed Brown if (i == 0) { /* west boundary */ 603c4762a1bSJed Brown *val = aeb*StokesExactVelocityX(y); 604c4762a1bSJed Brown } else { 605c4762a1bSJed Brown *val = 0.0; 606c4762a1bSJed Brown } 607c4762a1bSJed Brown PetscFunctionReturn(0); 608c4762a1bSJed Brown } 609c4762a1bSJed Brown 610c4762a1bSJed Brown PetscErrorCode StokesCalcResidual(Stokes *s) 611c4762a1bSJed Brown { 612c4762a1bSJed Brown PetscReal val; 613c4762a1bSJed Brown Vec b0, b1; 614c4762a1bSJed Brown PetscErrorCode ierr; 615c4762a1bSJed Brown 616c4762a1bSJed Brown PetscFunctionBeginUser; 617c4762a1bSJed Brown /* residual Ax-b (*warning* overwrites b) */ 618c4762a1bSJed Brown ierr = VecScale(s->b, -1.0);CHKERRQ(ierr); 619c4762a1bSJed Brown ierr = MatMultAdd(s->A, s->x, s->b, s->b);CHKERRQ(ierr); 620c4762a1bSJed Brown /* ierr = VecView(s->b, (PetscViewer)PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); */ 621c4762a1bSJed Brown 622c4762a1bSJed Brown /* residual velocity */ 623c4762a1bSJed Brown ierr = VecGetSubVector(s->b, s->isg[0], &b0);CHKERRQ(ierr); 624c4762a1bSJed Brown ierr = VecNorm(b0, NORM_2, &val);CHKERRQ(ierr); 625c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," residual u = %g\n",(double)val);CHKERRQ(ierr); 626c4762a1bSJed Brown ierr = VecRestoreSubVector(s->b, s->isg[0], &b0);CHKERRQ(ierr); 627c4762a1bSJed Brown 628c4762a1bSJed Brown /* residual pressure */ 629c4762a1bSJed Brown ierr = VecGetSubVector(s->b, s->isg[1], &b1);CHKERRQ(ierr); 630c4762a1bSJed Brown ierr = VecNorm(b1, NORM_2, &val);CHKERRQ(ierr); 631c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," residual p = %g\n",(double)val);CHKERRQ(ierr); 632c4762a1bSJed Brown ierr = VecRestoreSubVector(s->b, s->isg[1], &b1);CHKERRQ(ierr); 633c4762a1bSJed Brown 634c4762a1bSJed Brown /* total residual */ 635c4762a1bSJed Brown ierr = VecNorm(s->b, NORM_2, &val);CHKERRQ(ierr); 636c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," residual [u,p] = %g\n", (double)val);CHKERRQ(ierr); 637c4762a1bSJed Brown PetscFunctionReturn(0); 638c4762a1bSJed Brown } 639c4762a1bSJed Brown 640c4762a1bSJed Brown PetscErrorCode StokesCalcError(Stokes *s) 641c4762a1bSJed Brown { 642c4762a1bSJed Brown PetscScalar scale = PetscSqrtReal((double)s->nx*s->ny); 643c4762a1bSJed Brown PetscReal val; 644c4762a1bSJed Brown Vec y0, y1; 645c4762a1bSJed Brown PetscErrorCode ierr; 646c4762a1bSJed Brown 647c4762a1bSJed Brown PetscFunctionBeginUser; 648c4762a1bSJed Brown /* error y-x */ 649c4762a1bSJed Brown ierr = VecAXPY(s->y, -1.0, s->x);CHKERRQ(ierr); 650c4762a1bSJed Brown /* ierr = VecView(s->y, (PetscViewer)PETSC_VIEWER_DEFAULT);CHKERRQ(ierr); */ 651c4762a1bSJed Brown 652c4762a1bSJed Brown /* error in velocity */ 653c4762a1bSJed Brown ierr = VecGetSubVector(s->y, s->isg[0], &y0);CHKERRQ(ierr); 654c4762a1bSJed Brown ierr = VecNorm(y0, NORM_2, &val);CHKERRQ(ierr); 655c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," discretization error u = %g\n",(double)(PetscRealPart(val/scale)));CHKERRQ(ierr); 656c4762a1bSJed Brown ierr = VecRestoreSubVector(s->y, s->isg[0], &y0);CHKERRQ(ierr); 657c4762a1bSJed Brown 658c4762a1bSJed Brown /* error in pressure */ 659c4762a1bSJed Brown ierr = VecGetSubVector(s->y, s->isg[1], &y1);CHKERRQ(ierr); 660c4762a1bSJed Brown ierr = VecNorm(y1, NORM_2, &val);CHKERRQ(ierr); 661c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," discretization error p = %g\n",(double)(PetscRealPart(val/scale)));CHKERRQ(ierr); 662c4762a1bSJed Brown ierr = VecRestoreSubVector(s->y, s->isg[1], &y1);CHKERRQ(ierr); 663c4762a1bSJed Brown 664c4762a1bSJed Brown /* total error */ 665c4762a1bSJed Brown ierr = VecNorm(s->y, NORM_2, &val);CHKERRQ(ierr); 666c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," discretization error [u,p] = %g\n", (double)PetscRealPart((val/scale)));CHKERRQ(ierr); 667c4762a1bSJed Brown PetscFunctionReturn(0); 668c4762a1bSJed Brown } 669c4762a1bSJed Brown 670c4762a1bSJed Brown int main(int argc, char **argv) 671c4762a1bSJed Brown { 672c4762a1bSJed Brown Stokes s; 673c4762a1bSJed Brown KSP ksp; 674c4762a1bSJed Brown PetscErrorCode ierr; 675c4762a1bSJed Brown 676c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 677c4762a1bSJed Brown s.nx = 4; 678c4762a1bSJed Brown s.ny = 6; 679c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL, "-nx", &s.nx, NULL);CHKERRQ(ierr); 680c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL, "-ny", &s.ny, NULL);CHKERRQ(ierr); 681c4762a1bSJed Brown s.hx = 2.0/s.nx; 682c4762a1bSJed Brown s.hy = 1.0/s.ny; 683c4762a1bSJed Brown s.matsymmetric = PETSC_FALSE; 684c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL, "-mat_set_symmetric", &s.matsymmetric,NULL);CHKERRQ(ierr); 685c4762a1bSJed Brown s.userPC = s.userKSP = PETSC_FALSE; 686c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-user_pc", &s.userPC);CHKERRQ(ierr); 687c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-user_ksp", &s.userKSP);CHKERRQ(ierr); 688c4762a1bSJed Brown 689c4762a1bSJed Brown ierr = StokesSetupMatrix(&s);CHKERRQ(ierr); 690c4762a1bSJed Brown ierr = StokesSetupIndexSets(&s);CHKERRQ(ierr); 691c4762a1bSJed Brown ierr = StokesSetupVectors(&s);CHKERRQ(ierr); 692c4762a1bSJed Brown 693c4762a1bSJed Brown ierr = KSPCreate(PETSC_COMM_WORLD, &ksp);CHKERRQ(ierr); 694c4762a1bSJed Brown ierr = KSPSetOperators(ksp, s.A, s.A);CHKERRQ(ierr); 695c4762a1bSJed Brown ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 696c4762a1bSJed Brown ierr = StokesSetupPC(&s, ksp);CHKERRQ(ierr); 697c4762a1bSJed Brown ierr = KSPSolve(ksp, s.b, s.x);CHKERRQ(ierr); 698c4762a1bSJed Brown 699c4762a1bSJed Brown /* don't trust, verify! */ 700c4762a1bSJed Brown ierr = StokesCalcResidual(&s);CHKERRQ(ierr); 701c4762a1bSJed Brown ierr = StokesCalcError(&s);CHKERRQ(ierr); 702c4762a1bSJed Brown ierr = StokesWriteSolution(&s);CHKERRQ(ierr); 703c4762a1bSJed Brown 704c4762a1bSJed Brown ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 705c4762a1bSJed Brown ierr = MatDestroy(&s.subA[0]);CHKERRQ(ierr); 706c4762a1bSJed Brown ierr = MatDestroy(&s.subA[1]);CHKERRQ(ierr); 707c4762a1bSJed Brown ierr = MatDestroy(&s.subA[2]);CHKERRQ(ierr); 708c4762a1bSJed Brown ierr = MatDestroy(&s.subA[3]);CHKERRQ(ierr); 709c4762a1bSJed Brown ierr = MatDestroy(&s.A);CHKERRQ(ierr); 710c4762a1bSJed Brown ierr = VecDestroy(&s.x);CHKERRQ(ierr); 711c4762a1bSJed Brown ierr = VecDestroy(&s.b);CHKERRQ(ierr); 712c4762a1bSJed Brown ierr = VecDestroy(&s.y);CHKERRQ(ierr); 713c4762a1bSJed Brown ierr = MatDestroy(&s.myS);CHKERRQ(ierr); 714c4762a1bSJed Brown ierr = PetscFinalize(); 715c4762a1bSJed Brown return ierr; 716c4762a1bSJed Brown } 717c4762a1bSJed Brown 718c4762a1bSJed Brown 719c4762a1bSJed Brown /*TEST 720c4762a1bSJed Brown 721c4762a1bSJed Brown test: 722c4762a1bSJed Brown nsize: 2 723c4762a1bSJed 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 724c4762a1bSJed Brown 725c4762a1bSJed Brown test: 726c4762a1bSJed Brown suffix: 2 727c4762a1bSJed Brown nsize: 2 728c4762a1bSJed Brown args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc 729c4762a1bSJed Brown 730c4762a1bSJed Brown test: 731c4762a1bSJed Brown suffix: 3 732c4762a1bSJed Brown nsize: 2 733c4762a1bSJed Brown args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc 734c4762a1bSJed Brown 735c4762a1bSJed Brown test: 736c4762a1bSJed Brown suffix: 4 737c4762a1bSJed Brown nsize: 2 738c4762a1bSJed 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 739c4762a1bSJed Brown 740c4762a1bSJed Brown test: 741c4762a1bSJed Brown suffix: 4_pcksp 742c4762a1bSJed Brown nsize: 2 743c4762a1bSJed 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 744c4762a1bSJed Brown 745c4762a1bSJed Brown test: 746c4762a1bSJed Brown suffix: 5 747c4762a1bSJed Brown nsize: 2 748c4762a1bSJed 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 749c4762a1bSJed Brown 750c4762a1bSJed Brown test: 751c4762a1bSJed Brown suffix: 6 752c4762a1bSJed Brown nsize: 2 753c4762a1bSJed 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 754c4762a1bSJed Brown 755c4762a1bSJed Brown test: 756c4762a1bSJed Brown suffix: 7 757c4762a1bSJed Brown nsize: 2 758c4762a1bSJed 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 759c4762a1bSJed Brown 760c4762a1bSJed Brown 761c4762a1bSJed Brown TEST*/ 762