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 $a$ and $u$, given by 7c4762a1bSJed Brown \begin{align} 8c4762a1bSJed Brown L(u, a, \lambda) = \frac{1}{2} || Qu - d ||^2 + \frac{1}{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 exact control for the reference $a_r$. 14c4762a1bSJed Brown 15c4762a1bSJed Brown The PDE will be the Laplace equation with homogeneous boundary conditions 16c4762a1bSJed Brown \begin{align} 17c4762a1bSJed Brown -nabla \cdot a \nabla u = f 18c4762a1bSJed Brown \end{align} 19c4762a1bSJed Brown 20c4762a1bSJed Brown F*/ 21c4762a1bSJed Brown 22c4762a1bSJed Brown #include <petsc.h> 23c4762a1bSJed Brown #include <petscfe.h> 24c4762a1bSJed Brown 25*9371c9d4SSatish Balay typedef enum { 26*9371c9d4SSatish Balay RUN_FULL, 27*9371c9d4SSatish Balay RUN_TEST 28*9371c9d4SSatish Balay } RunType; 29c4762a1bSJed Brown 30c4762a1bSJed Brown typedef struct { 31c4762a1bSJed Brown RunType runType; /* Whether to run tests, or solve the full problem */ 32c4762a1bSJed Brown } AppCtx; 33c4762a1bSJed Brown 34*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 35c4762a1bSJed Brown const char *runTypes[2] = {"full", "test"}; 36c4762a1bSJed Brown PetscInt run; 37c4762a1bSJed Brown 38c4762a1bSJed Brown PetscFunctionBeginUser; 39c4762a1bSJed Brown options->runType = RUN_FULL; 40d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX"); 41c4762a1bSJed Brown run = options->runType; 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-run_type", "The run type", "ex1.c", runTypes, 2, runTypes[options->runType], &run, NULL)); 43c4762a1bSJed Brown options->runType = (RunType)run; 44d0609cedSBarry Smith PetscOptionsEnd(); 45c4762a1bSJed Brown PetscFunctionReturn(0); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 49c4762a1bSJed Brown PetscFunctionBeginUser; 509566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 519566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 529566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 539566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 54c4762a1bSJed Brown PetscFunctionReturn(0); 55c4762a1bSJed Brown } 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* u - (x^2 + y^2) */ 58*9371c9d4SSatish Balay 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[]) { 59c4762a1bSJed Brown f0[0] = u[0] - (x[0] * x[0] + x[1] * x[1]); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown /* a \nabla\lambda */ 62*9371c9d4SSatish Balay 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[]) { 63c4762a1bSJed Brown PetscInt d; 64c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u[1] * u_x[dim * 2 + d]; 65c4762a1bSJed Brown } 66c4762a1bSJed Brown /* I */ 67*9371c9d4SSatish Balay 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[]) { 68c4762a1bSJed Brown g0[0] = 1.0; 69c4762a1bSJed Brown } 70c4762a1bSJed Brown /* \nabla */ 71*9371c9d4SSatish Balay void g2_ua(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 g2[]) { 72c4762a1bSJed Brown PetscInt d; 73c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d] = u_x[dim * 2 + d]; 74c4762a1bSJed Brown } 75c4762a1bSJed Brown /* a */ 76*9371c9d4SSatish Balay 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[]) { 77c4762a1bSJed Brown PetscInt d; 78c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d * dim + d] = u[1]; 79c4762a1bSJed Brown } 80c4762a1bSJed Brown /* a - (x + y) */ 81*9371c9d4SSatish Balay 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[]) { 82c4762a1bSJed Brown f0[0] = u[1] - (x[0] + x[1]); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown /* \lambda \nabla u */ 85*9371c9d4SSatish Balay void f1_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 f1[]) { 86c4762a1bSJed Brown PetscInt d; 87c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u[2] * u_x[d]; 88c4762a1bSJed Brown } 89c4762a1bSJed Brown /* I */ 90*9371c9d4SSatish Balay 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[]) { 91c4762a1bSJed Brown g0[0] = 1.0; 92c4762a1bSJed Brown } 93c4762a1bSJed Brown /* 6 (x + y) */ 94*9371c9d4SSatish Balay 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[]) { 95c4762a1bSJed Brown f0[0] = 6.0 * (x[0] + x[1]); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown /* a \nabla u */ 98*9371c9d4SSatish Balay 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[]) { 99c4762a1bSJed Brown PetscInt d; 100c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u[1] * u_x[d]; 101c4762a1bSJed Brown } 102c4762a1bSJed Brown /* \nabla u */ 103*9371c9d4SSatish Balay void g2_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 g2[]) { 104c4762a1bSJed Brown PetscInt d; 105c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d] = u_x[d]; 106c4762a1bSJed Brown } 107c4762a1bSJed Brown /* a */ 108*9371c9d4SSatish Balay 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[]) { 109c4762a1bSJed Brown PetscInt d; 110c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d * dim + d] = u[1]; 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* 114c4762a1bSJed Brown In 2D for Dirichlet conditions with a variable coefficient, we use exact solution: 115c4762a1bSJed Brown 116c4762a1bSJed Brown u = x^2 + y^2 117c4762a1bSJed Brown f = 6 (x + y) 118c4762a1bSJed Brown kappa(a) = a = (x + y) 119c4762a1bSJed Brown 120c4762a1bSJed Brown so that 121c4762a1bSJed Brown 122c4762a1bSJed Brown -\div \kappa(a) \grad u + f = -6 (x + y) + 6 (x + y) = 0 123c4762a1bSJed Brown */ 124*9371c9d4SSatish Balay PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 125c4762a1bSJed Brown *u = x[0] * x[0] + x[1] * x[1]; 126c4762a1bSJed Brown return 0; 127c4762a1bSJed Brown } 128*9371c9d4SSatish Balay PetscErrorCode linear_a_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx) { 129c4762a1bSJed Brown *a = x[0] + x[1]; 130c4762a1bSJed Brown return 0; 131c4762a1bSJed Brown } 132*9371c9d4SSatish Balay PetscErrorCode zero(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx) { 133c4762a1bSJed Brown *l = 0.0; 134c4762a1bSJed Brown return 0; 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137*9371c9d4SSatish Balay PetscErrorCode SetupProblem(DM dm, AppCtx *user) { 13845480ffeSMatthew G. Knepley PetscDS ds; 13945480ffeSMatthew G. Knepley DMLabel label; 140c4762a1bSJed Brown const PetscInt id = 1; 141c4762a1bSJed Brown 142c4762a1bSJed Brown PetscFunctionBeginUser; 1439566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1449566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_u, f1_u)); 1459566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_a, f1_a)); 1469566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 2, f0_l, f1_l)); 1479566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, NULL)); 1489566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_ua, NULL)); 1499566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul)); 1509566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL)); 1519566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 2, 1, NULL, NULL, g2_la, NULL)); 1529566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu)); 153c4762a1bSJed Brown 1549566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL)); 1559566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, linear_a_2d, NULL)); 1569566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 2, zero, NULL)); 1579566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 1589566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))quadratic_u_2d, NULL, user, NULL)); 1599566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)(void))linear_a_2d, NULL, user, NULL)); 1609566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)(void))zero, NULL, user, NULL)); 161c4762a1bSJed Brown PetscFunctionReturn(0); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown 164*9371c9d4SSatish Balay PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) { 165c4762a1bSJed Brown DM cdm = dm; 166c4762a1bSJed Brown const PetscInt dim = 2; 167c4762a1bSJed Brown PetscFE fe[3]; 168c4762a1bSJed Brown PetscInt f; 169c4762a1bSJed Brown MPI_Comm comm; 170c4762a1bSJed Brown 171c4762a1bSJed Brown PetscFunctionBeginUser; 172c4762a1bSJed Brown /* Create finite element */ 1739566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 1749566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0])); 1759566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "potential")); 1769566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "conductivity_", -1, &fe[1])); 1779566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "conductivity")); 1789566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 1799566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2])); 1809566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[2], "multiplier")); 1819566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[2])); 182c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 1839566063dSJacob Faibussowitsch for (f = 0; f < 3; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe[f])); 1849566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 1859566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user)); 186c4762a1bSJed Brown while (cdm) { 1879566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 1889566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 189c4762a1bSJed Brown } 1909566063dSJacob Faibussowitsch for (f = 0; f < 3; ++f) PetscCall(PetscFEDestroy(&fe[f])); 191c4762a1bSJed Brown PetscFunctionReturn(0); 192c4762a1bSJed Brown } 193c4762a1bSJed Brown 194*9371c9d4SSatish Balay int main(int argc, char **argv) { 195c4762a1bSJed Brown DM dm; 196c4762a1bSJed Brown SNES snes; 197c4762a1bSJed Brown Vec u, r; 198c4762a1bSJed Brown AppCtx user; 199c4762a1bSJed Brown 200327415f7SBarry Smith PetscFunctionBeginUser; 2019566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2029566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 2039566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 2049566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 2059566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 2069566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &user)); 207c4762a1bSJed Brown 2089566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 2099566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "solution")); 2109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 2119566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 2129566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 213c4762a1bSJed Brown 2149566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 215c4762a1bSJed Brown if (user.runType == RUN_FULL) { 216348a1646SMatthew G. Knepley PetscDS ds; 217348a1646SMatthew G. Knepley PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 218c4762a1bSJed Brown PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx); 219c4762a1bSJed Brown PetscReal error; 220c4762a1bSJed Brown 2219566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 2229566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL)); 2239566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL)); 2249566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL)); 225c4762a1bSJed Brown initialGuess[0] = zero; 226c4762a1bSJed Brown initialGuess[1] = zero; 227c4762a1bSJed Brown initialGuess[2] = zero; 2289566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u)); 2299566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-initial_vec_view")); 2309566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 2319566063dSJacob Faibussowitsch if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n")); 23263a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", (double)error)); 2339566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 2349566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error)); 2359566063dSJacob Faibussowitsch if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n")); 23663a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", (double)error)); 237c4762a1bSJed Brown } 2389566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 239c4762a1bSJed Brown 2409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 2419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 2429566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 2439566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2449566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 245b122ec5aSJacob Faibussowitsch return 0; 246c4762a1bSJed Brown } 247c4762a1bSJed Brown 248c4762a1bSJed Brown /*TEST 249c4762a1bSJed Brown 250c4762a1bSJed Brown build: 251c4762a1bSJed Brown requires: !complex 252c4762a1bSJed Brown 253c4762a1bSJed Brown test: 254c4762a1bSJed Brown suffix: 0 255c4762a1bSJed Brown requires: triangle 256c4762a1bSJed Brown args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -conductivity_petscspace_degree 1 -multiplier_petscspace_degree 2 257c4762a1bSJed Brown 258c4762a1bSJed Brown test: 259c4762a1bSJed Brown suffix: 1 260c4762a1bSJed Brown requires: triangle 261c4762a1bSJed Brown args: -potential_petscspace_degree 2 -conductivity_petscspace_degree 1 -multiplier_petscspace_degree 2 -snes_monitor -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 -fieldsplit_multiplier_ksp_rtol 1.0e-10 -fieldsplit_multiplier_pc_type lu -sol_vec_view 262c4762a1bSJed Brown 263c4762a1bSJed Brown TEST*/ 264