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