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; 445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-run_type", "The run type", "ex2.c", runTypes, 2, runTypes[options->runType], &run, NULL)); 45c4762a1bSJed Brown options->runType = (RunType) run; 465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-use_dual_penalty", "Penalize deviation from both goals", "ex2.c", options->useDualPenalty, &options->useDualPenalty, NULL)); 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 PetscFunctionBeginUser; 545f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 58c4762a1bSJed Brown PetscFunctionReturn(0); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61c4762a1bSJed Brown void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 62c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 63c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 64c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 65c4762a1bSJed Brown { 66c4762a1bSJed Brown f0[0] = (u[0] - (x[0]*x[0] + x[1]*x[1])); 67c4762a1bSJed Brown } 68c4762a1bSJed Brown void f0_u_full(PetscInt dim, PetscInt Nf, PetscInt NfAux, 69c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 70c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 71c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 72c4762a1bSJed Brown { 73c4762a1bSJed 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]))) + 74c4762a1bSJed 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]))); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown void f1_u(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 f1[]) 80c4762a1bSJed Brown { 81c4762a1bSJed Brown PetscInt d; 82c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[dim*2+d]; 83c4762a1bSJed Brown } 84c4762a1bSJed Brown void g0_uu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 88c4762a1bSJed Brown { 89c4762a1bSJed Brown g0[0] = 1.0; 90c4762a1bSJed Brown } 91c4762a1bSJed Brown void g0_uu_full(PetscInt dim, PetscInt Nf, PetscInt NfAux, 92c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 93c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 94c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 95c4762a1bSJed Brown { 96c4762a1bSJed Brown g0[0] = PetscSqr(u[0] - sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1])) 97c4762a1bSJed Brown + PetscSqr(u[0] - (x[0]*x[0] + x[1]*x[1])) 98c4762a1bSJed 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] 99c4762a1bSJed 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])); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown void g3_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux, 102c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 103c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 104c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 105c4762a1bSJed Brown { 106c4762a1bSJed Brown PetscInt d; 107c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 110c4762a1bSJed Brown void f0_a(PetscInt dim, PetscInt Nf, PetscInt NfAux, 111c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 112c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 113c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 114c4762a1bSJed Brown { 115c4762a1bSJed Brown f0[0] = u[1] - 4.0 /* 0.0 */ + u[2]; 116c4762a1bSJed Brown } 117c4762a1bSJed Brown void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux, 118c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 119c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 120c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 121c4762a1bSJed Brown { 122c4762a1bSJed Brown g0[0] = 1.0; 123c4762a1bSJed Brown } 124c4762a1bSJed Brown void g0_al(PetscInt dim, PetscInt Nf, PetscInt NfAux, 125c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 126c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 127c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 128c4762a1bSJed Brown { 129c4762a1bSJed Brown g0[0] = 1.0; 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132c4762a1bSJed Brown void f0_l(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, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 136c4762a1bSJed Brown { 137c4762a1bSJed Brown f0[0] = u[1]; 138c4762a1bSJed Brown } 139c4762a1bSJed Brown void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, 140c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 141c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 142c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 143c4762a1bSJed Brown { 144c4762a1bSJed Brown PetscInt d; 145c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown void g0_la(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 151c4762a1bSJed Brown { 152c4762a1bSJed Brown g0[0] = 1.0; 153c4762a1bSJed Brown } 154c4762a1bSJed Brown void g3_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 155c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 156c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 157c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 158c4762a1bSJed Brown { 159c4762a1bSJed Brown PetscInt d; 160c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163c4762a1bSJed Brown /* 164c4762a1bSJed Brown In 2D for Dirichlet conditions with a variable coefficient, we use exact solution: 165c4762a1bSJed Brown 166c4762a1bSJed Brown u = x^2 + y^2 167c4762a1bSJed Brown a = 4 168c4762a1bSJed Brown d_A = 4 169c4762a1bSJed Brown d_B = sin(2*pi*x[0]) * sin(2*pi*x[1]) 170c4762a1bSJed Brown 171c4762a1bSJed Brown so that 172c4762a1bSJed Brown 173c4762a1bSJed Brown -\Delta u + a = -4 + 4 = 0 174c4762a1bSJed Brown */ 175c4762a1bSJed Brown PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 176c4762a1bSJed Brown { 177c4762a1bSJed Brown *u = x[0]*x[0] + x[1]*x[1]; 178c4762a1bSJed Brown return 0; 179c4762a1bSJed Brown } 180c4762a1bSJed Brown PetscErrorCode constant_a_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx) 181c4762a1bSJed Brown { 182c4762a1bSJed Brown *a = 4; 183c4762a1bSJed Brown return 0; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx) 186c4762a1bSJed Brown { 187c4762a1bSJed Brown *l = 0.0; 188c4762a1bSJed Brown return 0; 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191c4762a1bSJed Brown PetscErrorCode SetupProblem(DM dm, AppCtx *user) 192c4762a1bSJed Brown { 19345480ffeSMatthew G. Knepley PetscDS ds; 19445480ffeSMatthew G. Knepley DMLabel label; 195c4762a1bSJed Brown const PetscInt id = 1; 196c4762a1bSJed Brown 197c4762a1bSJed Brown PetscFunctionBeginUser; 1985f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 1995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u)); 2005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 1, f0_a, NULL)); 2015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 2, f0_l, f1_l)); 2025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL)); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul)); 2045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 1, 2, g0_al, NULL, NULL, NULL)); 2065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 2, 1, g0_la, NULL, NULL, NULL)); 2075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu)); 208c4762a1bSJed Brown 2095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(ds, 1, constant_a_2d, NULL)); 2115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(ds, 2, zero, NULL)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 2135f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)()) quadratic_u_2d, NULL, user, NULL)); 2145f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)()) constant_a_2d, NULL, user, NULL)); 2155f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)()) zero, NULL, user, NULL)); 216c4762a1bSJed Brown PetscFunctionReturn(0); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 219c4762a1bSJed Brown PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 220c4762a1bSJed Brown { 221c4762a1bSJed Brown DM cdm = dm; 222c4762a1bSJed Brown const PetscInt dim = 2; 223c4762a1bSJed Brown PetscFE fe[3]; 224c4762a1bSJed Brown PetscInt f; 225c4762a1bSJed Brown MPI_Comm comm; 226c4762a1bSJed Brown 227c4762a1bSJed Brown PetscFunctionBeginUser; 228c4762a1bSJed Brown /* Create finite element */ 2295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm)); 2305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0])); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe[0], "potential")); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "charge_", -1, &fe[1])); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe[1], "charge")); 2345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECopyQuadrature(fe[0], fe[1])); 2355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2])); 2365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe[2], "multiplier")); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECopyQuadrature(fe[0], fe[2])); 238c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 2395f80ce2aSJacob Faibussowitsch for (f = 0; f < 3; ++f) CHKERRQ(DMSetField(dm, f, NULL, (PetscObject) fe[f])); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(cdm)); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(SetupProblem(dm, user)); 242c4762a1bSJed Brown while (cdm) { 2435f80ce2aSJacob Faibussowitsch CHKERRQ(DMCopyDisc(dm, cdm)); 2445f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 245c4762a1bSJed Brown } 2465f80ce2aSJacob Faibussowitsch for (f = 0; f < 3; ++f) CHKERRQ(PetscFEDestroy(&fe[f])); 247c4762a1bSJed Brown PetscFunctionReturn(0); 248c4762a1bSJed Brown } 249c4762a1bSJed Brown 250c4762a1bSJed Brown int main(int argc, char **argv) 251c4762a1bSJed Brown { 252c4762a1bSJed Brown DM dm; 253c4762a1bSJed Brown SNES snes; 254c4762a1bSJed Brown Vec u, r; 255c4762a1bSJed Brown AppCtx user; 256c4762a1bSJed Brown 257*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL,help)); 2585f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 2595f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD, &snes)); 2605f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2615f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetDM(snes, dm)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(SetupDiscretization(dm, &user)); 263c4762a1bSJed Brown 2645f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &u)); 2655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) u, "solution")); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u, &r)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetSNESLocalFEM(dm,&user,&user,&user)); 2685f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 269c4762a1bSJed Brown 2705f80ce2aSJacob Faibussowitsch CHKERRQ(DMSNESCheckFromOptions(snes, u)); 271c4762a1bSJed Brown if (user.runType == RUN_FULL) { 272348a1646SMatthew G. Knepley PetscDS ds; 273348a1646SMatthew G. Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 274c4762a1bSJed Brown PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx); 275c4762a1bSJed Brown PetscReal error; 276c4762a1bSJed Brown 2775f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 2785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL)); 2795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL)); 2805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL)); 281c4762a1bSJed Brown initialGuess[0] = zero; 282c4762a1bSJed Brown initialGuess[1] = zero; 283c4762a1bSJed Brown initialGuess[2] = zero; 2845f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u)); 2855f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-initial_vec_view")); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 2875f80ce2aSJacob Faibussowitsch if (error < 1.0e-11) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n")); 2885f80ce2aSJacob Faibussowitsch else CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", error)); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes, NULL, u)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 2915f80ce2aSJacob Faibussowitsch if (error < 1.0e-11) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n")); 2925f80ce2aSJacob Faibussowitsch else CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", error)); 293c4762a1bSJed Brown } 2945f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-sol_vec_view")); 295c4762a1bSJed Brown 2965f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 2975f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 2985f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 2995f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 300*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 301*b122ec5aSJacob Faibussowitsch return 0; 302c4762a1bSJed Brown } 303c4762a1bSJed Brown 304c4762a1bSJed Brown /*TEST 305c4762a1bSJed Brown 306c4762a1bSJed Brown build: 307c4762a1bSJed Brown requires: !complex triangle 308c4762a1bSJed Brown 309c4762a1bSJed Brown test: 310c4762a1bSJed Brown suffix: 0 311c4762a1bSJed Brown args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1 312c4762a1bSJed Brown 313c4762a1bSJed Brown test: 314c4762a1bSJed Brown suffix: 1 315c4762a1bSJed 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 316c4762a1bSJed Brown 317c4762a1bSJed Brown test: 318c4762a1bSJed Brown suffix: 2 319c4762a1bSJed 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 320c4762a1bSJed Brown 321c4762a1bSJed Brown TEST*/ 322