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 259371c9d4SSatish Balay typedef enum { 269371c9d4SSatish Balay RUN_FULL, 279371c9d4SSatish Balay RUN_TEST 289371c9d4SSatish Balay } RunType; 29c4762a1bSJed Brown 30c4762a1bSJed Brown typedef struct { 31c4762a1bSJed Brown RunType runType; /* Whether to run tests, or solve the full problem */ 32c4762a1bSJed Brown PetscBool useDualPenalty; /* Penalize deviation from both goals */ 33c4762a1bSJed Brown } AppCtx; 34c4762a1bSJed Brown 35d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 36d71ae5a4SJacob Faibussowitsch { 37c4762a1bSJed Brown const char *runTypes[2] = {"full", "test"}; 38c4762a1bSJed Brown PetscInt run; 39c4762a1bSJed Brown 40c4762a1bSJed Brown PetscFunctionBeginUser; 41c4762a1bSJed Brown options->runType = RUN_FULL; 42c4762a1bSJed Brown options->useDualPenalty = PETSC_FALSE; 43d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX"); 44c4762a1bSJed Brown run = options->runType; 459566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-run_type", "The run type", "ex2.c", runTypes, 2, runTypes[options->runType], &run, NULL)); 46c4762a1bSJed Brown options->runType = (RunType)run; 479566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_dual_penalty", "Penalize deviation from both goals", "ex2.c", options->useDualPenalty, &options->useDualPenalty, NULL)); 48d0609cedSBarry Smith PetscOptionsEnd(); 493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown 52d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 53d71ae5a4SJacob Faibussowitsch { 54c4762a1bSJed Brown PetscFunctionBeginUser; 559566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 569566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 579566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 589566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown 62d71ae5a4SJacob Faibussowitsch void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 63d71ae5a4SJacob Faibussowitsch { 64c4762a1bSJed Brown f0[0] = (u[0] - (x[0] * x[0] + x[1] * x[1])); 65c4762a1bSJed Brown } 66d71ae5a4SJacob Faibussowitsch void f0_u_full(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 67d71ae5a4SJacob Faibussowitsch { 689371c9d4SSatish Balay 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]))) + 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]))); 69c4762a1bSJed Brown } 70d71ae5a4SJacob Faibussowitsch void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 71d71ae5a4SJacob Faibussowitsch { 72c4762a1bSJed Brown PetscInt d; 73c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[dim * 2 + d]; 74c4762a1bSJed Brown } 75d71ae5a4SJacob Faibussowitsch void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 76d71ae5a4SJacob Faibussowitsch { 77c4762a1bSJed Brown g0[0] = 1.0; 78c4762a1bSJed Brown } 79d71ae5a4SJacob Faibussowitsch void g0_uu_full(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 80d71ae5a4SJacob Faibussowitsch { 819371c9d4SSatish Balay g0[0] = PetscSqr(u[0] - sin(2.0 * PETSC_PI * x[0]) * sin(2.0 * PETSC_PI * x[1])) + PetscSqr(u[0] - (x[0] * x[0] + x[1] * x[1])) - 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] + 4.0 * (x[0] * x[0] + x[1] * x[1]) * (sin(2.0 * PETSC_PI * x[0]) * sin(2.0 * PETSC_PI * x[1])); 82c4762a1bSJed Brown } 83d71ae5a4SJacob Faibussowitsch void g3_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 84d71ae5a4SJacob Faibussowitsch { 85c4762a1bSJed Brown PetscInt d; 86c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 89d71ae5a4SJacob Faibussowitsch void f0_a(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 90d71ae5a4SJacob Faibussowitsch { 91c4762a1bSJed Brown f0[0] = u[1] - 4.0 /* 0.0 */ + u[2]; 92c4762a1bSJed Brown } 93d71ae5a4SJacob Faibussowitsch void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 94d71ae5a4SJacob Faibussowitsch { 95c4762a1bSJed Brown g0[0] = 1.0; 96c4762a1bSJed Brown } 97d71ae5a4SJacob Faibussowitsch void g0_al(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 98d71ae5a4SJacob Faibussowitsch { 99c4762a1bSJed Brown g0[0] = 1.0; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102d71ae5a4SJacob Faibussowitsch void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 103d71ae5a4SJacob Faibussowitsch { 104c4762a1bSJed Brown f0[0] = u[1]; 105c4762a1bSJed Brown } 106d71ae5a4SJacob Faibussowitsch void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 107d71ae5a4SJacob Faibussowitsch { 108c4762a1bSJed Brown PetscInt d; 109c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 110c4762a1bSJed Brown } 111d71ae5a4SJacob Faibussowitsch void g0_la(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 112d71ae5a4SJacob Faibussowitsch { 113c4762a1bSJed Brown g0[0] = 1.0; 114c4762a1bSJed Brown } 115d71ae5a4SJacob Faibussowitsch void g3_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 116d71ae5a4SJacob Faibussowitsch { 117c4762a1bSJed Brown PetscInt d; 118c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* 122c4762a1bSJed Brown In 2D for Dirichlet conditions with a variable coefficient, we use exact solution: 123c4762a1bSJed Brown 124c4762a1bSJed Brown u = x^2 + y^2 125c4762a1bSJed Brown a = 4 126c4762a1bSJed Brown d_A = 4 127c4762a1bSJed Brown d_B = sin(2*pi*x[0]) * sin(2*pi*x[1]) 128c4762a1bSJed Brown 129c4762a1bSJed Brown so that 130c4762a1bSJed Brown 131c4762a1bSJed Brown -\Delta u + a = -4 + 4 = 0 132c4762a1bSJed Brown */ 133d71ae5a4SJacob Faibussowitsch PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 134d71ae5a4SJacob Faibussowitsch { 135c4762a1bSJed Brown *u = x[0] * x[0] + x[1] * x[1]; 1363ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 137c4762a1bSJed Brown } 138d71ae5a4SJacob Faibussowitsch PetscErrorCode constant_a_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx) 139d71ae5a4SJacob Faibussowitsch { 140c4762a1bSJed Brown *a = 4; 1413ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 142c4762a1bSJed Brown } 143d71ae5a4SJacob Faibussowitsch PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx) 144d71ae5a4SJacob Faibussowitsch { 145c4762a1bSJed Brown *l = 0.0; 1463ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149d71ae5a4SJacob Faibussowitsch PetscErrorCode SetupProblem(DM dm, AppCtx *user) 150d71ae5a4SJacob Faibussowitsch { 15145480ffeSMatthew G. Knepley PetscDS ds; 15245480ffeSMatthew G. Knepley DMLabel label; 153c4762a1bSJed Brown const PetscInt id = 1; 154c4762a1bSJed Brown 155c4762a1bSJed Brown PetscFunctionBeginUser; 1569566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1579566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u)); 1589566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_a, NULL)); 1599566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 2, f0_l, f1_l)); 1609566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL)); 1619566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul)); 1629566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL)); 1639566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 2, g0_al, NULL, NULL, NULL)); 1649566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 2, 1, g0_la, NULL, NULL, NULL)); 1659566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu)); 166c4762a1bSJed Brown 1679566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL)); 1689566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, constant_a_2d, NULL)); 1699566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 2, zero, NULL)); 1709566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 171*2b2f8cc6SPierre Jolivet PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)quadratic_u_2d, NULL, user, NULL)); 172*2b2f8cc6SPierre Jolivet PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (PetscVoidFn *)constant_a_2d, NULL, user, NULL)); 173*2b2f8cc6SPierre Jolivet PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (PetscVoidFn *)zero, NULL, user, NULL)); 1743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 177d71ae5a4SJacob Faibussowitsch PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 178d71ae5a4SJacob Faibussowitsch { 179c4762a1bSJed Brown DM cdm = dm; 180c4762a1bSJed Brown const PetscInt dim = 2; 181c4762a1bSJed Brown PetscFE fe[3]; 182c4762a1bSJed Brown PetscInt f; 183c4762a1bSJed Brown MPI_Comm comm; 184c4762a1bSJed Brown 185c4762a1bSJed Brown PetscFunctionBeginUser; 186c4762a1bSJed Brown /* Create finite element */ 1879566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 1889566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0])); 1899566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "potential")); 1909566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "charge_", -1, &fe[1])); 1919566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "charge")); 1929566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 1939566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2])); 1949566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[2], "multiplier")); 1959566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[2])); 196c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 1979566063dSJacob Faibussowitsch for (f = 0; f < 3; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe[f])); 1989566063dSJacob Faibussowitsch PetscCall(DMCreateDS(cdm)); 1999566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user)); 200c4762a1bSJed Brown while (cdm) { 2019566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 2029566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 203c4762a1bSJed Brown } 2049566063dSJacob Faibussowitsch for (f = 0; f < 3; ++f) PetscCall(PetscFEDestroy(&fe[f])); 2053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 208d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 209d71ae5a4SJacob Faibussowitsch { 210c4762a1bSJed Brown DM dm; 211c4762a1bSJed Brown SNES snes; 212c4762a1bSJed Brown Vec u, r; 213c4762a1bSJed Brown AppCtx user; 214c4762a1bSJed Brown 215327415f7SBarry Smith PetscFunctionBeginUser; 2169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2179566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 2189566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 2199566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2209566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 2219566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &user)); 222c4762a1bSJed Brown 2239566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 2249566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "solution")); 2259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 2266493148fSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user)); 2279566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 228c4762a1bSJed Brown 2299566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 230c4762a1bSJed Brown if (user.runType == RUN_FULL) { 231348a1646SMatthew G. Knepley PetscDS ds; 232348a1646SMatthew G. Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 233c4762a1bSJed Brown PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx); 234c4762a1bSJed Brown PetscReal error; 235c4762a1bSJed Brown 2369566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2379566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL)); 2389566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL)); 2399566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL)); 240c4762a1bSJed Brown initialGuess[0] = zero; 241c4762a1bSJed Brown initialGuess[1] = zero; 242c4762a1bSJed Brown initialGuess[2] = zero; 2439566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u)); 2449566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-initial_vec_view")); 2459566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 2469566063dSJacob Faibussowitsch if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n")); 24763a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", (double)error)); 2489566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 2499566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 2509566063dSJacob Faibussowitsch if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n")); 25163a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", (double)error)); 252c4762a1bSJed Brown } 2539566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 254c4762a1bSJed Brown 2559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 2569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 2579566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 2589566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2599566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 260b122ec5aSJacob Faibussowitsch return 0; 261c4762a1bSJed Brown } 262c4762a1bSJed Brown 263c4762a1bSJed Brown /*TEST 264c4762a1bSJed Brown 265c4762a1bSJed Brown build: 266c4762a1bSJed Brown requires: !complex triangle 267c4762a1bSJed Brown 268c4762a1bSJed Brown test: 269c4762a1bSJed Brown suffix: 0 270c4762a1bSJed Brown args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1 271c4762a1bSJed Brown 272c4762a1bSJed Brown test: 273c4762a1bSJed Brown suffix: 1 274c4762a1bSJed 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 275c4762a1bSJed Brown 276c4762a1bSJed Brown test: 277c4762a1bSJed Brown suffix: 2 278c4762a1bSJed 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 279c4762a1bSJed Brown 280c4762a1bSJed Brown TEST*/ 281