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 7c4762a1bSJed Brown #include <petscdmplex.h> 8c4762a1bSJed Brown #include <petscsnes.h> 9c4762a1bSJed Brown #include <petscds.h> 10c4762a1bSJed Brown #include <petscconvest.h> 11c4762a1bSJed Brown 12c4762a1bSJed Brown typedef enum {SOL_LINEAR, SOL_QUADRATIC, SOL_QUARTIC, SOL_UNKNOWN, NUM_SOLTYPE} SolType; 13c4762a1bSJed Brown const char *SolTypeNames[NUM_SOLTYPE+3] = {"linear", "quadratic", "quartic", "unknown", "SolType", "SOL_", NULL}; 14c4762a1bSJed Brown 15c4762a1bSJed Brown typedef struct { 16c4762a1bSJed Brown /* Domain and mesh definition */ 17c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 18c4762a1bSJed Brown PetscBool simplex; /* Simplicial mesh */ 19c4762a1bSJed Brown SolType solType; /* The type of exact solution */ 20c4762a1bSJed Brown } AppCtx; 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* 2D Dirichlet potential example 23c4762a1bSJed Brown 24c4762a1bSJed Brown u = x 25c4762a1bSJed Brown q = <1, 0> 26c4762a1bSJed Brown f = 0 27c4762a1bSJed Brown 28c4762a1bSJed Brown We will need a boundary integral of u over \Gamma. 29c4762a1bSJed Brown */ 30c4762a1bSJed Brown static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 31c4762a1bSJed Brown { 32c4762a1bSJed Brown u[0] = x[0]; 33c4762a1bSJed Brown return 0; 34c4762a1bSJed Brown } 35c4762a1bSJed Brown static PetscErrorCode linear_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 36c4762a1bSJed Brown { 37c4762a1bSJed Brown PetscInt c; 38c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = c ? 0.0 : 1.0; 39c4762a1bSJed Brown return 0; 40c4762a1bSJed Brown } 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* 2D Dirichlet potential example 43c4762a1bSJed Brown 44c4762a1bSJed Brown u = x^2 + y^2 45c4762a1bSJed Brown q = <2x, 2y> 46c4762a1bSJed Brown f = 4 47c4762a1bSJed Brown 48c4762a1bSJed Brown We will need a boundary integral of u over \Gamma. 49c4762a1bSJed Brown */ 50c4762a1bSJed Brown static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 51c4762a1bSJed Brown { 52c4762a1bSJed Brown PetscInt d; 53c4762a1bSJed Brown 54c4762a1bSJed Brown u[0] = 0.0; 55c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[0] += x[d]*x[d]; 56c4762a1bSJed Brown return 0; 57c4762a1bSJed Brown } 58c4762a1bSJed Brown static PetscErrorCode quadratic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 59c4762a1bSJed Brown { 60c4762a1bSJed Brown PetscInt c; 61c4762a1bSJed Brown for (c = 0; c < Nc; ++c) u[c] = 2.0*x[c]; 62c4762a1bSJed Brown return 0; 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* 2D Dirichlet potential example 66c4762a1bSJed Brown 67c4762a1bSJed Brown u = x (1-x) y (1-y) 68c4762a1bSJed Brown q = <(1-2x) y (1-y), x (1-x) (1-2y)> 69c4762a1bSJed Brown f = -y (1-y) - x (1-x) 70c4762a1bSJed Brown 71c4762a1bSJed Brown u|_\Gamma = 0 so that the boundary integral vanishes 72c4762a1bSJed Brown */ 73c4762a1bSJed Brown static PetscErrorCode quartic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 74c4762a1bSJed Brown { 75c4762a1bSJed Brown PetscInt d; 76c4762a1bSJed Brown 77c4762a1bSJed Brown u[0] = 1.0; 78c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[0] *= x[d]*(1.0 - x[d]); 79c4762a1bSJed Brown return 0; 80c4762a1bSJed Brown } 81c4762a1bSJed Brown static PetscErrorCode quartic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 82c4762a1bSJed Brown { 83c4762a1bSJed Brown PetscInt c, d; 84c4762a1bSJed Brown 85c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 86c4762a1bSJed Brown u[c] = 1.0; 87c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 88c4762a1bSJed Brown if (c == d) u[c] *= 1 - 2.0*x[d]; 89c4762a1bSJed Brown else u[c] *= x[d]*(1.0 - x[d]); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown } 92c4762a1bSJed Brown return 0; 93c4762a1bSJed Brown } 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* <v, -\nabla\cdot q> + <v, f> */ 96c4762a1bSJed Brown static void f0_linear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 97c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 98c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 99c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 100c4762a1bSJed Brown { 101c4762a1bSJed Brown f0[0] = 0.0; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown static void f0_bd_linear_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 104c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 105c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 106c4762a1bSJed Brown PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 107c4762a1bSJed Brown { 108c4762a1bSJed Brown PetscScalar potential; 109c4762a1bSJed Brown PetscInt d; 110c4762a1bSJed Brown 111c4762a1bSJed Brown linear_u(dim, t, x, dim, &potential, NULL); 112c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] = -potential*n[d]; 113c4762a1bSJed Brown } 114c4762a1bSJed Brown static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 115c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 116c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 117c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 118c4762a1bSJed Brown { 119c4762a1bSJed Brown PetscInt d; 120c4762a1bSJed Brown f0[0] = 0.0; 121c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 122c4762a1bSJed Brown f0[0] -= u_x[uOff_x[0]+d*dim+d]; 123c4762a1bSJed Brown } 124c4762a1bSJed Brown f0[0] += 4.0; 125c4762a1bSJed Brown } 126c4762a1bSJed Brown static void f0_bd_quadratic_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 127c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 128c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 129c4762a1bSJed Brown PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 130c4762a1bSJed Brown { 131c4762a1bSJed Brown PetscScalar potential; 132c4762a1bSJed Brown PetscInt d; 133c4762a1bSJed Brown 134c4762a1bSJed Brown quadratic_u(dim, t, x, dim, &potential, NULL); 135c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] = -potential*n[d]; 136c4762a1bSJed Brown } 137c4762a1bSJed Brown static void f0_quartic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 138c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 139c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 140c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 141c4762a1bSJed Brown { 142c4762a1bSJed Brown PetscInt d; 143c4762a1bSJed Brown f0[0] = 0.0; 144c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 145c4762a1bSJed Brown f0[0] -= u_x[uOff_x[0]+d*dim+d]; 146c4762a1bSJed Brown f0[0] += -2.0*x[d]*(1.0 - x[d]); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* <w, q> */ 151c4762a1bSJed Brown static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 152c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 153c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 154c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 155c4762a1bSJed Brown { 156c4762a1bSJed Brown PetscInt c; 157c4762a1bSJed Brown 158c4762a1bSJed Brown for (c = 0; c < dim; ++c) { 159c4762a1bSJed Brown f0[c] = u[uOff[0]+c]; 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163c4762a1bSJed Brown /* <\nabla\cdot w, u> = <\nabla w, Iu> */ 164c4762a1bSJed Brown static void f1_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, 165c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 166c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 167c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 168c4762a1bSJed Brown { 169c4762a1bSJed Brown PetscInt c, d; 170c4762a1bSJed Brown for (c = 0; c < dim; ++c) { 171c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 172c4762a1bSJed Brown if (c == d) f1[c*dim+d] = u[uOff[1]]; 173c4762a1bSJed Brown } 174c4762a1bSJed Brown } 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* <w, q> */ 178c4762a1bSJed Brown static void g0_qq(PetscInt dim, PetscInt Nf, PetscInt NfAux, 179c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 180c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 181c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 182c4762a1bSJed Brown { 183c4762a1bSJed Brown PetscInt c; 184c4762a1bSJed Brown for (c = 0; c < dim; ++c) g0[c*dim+c] = 1.0; 185c4762a1bSJed Brown } 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* <\nabla\cdot w, u> = <\nabla w, Iu> */ 188c4762a1bSJed Brown static void g2_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 189c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 190c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 191c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 192c4762a1bSJed Brown { 193c4762a1bSJed Brown PetscInt d; 194c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d*dim+d] = 1.0; 195c4762a1bSJed Brown } 196c4762a1bSJed Brown 197c4762a1bSJed Brown /* <v, -\nabla\cdot q> */ 198c4762a1bSJed Brown static void g1_uq(PetscInt dim, PetscInt Nf, PetscInt NfAux, 199c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 200c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 201c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 202c4762a1bSJed Brown { 203c4762a1bSJed Brown PetscInt d; 204c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d*dim+d] = -1.0; 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 208c4762a1bSJed Brown { 209c4762a1bSJed Brown PetscErrorCode ierr; 210c4762a1bSJed Brown 211c4762a1bSJed Brown PetscFunctionBeginUser; 212c4762a1bSJed Brown options->dim = 2; 213c4762a1bSJed Brown options->simplex = PETSC_TRUE; 214c4762a1bSJed Brown options->solType = SOL_LINEAR; 215c4762a1bSJed Brown 216c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex24.c", options->dim, &options->dim, NULL);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = PetscOptionsBool("-simplex", "Simplicial (true) or tensor (false) mesh", "ex24.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = PetscOptionsEnum("-sol_type", "Type of exact solution", "ex24.c", SolTypeNames, (PetscEnum) options->solType, (PetscEnum *) &options->solType, NULL);CHKERRQ(ierr); 220c4762a1bSJed Brown 221c4762a1bSJed Brown ierr = PetscOptionsEnd(); 222c4762a1bSJed Brown PetscFunctionReturn(0); 223c4762a1bSJed Brown } 224c4762a1bSJed Brown 225c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 226c4762a1bSJed Brown { 227c4762a1bSJed Brown PetscErrorCode ierr; 228c4762a1bSJed Brown 229c4762a1bSJed Brown PetscFunctionBeginUser; 230c4762a1bSJed Brown if (0) { 231c4762a1bSJed Brown DMLabel label; 232c4762a1bSJed Brown const char *name = "marker"; 233c4762a1bSJed Brown 234c4762a1bSJed Brown ierr = DMPlexCreateReferenceCell(comm, user->dim, user->simplex, dm);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = DMCreateLabel(*dm, name);CHKERRQ(ierr); 236c4762a1bSJed Brown ierr = DMGetLabel(*dm, name, &label);CHKERRQ(ierr); 237c4762a1bSJed Brown ierr = DMPlexMarkBoundaryFaces(*dm, 1, label);CHKERRQ(ierr); 238c4762a1bSJed Brown ierr = DMPlexLabelComplete(*dm, label);CHKERRQ(ierr); 239c4762a1bSJed Brown } else { 240c4762a1bSJed Brown /* Create box mesh */ 241c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, user->dim, user->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 242c4762a1bSJed Brown } 243c4762a1bSJed Brown /* Distribute mesh over processes */ 244c4762a1bSJed Brown { 245c4762a1bSJed Brown DM dmDist = NULL; 246c4762a1bSJed Brown PetscPartitioner part; 247c4762a1bSJed Brown 248c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 249c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 250c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &dmDist);CHKERRQ(ierr); 251c4762a1bSJed Brown if (dmDist) { 252c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 253c4762a1bSJed Brown *dm = dmDist; 254c4762a1bSJed Brown } 255c4762a1bSJed Brown } 256c4762a1bSJed Brown /* TODO: This should be pulled into the library */ 257c4762a1bSJed Brown { 258c4762a1bSJed Brown char convType[256]; 259c4762a1bSJed Brown PetscBool flg; 260c4762a1bSJed Brown 261c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr); 262c4762a1bSJed Brown ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr); 263c4762a1bSJed Brown ierr = PetscOptionsEnd(); 264c4762a1bSJed Brown if (flg) { 265c4762a1bSJed Brown DM dmConv; 266c4762a1bSJed Brown 267c4762a1bSJed Brown ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr); 268c4762a1bSJed Brown if (dmConv) { 269c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 270c4762a1bSJed Brown *dm = dmConv; 271c4762a1bSJed Brown } 272c4762a1bSJed Brown } 273c4762a1bSJed Brown } 274c4762a1bSJed Brown /* TODO: This should be pulled into the library */ 275c4762a1bSJed Brown ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); 276c4762a1bSJed Brown 277c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 278c4762a1bSJed Brown ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr); 279c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 280c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 281c4762a1bSJed Brown PetscFunctionReturn(0); 282c4762a1bSJed Brown } 283c4762a1bSJed Brown 284c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 285c4762a1bSJed Brown { 286c4762a1bSJed Brown PetscDS prob; 287c4762a1bSJed Brown const PetscInt id = 1; 288c4762a1bSJed Brown PetscErrorCode ierr; 289c4762a1bSJed Brown 290c4762a1bSJed Brown PetscFunctionBeginUser; 291c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 292c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_q, f1_q);CHKERRQ(ierr); 293c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, g0_qq, NULL, NULL, NULL);CHKERRQ(ierr); 294c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_qu, NULL);CHKERRQ(ierr); 295c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_uq, NULL, NULL);CHKERRQ(ierr); 296c4762a1bSJed Brown switch (user->solType) 297c4762a1bSJed Brown { 298c4762a1bSJed Brown case SOL_LINEAR: 299c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 1, f0_linear_u, NULL);CHKERRQ(ierr); 300c4762a1bSJed Brown ierr = PetscDSSetBdResidual(prob, 0, f0_bd_linear_q, NULL);CHKERRQ(ierr); 301408cafa0SMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", "marker", 0, 0, NULL, (void (*)(void)) NULL, 1, &id, user);CHKERRQ(ierr); 302c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob, 0, linear_q, user);CHKERRQ(ierr); 303c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob, 1, linear_u, user);CHKERRQ(ierr); 304c4762a1bSJed Brown break; 305c4762a1bSJed Brown case SOL_QUADRATIC: 306c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 1, f0_quadratic_u, NULL);CHKERRQ(ierr); 307c4762a1bSJed Brown ierr = PetscDSSetBdResidual(prob, 0, f0_bd_quadratic_q, NULL);CHKERRQ(ierr); 308408cafa0SMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", "marker", 0, 0, NULL, (void (*)(void)) NULL, 1, &id, user);CHKERRQ(ierr); 309c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob, 0, quadratic_q, user);CHKERRQ(ierr); 310c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob, 1, quadratic_u, user);CHKERRQ(ierr); 311c4762a1bSJed Brown break; 312c4762a1bSJed Brown case SOL_QUARTIC: 313c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 1, f0_quartic_u, NULL);CHKERRQ(ierr); 314c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob, 0, quartic_q, user);CHKERRQ(ierr); 315c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob, 1, quartic_u, user);CHKERRQ(ierr); 316c4762a1bSJed Brown break; 317c4762a1bSJed Brown default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid exact solution type %s", SolTypeNames[PetscMin(user->solType, SOL_UNKNOWN)]); 318c4762a1bSJed Brown } 319c4762a1bSJed Brown PetscFunctionReturn(0); 320c4762a1bSJed Brown } 321c4762a1bSJed Brown 322c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user) 323c4762a1bSJed Brown { 324c4762a1bSJed Brown DM cdm = dm; 325c4762a1bSJed Brown PetscFE feq, feu; 326c4762a1bSJed Brown const PetscInt dim = user->dim; 327c4762a1bSJed Brown PetscErrorCode ierr; 328c4762a1bSJed Brown 329c4762a1bSJed Brown PetscFunctionBeginUser; 330c4762a1bSJed Brown /* Create finite element */ 331c4762a1bSJed Brown ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, dim, user->simplex, "field_", -1, &feq);CHKERRQ(ierr); 332c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) feq, "field");CHKERRQ(ierr); 333c4762a1bSJed Brown ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, user->simplex, "potential_", -1, &feu);CHKERRQ(ierr); 334c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) feu, "potential");CHKERRQ(ierr); 335c4762a1bSJed Brown ierr = PetscFECopyQuadrature(feq, feu);CHKERRQ(ierr); 336c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 337c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) feq);CHKERRQ(ierr); 338c4762a1bSJed Brown ierr = DMSetField(dm, 1, NULL, (PetscObject) feu);CHKERRQ(ierr); 339c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 340c4762a1bSJed Brown ierr = (*setup)(dm, user);CHKERRQ(ierr); 341c4762a1bSJed Brown while (cdm) { 342c4762a1bSJed Brown ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr); 343c4762a1bSJed Brown /* TODO: Check whether the boundary of coarse meshes is marked */ 344c4762a1bSJed Brown ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 345c4762a1bSJed Brown } 346c4762a1bSJed Brown ierr = PetscFEDestroy(&feq);CHKERRQ(ierr); 347c4762a1bSJed Brown ierr = PetscFEDestroy(&feu);CHKERRQ(ierr); 348c4762a1bSJed Brown PetscFunctionReturn(0); 349c4762a1bSJed Brown } 350c4762a1bSJed Brown 351c4762a1bSJed Brown int main(int argc, char **argv) 352c4762a1bSJed Brown { 353c4762a1bSJed Brown DM dm; /* Problem specification */ 354c4762a1bSJed Brown SNES snes; /* Nonlinear solver */ 355c4762a1bSJed Brown Vec u; /* Solutions */ 356c4762a1bSJed Brown AppCtx user; /* User-defined work context */ 357c4762a1bSJed Brown PetscErrorCode ierr; 358c4762a1bSJed Brown 359c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 360c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 361c4762a1bSJed Brown /* Primal system */ 362c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 363c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 364c4762a1bSJed Brown ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 365c4762a1bSJed Brown ierr = SetupDiscretization(dm, SetupPrimalProblem, &user);CHKERRQ(ierr); 366c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 367c4762a1bSJed Brown ierr = VecSet(u, 0.0);CHKERRQ(ierr); 368c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr); 369c4762a1bSJed Brown ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr); 370c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 371*348a1646SMatthew G. Knepley ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr); 372c4762a1bSJed Brown ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 373c4762a1bSJed Brown ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr); 374c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr); 375c4762a1bSJed Brown /* Cleanup */ 376c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 377c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 378c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 379c4762a1bSJed Brown ierr = PetscFinalize(); 380c4762a1bSJed Brown return ierr; 381c4762a1bSJed Brown } 382c4762a1bSJed Brown 383c4762a1bSJed Brown /*TEST 384c4762a1bSJed Brown 385c4762a1bSJed Brown test: 386c4762a1bSJed Brown suffix: 2d_bdm1_p0_0 387c4762a1bSJed Brown requires: triangle 388c4762a1bSJed Brown args: -sol_type linear \ 389c4762a1bSJed Brown -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 390c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 391c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 392c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 393c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 394c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 395c4762a1bSJed Brown test: 396c4762a1bSJed Brown suffix: 2d_bdm1_p0_1 397c4762a1bSJed Brown requires: triangle 398c4762a1bSJed Brown args: -sol_type quadratic \ 399c4762a1bSJed Brown -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 400c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 401c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 402c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 403c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 404c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 405c4762a1bSJed Brown test: 406c4762a1bSJed Brown suffix: 2d_bdm1_p0_2 407c4762a1bSJed Brown requires: triangle 408c4762a1bSJed Brown args: -sol_type quartic \ 409c4762a1bSJed Brown -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 410c4762a1bSJed Brown -snes_error_if_not_converged \ 411c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 412c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 413c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 414c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 415c4762a1bSJed Brown 416c4762a1bSJed Brown test: 417c4762a1bSJed Brown suffix: 2d_p2_p0_vtk 418c4762a1bSJed Brown requires: triangle 419c4762a1bSJed Brown args: -sol_type linear \ 420c4762a1bSJed Brown -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 421c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 422c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 423c4762a1bSJed Brown -potential_view vtk: -exact_vec_view vtk: \ 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 427c4762a1bSJed Brown 428c4762a1bSJed Brown test: 429c4762a1bSJed Brown suffix: 2d_p2_p0_vtu 430c4762a1bSJed Brown requires: triangle 431c4762a1bSJed Brown args: -sol_type linear \ 432c4762a1bSJed Brown -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \ 433c4762a1bSJed Brown -dmsnes_check .001 -snes_error_if_not_converged \ 434c4762a1bSJed Brown -ksp_rtol 1e-10 -ksp_error_if_not_converged \ 435c4762a1bSJed Brown -potential_view vtk:multifield.vtu -exact_vec_view vtk:exact.vtu \ 436c4762a1bSJed Brown -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \ 437c4762a1bSJed Brown -fieldsplit_field_pc_type lu \ 438c4762a1bSJed Brown -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu 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 446c4762a1bSJed Brown args: -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 456c4762a1bSJed Brown args: -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