1c4762a1bSJed Brown static char help[] = "One-Shot Multigrid for Parameter Estimation Problem for the Poisson Equation.\n\ 2c4762a1bSJed Brown Using the Interior Point Method.\n\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*F 5c4762a1bSJed Brown We are solving the parameter estimation problem for the Laplacian. We will ask to minimize a Lagrangian 6c4762a1bSJed Brown function over $y$ and $u$, given by 7c4762a1bSJed Brown \begin{align} 8c4762a1bSJed Brown L(u, a, \lambda) = \frac{1}{2} || Qu - d_A ||^2 || Qu - d_B ||^2 + \frac{\beta}{2} || L (a - a_r) ||^2 + \lambda F(u; a) 9c4762a1bSJed Brown \end{align} 10c4762a1bSJed Brown where $Q$ is a sampling operator, $L$ is a regularization operator, $F$ defines the PDE. 11c4762a1bSJed Brown 12c4762a1bSJed Brown Currently, we have perfect information, meaning $Q = I$, and then we need no regularization, $L = I$. We 13c4762a1bSJed Brown also give the null vector for the reference control $a_r$. Right now $\beta = 1$. 14c4762a1bSJed Brown 15c4762a1bSJed Brown The PDE will be the Laplace equation with homogeneous boundary conditions 16c4762a1bSJed Brown \begin{align} 17c4762a1bSJed Brown -Delta u = a 18c4762a1bSJed Brown \end{align} 19c4762a1bSJed Brown 20c4762a1bSJed Brown F*/ 21c4762a1bSJed Brown 22c4762a1bSJed Brown #include <petsc.h> 23c4762a1bSJed Brown #include <petscfe.h> 24c4762a1bSJed Brown 25c4762a1bSJed Brown typedef enum {RUN_FULL, RUN_TEST} RunType; 26c4762a1bSJed Brown 27c4762a1bSJed Brown typedef struct { 28c4762a1bSJed Brown RunType runType; /* Whether to run tests, or solve the full problem */ 29c4762a1bSJed Brown PetscBool useDualPenalty; /* Penalize deviation from both goals */ 30c4762a1bSJed Brown } AppCtx; 31c4762a1bSJed Brown 32c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 33c4762a1bSJed Brown { 34c4762a1bSJed Brown const char *runTypes[2] = {"full", "test"}; 35c4762a1bSJed Brown PetscInt run; 36c4762a1bSJed Brown PetscErrorCode ierr; 37c4762a1bSJed Brown 38c4762a1bSJed Brown PetscFunctionBeginUser; 39c4762a1bSJed Brown options->runType = RUN_FULL; 40c4762a1bSJed Brown options->useDualPenalty = PETSC_FALSE; 41c4762a1bSJed Brown 42c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");CHKERRQ(ierr); 43c4762a1bSJed Brown run = options->runType; 44c4762a1bSJed Brown ierr = PetscOptionsEList("-run_type", "The run type", "ex2.c", runTypes, 2, runTypes[options->runType], &run, NULL);CHKERRQ(ierr); 45c4762a1bSJed Brown options->runType = (RunType) run; 46c4762a1bSJed Brown ierr = PetscOptionsBool("-use_dual_penalty", "Penalize deviation from both goals", "ex2.c", options->useDualPenalty, &options->useDualPenalty, NULL);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 48c4762a1bSJed Brown PetscFunctionReturn(0); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 51c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 52c4762a1bSJed Brown { 53c4762a1bSJed Brown DM distributedMesh = NULL; 54c4762a1bSJed Brown PetscErrorCode ierr; 55c4762a1bSJed Brown 56c4762a1bSJed Brown PetscFunctionBeginUser; 57c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 59c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr); 60c4762a1bSJed Brown if (distributedMesh) { 61c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 62c4762a1bSJed Brown *dm = distributedMesh; 63c4762a1bSJed Brown } 64c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 65c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 66c4762a1bSJed Brown PetscFunctionReturn(0); 67c4762a1bSJed Brown } 68c4762a1bSJed Brown 69c4762a1bSJed Brown void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 70c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 71c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 72c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 73c4762a1bSJed Brown { 74c4762a1bSJed Brown f0[0] = (u[0] - (x[0]*x[0] + x[1]*x[1])); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown void f0_u_full(PetscInt dim, PetscInt Nf, PetscInt NfAux, 77c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 78c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 79c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 80c4762a1bSJed Brown { 81c4762a1bSJed Brown f0[0] = (u[0] - (x[0]*x[0] + x[1]*x[1]))*PetscSqr(u[0] - (sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1]))) + 82c4762a1bSJed Brown PetscSqr(u[0] - (x[0]*x[0] + x[1]*x[1]))*(u[0] - (sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1]))); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 85c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 86c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 87c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 88c4762a1bSJed Brown { 89c4762a1bSJed Brown PetscInt d; 90c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[dim*2+d]; 91c4762a1bSJed Brown } 92c4762a1bSJed Brown void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 93c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 94c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 95c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 96c4762a1bSJed Brown { 97c4762a1bSJed Brown g0[0] = 1.0; 98c4762a1bSJed Brown } 99c4762a1bSJed Brown void g0_uu_full(PetscInt dim, PetscInt Nf, PetscInt NfAux, 100c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 101c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 102c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 103c4762a1bSJed Brown { 104c4762a1bSJed Brown g0[0] = PetscSqr(u[0] - sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1])) 105c4762a1bSJed Brown + PetscSqr(u[0] - (x[0]*x[0] + x[1]*x[1])) 106c4762a1bSJed Brown - 2.0*((x[0]*x[0] + x[1]*x[1]) + (sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1])))*u[0] 107c4762a1bSJed Brown + 4.0*(x[0]*x[0] + x[1]*x[1])*(sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1])); 108c4762a1bSJed Brown } 109c4762a1bSJed Brown void g3_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux, 110c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 111c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 112c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 113c4762a1bSJed Brown { 114c4762a1bSJed Brown PetscInt d; 115c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 116c4762a1bSJed Brown } 117c4762a1bSJed Brown 118c4762a1bSJed Brown void f0_a(PetscInt dim, PetscInt Nf, PetscInt NfAux, 119c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 120c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 121c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 122c4762a1bSJed Brown { 123c4762a1bSJed Brown f0[0] = u[1] - 4.0 /* 0.0 */ + u[2]; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux, 126c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 127c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 128c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 129c4762a1bSJed Brown { 130c4762a1bSJed Brown g0[0] = 1.0; 131c4762a1bSJed Brown } 132c4762a1bSJed Brown void g0_al(PetscInt dim, PetscInt Nf, PetscInt NfAux, 133c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 134c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 135c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 136c4762a1bSJed Brown { 137c4762a1bSJed Brown g0[0] = 1.0; 138c4762a1bSJed Brown } 139c4762a1bSJed Brown 140c4762a1bSJed Brown void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, 141c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 142c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 143c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 144c4762a1bSJed Brown { 145c4762a1bSJed Brown f0[0] = u[1]; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, 148c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 149c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 150c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 151c4762a1bSJed Brown { 152c4762a1bSJed Brown PetscInt d; 153c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 154c4762a1bSJed Brown } 155c4762a1bSJed Brown void g0_la(PetscInt dim, PetscInt Nf, PetscInt NfAux, 156c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 157c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 158c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 159c4762a1bSJed Brown { 160c4762a1bSJed Brown g0[0] = 1.0; 161c4762a1bSJed Brown } 162c4762a1bSJed Brown void g3_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 163c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 164c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 165c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 166c4762a1bSJed Brown { 167c4762a1bSJed Brown PetscInt d; 168c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 169c4762a1bSJed Brown } 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* 172c4762a1bSJed Brown In 2D for Dirichlet conditions with a variable coefficient, we use exact solution: 173c4762a1bSJed Brown 174c4762a1bSJed Brown u = x^2 + y^2 175c4762a1bSJed Brown a = 4 176c4762a1bSJed Brown d_A = 4 177c4762a1bSJed Brown d_B = sin(2*pi*x[0]) * sin(2*pi*x[1]) 178c4762a1bSJed Brown 179c4762a1bSJed Brown so that 180c4762a1bSJed Brown 181c4762a1bSJed Brown -\Delta u + a = -4 + 4 = 0 182c4762a1bSJed Brown */ 183c4762a1bSJed Brown PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 184c4762a1bSJed Brown { 185c4762a1bSJed Brown *u = x[0]*x[0] + x[1]*x[1]; 186c4762a1bSJed Brown return 0; 187c4762a1bSJed Brown } 188c4762a1bSJed Brown PetscErrorCode constant_a_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx) 189c4762a1bSJed Brown { 190c4762a1bSJed Brown *a = 4; 191c4762a1bSJed Brown return 0; 192c4762a1bSJed Brown } 193c4762a1bSJed Brown PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx) 194c4762a1bSJed Brown { 195c4762a1bSJed Brown *l = 0.0; 196c4762a1bSJed Brown return 0; 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199c4762a1bSJed Brown PetscErrorCode SetupProblem(DM dm, AppCtx *user) 200c4762a1bSJed Brown { 201*45480ffeSMatthew G. Knepley PetscDS ds; 202*45480ffeSMatthew G. Knepley DMLabel label; 203c4762a1bSJed Brown const PetscInt id = 1; 204c4762a1bSJed Brown PetscErrorCode ierr; 205c4762a1bSJed Brown 206c4762a1bSJed Brown PetscFunctionBeginUser; 207*45480ffeSMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 208*45480ffeSMatthew G. Knepley ierr = PetscDSSetResidual(ds, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u);CHKERRQ(ierr); 209*45480ffeSMatthew G. Knepley ierr = PetscDSSetResidual(ds, 1, f0_a, NULL);CHKERRQ(ierr); 210*45480ffeSMatthew G. Knepley ierr = PetscDSSetResidual(ds, 2, f0_l, f1_l);CHKERRQ(ierr); 211*45480ffeSMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL);CHKERRQ(ierr); 212*45480ffeSMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul);CHKERRQ(ierr); 213*45480ffeSMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL);CHKERRQ(ierr); 214*45480ffeSMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 1, 2, g0_al, NULL, NULL, NULL);CHKERRQ(ierr); 215*45480ffeSMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 2, 1, g0_la, NULL, NULL, NULL);CHKERRQ(ierr); 216*45480ffeSMatthew G. Knepley ierr = PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu);CHKERRQ(ierr); 217c4762a1bSJed Brown 218*45480ffeSMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL);CHKERRQ(ierr); 219*45480ffeSMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 1, constant_a_2d, NULL);CHKERRQ(ierr); 220*45480ffeSMatthew G. Knepley ierr = PetscDSSetExactSolution(ds, 2, zero, NULL);CHKERRQ(ierr); 221*45480ffeSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 222*45480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)()) quadratic_u_2d, NULL, user, NULL);CHKERRQ(ierr); 223*45480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)()) constant_a_2d, NULL, user, NULL);CHKERRQ(ierr); 224*45480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)()) zero, NULL, user, NULL);CHKERRQ(ierr); 225c4762a1bSJed Brown PetscFunctionReturn(0); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown 228c4762a1bSJed Brown PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 229c4762a1bSJed Brown { 230c4762a1bSJed Brown DM cdm = dm; 231c4762a1bSJed Brown const PetscInt dim = 2; 232c4762a1bSJed Brown PetscFE fe[3]; 233c4762a1bSJed Brown PetscInt f; 234c4762a1bSJed Brown MPI_Comm comm; 235c4762a1bSJed Brown PetscErrorCode ierr; 236c4762a1bSJed Brown 237c4762a1bSJed Brown PetscFunctionBeginUser; 238c4762a1bSJed Brown /* Create finite element */ 239c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]);CHKERRQ(ierr); 241c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe[0], "potential");CHKERRQ(ierr); 242c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "charge_", -1, &fe[1]);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe[1], "charge");CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 245c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]);CHKERRQ(ierr); 246c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe[2], "multiplier");CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr); 248c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 249c4762a1bSJed Brown for (f = 0; f < 3; ++f) {ierr = DMSetField(dm, f, NULL, (PetscObject) fe[f]);CHKERRQ(ierr);} 250c4762a1bSJed Brown ierr = DMCreateDS(cdm);CHKERRQ(ierr); 251c4762a1bSJed Brown ierr = SetupProblem(dm, user);CHKERRQ(ierr); 252c4762a1bSJed Brown while (cdm) { 253c4762a1bSJed Brown ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 254c4762a1bSJed Brown ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 255c4762a1bSJed Brown } 256c4762a1bSJed Brown for (f = 0; f < 3; ++f) {ierr = PetscFEDestroy(&fe[f]);CHKERRQ(ierr);} 257c4762a1bSJed Brown PetscFunctionReturn(0); 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 260c4762a1bSJed Brown int main(int argc, char **argv) 261c4762a1bSJed Brown { 262c4762a1bSJed Brown DM dm; 263c4762a1bSJed Brown SNES snes; 264c4762a1bSJed Brown Vec u, r; 265c4762a1bSJed Brown AppCtx user; 266c4762a1bSJed Brown PetscErrorCode ierr; 267c4762a1bSJed Brown 268c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 269c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 270c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 272c4762a1bSJed Brown ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 273c4762a1bSJed Brown ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr); 274c4762a1bSJed Brown 275c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 276c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "solution");CHKERRQ(ierr); 277c4762a1bSJed Brown ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 278c4762a1bSJed Brown ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr); 279c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 280c4762a1bSJed Brown 281348a1646SMatthew G. Knepley ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr); 282c4762a1bSJed Brown if (user.runType == RUN_FULL) { 283348a1646SMatthew G. Knepley PetscDS ds; 284348a1646SMatthew G. Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 285c4762a1bSJed Brown PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx); 286c4762a1bSJed Brown PetscReal error; 287c4762a1bSJed Brown 288348a1646SMatthew G. Knepley ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 289348a1646SMatthew G. Knepley ierr = PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL);CHKERRQ(ierr); 290348a1646SMatthew G. Knepley ierr = PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL);CHKERRQ(ierr); 291348a1646SMatthew G. Knepley ierr = PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL);CHKERRQ(ierr); 292c4762a1bSJed Brown initialGuess[0] = zero; 293c4762a1bSJed Brown initialGuess[1] = zero; 294c4762a1bSJed Brown initialGuess[2] = zero; 295c4762a1bSJed Brown ierr = DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr); 296c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-initial_vec_view");CHKERRQ(ierr); 297348a1646SMatthew G. Knepley ierr = DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);CHKERRQ(ierr); 298c4762a1bSJed Brown if (error < 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);} 299c4762a1bSJed Brown else {ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", error);CHKERRQ(ierr);} 300c4762a1bSJed Brown ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 301348a1646SMatthew G. Knepley ierr = DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);CHKERRQ(ierr); 302c4762a1bSJed Brown if (error < 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);} 303c4762a1bSJed Brown else {ierr = PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", error);CHKERRQ(ierr);} 304c4762a1bSJed Brown } 305c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 306c4762a1bSJed Brown 307c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 308c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 309c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 310c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 311c4762a1bSJed Brown ierr = PetscFinalize(); 312c4762a1bSJed Brown return ierr; 313c4762a1bSJed Brown } 314c4762a1bSJed Brown 315c4762a1bSJed Brown /*TEST 316c4762a1bSJed Brown 317c4762a1bSJed Brown build: 318c4762a1bSJed Brown requires: !complex triangle 319c4762a1bSJed Brown 320c4762a1bSJed Brown test: 321c4762a1bSJed Brown suffix: 0 322c4762a1bSJed Brown args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1 323c4762a1bSJed Brown 324c4762a1bSJed Brown test: 325c4762a1bSJed Brown suffix: 1 326c4762a1bSJed Brown args: -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1 -snes_monitor -snes_converged_reason -pc_type fieldsplit -pc_fieldsplit_0_fields 0,1 -pc_fieldsplit_1_fields 2 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition selfp -fieldsplit_0_pc_type lu -sol_vec_view 327c4762a1bSJed Brown 328c4762a1bSJed Brown test: 329c4762a1bSJed Brown suffix: 2 330c4762a1bSJed Brown args: -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1 -snes_monitor -snes_converged_reason -snes_fd -pc_type fieldsplit -pc_fieldsplit_0_fields 0,1 -pc_fieldsplit_1_fields 2 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition selfp -fieldsplit_0_pc_type lu -sol_vec_view 331c4762a1bSJed Brown 332c4762a1bSJed Brown TEST*/ 333