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