1c4762a1bSJed Brown static char help[] = "Poisson Problem in mixed form with 2d and 3d with finite elements.\n\ 2c4762a1bSJed Brown We solve the Poisson problem in a rectangular\n\ 3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4c4762a1bSJed Brown This example supports automatic convergence estimation\n\ 5c4762a1bSJed Brown and Hdiv elements.\n\n\n"; 6c4762a1bSJed Brown 769b95281SMatthew G. Knepley /* 869b95281SMatthew G. Knepley 969b95281SMatthew G. Knepley The mixed form of Poisson's equation, e.g. https://firedrakeproject.org/demos/poisson_mixed.py.html, is given 1069b95281SMatthew G. Knepley in the strong form by 1169b95281SMatthew G. Knepley \begin{align} 1269b95281SMatthew G. Knepley \vb{q} - \nabla u &= 0 \\ 1369b95281SMatthew G. Knepley \nabla \cdot \vb{q} &= f 1469b95281SMatthew G. Knepley \end{align} 1569b95281SMatthew G. Knepley where $u$ is the potential, as in the original problem, but we also discretize the gradient of potential $\vb{q}$, 1669b95281SMatthew G. Knepley or flux, directly. The weak form is then 1769b95281SMatthew G. Knepley \begin{align} 1869b95281SMatthew G. Knepley <t, \vb{q}> + <\nabla \cdot t, u> - <t_n, u>_\Gamma &= 0 \\ 1969b95281SMatthew G. Knepley <v, \nabla \cdot \vb{q}> &= <v, f> 2069b95281SMatthew G. Knepley \end{align} 2169b95281SMatthew G. Knepley For the original Poisson problem, the Dirichlet boundary forces an essential boundary condition on the potential space, 2269b95281SMatthew G. Knepley and the Neumann boundary gives a natural boundary condition in the weak form. In the mixed formulation, the Neumann 2369b95281SMatthew G. Knepley boundary gives an essential boundary condition on the flux space, $\vb{q} \cdot \vu{n} = h$, and the Dirichlet condition 2469b95281SMatthew G. Knepley becomes a natural condition in the weak form, <t_n, g>_\Gamma. 2569b95281SMatthew G. Knepley */ 2669b95281SMatthew G. Knepley 27c4762a1bSJed Brown #include <petscdmplex.h> 28c4762a1bSJed Brown #include <petscsnes.h> 29c4762a1bSJed Brown #include <petscds.h> 30c4762a1bSJed Brown #include <petscconvest.h> 31c4762a1bSJed Brown 32*9371c9d4SSatish Balay typedef enum { 33*9371c9d4SSatish Balay SOL_LINEAR, 34*9371c9d4SSatish Balay SOL_QUADRATIC, 35*9371c9d4SSatish Balay SOL_QUARTIC, 36*9371c9d4SSatish Balay SOL_QUARTIC_NEUMANN, 37*9371c9d4SSatish Balay SOL_UNKNOWN, 38*9371c9d4SSatish Balay NUM_SOLTYPE 39*9371c9d4SSatish Balay } SolType; 4069b95281SMatthew G. Knepley const char *SolTypeNames[NUM_SOLTYPE + 3] = {"linear", "quadratic", "quartic", "quartic_neumann", "unknown", "SolType", "SOL_", NULL}; 41c4762a1bSJed Brown 42c4762a1bSJed Brown typedef struct { 43c4762a1bSJed Brown SolType solType; /* The type of exact solution */ 44c4762a1bSJed Brown } AppCtx; 45c4762a1bSJed Brown 46*9371c9d4SSatish Balay static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 4769b95281SMatthew G. Knepley PetscInt d; 4869b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += u_x[uOff_x[0] + d * dim + d]; 4969b95281SMatthew G. Knepley } 5069b95281SMatthew G. Knepley 51c4762a1bSJed Brown /* 2D Dirichlet potential example 52c4762a1bSJed Brown 53c4762a1bSJed Brown u = x 54c4762a1bSJed Brown q = <1, 0> 55c4762a1bSJed Brown f = 0 56c4762a1bSJed Brown 57c4762a1bSJed Brown We will need a boundary integral of u over \Gamma. 58c4762a1bSJed Brown */ 59*9371c9d4SSatish Balay static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 60c4762a1bSJed Brown u[0] = x[0]; 61c4762a1bSJed Brown return 0; 62c4762a1bSJed Brown } 63*9371c9d4SSatish Balay static PetscErrorCode linear_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 64c4762a1bSJed Brown PetscInt c; 65c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = c ? 0.0 : 1.0; 66c4762a1bSJed Brown return 0; 67c4762a1bSJed Brown } 68c4762a1bSJed Brown 69*9371c9d4SSatish Balay static void f0_linear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 7069b95281SMatthew G. Knepley f0[0] = 0.0; 7169b95281SMatthew G. Knepley 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); 7269b95281SMatthew G. Knepley } 73*9371c9d4SSatish Balay static void f0_bd_linear_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 7469b95281SMatthew G. Knepley PetscScalar potential; 7569b95281SMatthew G. Knepley PetscInt d; 7669b95281SMatthew G. Knepley 7769b95281SMatthew G. Knepley linear_u(dim, t, x, dim, &potential, NULL); 7869b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[d] = -potential * n[d]; 7969b95281SMatthew G. Knepley } 8069b95281SMatthew G. Knepley 81c4762a1bSJed Brown /* 2D Dirichlet potential example 82c4762a1bSJed Brown 83c4762a1bSJed Brown u = x^2 + y^2 84c4762a1bSJed Brown q = <2x, 2y> 85c4762a1bSJed Brown f = 4 86c4762a1bSJed Brown 87c4762a1bSJed Brown We will need a boundary integral of u over \Gamma. 88c4762a1bSJed Brown */ 89*9371c9d4SSatish Balay static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 90c4762a1bSJed Brown PetscInt d; 91c4762a1bSJed Brown 92c4762a1bSJed Brown u[0] = 0.0; 93c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[0] += x[d] * x[d]; 94c4762a1bSJed Brown return 0; 95c4762a1bSJed Brown } 96*9371c9d4SSatish Balay static PetscErrorCode quadratic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 97c4762a1bSJed Brown PetscInt c; 98c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = 2.0 * x[c]; 99c4762a1bSJed Brown return 0; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102*9371c9d4SSatish Balay static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 10369b95281SMatthew G. Knepley f0[0] = -4.0; 10469b95281SMatthew G. Knepley 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); 10569b95281SMatthew G. Knepley } 106*9371c9d4SSatish Balay static void f0_bd_quadratic_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 10769b95281SMatthew G. Knepley PetscScalar potential; 10869b95281SMatthew G. Knepley PetscInt d; 10969b95281SMatthew G. Knepley 11069b95281SMatthew G. Knepley quadratic_u(dim, t, x, dim, &potential, NULL); 11169b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[d] = -potential * n[d]; 11269b95281SMatthew G. Knepley } 11369b95281SMatthew G. Knepley 114c4762a1bSJed Brown /* 2D Dirichlet potential example 115c4762a1bSJed Brown 116c4762a1bSJed Brown u = x (1-x) y (1-y) 117c4762a1bSJed Brown q = <(1-2x) y (1-y), x (1-x) (1-2y)> 118c4762a1bSJed Brown f = -y (1-y) - x (1-x) 119c4762a1bSJed Brown 120c4762a1bSJed Brown u|_\Gamma = 0 so that the boundary integral vanishes 121c4762a1bSJed Brown */ 122*9371c9d4SSatish Balay static PetscErrorCode quartic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 123c4762a1bSJed Brown PetscInt d; 124c4762a1bSJed Brown 125c4762a1bSJed Brown u[0] = 1.0; 126c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[0] *= x[d] * (1.0 - x[d]); 127c4762a1bSJed Brown return 0; 128c4762a1bSJed Brown } 129*9371c9d4SSatish Balay static PetscErrorCode quartic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 130c4762a1bSJed Brown PetscInt c, d; 131c4762a1bSJed Brown 132c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 133c4762a1bSJed Brown u[c] = 1.0; 134c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 135c4762a1bSJed Brown if (c == d) u[c] *= 1 - 2.0 * x[d]; 136c4762a1bSJed Brown else u[c] *= x[d] * (1.0 - x[d]); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown } 139c4762a1bSJed Brown return 0; 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142*9371c9d4SSatish Balay static void f0_quartic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 143c4762a1bSJed Brown PetscInt d; 144c4762a1bSJed Brown f0[0] = 0.0; 14569b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += 2.0 * x[d] * (1.0 - x[d]); 14669b95281SMatthew G. Knepley 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); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 14969b95281SMatthew G. Knepley /* 2D Dirichlet potential example 15069b95281SMatthew G. Knepley 15169b95281SMatthew G. Knepley u = x (1-x) y (1-y) 15269b95281SMatthew G. Knepley q = <(1-2x) y (1-y), x (1-x) (1-2y)> 15369b95281SMatthew G. Knepley f = -y (1-y) - x (1-x) 15469b95281SMatthew G. Knepley 15569b95281SMatthew G. Knepley du/dn_\Gamma = {(1-2x) y (1-y), x (1-x) (1-2y)} . n produces an essential condition on q 15669b95281SMatthew G. Knepley */ 15769b95281SMatthew G. Knepley 158*9371c9d4SSatish Balay static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 159c4762a1bSJed Brown PetscInt c; 16069b95281SMatthew G. Knepley for (c = 0; c < dim; ++c) f0[c] = u[uOff[0] + c]; 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163c4762a1bSJed Brown /* <\nabla\cdot w, u> = <\nabla w, Iu> */ 164*9371c9d4SSatish Balay static void f1_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) { 16569b95281SMatthew G. Knepley PetscInt c; 16669b95281SMatthew G. Knepley for (c = 0; c < dim; ++c) f1[c * dim + c] = u[uOff[1]]; 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 16969b95281SMatthew G. Knepley /* <t, q> */ 170*9371c9d4SSatish Balay static void g0_qq(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) { 171c4762a1bSJed Brown PetscInt c; 172c4762a1bSJed Brown for (c = 0; c < dim; ++c) g0[c * dim + c] = 1.0; 173c4762a1bSJed Brown } 174c4762a1bSJed Brown 17569b95281SMatthew G. Knepley /* <\nabla\cdot t, u> = <\nabla t, Iu> */ 176*9371c9d4SSatish Balay static void g2_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) { 177c4762a1bSJed Brown PetscInt d; 178c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d * dim + d] = 1.0; 179c4762a1bSJed Brown } 180c4762a1bSJed Brown 18169b95281SMatthew G. Knepley /* <v, \nabla\cdot q> */ 182*9371c9d4SSatish Balay static void g1_uq(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) { 183c4762a1bSJed Brown PetscInt d; 18469b95281SMatthew G. Knepley for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; 185c4762a1bSJed Brown } 186c4762a1bSJed Brown 187*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 188c4762a1bSJed Brown PetscFunctionBeginUser; 189c4762a1bSJed Brown options->solType = SOL_LINEAR; 190d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX"); 1919566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-sol_type", "Type of exact solution", "ex24.c", SolTypeNames, (PetscEnum)options->solType, (PetscEnum *)&options->solType, NULL)); 192d0609cedSBarry Smith PetscOptionsEnd(); 193c4762a1bSJed Brown PetscFunctionReturn(0); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 197c4762a1bSJed Brown PetscFunctionBeginUser; 1989566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1999566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 2009566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Example Mesh")); 2019566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 2029566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2039566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 204c4762a1bSJed Brown PetscFunctionReturn(0); 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207*9371c9d4SSatish Balay static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) { 20845480ffeSMatthew G. Knepley PetscDS ds; 20945480ffeSMatthew G. Knepley DMLabel label; 21045480ffeSMatthew G. Knepley PetscWeakForm wf; 211c4762a1bSJed Brown const PetscInt id = 1; 21245480ffeSMatthew G. Knepley PetscInt bd; 213c4762a1bSJed Brown 214c4762a1bSJed Brown PetscFunctionBeginUser; 2159566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 2169566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2179566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_q, f1_q)); 2189566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_qq, NULL, NULL, NULL)); 2199566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_qu, NULL)); 2209566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_uq, NULL, NULL)); 221*9371c9d4SSatish Balay switch (user->solType) { 222c4762a1bSJed Brown case SOL_LINEAR: 2239566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_linear_u, NULL)); 2249566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 2259566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2269566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_linear_q, 0, NULL)); 2279566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, linear_q, user)); 2289566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, linear_u, user)); 229c4762a1bSJed Brown break; 230c4762a1bSJed Brown case SOL_QUADRATIC: 2319566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_quadratic_u, NULL)); 2329566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 2339566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2349566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_quadratic_q, 0, NULL)); 2359566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_q, user)); 2369566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, quadratic_u, user)); 237c4762a1bSJed Brown break; 238c4762a1bSJed Brown case SOL_QUARTIC: 2399566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_quartic_u, NULL)); 2409566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quartic_q, user)); 2419566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, quartic_u, user)); 242c4762a1bSJed Brown break; 24369b95281SMatthew G. Knepley case SOL_QUARTIC_NEUMANN: 2449566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_quartic_u, NULL)); 2459566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quartic_q, user)); 2469566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, quartic_u, user)); 2479566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "Flux condition", label, 1, &id, 0, 0, NULL, (void (*)(void))quartic_q, NULL, user, NULL)); 24869b95281SMatthew G. Knepley break; 24998921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid exact solution type %s", SolTypeNames[PetscMin(user->solType, SOL_UNKNOWN)]); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown PetscFunctionReturn(0); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown 254*9371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) { 255c4762a1bSJed Brown DM cdm = dm; 256c4762a1bSJed Brown PetscFE feq, feu; 25769b95281SMatthew G. Knepley DMPolytopeType ct; 25869b95281SMatthew G. Knepley PetscBool simplex; 25969b95281SMatthew G. Knepley PetscInt dim, cStart; 260c4762a1bSJed Brown 261c4762a1bSJed Brown PetscFunctionBeginUser; 2629566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 2639566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 2649566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 26569b95281SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 266c4762a1bSJed Brown /* Create finite element */ 2679566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "field_", -1, &feq)); 2689566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)feq, "field")); 2699566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "potential_", -1, &feu)); 2709566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)feu, "potential")); 2719566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(feq, feu)); 272c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 2739566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)feq)); 2749566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject)feu)); 2759566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2769566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 277c4762a1bSJed Brown while (cdm) { 2789566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 2799566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 280c4762a1bSJed Brown } 2819566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feq)); 2829566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feu)); 283c4762a1bSJed Brown PetscFunctionReturn(0); 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286*9371c9d4SSatish Balay int main(int argc, char **argv) { 287c4762a1bSJed Brown DM dm; /* Problem specification */ 288c4762a1bSJed Brown SNES snes; /* Nonlinear solver */ 289c4762a1bSJed Brown Vec u; /* Solutions */ 290c4762a1bSJed Brown AppCtx user; /* User-defined work context */ 291c4762a1bSJed Brown 292327415f7SBarry Smith PetscFunctionBeginUser; 2939566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2949566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 295c4762a1bSJed Brown /* Primal system */ 2969566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 2979566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2989566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 2999566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, SetupPrimalProblem, &user)); 3009566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 3019566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 3029566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "potential")); 3039566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 3049566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 3059566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 3069566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 3079566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 3089566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-potential_view")); 309c4762a1bSJed Brown /* Cleanup */ 3109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 3119566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 3129566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3139566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 314b122ec5aSJacob Faibussowitsch return 0; 315c4762a1bSJed Brown } 316c4762a1bSJed Brown 317c4762a1bSJed Brown /*TEST 318c4762a1bSJed Brown 319c4762a1bSJed Brown test: 32069b95281SMatthew G. Knepley suffix: 2d_bdm1_p0 321c4762a1bSJed Brown requires: triangle 322c4762a1bSJed Brown args: -sol_type linear \ 32369b95281SMatthew G. Knepley -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 1 \ 324c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 325c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 326c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 327c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 328c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 329c4762a1bSJed Brown test: 33069b95281SMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: [2.0, 1.0] 33169b95281SMatthew G. Knepley # Using -sol_type quadratic -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: [2.9, 1.0] 33269b95281SMatthew G. Knepley suffix: 2d_bdm1_p0_conv 333c4762a1bSJed Brown requires: triangle 334c4762a1bSJed Brown args: -sol_type quartic \ 335c4762a1bSJed Brown -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 336c4762a1bSJed Brown -snes_error_if_not_converged \ 337c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 338c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 339c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 340c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 341c4762a1bSJed Brown 342c4762a1bSJed Brown test: 34369b95281SMatthew G. Knepley # 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 34469b95281SMatthew G. Knepley # VTK output: -potential_view vtk: -exact_vec_view vtk: 34569b95281SMatthew G. Knepley # VTU output: -potential_view vtk:multifield.vtu -exact_vec_view vtk:exact.vtu 34669b95281SMatthew G. Knepley suffix: 2d_q2_p0 347c4762a1bSJed Brown requires: triangle 34830602db0SMatthew G. Knepley args: -sol_type linear -dm_plex_simplex 0 \ 34969b95281SMatthew G. Knepley -field_petscspace_degree 2 -dm_refine 0 \ 350c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 351c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 352c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 353c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 354c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 355c4762a1bSJed Brown test: 35669b95281SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: [0.0057, 1.0] 35769b95281SMatthew G. Knepley suffix: 2d_q2_p0_conv 358c4762a1bSJed Brown requires: triangle 35930602db0SMatthew G. Knepley args: -sol_type linear -dm_plex_simplex 0 \ 360c4762a1bSJed Brown -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 36169b95281SMatthew G. Knepley -snes_error_if_not_converged \ 362c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 363c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 364c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 365c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 36669b95281SMatthew G. Knepley test: 36769b95281SMatthew G. Knepley # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: [-0.014, 0.0066] 36869b95281SMatthew G. Knepley suffix: 2d_q2_p0_neumann_conv 36969b95281SMatthew G. Knepley requires: triangle 37030602db0SMatthew G. Knepley args: -sol_type quartic_neumann -dm_plex_simplex 0 \ 37169b95281SMatthew G. Knepley -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 37269b95281SMatthew G. Knepley -snes_error_if_not_converged \ 37369b95281SMatthew G. Knepley -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 37469b95281SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 37569b95281SMatthew G. Knepley -fieldsplit_field_pc_type lu \ 37669b95281SMatthew G. Knepley -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type svd 37769b95281SMatthew G. Knepley 378c4762a1bSJed Brown TEST*/ 379c4762a1bSJed Brown 380c4762a1bSJed Brown /* These tests will be active once tensor P^- is working 381c4762a1bSJed Brown 382c4762a1bSJed Brown test: 383c4762a1bSJed Brown suffix: 2d_bdmq1_p0_0 384c4762a1bSJed Brown requires: triangle 38530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -sol_type linear \ 386c4762a1bSJed Brown -field_petscspace_poly_type pminus_hdiv -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate \ 387c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 388c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 389c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 390c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 391c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 392c4762a1bSJed Brown test: 393c4762a1bSJed Brown suffix: 2d_bdmq1_p0_2 394c4762a1bSJed Brown requires: triangle 39530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -sol_type quartic \ 396c4762a1bSJed Brown -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 \ 397c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 398c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 399c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 400c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 401c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 402c4762a1bSJed Brown 403c4762a1bSJed Brown */ 404