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