1 static char help[] = "Stokes Problem discretized with finite elements,\n\ 2 using a parallel unstructured mesh (DMPLEX) to represent the domain.\n\n\n"; 3 4 /* 5 For the isoviscous Stokes problem, which we discretize using the finite 6 element method on an unstructured mesh, the weak form equations are 7 8 < \nabla v, \nabla u + {\nabla u}^T > - < \nabla\cdot v, p > - < v, f > = 0 9 < q, -\nabla\cdot u > = 0 10 11 Viewing: 12 13 To produce nice output, use 14 15 -dm_refine 3 -dm_view hdf5:sol1.h5 -error_vec_view hdf5:sol1.h5::append -snes_view_solution hdf5:sol1.h5::append -exact_vec_view hdf5:sol1.h5::append 16 17 You can get a LaTeX view of the mesh, with point numbering using 18 19 -dm_view :mesh.tex:ascii_latex -dm_plex_view_scale 8.0 20 21 The data layout can be viewed using 22 23 -dm_petscsection_view 24 25 Lots of information about the FEM assembly can be printed using 26 27 -dm_plex_print_fem 3 28 */ 29 30 #include <petscdmplex.h> 31 #include <petscsnes.h> 32 #include <petscds.h> 33 #include <petscbag.h> 34 35 // TODO: Plot residual by fields after each smoother iterate 36 37 typedef enum {SOL_QUADRATIC, SOL_TRIG, SOL_UNKNOWN} SolType; 38 const char *SolTypes[] = {"quadratic", "trig", "unknown", "SolType", "SOL_", 0}; 39 40 typedef struct { 41 PetscScalar mu; /* dynamic shear viscosity */ 42 } Parameter; 43 44 typedef struct { 45 PetscBag bag; /* Problem parameters */ 46 SolType sol; /* MMS solution */ 47 } AppCtx; 48 49 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 50 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 51 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 52 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 53 { 54 const PetscReal mu = PetscRealPart(constants[0]); 55 const PetscInt Nc = uOff[1]-uOff[0]; 56 PetscInt c, d; 57 58 for (c = 0; c < Nc; ++c) { 59 for (d = 0; d < dim; ++d) { 60 f1[c*dim+d] = mu * (u_x[c*dim+d] + u_x[d*dim+c]); 61 } 62 f1[c*dim+c] -= u[uOff[1]]; 63 } 64 } 65 66 static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, 67 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 68 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 69 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 70 { 71 PetscInt d; 72 for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] -= u_x[d*dim+d]; 73 } 74 75 static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 76 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 77 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 78 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 79 { 80 PetscInt d; 81 for (d = 0; d < dim; ++d) g1[d*dim+d] = -1.0; /* < q, -\nabla\cdot u > */ 82 } 83 84 static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux, 85 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 86 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 87 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 88 { 89 PetscInt d; 90 for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* -< \nabla\cdot v, p > */ 91 } 92 93 static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 94 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 95 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 96 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 97 { 98 const PetscReal mu = PetscRealPart(constants[0]); 99 const PetscInt Nc = uOff[1]-uOff[0]; 100 PetscInt c, d; 101 102 for (c = 0; c < Nc; ++c) { 103 for (d = 0; d < dim; ++d) { 104 g3[((c*Nc+c)*dim+d)*dim+d] += mu; /* < \nabla v, \nabla u > */ 105 g3[((c*Nc+d)*dim+d)*dim+c] += mu; /* < \nabla v, {\nabla u}^T > */ 106 } 107 } 108 } 109 110 static void g0_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux, 111 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 112 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 113 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 114 { 115 const PetscReal mu = PetscRealPart(constants[0]); 116 117 g0[0] = 1.0/mu; 118 } 119 120 /* Quadratic MMS Solution 121 2D: 122 123 u = x^2 + y^2 124 v = 2 x^2 - 2xy 125 p = x + y - 1 126 f = <1 - 4 mu, 1 - 4 mu> 127 128 so that 129 130 e(u) = (grad u + grad u^T) = / 4x 4x \ 131 \ 4x -4x / 132 div mu e(u) - \nabla p + f = mu <4, 4> - <1, 1> + <1 - 4 mu, 1 - 4 mu> = 0 133 \nabla \cdot u = 2x - 2x = 0 134 135 3D: 136 137 u = 2 x^2 + y^2 + z^2 138 v = 2 x^2 - 2xy 139 w = 2 x^2 - 2xz 140 p = x + y + z - 3/2 141 f = <1 - 8 mu, 1 - 4 mu, 1 - 4 mu> 142 143 so that 144 145 e(u) = (grad u + grad u^T) = / 8x 4x 4x \ 146 | 4x -4x 0 | 147 \ 4x 0 -4x / 148 div mu e(u) - \nabla p + f = mu <8, 4, 4> - <1, 1, 1> + <1 - 8 mu, 1 - 4 mu, 1 - 4 mu> = 0 149 \nabla \cdot u = 4x - 2x - 2x = 0 150 */ 151 static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 152 { 153 PetscInt c; 154 155 u[0] = (dim-1)*PetscSqr(x[0]); 156 for (c = 1; c < Nc; ++c) { 157 u[0] += PetscSqr(x[c]); 158 u[c] = 2.0*PetscSqr(x[0]) - 2.0*x[0]*x[c]; 159 } 160 return 0; 161 } 162 163 static PetscErrorCode quadratic_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 164 { 165 PetscInt d; 166 167 u[0] = -0.5*dim; 168 for (d = 0; d < dim; ++d) u[0] += x[d]; 169 return 0; 170 } 171 172 static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 173 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 174 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 175 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 176 { 177 const PetscReal mu = PetscRealPart(constants[0]); 178 PetscInt d; 179 180 f0[0] = (dim-1)*4.0*mu - 1.0; 181 for (d = 1; d < dim; ++d) f0[d] = 4.0*mu - 1.0; 182 } 183 184 /* Trigonometric MMS Solution 185 2D: 186 187 u = sin(pi x) + sin(pi y) 188 v = -pi cos(pi x) y 189 p = sin(2 pi x) + sin(2 pi y) 190 f = <2pi cos(2 pi x) + mu pi^2 sin(pi x) + mu pi^2 sin(pi y), 2pi cos(2 pi y) - mu pi^3 cos(pi x) y> 191 192 so that 193 194 e(u) = (grad u + grad u^T) = / 2pi cos(pi x) pi cos(pi y) + pi^2 sin(pi x) y \ 195 \ pi cos(pi y) + pi^2 sin(pi x) y -2pi cos(pi x) / 196 div mu e(u) - \nabla p + f = mu <-pi^2 sin(pi x) - pi^2 sin(pi y), pi^3 cos(pi x) y> - <2pi cos(2 pi x), 2pi cos(2 pi y)> + <f_x, f_y> = 0 197 \nabla \cdot u = pi cos(pi x) - pi cos(pi x) = 0 198 199 3D: 200 201 u = 2 sin(pi x) + sin(pi y) + sin(pi z) 202 v = -pi cos(pi x) y 203 w = -pi cos(pi x) z 204 p = sin(2 pi x) + sin(2 pi y) + sin(2 pi z) 205 f = <2pi cos(2 pi x) + mu 2pi^2 sin(pi x) + mu pi^2 sin(pi y) + mu pi^2 sin(pi z), 2pi cos(2 pi y) - mu pi^3 cos(pi x) y, 2pi cos(2 pi z) - mu pi^3 cos(pi x) z> 206 207 so that 208 209 e(u) = (grad u + grad u^T) = / 4pi cos(pi x) pi cos(pi y) + pi^2 sin(pi x) y pi cos(pi z) + pi^2 sin(pi x) z \ 210 | pi cos(pi y) + pi^2 sin(pi x) y -2pi cos(pi x) 0 | 211 \ pi cos(pi z) + pi^2 sin(pi x) z 0 -2pi cos(pi x) / 212 div mu e(u) - \nabla p + f = mu <-2pi^2 sin(pi x) - pi^2 sin(pi y) - pi^2 sin(pi z), pi^3 cos(pi x) y, pi^3 cos(pi x) z> - <2pi cos(2 pi x), 2pi cos(2 pi y), 2pi cos(2 pi z)> + <f_x, f_y, f_z> = 0 213 \nabla \cdot u = 2 pi cos(pi x) - pi cos(pi x) - pi cos(pi x) = 0 214 */ 215 static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 216 { 217 PetscInt c; 218 219 u[0] = (dim-1)*PetscSinReal(PETSC_PI*x[0]); 220 for (c = 1; c < Nc; ++c) { 221 u[0] += PetscSinReal(PETSC_PI*x[c]); 222 u[c] = -PETSC_PI*PetscCosReal(PETSC_PI*x[0]) * x[c]; 223 } 224 return 0; 225 } 226 227 static PetscErrorCode trig_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 228 { 229 PetscInt d; 230 231 for (d = 0, u[0] = 0.0; d < dim; ++d) u[0] += PetscSinReal(2.0*PETSC_PI*x[d]); 232 return 0; 233 } 234 235 static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 236 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 237 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 238 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 239 { 240 const PetscReal mu = PetscRealPart(constants[0]); 241 PetscInt d; 242 243 f0[0] = -2.0*PETSC_PI*PetscCosReal(2.0*PETSC_PI*x[0]) - (dim-1)*mu*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x[0]); 244 for (d = 1; d < dim; ++d) { 245 f0[0] -= mu*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x[d]); 246 f0[d] = -2.0*PETSC_PI*PetscCosReal(2.0*PETSC_PI*x[d]) + mu*PetscPowRealInt(PETSC_PI, 3)*PetscCosReal(PETSC_PI*x[0])*x[d]; 247 } 248 } 249 250 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 251 { 252 PetscInt sol; 253 PetscErrorCode ierr; 254 255 PetscFunctionBeginUser; 256 options->sol = SOL_QUADRATIC; 257 258 ierr = PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");CHKERRQ(ierr); 259 sol = options->sol; 260 CHKERRQ(PetscOptionsEList("-sol", "The MMS solution", "ex62.c", SolTypes, (sizeof(SolTypes)/sizeof(SolTypes[0]))-3, SolTypes[options->sol], &sol, NULL)); 261 options->sol = (SolType) sol; 262 ierr = PetscOptionsEnd();CHKERRQ(ierr); 263 PetscFunctionReturn(0); 264 } 265 266 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 267 { 268 PetscFunctionBeginUser; 269 CHKERRQ(DMCreate(comm, dm)); 270 CHKERRQ(DMSetType(*dm, DMPLEX)); 271 CHKERRQ(DMSetFromOptions(*dm)); 272 CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 273 PetscFunctionReturn(0); 274 } 275 276 static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx) 277 { 278 Parameter *p; 279 280 PetscFunctionBeginUser; 281 /* setup PETSc parameter bag */ 282 CHKERRQ(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &ctx->bag)); 283 CHKERRQ(PetscBagGetData(ctx->bag, (void **) &p)); 284 CHKERRQ(PetscBagSetName(ctx->bag, "par", "Stokes Parameters")); 285 CHKERRQ(PetscBagRegisterScalar(ctx->bag, &p->mu, 1.0, "mu", "Dynamic Shear Viscosity, Pa s")); 286 CHKERRQ(PetscBagSetFromOptions(ctx->bag)); 287 { 288 PetscViewer viewer; 289 PetscViewerFormat format; 290 PetscBool flg; 291 292 CHKERRQ(PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg)); 293 if (flg) { 294 CHKERRQ(PetscViewerPushFormat(viewer, format)); 295 CHKERRQ(PetscBagView(ctx->bag, viewer)); 296 CHKERRQ(PetscViewerFlush(viewer)); 297 CHKERRQ(PetscViewerPopFormat(viewer)); 298 CHKERRQ(PetscViewerDestroy(&viewer)); 299 } 300 } 301 PetscFunctionReturn(0); 302 } 303 304 static PetscErrorCode SetupEqn(DM dm, AppCtx *user) 305 { 306 PetscErrorCode (*exactFuncs[2])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 307 PetscDS ds; 308 DMLabel label; 309 const PetscInt id = 1; 310 311 PetscFunctionBeginUser; 312 CHKERRQ(DMGetDS(dm, &ds)); 313 switch (user->sol) { 314 case SOL_QUADRATIC: 315 CHKERRQ(PetscDSSetResidual(ds, 0, f0_quadratic_u, f1_u)); 316 exactFuncs[0] = quadratic_u; 317 exactFuncs[1] = quadratic_p; 318 break; 319 case SOL_TRIG: 320 CHKERRQ(PetscDSSetResidual(ds, 0, f0_trig_u, f1_u)); 321 exactFuncs[0] = trig_u; 322 exactFuncs[1] = trig_p; 323 break; 324 default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", SolTypes[PetscMin(user->sol, SOL_UNKNOWN)], user->sol); 325 } 326 CHKERRQ(PetscDSSetResidual(ds, 1, f0_p, NULL)); 327 CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 328 CHKERRQ(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL)); 329 CHKERRQ(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL)); 330 CHKERRQ(PetscDSSetJacobianPreconditioner(ds, 0, 0, NULL, NULL, NULL, g3_uu)); 331 CHKERRQ(PetscDSSetJacobianPreconditioner(ds, 1, 1, g0_pp, NULL, NULL, NULL)); 332 333 CHKERRQ(PetscDSSetExactSolution(ds, 0, exactFuncs[0], user)); 334 CHKERRQ(PetscDSSetExactSolution(ds, 1, exactFuncs[1], user)); 335 336 CHKERRQ(DMGetLabel(dm, "marker", &label)); 337 CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, user, NULL)); 338 339 /* Make constant values available to pointwise functions */ 340 { 341 Parameter *param; 342 PetscScalar constants[1]; 343 344 CHKERRQ(PetscBagGetData(user->bag, (void **) ¶m)); 345 constants[0] = param->mu; /* dynamic shear viscosity, Pa s */ 346 CHKERRQ(PetscDSSetConstants(ds, 1, constants)); 347 } 348 PetscFunctionReturn(0); 349 } 350 351 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 352 { 353 PetscInt c; 354 for (c = 0; c < Nc; ++c) u[c] = 0.0; 355 return 0; 356 } 357 static PetscErrorCode one(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 358 { 359 PetscInt c; 360 for (c = 0; c < Nc; ++c) u[c] = 1.0; 361 return 0; 362 } 363 364 static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace) 365 { 366 Vec vec; 367 PetscErrorCode (*funcs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx) = {zero, one}; 368 369 PetscFunctionBeginUser; 370 PetscCheckFalse(origField != 1,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Field %D should be 1 for pressure", origField); 371 funcs[field] = one; 372 { 373 PetscDS ds; 374 CHKERRQ(DMGetDS(dm, &ds)); 375 CHKERRQ(PetscObjectViewFromOptions((PetscObject) ds, NULL, "-ds_view")); 376 } 377 CHKERRQ(DMCreateGlobalVector(dm, &vec)); 378 CHKERRQ(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec)); 379 CHKERRQ(VecNormalize(vec, NULL)); 380 CHKERRQ(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullspace)); 381 CHKERRQ(VecDestroy(&vec)); 382 /* New style for field null spaces */ 383 { 384 PetscObject pressure; 385 MatNullSpace nullspacePres; 386 387 CHKERRQ(DMGetField(dm, field, NULL, &pressure)); 388 CHKERRQ(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres)); 389 CHKERRQ(PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres)); 390 CHKERRQ(MatNullSpaceDestroy(&nullspacePres)); 391 } 392 PetscFunctionReturn(0); 393 } 394 395 static PetscErrorCode SetupProblem(DM dm, PetscErrorCode (*setupEqn)(DM, AppCtx *), AppCtx *user) 396 { 397 DM cdm = dm; 398 PetscQuadrature q = NULL; 399 PetscBool simplex; 400 PetscInt dim, Nf = 2, f, Nc[2]; 401 const char *name[2] = {"velocity", "pressure"}; 402 const char *prefix[2] = {"vel_", "pres_"}; 403 404 PetscFunctionBegin; 405 CHKERRQ(DMGetDimension(dm, &dim)); 406 CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 407 Nc[0] = dim; 408 Nc[1] = 1; 409 for (f = 0; f < Nf; ++f) { 410 PetscFE fe; 411 412 CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc[f], simplex, prefix[f], -1, &fe)); 413 CHKERRQ(PetscObjectSetName((PetscObject) fe, name[f])); 414 if (!q) CHKERRQ(PetscFEGetQuadrature(fe, &q)); 415 CHKERRQ(PetscFESetQuadrature(fe, q)); 416 CHKERRQ(DMSetField(dm, f, NULL, (PetscObject) fe)); 417 CHKERRQ(PetscFEDestroy(&fe)); 418 } 419 CHKERRQ(DMCreateDS(dm)); 420 CHKERRQ((*setupEqn)(dm, user)); 421 while (cdm) { 422 CHKERRQ(DMCopyDisc(dm, cdm)); 423 CHKERRQ(DMSetNullSpaceConstructor(cdm, 1, CreatePressureNullSpace)); 424 CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 425 } 426 PetscFunctionReturn(0); 427 } 428 429 int main(int argc, char **argv) 430 { 431 SNES snes; 432 DM dm; 433 Vec u; 434 AppCtx user; 435 PetscErrorCode ierr; 436 437 ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 438 CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 439 CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 440 CHKERRQ(SNESCreate(PetscObjectComm((PetscObject) dm), &snes)); 441 CHKERRQ(SNESSetDM(snes, dm)); 442 CHKERRQ(DMSetApplicationContext(dm, &user)); 443 444 CHKERRQ(SetupParameters(PETSC_COMM_WORLD, &user)); 445 CHKERRQ(SetupProblem(dm, SetupEqn, &user)); 446 CHKERRQ(DMPlexCreateClosureIndex(dm, NULL)); 447 448 CHKERRQ(DMCreateGlobalVector(dm, &u)); 449 CHKERRQ(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 450 CHKERRQ(SNESSetFromOptions(snes)); 451 CHKERRQ(DMSNESCheckFromOptions(snes, u)); 452 CHKERRQ(PetscObjectSetName((PetscObject) u, "Solution")); 453 { 454 Mat J; 455 MatNullSpace sp; 456 457 CHKERRQ(SNESSetUp(snes)); 458 CHKERRQ(CreatePressureNullSpace(dm, 1, 1, &sp)); 459 CHKERRQ(SNESGetJacobian(snes, &J, NULL, NULL, NULL)); 460 CHKERRQ(MatSetNullSpace(J, sp)); 461 CHKERRQ(MatNullSpaceDestroy(&sp)); 462 CHKERRQ(PetscObjectSetName((PetscObject) J, "Jacobian")); 463 CHKERRQ(MatViewFromOptions(J, NULL, "-J_view")); 464 } 465 CHKERRQ(SNESSolve(snes, NULL, u)); 466 467 CHKERRQ(VecDestroy(&u)); 468 CHKERRQ(SNESDestroy(&snes)); 469 CHKERRQ(DMDestroy(&dm)); 470 CHKERRQ(PetscBagDestroy(&user.bag)); 471 ierr = PetscFinalize(); 472 return ierr; 473 } 474 /*TEST 475 476 test: 477 suffix: 2d_p2_p1_check 478 requires: triangle 479 args: -sol quadratic -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 480 481 test: 482 suffix: 2d_p2_p1_check_parallel 483 nsize: {{2 3 5}} 484 requires: triangle 485 args: -sol quadratic -dm_refine 2 -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 486 487 test: 488 suffix: 3d_p2_p1_check 489 requires: ctetgen 490 args: -sol quadratic -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 491 492 test: 493 suffix: 3d_p2_p1_check_parallel 494 nsize: {{2 3 5}} 495 requires: ctetgen 496 args: -sol quadratic -dm_refine 2 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 497 498 test: 499 suffix: 2d_p2_p1_conv 500 requires: triangle 501 # Using -dm_refine 3 gives L_2 convergence rate: [3.0, 2.1] 502 args: -sol trig -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 -ksp_error_if_not_converged \ 503 -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 504 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \ 505 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 506 507 test: 508 suffix: 2d_p2_p1_conv_gamg 509 requires: triangle 510 args: -sol trig -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 \ 511 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \ 512 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_coarse_pc_type svd 513 514 test: 515 suffix: 3d_p2_p1_conv 516 requires: ctetgen !single 517 # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.8] 518 args: -sol trig -dm_plex_dim 3 -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 \ 519 -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 520 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \ 521 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 522 523 test: 524 suffix: 2d_q2_q1_check 525 args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 526 527 test: 528 suffix: 3d_q2_q1_check 529 args: -sol quadratic -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 530 531 test: 532 suffix: 2d_q2_q1_conv 533 # Using -dm_refine 3 -convest_num_refine 1 gives L_2 convergence rate: [3.0, 2.1] 534 args: -sol trig -dm_plex_simplex 0 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -ksp_error_if_not_converged \ 535 -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 536 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \ 537 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 538 539 test: 540 suffix: 3d_q2_q1_conv 541 requires: !single 542 # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.4] 543 args: -sol trig -dm_plex_simplex 0 -dm_plex_dim 3 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 \ 544 -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 545 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \ 546 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 547 548 test: 549 suffix: 2d_p3_p2_check 550 requires: triangle 551 args: -sol quadratic -vel_petscspace_degree 3 -pres_petscspace_degree 2 -dmsnes_check 0.0001 552 553 test: 554 suffix: 3d_p3_p2_check 555 requires: ctetgen !single 556 args: -sol quadratic -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -vel_petscspace_degree 3 -pres_petscspace_degree 2 -dmsnes_check 0.0001 557 558 test: 559 suffix: 2d_p3_p2_conv 560 requires: triangle 561 # Using -dm_refine 2 gives L_2 convergence rate: [3.8, 3.0] 562 args: -sol trig -vel_petscspace_degree 3 -pres_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 -ksp_error_if_not_converged \ 563 -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 564 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \ 565 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 566 567 test: 568 suffix: 3d_p3_p2_conv 569 requires: ctetgen long_runtime 570 # Using -dm_refine 1 -convest_num_refine 2 gives L_2 convergence rate: [3.6, 3.9] 571 args: -sol trig -dm_plex_dim 3 -dm_refine 1 -vel_petscspace_degree 3 -pres_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 \ 572 -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \ 573 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \ 574 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu 575 576 test: 577 suffix: 2d_q1_p0_conv 578 requires: !single 579 # Using -dm_refine 3 gives L_2 convergence rate: [1.9, 1.0] 580 args: -sol quadratic -dm_plex_simplex 0 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -snes_convergence_estimate -convest_num_refine 2 \ 581 -ksp_atol 1e-10 -petscds_jac_pre 0 \ 582 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \ 583 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_levels_pc_type jacobi -fieldsplit_pressure_mg_coarse_pc_type svd 584 585 test: 586 suffix: 3d_q1_p0_conv 587 requires: !single 588 # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [1.7, 1.0] 589 args: -sol quadratic -dm_plex_simplex 0 -dm_plex_dim 3 -dm_refine 1 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -snes_convergence_estimate -convest_num_refine 1 \ 590 -ksp_atol 1e-10 -petscds_jac_pre 0 \ 591 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \ 592 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_levels_pc_type jacobi -fieldsplit_pressure_mg_coarse_pc_type svd 593 594 # Stokes preconditioners 595 # Block diagonal \begin{pmatrix} A & 0 \\ 0 & I \end{pmatrix} 596 test: 597 suffix: 2d_p2_p1_block_diagonal 598 requires: triangle 599 args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 600 -snes_error_if_not_converged \ 601 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-4 -ksp_error_if_not_converged \ 602 -pc_type fieldsplit -pc_fieldsplit_type additive -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi 603 # Block triangular \begin{pmatrix} A & B \\ 0 & I \end{pmatrix} 604 test: 605 suffix: 2d_p2_p1_block_triangular 606 requires: triangle 607 args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 608 -snes_error_if_not_converged \ 609 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 610 -pc_type fieldsplit -pc_fieldsplit_type multiplicative -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi 611 # Diagonal Schur complement \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix} 612 test: 613 suffix: 2d_p2_p1_schur_diagonal 614 requires: triangle 615 args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 616 -snes_error_if_not_converged \ 617 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \ 618 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type diag -pc_fieldsplit_off_diag_use_amat \ 619 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 620 # Upper triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix} 621 test: 622 suffix: 2d_p2_p1_schur_upper 623 requires: triangle 624 args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 \ 625 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \ 626 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -pc_fieldsplit_off_diag_use_amat \ 627 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 628 # Lower triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix} 629 test: 630 suffix: 2d_p2_p1_schur_lower 631 requires: triangle 632 args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 633 -snes_error_if_not_converged \ 634 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \ 635 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type lower -pc_fieldsplit_off_diag_use_amat \ 636 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 637 # Full Schur complement \begin{pmatrix} I & 0 \\ B^T A^{-1} & I \end{pmatrix} \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix} \begin{pmatrix} I & A^{-1} B \\ 0 & I \end{pmatrix} 638 test: 639 suffix: 2d_p2_p1_schur_full 640 requires: triangle 641 args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 642 -snes_error_if_not_converged \ 643 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \ 644 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_off_diag_use_amat \ 645 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 646 # Full Schur + Velocity GMG 647 test: 648 suffix: 2d_p2_p1_gmg_vcycle 649 requires: triangle 650 args: -sol quadratic -dm_refine_hierarchy 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \ 651 -ksp_type fgmres -ksp_atol 1e-9 -snes_error_if_not_converged -pc_use_amat \ 652 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_off_diag_use_amat \ 653 -fieldsplit_velocity_pc_type mg -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_pc_gamg_esteig_ksp_max_it 10 -fieldsplit_pressure_mg_levels_pc_type sor -fieldsplit_pressure_mg_coarse_pc_type svd 654 # SIMPLE \begin{pmatrix} I & 0 \\ B^T A^{-1} & I \end{pmatrix} \begin{pmatrix} A & 0 \\ 0 & B^T diag(A)^{-1} B \end{pmatrix} \begin{pmatrix} I & diag(A)^{-1} B \\ 0 & I \end{pmatrix} 655 test: 656 suffix: 2d_p2_p1_simple 657 requires: triangle 658 args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 659 -snes_error_if_not_converged \ 660 -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \ 661 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \ 662 -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \ 663 -fieldsplit_pressure_inner_ksp_type preonly -fieldsplit_pressure_inner_pc_type jacobi -fieldsplit_pressure_upper_ksp_type preonly -fieldsplit_pressure_upper_pc_type jacobi 664 # FETI-DP solvers (these solvers are quite inefficient, they are here to exercise the code) 665 test: 666 suffix: 2d_p2_p1_fetidp 667 requires: triangle mumps 668 nsize: 5 669 args: -sol quadratic -dm_refine 2 -dm_mat_type is -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 670 -snes_error_if_not_converged \ 671 -ksp_type fetidp -ksp_rtol 1.0e-8 \ 672 -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \ 673 -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 200 -fetidp_fieldsplit_p_pc_type none \ 674 -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type mumps -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type mumps -fetidp_fieldsplit_lag_ksp_type preonly 675 test: 676 suffix: 2d_q2_q1_fetidp 677 requires: mumps 678 nsize: 5 679 args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -dm_mat_type is -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 680 -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \ 681 -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \ 682 -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 200 -fetidp_fieldsplit_p_pc_type none \ 683 -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type mumps -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type mumps -fetidp_fieldsplit_lag_ksp_type preonly 684 test: 685 suffix: 3d_p2_p1_fetidp 686 requires: ctetgen mumps suitesparse 687 nsize: 5 688 args: -sol quadratic -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_refine 1 -dm_mat_type is -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 689 -snes_error_if_not_converged \ 690 -ksp_type fetidp -ksp_rtol 1.0e-9 \ 691 -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \ 692 -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 1000 -fetidp_fieldsplit_p_pc_type none \ 693 -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_benign_trick -fetidp_bddc_pc_bddc_deluxe_singlemat \ 694 -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \ 695 -fetidp_bddelta_pc_factor_mat_solver_type umfpack -fetidp_fieldsplit_lag_ksp_type preonly -fetidp_bddc_sub_schurs_mat_solver_type mumps -fetidp_bddc_sub_schurs_mat_mumps_icntl_14 100000 \ 696 -fetidp_bddelta_pc_factor_mat_ordering_type external \ 697 -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \ 698 -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external 699 test: 700 suffix: 3d_q2_q1_fetidp 701 requires: suitesparse 702 nsize: 5 703 args: -sol quadratic -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_refine 1 -dm_mat_type is -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 704 -snes_error_if_not_converged \ 705 -ksp_type fetidp -ksp_rtol 1.0e-8 \ 706 -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \ 707 -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 2000 -fetidp_fieldsplit_p_pc_type none \ 708 -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \ 709 -fetidp_bddc_pc_bddc_symmetric -fetidp_fieldsplit_lag_ksp_type preonly \ 710 -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \ 711 -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external 712 # BDDC solvers (these solvers are quite inefficient, they are here to exercise the code) 713 test: 714 suffix: 2d_p2_p1_bddc 715 nsize: 2 716 requires: triangle !single 717 args: -sol quadratic -dm_plex_box_faces 2,2,2 -dm_refine 1 -dm_mat_type is -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \ 718 -snes_error_if_not_converged \ 719 -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \ 720 -pc_type bddc -pc_bddc_corner_selection -pc_bddc_dirichlet_pc_type svd -pc_bddc_neumann_pc_type svd -pc_bddc_coarse_redundant_pc_type svd 721 # Vanka 722 test: 723 suffix: 2d_q1_p0_vanka 724 requires: double !complex 725 args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \ 726 -snes_rtol 1.0e-4 \ 727 -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \ 728 -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \ 729 -sub_ksp_type preonly -sub_pc_type lu 730 test: 731 suffix: 2d_q1_p0_vanka_denseinv 732 requires: double !complex 733 args: -sol quadratic -dm_plex_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \ 734 -snes_rtol 1.0e-4 \ 735 -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \ 736 -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \ 737 -pc_patch_dense_inverse -pc_patch_sub_mat_type seqdense 738 # Vanka smoother 739 test: 740 suffix: 2d_q1_p0_gmg_vanka 741 requires: double !complex 742 args: -sol quadratic -dm_plex_simplex 0 -dm_refine_hierarchy 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \ 743 -snes_rtol 1.0e-4 \ 744 -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \ 745 -pc_type mg \ 746 -mg_levels_ksp_type gmres -mg_levels_ksp_max_it 30 \ 747 -mg_levels_pc_type patch -mg_levels_pc_patch_partition_of_unity 0 -mg_levels_pc_patch_construct_codim 0 -mg_levels_pc_patch_construct_type vanka \ 748 -mg_levels_sub_ksp_type preonly -mg_levels_sub_pc_type lu \ 749 -mg_coarse_pc_type svd 750 751 TEST*/ 752