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