xref: /petsc/src/tao/tutorials/ex1.c (revision 348a1646b1a7771dc468a628d08efb39f735eb8e)
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 
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 } AppCtx;
30c4762a1bSJed Brown 
31c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
32c4762a1bSJed Brown {
33c4762a1bSJed Brown   const char    *runTypes[2] = {"full", "test"};
34c4762a1bSJed Brown   PetscInt       run;
35c4762a1bSJed Brown   PetscErrorCode ierr;
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   PetscFunctionBeginUser;
38c4762a1bSJed Brown   options->runType = RUN_FULL;
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");CHKERRQ(ierr);
41c4762a1bSJed Brown   run  = options->runType;
42c4762a1bSJed Brown   ierr = PetscOptionsEList("-run_type", "The run type", "ex1.c", runTypes, 2, runTypes[options->runType], &run, NULL);CHKERRQ(ierr);
43c4762a1bSJed Brown   options->runType = (RunType) run;
44c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
45c4762a1bSJed Brown   PetscFunctionReturn(0);
46c4762a1bSJed Brown }
47c4762a1bSJed Brown 
48c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
49c4762a1bSJed Brown {
50c4762a1bSJed Brown   DM             distributedMesh = NULL;
51c4762a1bSJed Brown   PetscErrorCode ierr;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   PetscFunctionBeginUser;
54c4762a1bSJed Brown   ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
55c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
56c4762a1bSJed Brown   ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr);
57c4762a1bSJed Brown   if (distributedMesh) {
58c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
59c4762a1bSJed Brown     *dm  = distributedMesh;
60c4762a1bSJed Brown   }
61c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
62c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
63c4762a1bSJed Brown   PetscFunctionReturn(0);
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown /* u - (x^2 + y^2) */
67c4762a1bSJed Brown void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
68c4762a1bSJed Brown           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
69c4762a1bSJed Brown           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
70c4762a1bSJed Brown           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
71c4762a1bSJed Brown {
72c4762a1bSJed Brown   f0[0] = u[0] - (x[0]*x[0] + x[1]*x[1]);
73c4762a1bSJed Brown }
74c4762a1bSJed Brown /* a \nabla\lambda */
75c4762a1bSJed Brown void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
76c4762a1bSJed Brown           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
77c4762a1bSJed Brown           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
78c4762a1bSJed Brown           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
79c4762a1bSJed Brown {
80c4762a1bSJed Brown   PetscInt d;
81c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[dim*2+d];
82c4762a1bSJed Brown }
83c4762a1bSJed Brown /* I */
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 /* \nabla */
92c4762a1bSJed Brown void g2_ua(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 g2[])
96c4762a1bSJed Brown {
97c4762a1bSJed Brown   PetscInt d;
98c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g2[d] = u_x[dim*2+d];
99c4762a1bSJed Brown }
100c4762a1bSJed Brown /* a */
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] = u[1];
108c4762a1bSJed Brown }
109c4762a1bSJed Brown /* a - (x + y) */
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] - (x[0] + x[1]);
116c4762a1bSJed Brown }
117c4762a1bSJed Brown /* \lambda \nabla u */
118c4762a1bSJed Brown void f1_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 f1[])
122c4762a1bSJed Brown {
123c4762a1bSJed Brown   PetscInt d;
124c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = u[2]*u_x[d];
125c4762a1bSJed Brown }
126c4762a1bSJed Brown /* I */
127c4762a1bSJed Brown void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux,
128c4762a1bSJed Brown            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
129c4762a1bSJed Brown            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
130c4762a1bSJed Brown            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
131c4762a1bSJed Brown {
132c4762a1bSJed Brown   g0[0] = 1.0;
133c4762a1bSJed Brown }
134c4762a1bSJed Brown /* 6 (x + y) */
135c4762a1bSJed Brown void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
136c4762a1bSJed Brown           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
137c4762a1bSJed Brown           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
138c4762a1bSJed Brown           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
139c4762a1bSJed Brown {
140c4762a1bSJed Brown   f0[0] = 6.0*(x[0] + x[1]);
141c4762a1bSJed Brown }
142c4762a1bSJed Brown /* a \nabla u */
143c4762a1bSJed Brown void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
144c4762a1bSJed Brown           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
145c4762a1bSJed Brown           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
146c4762a1bSJed Brown           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
147c4762a1bSJed Brown {
148c4762a1bSJed Brown   PetscInt d;
149c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[d];
150c4762a1bSJed Brown }
151c4762a1bSJed Brown /* \nabla u */
152c4762a1bSJed Brown void g2_la(PetscInt dim, PetscInt Nf, PetscInt NfAux,
153c4762a1bSJed Brown            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
154c4762a1bSJed Brown            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
155c4762a1bSJed Brown            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
156c4762a1bSJed Brown {
157c4762a1bSJed Brown   PetscInt d;
158c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g2[d] = u_x[d];
159c4762a1bSJed Brown }
160c4762a1bSJed Brown /* a */
161c4762a1bSJed Brown void g3_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
162c4762a1bSJed Brown            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
163c4762a1bSJed Brown            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
164c4762a1bSJed Brown            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
165c4762a1bSJed Brown {
166c4762a1bSJed Brown   PetscInt d;
167c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g3[d*dim+d] = u[1];
168c4762a1bSJed Brown }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown /*
171c4762a1bSJed Brown   In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:
172c4762a1bSJed Brown 
173c4762a1bSJed Brown     u  = x^2 + y^2
174c4762a1bSJed Brown     f  = 6 (x + y)
175c4762a1bSJed Brown     kappa(a) = a = (x + y)
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   so that
178c4762a1bSJed Brown 
179c4762a1bSJed Brown     -\div \kappa(a) \grad u + f = -6 (x + y) + 6 (x + y) = 0
180c4762a1bSJed Brown */
181c4762a1bSJed Brown PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
182c4762a1bSJed Brown {
183c4762a1bSJed Brown   *u = x[0]*x[0] + x[1]*x[1];
184c4762a1bSJed Brown   return 0;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown PetscErrorCode linear_a_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx)
187c4762a1bSJed Brown {
188c4762a1bSJed Brown   *a = x[0] + x[1];
189c4762a1bSJed Brown   return 0;
190c4762a1bSJed Brown }
191c4762a1bSJed Brown PetscErrorCode zero(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx)
192c4762a1bSJed Brown {
193c4762a1bSJed Brown   *l = 0.0;
194c4762a1bSJed Brown   return 0;
195c4762a1bSJed Brown }
196c4762a1bSJed Brown 
197c4762a1bSJed Brown PetscErrorCode SetupProblem(DM dm, AppCtx *user)
198c4762a1bSJed Brown {
199c4762a1bSJed Brown   PetscDS        prob;
200c4762a1bSJed Brown   const PetscInt id = 1;
201c4762a1bSJed Brown   PetscErrorCode ierr;
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   PetscFunctionBeginUser;
204c4762a1bSJed Brown   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
205c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 0, f0_u, f1_u);CHKERRQ(ierr);
206c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 1, f0_a, f1_a);CHKERRQ(ierr);
207c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 2, f0_l, f1_l);CHKERRQ(ierr);
208c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 0, g0_uu, NULL, NULL, NULL);CHKERRQ(ierr);
209c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_ua, NULL);CHKERRQ(ierr);
210c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 2, NULL, NULL, NULL, g3_ul);CHKERRQ(ierr);
211c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 1, 1, g0_aa, NULL, NULL, NULL);CHKERRQ(ierr);
212c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 2, 1, NULL, NULL, g2_la, NULL);CHKERRQ(ierr);
213c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 2, 0, NULL, NULL, NULL, g3_lu);CHKERRQ(ierr);
214c4762a1bSJed Brown 
215*348a1646SMatthew G. Knepley   ierr = PetscDSSetExactSolution(prob, 0, quadratic_u_2d, NULL);CHKERRQ(ierr);
216*348a1646SMatthew G. Knepley   ierr = PetscDSSetExactSolution(prob, 1, linear_a_2d, NULL);CHKERRQ(ierr);
217*348a1646SMatthew G. Knepley   ierr = PetscDSSetExactSolution(prob, 2, zero, NULL);CHKERRQ(ierr);
218*348a1646SMatthew G. Knepley   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) quadratic_u_2d, 1, &id, user);CHKERRQ(ierr);
219*348a1646SMatthew G. Knepley   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 1, 0, NULL, (void (*)(void)) linear_a_2d, 1, &id, user);CHKERRQ(ierr);
220*348a1646SMatthew G. Knepley   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 2, 0, NULL, (void (*)(void)) zero, 1, &id, user);CHKERRQ(ierr);
221c4762a1bSJed Brown   PetscFunctionReturn(0);
222c4762a1bSJed Brown }
223c4762a1bSJed Brown 
224c4762a1bSJed Brown PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
225c4762a1bSJed Brown {
226c4762a1bSJed Brown   DM              cdm = dm;
227c4762a1bSJed Brown   const PetscInt  dim = 2;
228c4762a1bSJed Brown   PetscFE         fe[3];
229c4762a1bSJed Brown   PetscInt        f;
230c4762a1bSJed Brown   MPI_Comm        comm;
231c4762a1bSJed Brown   PetscErrorCode  ierr;
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   PetscFunctionBeginUser;
234c4762a1bSJed Brown   /* Create finite element */
235c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
236c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]);CHKERRQ(ierr);
237c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[0], "potential");CHKERRQ(ierr);
238c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "conductivity_", -1, &fe[1]);CHKERRQ(ierr);
239c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[1], "conductivity");CHKERRQ(ierr);
240c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
241c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]);CHKERRQ(ierr);
242c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[2], "multiplier");CHKERRQ(ierr);
243c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr);
244c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
245c4762a1bSJed Brown   for (f = 0; f < 3; ++f) {ierr = DMSetField(dm, f, NULL, (PetscObject) fe[f]);CHKERRQ(ierr);}
246c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
247c4762a1bSJed Brown   ierr = SetupProblem(dm, user);CHKERRQ(ierr);
248c4762a1bSJed Brown   while (cdm) {
249c4762a1bSJed Brown     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
250c4762a1bSJed Brown     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
251c4762a1bSJed Brown   }
252c4762a1bSJed Brown   for (f = 0; f < 3; ++f) {ierr = PetscFEDestroy(&fe[f]);CHKERRQ(ierr);}
253c4762a1bSJed Brown   PetscFunctionReturn(0);
254c4762a1bSJed Brown }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown int main(int argc, char **argv)
257c4762a1bSJed Brown {
258c4762a1bSJed Brown   DM             dm;
259c4762a1bSJed Brown   SNES           snes;
260c4762a1bSJed Brown   Vec            u, r;
261c4762a1bSJed Brown   AppCtx         user;
262c4762a1bSJed Brown   PetscErrorCode ierr;
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
265c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
266c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
267c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
268c4762a1bSJed Brown   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
269c4762a1bSJed Brown   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
272c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) u, "solution");CHKERRQ(ierr);
273c4762a1bSJed Brown   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
274c4762a1bSJed Brown   ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr);
275c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
276c4762a1bSJed Brown 
277*348a1646SMatthew G. Knepley   ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr);
278c4762a1bSJed Brown   if (user.runType == RUN_FULL) {
279*348a1646SMatthew G. Knepley     PetscDS          ds;
280*348a1646SMatthew G. Knepley     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
281c4762a1bSJed Brown     PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
282c4762a1bSJed Brown     PetscReal        error;
283c4762a1bSJed Brown 
284*348a1646SMatthew G. Knepley     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
285*348a1646SMatthew G. Knepley     ierr = PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL);CHKERRQ(ierr);
286*348a1646SMatthew G. Knepley     ierr = PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL);CHKERRQ(ierr);
287*348a1646SMatthew G. Knepley     ierr = PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL);CHKERRQ(ierr);
288c4762a1bSJed Brown     initialGuess[0] = zero;
289c4762a1bSJed Brown     initialGuess[1] = zero;
290c4762a1bSJed Brown     initialGuess[2] = zero;
291c4762a1bSJed Brown     ierr = DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr);
292c4762a1bSJed Brown     ierr = VecViewFromOptions(u, NULL, "-initial_vec_view");CHKERRQ(ierr);
293*348a1646SMatthew G. Knepley     ierr = DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);CHKERRQ(ierr);
294c4762a1bSJed Brown     if (error < 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
295c4762a1bSJed Brown     else                 {ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", error);CHKERRQ(ierr);}
296c4762a1bSJed Brown     ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
297*348a1646SMatthew 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, "Final L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
299c4762a1bSJed Brown     else                 {ierr = PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", error);CHKERRQ(ierr);}
300c4762a1bSJed Brown   }
301c4762a1bSJed Brown   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
304c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
305c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
306c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
307c4762a1bSJed Brown   ierr = PetscFinalize();
308c4762a1bSJed Brown   return ierr;
309c4762a1bSJed Brown }
310c4762a1bSJed Brown 
311c4762a1bSJed Brown /*TEST
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   build:
314c4762a1bSJed Brown     requires: !complex
315c4762a1bSJed Brown 
316c4762a1bSJed Brown   test:
317c4762a1bSJed Brown     suffix: 0
318c4762a1bSJed Brown     requires: triangle
319c4762a1bSJed Brown     args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -conductivity_petscspace_degree 1 -multiplier_petscspace_degree 2
320c4762a1bSJed Brown 
321c4762a1bSJed Brown   test:
322c4762a1bSJed Brown     suffix: 1
323c4762a1bSJed Brown     requires: triangle
324c4762a1bSJed 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
325c4762a1bSJed Brown 
326c4762a1bSJed Brown TEST*/
327