xref: /petsc/src/tao/tutorials/ex2.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
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   PetscBool useDualPenalty; /* Penalize deviation from both goals */
33c4762a1bSJed Brown } AppCtx;
34c4762a1bSJed Brown 
35*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
36c4762a1bSJed Brown   const char *runTypes[2] = {"full", "test"};
37c4762a1bSJed Brown   PetscInt    run;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   PetscFunctionBeginUser;
40c4762a1bSJed Brown   options->runType        = RUN_FULL;
41c4762a1bSJed Brown   options->useDualPenalty = PETSC_FALSE;
42d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");
43c4762a1bSJed Brown   run = options->runType;
449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-run_type", "The run type", "ex2.c", runTypes, 2, runTypes[options->runType], &run, NULL));
45c4762a1bSJed Brown   options->runType = (RunType)run;
469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-use_dual_penalty", "Penalize deviation from both goals", "ex2.c", options->useDualPenalty, &options->useDualPenalty, NULL));
47d0609cedSBarry Smith   PetscOptionsEnd();
48c4762a1bSJed Brown   PetscFunctionReturn(0);
49c4762a1bSJed Brown }
50c4762a1bSJed Brown 
51*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
52c4762a1bSJed Brown   PetscFunctionBeginUser;
539566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
549566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
559566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
569566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
57c4762a1bSJed Brown   PetscFunctionReturn(0);
58c4762a1bSJed Brown }
59c4762a1bSJed Brown 
60*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[]) {
61c4762a1bSJed Brown   f0[0] = (u[0] - (x[0] * x[0] + x[1] * x[1]));
62c4762a1bSJed Brown }
63*9371c9d4SSatish Balay 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[]) {
64*9371c9d4SSatish 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])));
65c4762a1bSJed Brown }
66*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[]) {
67c4762a1bSJed Brown   PetscInt d;
68c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = u_x[dim * 2 + d];
69c4762a1bSJed Brown }
70*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[]) {
71c4762a1bSJed Brown   g0[0] = 1.0;
72c4762a1bSJed Brown }
73*9371c9d4SSatish Balay 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[]) {
74*9371c9d4SSatish 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]));
75c4762a1bSJed Brown }
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] = 1.0;
79c4762a1bSJed Brown }
80c4762a1bSJed Brown 
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] - 4.0 /* 0.0 */ + u[2];
83c4762a1bSJed Brown }
84*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[]) {
85c4762a1bSJed Brown   g0[0] = 1.0;
86c4762a1bSJed Brown }
87*9371c9d4SSatish Balay 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[]) {
88c4762a1bSJed Brown   g0[0] = 1.0;
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
91*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[]) {
92c4762a1bSJed Brown   f0[0] = u[1];
93c4762a1bSJed Brown }
94*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[]) {
95c4762a1bSJed Brown   PetscInt d;
96c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
97c4762a1bSJed Brown }
98*9371c9d4SSatish Balay 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[]) {
99c4762a1bSJed Brown   g0[0] = 1.0;
100c4762a1bSJed Brown }
101*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[]) {
102c4762a1bSJed Brown   PetscInt d;
103c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
104c4762a1bSJed Brown }
105c4762a1bSJed Brown 
106c4762a1bSJed Brown /*
107c4762a1bSJed Brown   In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:
108c4762a1bSJed Brown 
109c4762a1bSJed Brown     u   = x^2 + y^2
110c4762a1bSJed Brown     a   = 4
111c4762a1bSJed Brown     d_A = 4
112c4762a1bSJed Brown     d_B = sin(2*pi*x[0]) * sin(2*pi*x[1])
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   so that
115c4762a1bSJed Brown 
116c4762a1bSJed Brown     -\Delta u + a = -4 + 4 = 0
117c4762a1bSJed Brown */
118*9371c9d4SSatish Balay PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) {
119c4762a1bSJed Brown   *u = x[0] * x[0] + x[1] * x[1];
120c4762a1bSJed Brown   return 0;
121c4762a1bSJed Brown }
122*9371c9d4SSatish Balay PetscErrorCode constant_a_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx) {
123c4762a1bSJed Brown   *a = 4;
124c4762a1bSJed Brown   return 0;
125c4762a1bSJed Brown }
126*9371c9d4SSatish Balay PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx) {
127c4762a1bSJed Brown   *l = 0.0;
128c4762a1bSJed Brown   return 0;
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131*9371c9d4SSatish Balay PetscErrorCode SetupProblem(DM dm, AppCtx *user) {
13245480ffeSMatthew G. Knepley   PetscDS        ds;
13345480ffeSMatthew G. Knepley   DMLabel        label;
134c4762a1bSJed Brown   const PetscInt id = 1;
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   PetscFunctionBeginUser;
1379566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1389566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u));
1399566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds, 1, f0_a, NULL));
1409566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds, 2, f0_l, f1_l));
1419566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL));
1429566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul));
1439566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL));
1449566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 1, 2, g0_al, NULL, NULL, NULL));
1459566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 2, 1, g0_la, NULL, NULL, NULL));
1469566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu));
147c4762a1bSJed Brown 
1489566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL));
1499566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 1, constant_a_2d, NULL));
1509566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 2, zero, NULL));
1519566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
1529566063dSJacob Faibussowitsch   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)())quadratic_u_2d, NULL, user, NULL));
1539566063dSJacob Faibussowitsch   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)())constant_a_2d, NULL, user, NULL));
1549566063dSJacob Faibussowitsch   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)())zero, NULL, user, NULL));
155c4762a1bSJed Brown   PetscFunctionReturn(0);
156c4762a1bSJed Brown }
157c4762a1bSJed Brown 
158*9371c9d4SSatish Balay PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) {
159c4762a1bSJed Brown   DM             cdm = dm;
160c4762a1bSJed Brown   const PetscInt dim = 2;
161c4762a1bSJed Brown   PetscFE        fe[3];
162c4762a1bSJed Brown   PetscInt       f;
163c4762a1bSJed Brown   MPI_Comm       comm;
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   PetscFunctionBeginUser;
166c4762a1bSJed Brown   /* Create finite element */
1679566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1689566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]));
1699566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe[0], "potential"));
1709566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "charge_", -1, &fe[1]));
1719566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe[1], "charge"));
1729566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
1739566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]));
1749566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe[2], "multiplier"));
1759566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fe[0], fe[2]));
176c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
1779566063dSJacob Faibussowitsch   for (f = 0; f < 3; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe[f]));
1789566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(cdm));
1799566063dSJacob Faibussowitsch   PetscCall(SetupProblem(dm, user));
180c4762a1bSJed Brown   while (cdm) {
1819566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
1829566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
183c4762a1bSJed Brown   }
1849566063dSJacob Faibussowitsch   for (f = 0; f < 3; ++f) PetscCall(PetscFEDestroy(&fe[f]));
185c4762a1bSJed Brown   PetscFunctionReturn(0);
186c4762a1bSJed Brown }
187c4762a1bSJed Brown 
188*9371c9d4SSatish Balay int main(int argc, char **argv) {
189c4762a1bSJed Brown   DM     dm;
190c4762a1bSJed Brown   SNES   snes;
191c4762a1bSJed Brown   Vec    u, r;
192c4762a1bSJed Brown   AppCtx user;
193c4762a1bSJed Brown 
194327415f7SBarry Smith   PetscFunctionBeginUser;
1959566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1969566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1979566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
1989566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1999566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, dm));
2009566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, &user));
201c4762a1bSJed Brown 
2029566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
2039566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "solution"));
2049566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
2059566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
2069566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
207c4762a1bSJed Brown 
2089566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckFromOptions(snes, u));
209c4762a1bSJed Brown   if (user.runType == RUN_FULL) {
210348a1646SMatthew G. Knepley     PetscDS ds;
211348a1646SMatthew G. Knepley     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
212c4762a1bSJed Brown     PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
213c4762a1bSJed Brown     PetscReal error;
214c4762a1bSJed Brown 
2159566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &ds));
2169566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL));
2179566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL));
2189566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL));
219c4762a1bSJed Brown     initialGuess[0] = zero;
220c4762a1bSJed Brown     initialGuess[1] = zero;
221c4762a1bSJed Brown     initialGuess[2] = zero;
2229566063dSJacob Faibussowitsch     PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u));
2239566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(u, NULL, "-initial_vec_view"));
2249566063dSJacob Faibussowitsch     PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error));
2259566063dSJacob Faibussowitsch     if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n"));
22663a3b9bcSJacob Faibussowitsch     else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", (double)error));
2279566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, NULL, u));
2289566063dSJacob Faibussowitsch     PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error));
2299566063dSJacob Faibussowitsch     if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n"));
23063a3b9bcSJacob Faibussowitsch     else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", (double)error));
231c4762a1bSJed Brown   }
2329566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
233c4762a1bSJed Brown 
2349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
2359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
2369566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
2379566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
2389566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
239b122ec5aSJacob Faibussowitsch   return 0;
240c4762a1bSJed Brown }
241c4762a1bSJed Brown 
242c4762a1bSJed Brown /*TEST
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   build:
245c4762a1bSJed Brown     requires: !complex triangle
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   test:
248c4762a1bSJed Brown     suffix: 0
249c4762a1bSJed Brown     args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   test:
252c4762a1bSJed Brown     suffix: 1
253c4762a1bSJed 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
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   test:
256c4762a1bSJed Brown     suffix: 2
257c4762a1bSJed 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
258c4762a1bSJed Brown 
259c4762a1bSJed Brown TEST*/
260