1 static char help[] = "Poisson Problem in mixed form with 2d and 3d with finite elements.\n\ 2 We solve the Poisson problem in a rectangular\n\ 3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4 This example supports automatic convergence estimation\n\ 5 and Hdiv elements.\n\n\n"; 6 7 /* 8 9 The mixed form of Poisson's equation, e.g. https://firedrakeproject.org/demos/poisson_mixed.py.html, is given 10 in the strong form by 11 \begin{align} 12 \vb{q} - \nabla u &= 0 \\ 13 \nabla \cdot \vb{q} &= f 14 \end{align} 15 where $u$ is the potential, as in the original problem, but we also discretize the gradient of potential $\vb{q}$, 16 or flux, directly. The weak form is then 17 \begin{align} 18 <t, \vb{q}> + <\nabla \cdot t, u> - <t_n, u>_\Gamma &= 0 \\ 19 <v, \nabla \cdot \vb{q}> &= <v, f> 20 \end{align} 21 For the original Poisson problem, the Dirichlet boundary forces an essential boundary condition on the potential space, 22 and the Neumann boundary gives a natural boundary condition in the weak form. In the mixed formulation, the Neumann 23 boundary gives an essential boundary condition on the flux space, $\vb{q} \cdot \vu{n} = h$, and the Dirichlet condition 24 becomes a natural condition in the weak form, <t_n, g>_\Gamma. 25 */ 26 27 #include <petscdmplex.h> 28 #include <petscsnes.h> 29 #include <petscds.h> 30 #include <petscconvest.h> 31 32 typedef enum {SOL_LINEAR, SOL_QUADRATIC, SOL_QUARTIC, SOL_QUARTIC_NEUMANN, SOL_UNKNOWN, NUM_SOLTYPE} SolType; 33 const char *SolTypeNames[NUM_SOLTYPE+3] = {"linear", "quadratic", "quartic", "quartic_neumann", "unknown", "SolType", "SOL_", NULL}; 34 35 typedef struct { 36 SolType solType; /* The type of exact solution */ 37 } AppCtx; 38 39 static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 40 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 41 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 42 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 43 { 44 PetscInt d; 45 for (d = 0; d < dim; ++d) f0[0] += u_x[uOff_x[0]+d*dim+d]; 46 } 47 48 /* 2D Dirichlet potential example 49 50 u = x 51 q = <1, 0> 52 f = 0 53 54 We will need a boundary integral of u over \Gamma. 55 */ 56 static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 57 { 58 u[0] = x[0]; 59 return 0; 60 } 61 static PetscErrorCode linear_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 62 { 63 PetscInt c; 64 for (c = 0; c < Nc; ++c) u[c] = c ? 0.0 : 1.0; 65 return 0; 66 } 67 68 static void f0_linear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 69 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 70 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 71 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 72 { 73 f0[0] = 0.0; 74 f0_u(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, f0); 75 } 76 static void f0_bd_linear_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 77 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 78 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 79 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 80 { 81 PetscScalar potential; 82 PetscInt d; 83 84 linear_u(dim, t, x, dim, &potential, NULL); 85 for (d = 0; d < dim; ++d) f0[d] = -potential*n[d]; 86 } 87 88 /* 2D Dirichlet potential example 89 90 u = x^2 + y^2 91 q = <2x, 2y> 92 f = 4 93 94 We will need a boundary integral of u over \Gamma. 95 */ 96 static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 97 { 98 PetscInt d; 99 100 u[0] = 0.0; 101 for (d = 0; d < dim; ++d) u[0] += x[d]*x[d]; 102 return 0; 103 } 104 static PetscErrorCode quadratic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 105 { 106 PetscInt c; 107 for (c = 0; c < Nc; ++c) u[c] = 2.0*x[c]; 108 return 0; 109 } 110 111 static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 112 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 113 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 114 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 115 { 116 f0[0] = -4.0; 117 f0_u(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, f0); 118 } 119 static void f0_bd_quadratic_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 120 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 121 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 122 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 123 { 124 PetscScalar potential; 125 PetscInt d; 126 127 quadratic_u(dim, t, x, dim, &potential, NULL); 128 for (d = 0; d < dim; ++d) f0[d] = -potential*n[d]; 129 } 130 131 /* 2D Dirichlet potential example 132 133 u = x (1-x) y (1-y) 134 q = <(1-2x) y (1-y), x (1-x) (1-2y)> 135 f = -y (1-y) - x (1-x) 136 137 u|_\Gamma = 0 so that the boundary integral vanishes 138 */ 139 static PetscErrorCode quartic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 140 { 141 PetscInt d; 142 143 u[0] = 1.0; 144 for (d = 0; d < dim; ++d) u[0] *= x[d]*(1.0 - x[d]); 145 return 0; 146 } 147 static PetscErrorCode quartic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 148 { 149 PetscInt c, d; 150 151 for (c = 0; c < Nc; ++c) { 152 u[c] = 1.0; 153 for (d = 0; d < dim; ++d) { 154 if (c == d) u[c] *= 1 - 2.0*x[d]; 155 else u[c] *= x[d]*(1.0 - x[d]); 156 } 157 } 158 return 0; 159 } 160 161 static void f0_quartic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 162 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 163 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 164 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 165 { 166 PetscInt d; 167 f0[0] = 0.0; 168 for (d = 0; d < dim; ++d) f0[0] += 2.0*x[d]*(1.0 - x[d]); 169 f0_u(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, f0); 170 } 171 172 /* 2D Dirichlet potential example 173 174 u = x (1-x) y (1-y) 175 q = <(1-2x) y (1-y), x (1-x) (1-2y)> 176 f = -y (1-y) - x (1-x) 177 178 du/dn_\Gamma = {(1-2x) y (1-y), x (1-x) (1-2y)} . n produces an essential condition on q 179 */ 180 181 static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 182 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 183 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 184 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 185 { 186 PetscInt c; 187 for (c = 0; c < dim; ++c) f0[c] = u[uOff[0]+c]; 188 } 189 190 /* <\nabla\cdot w, u> = <\nabla w, Iu> */ 191 static void f1_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 192 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 193 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 194 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 195 { 196 PetscInt c; 197 for (c = 0; c < dim; ++c) f1[c*dim+c] = u[uOff[1]]; 198 } 199 200 /* <t, q> */ 201 static void g0_qq(PetscInt dim, PetscInt Nf, PetscInt NfAux, 202 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 203 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 204 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 205 { 206 PetscInt c; 207 for (c = 0; c < dim; ++c) g0[c*dim+c] = 1.0; 208 } 209 210 /* <\nabla\cdot t, u> = <\nabla t, Iu> */ 211 static void g2_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 212 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 213 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 214 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 215 { 216 PetscInt d; 217 for (d = 0; d < dim; ++d) g2[d*dim+d] = 1.0; 218 } 219 220 /* <v, \nabla\cdot q> */ 221 static void g1_uq(PetscInt dim, PetscInt Nf, PetscInt NfAux, 222 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 223 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 224 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 225 { 226 PetscInt d; 227 for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; 228 } 229 230 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 231 { 232 PetscFunctionBeginUser; 233 options->solType = SOL_LINEAR; 234 PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX"); 235 PetscCall(PetscOptionsEnum("-sol_type", "Type of exact solution", "ex24.c", SolTypeNames, (PetscEnum) options->solType, (PetscEnum *) &options->solType, NULL)); 236 PetscOptionsEnd(); 237 PetscFunctionReturn(0); 238 } 239 240 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 241 { 242 PetscFunctionBeginUser; 243 PetscCall(DMCreate(comm, dm)); 244 PetscCall(DMSetType(*dm, DMPLEX)); 245 PetscCall(PetscObjectSetName((PetscObject)*dm, "Example Mesh")); 246 PetscCall(DMSetApplicationContext(*dm, user)); 247 PetscCall(DMSetFromOptions(*dm)); 248 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 249 PetscFunctionReturn(0); 250 } 251 252 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 253 { 254 PetscDS ds; 255 DMLabel label; 256 PetscWeakForm wf; 257 const PetscInt id = 1; 258 PetscInt bd; 259 260 PetscFunctionBeginUser; 261 PetscCall(DMGetLabel(dm, "marker", &label)); 262 PetscCall(DMGetDS(dm, &ds)); 263 PetscCall(PetscDSSetResidual(ds, 0, f0_q, f1_q)); 264 PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_qq, NULL, NULL, NULL)); 265 PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_qu, NULL)); 266 PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_uq, NULL, NULL)); 267 switch (user->solType) 268 { 269 case SOL_LINEAR: 270 PetscCall(PetscDSSetResidual(ds, 1, f0_linear_u, NULL)); 271 PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 272 PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 273 PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_linear_q, 0, NULL)); 274 PetscCall(PetscDSSetExactSolution(ds, 0, linear_q, user)); 275 PetscCall(PetscDSSetExactSolution(ds, 1, linear_u, user)); 276 break; 277 case SOL_QUADRATIC: 278 PetscCall(PetscDSSetResidual(ds, 1, f0_quadratic_u, NULL)); 279 PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 280 PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 281 PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_quadratic_q, 0, NULL)); 282 PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_q, user)); 283 PetscCall(PetscDSSetExactSolution(ds, 1, quadratic_u, user)); 284 break; 285 case SOL_QUARTIC: 286 PetscCall(PetscDSSetResidual(ds, 1, f0_quartic_u, NULL)); 287 PetscCall(PetscDSSetExactSolution(ds, 0, quartic_q, user)); 288 PetscCall(PetscDSSetExactSolution(ds, 1, quartic_u, user)); 289 break; 290 case SOL_QUARTIC_NEUMANN: 291 PetscCall(PetscDSSetResidual(ds, 1, f0_quartic_u, NULL)); 292 PetscCall(PetscDSSetExactSolution(ds, 0, quartic_q, user)); 293 PetscCall(PetscDSSetExactSolution(ds, 1, quartic_u, user)); 294 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "Flux condition", label, 1, &id, 0, 0, NULL, (void (*)(void)) quartic_q, NULL, user, NULL)); 295 break; 296 default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid exact solution type %s", SolTypeNames[PetscMin(user->solType, SOL_UNKNOWN)]); 297 } 298 PetscFunctionReturn(0); 299 } 300 301 static PetscErrorCode SetupDiscretization(DM dm, PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 302 { 303 DM cdm = dm; 304 PetscFE feq, feu; 305 DMPolytopeType ct; 306 PetscBool simplex; 307 PetscInt dim, cStart; 308 309 PetscFunctionBeginUser; 310 PetscCall(DMGetDimension(dm, &dim)); 311 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 312 PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 313 simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 314 /* Create finite element */ 315 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "field_", -1, &feq)); 316 PetscCall(PetscObjectSetName((PetscObject) feq, "field")); 317 PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "potential_", -1, &feu)); 318 PetscCall(PetscObjectSetName((PetscObject) feu, "potential")); 319 PetscCall(PetscFECopyQuadrature(feq, feu)); 320 /* Set discretization and boundary conditions for each mesh */ 321 PetscCall(DMSetField(dm, 0, NULL, (PetscObject) feq)); 322 PetscCall(DMSetField(dm, 1, NULL, (PetscObject) feu)); 323 PetscCall(DMCreateDS(dm)); 324 PetscCall((*setup)(dm, user)); 325 while (cdm) { 326 PetscCall(DMCopyDisc(dm,cdm)); 327 PetscCall(DMGetCoarseDM(cdm, &cdm)); 328 } 329 PetscCall(PetscFEDestroy(&feq)); 330 PetscCall(PetscFEDestroy(&feu)); 331 PetscFunctionReturn(0); 332 } 333 334 int main(int argc, char **argv) 335 { 336 DM dm; /* Problem specification */ 337 SNES snes; /* Nonlinear solver */ 338 Vec u; /* Solutions */ 339 AppCtx user; /* User-defined work context */ 340 341 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 342 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 343 /* Primal system */ 344 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 345 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 346 PetscCall(SNESSetDM(snes, dm)); 347 PetscCall(SetupDiscretization(dm, SetupPrimalProblem, &user)); 348 PetscCall(DMCreateGlobalVector(dm, &u)); 349 PetscCall(VecSet(u, 0.0)); 350 PetscCall(PetscObjectSetName((PetscObject) u, "potential")); 351 PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 352 PetscCall(SNESSetFromOptions(snes)); 353 PetscCall(DMSNESCheckFromOptions(snes, u)); 354 PetscCall(SNESSolve(snes, NULL, u)); 355 PetscCall(SNESGetSolution(snes, &u)); 356 PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 357 /* Cleanup */ 358 PetscCall(VecDestroy(&u)); 359 PetscCall(SNESDestroy(&snes)); 360 PetscCall(DMDestroy(&dm)); 361 PetscCall(PetscFinalize()); 362 return 0; 363 } 364 365 /*TEST 366 367 test: 368 suffix: 2d_bdm1_p0 369 requires: triangle 370 args: -sol_type linear \ 371 -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 1 \ 372 -dmsnes_check .001 -snes_error_if_not_converged \ 373 -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 374 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 375 -fieldsplit_field_pc_type lu \ 376 -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 377 test: 378 # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: [2.0, 1.0] 379 # Using -sol_type quadratic -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: [2.9, 1.0] 380 suffix: 2d_bdm1_p0_conv 381 requires: triangle 382 args: -sol_type quartic \ 383 -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 384 -snes_error_if_not_converged \ 385 -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 386 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 387 -fieldsplit_field_pc_type lu \ 388 -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 389 390 test: 391 # HDF5 output: -dm_view hdf5:${PETSC_DIR}/sol.h5 -potential_view hdf5:${PETSC_DIR}/sol.h5::append -exact_vec_view hdf5:${PETSC_DIR}/sol.h5::append 392 # VTK output: -potential_view vtk: -exact_vec_view vtk: 393 # VTU output: -potential_view vtk:multifield.vtu -exact_vec_view vtk:exact.vtu 394 suffix: 2d_q2_p0 395 requires: triangle 396 args: -sol_type linear -dm_plex_simplex 0 \ 397 -field_petscspace_degree 2 -dm_refine 0 \ 398 -dmsnes_check .001 -snes_error_if_not_converged \ 399 -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 400 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 401 -fieldsplit_field_pc_type lu \ 402 -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 403 test: 404 # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: [0.0057, 1.0] 405 suffix: 2d_q2_p0_conv 406 requires: triangle 407 args: -sol_type linear -dm_plex_simplex 0 \ 408 -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 409 -snes_error_if_not_converged \ 410 -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 411 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 412 -fieldsplit_field_pc_type lu \ 413 -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 414 test: 415 # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: [-0.014, 0.0066] 416 suffix: 2d_q2_p0_neumann_conv 417 requires: triangle 418 args: -sol_type quartic_neumann -dm_plex_simplex 0 \ 419 -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 420 -snes_error_if_not_converged \ 421 -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 422 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 423 -fieldsplit_field_pc_type lu \ 424 -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type svd 425 426 TEST*/ 427 428 /* These tests will be active once tensor P^- is working 429 430 test: 431 suffix: 2d_bdmq1_p0_0 432 requires: triangle 433 args: -dm_plex_simplex 0 -sol_type linear \ 434 -field_petscspace_poly_type pminus_hdiv -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate \ 435 -dmsnes_check .001 -snes_error_if_not_converged \ 436 -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 437 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 438 -fieldsplit_field_pc_type lu \ 439 -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 440 test: 441 suffix: 2d_bdmq1_p0_2 442 requires: triangle 443 args: -dm_plex_simplex 0 -sol_type quartic \ 444 -field_petscspace_poly_type_no pminus_hdiv -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate \ 445 -dmsnes_check .001 -snes_error_if_not_converged \ 446 -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 447 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 448 -fieldsplit_field_pc_type lu \ 449 -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 450 451 */ 452