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