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
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 } AppCtx;
33c4762a1bSJed Brown
ProcessOptions(MPI_Comm comm,AppCtx * options)34d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
35d71ae5a4SJacob Faibussowitsch {
36c4762a1bSJed Brown const char *runTypes[2] = {"full", "test"};
37c4762a1bSJed Brown PetscInt run;
38c4762a1bSJed Brown
39c4762a1bSJed Brown PetscFunctionBeginUser;
40c4762a1bSJed Brown options->runType = RUN_FULL;
41d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");
42c4762a1bSJed Brown run = options->runType;
439566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-run_type", "The run type", "ex1.c", runTypes, 2, runTypes[options->runType], &run, NULL));
44c4762a1bSJed Brown options->runType = (RunType)run;
45d0609cedSBarry Smith PetscOptionsEnd();
463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
47c4762a1bSJed Brown }
48c4762a1bSJed Brown
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)49d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
50d71ae5a4SJacob Faibussowitsch {
51c4762a1bSJed Brown PetscFunctionBeginUser;
529566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
539566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
549566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
559566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
57c4762a1bSJed Brown }
58c4762a1bSJed Brown
59c4762a1bSJed Brown /* u - (x^2 + y^2) */
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[])60d71ae5a4SJacob 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[])
61d71ae5a4SJacob Faibussowitsch {
62c4762a1bSJed Brown f0[0] = u[0] - (x[0] * x[0] + x[1] * x[1]);
63c4762a1bSJed Brown }
64c4762a1bSJed Brown /* a \nabla\lambda */
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[])65d71ae5a4SJacob 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[])
66d71ae5a4SJacob Faibussowitsch {
67c4762a1bSJed Brown PetscInt d;
68c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u[1] * u_x[dim * 2 + d];
69c4762a1bSJed Brown }
70c4762a1bSJed Brown /* I */
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[])71d71ae5a4SJacob 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[])
72d71ae5a4SJacob Faibussowitsch {
73c4762a1bSJed Brown g0[0] = 1.0;
74c4762a1bSJed Brown }
75c4762a1bSJed Brown /* \nabla */
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[])76d71ae5a4SJacob Faibussowitsch 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[])
77d71ae5a4SJacob Faibussowitsch {
78c4762a1bSJed Brown PetscInt d;
79c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d] = u_x[dim * 2 + d];
80c4762a1bSJed Brown }
81c4762a1bSJed Brown /* a */
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[])82d71ae5a4SJacob 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[])
83d71ae5a4SJacob Faibussowitsch {
84c4762a1bSJed Brown PetscInt d;
85c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d * dim + d] = u[1];
86c4762a1bSJed Brown }
87c4762a1bSJed Brown /* a - (x + y) */
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[])88d71ae5a4SJacob 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[])
89d71ae5a4SJacob Faibussowitsch {
90c4762a1bSJed Brown f0[0] = u[1] - (x[0] + x[1]);
91c4762a1bSJed Brown }
92c4762a1bSJed Brown /* \lambda \nabla u */
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[])93d71ae5a4SJacob Faibussowitsch 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[])
94d71ae5a4SJacob Faibussowitsch {
95c4762a1bSJed Brown PetscInt d;
96c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u[2] * u_x[d];
97c4762a1bSJed Brown }
98c4762a1bSJed Brown /* I */
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[])99d71ae5a4SJacob 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[])
100d71ae5a4SJacob Faibussowitsch {
101c4762a1bSJed Brown g0[0] = 1.0;
102c4762a1bSJed Brown }
103c4762a1bSJed Brown /* 6 (x + y) */
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[])104d71ae5a4SJacob 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[])
105d71ae5a4SJacob Faibussowitsch {
106c4762a1bSJed Brown f0[0] = 6.0 * (x[0] + x[1]);
107c4762a1bSJed Brown }
108c4762a1bSJed Brown /* a \nabla u */
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[])109d71ae5a4SJacob 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[])
110d71ae5a4SJacob Faibussowitsch {
111c4762a1bSJed Brown PetscInt d;
112c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = u[1] * u_x[d];
113c4762a1bSJed Brown }
114c4762a1bSJed Brown /* \nabla u */
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[])115d71ae5a4SJacob Faibussowitsch 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[])
116d71ae5a4SJacob Faibussowitsch {
117c4762a1bSJed Brown PetscInt d;
118c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d] = u_x[d];
119c4762a1bSJed Brown }
120c4762a1bSJed Brown /* a */
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[])121d71ae5a4SJacob 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[])
122d71ae5a4SJacob Faibussowitsch {
123c4762a1bSJed Brown PetscInt d;
124c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d * dim + d] = u[1];
125c4762a1bSJed Brown }
126c4762a1bSJed Brown
127c4762a1bSJed Brown /*
128c4762a1bSJed Brown In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:
129c4762a1bSJed Brown
130c4762a1bSJed Brown u = x^2 + y^2
131c4762a1bSJed Brown f = 6 (x + y)
132c4762a1bSJed Brown kappa(a) = a = (x + y)
133c4762a1bSJed Brown
134c4762a1bSJed Brown so that
135c4762a1bSJed Brown
136c4762a1bSJed Brown -\div \kappa(a) \grad u + f = -6 (x + y) + 6 (x + y) = 0
137c4762a1bSJed Brown */
quadratic_u_2d(PetscInt dim,PetscReal t,const PetscReal x[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)138*2a8381b2SBarry Smith PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
139d71ae5a4SJacob Faibussowitsch {
140c4762a1bSJed Brown *u = x[0] * x[0] + x[1] * x[1];
1413ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
142c4762a1bSJed Brown }
linear_a_2d(PetscInt dim,PetscReal t,const PetscReal x[],PetscInt Nf,PetscScalar * a,PetscCtx ctx)143*2a8381b2SBarry Smith PetscErrorCode linear_a_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *a, PetscCtx ctx)
144d71ae5a4SJacob Faibussowitsch {
145c4762a1bSJed Brown *a = x[0] + x[1];
1463ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
147c4762a1bSJed Brown }
zero(PetscInt dim,PetscReal t,const PetscReal x[],PetscInt Nf,PetscScalar * l,PetscCtx ctx)148*2a8381b2SBarry Smith PetscErrorCode zero(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *l, PetscCtx ctx)
149d71ae5a4SJacob Faibussowitsch {
150c4762a1bSJed Brown *l = 0.0;
1513ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
152c4762a1bSJed Brown }
153c4762a1bSJed Brown
SetupProblem(DM dm,AppCtx * user)154d71ae5a4SJacob Faibussowitsch PetscErrorCode SetupProblem(DM dm, AppCtx *user)
155d71ae5a4SJacob Faibussowitsch {
15645480ffeSMatthew G. Knepley PetscDS ds;
15745480ffeSMatthew G. Knepley DMLabel label;
158c4762a1bSJed Brown const PetscInt id = 1;
159c4762a1bSJed Brown
160c4762a1bSJed Brown PetscFunctionBeginUser;
1619566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
1629566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_u, f1_u));
1639566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_a, f1_a));
1649566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 2, f0_l, f1_l));
1659566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, NULL));
1669566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_ua, NULL));
1679566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul));
1689566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL));
1699566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 2, 1, NULL, NULL, g2_la, NULL));
1709566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu));
171c4762a1bSJed Brown
1729566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL));
1739566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, linear_a_2d, NULL));
1749566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 2, zero, NULL));
1759566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label));
17657d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)quadratic_u_2d, NULL, user, NULL));
17757d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (PetscVoidFn *)linear_a_2d, NULL, user, NULL));
17857d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (PetscVoidFn *)zero, NULL, user, NULL));
1793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
180c4762a1bSJed Brown }
181c4762a1bSJed Brown
SetupDiscretization(DM dm,AppCtx * user)182d71ae5a4SJacob Faibussowitsch PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
183d71ae5a4SJacob Faibussowitsch {
184c4762a1bSJed Brown DM cdm = dm;
185c4762a1bSJed Brown const PetscInt dim = 2;
186c4762a1bSJed Brown PetscFE fe[3];
187c4762a1bSJed Brown PetscInt f;
188c4762a1bSJed Brown MPI_Comm comm;
189c4762a1bSJed Brown
190c4762a1bSJed Brown PetscFunctionBeginUser;
191c4762a1bSJed Brown /* Create finite element */
1929566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1939566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]));
1949566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "potential"));
1959566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "conductivity_", -1, &fe[1]));
1969566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "conductivity"));
1979566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
1989566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]));
1999566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[2], "multiplier"));
2009566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[2]));
201c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */
2029566063dSJacob Faibussowitsch for (f = 0; f < 3; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe[f]));
2039566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm));
2049566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user));
205c4762a1bSJed Brown while (cdm) {
2069566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm));
2079566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm));
208c4762a1bSJed Brown }
2099566063dSJacob Faibussowitsch for (f = 0; f < 3; ++f) PetscCall(PetscFEDestroy(&fe[f]));
2103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
211c4762a1bSJed Brown }
212c4762a1bSJed Brown
main(int argc,char ** argv)213d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
214d71ae5a4SJacob Faibussowitsch {
215c4762a1bSJed Brown DM dm;
216c4762a1bSJed Brown SNES snes;
217c4762a1bSJed Brown Vec u, r;
218c4762a1bSJed Brown AppCtx user;
219c4762a1bSJed Brown
220327415f7SBarry Smith PetscFunctionBeginUser;
2219566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2229566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
2239566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
2249566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
2259566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm));
2269566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &user));
227c4762a1bSJed Brown
2289566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u));
2299566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "solution"));
2309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r));
2316493148fSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
2329566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
233c4762a1bSJed Brown
2349566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u));
235c4762a1bSJed Brown if (user.runType == RUN_FULL) {
236348a1646SMatthew G. Knepley PetscDS ds;
237*2a8381b2SBarry Smith PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
238*2a8381b2SBarry Smith PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar u[], PetscCtx ctx);
239c4762a1bSJed Brown PetscReal error;
240c4762a1bSJed Brown
2419566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
2429566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL));
2439566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL));
2449566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL));
245c4762a1bSJed Brown initialGuess[0] = zero;
246c4762a1bSJed Brown initialGuess[1] = zero;
247c4762a1bSJed Brown initialGuess[2] = zero;
2489566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u));
2499566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-initial_vec_view"));
2509566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error));
2519566063dSJacob Faibussowitsch if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n"));
25263a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", (double)error));
2539566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u));
2549566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error));
2559566063dSJacob Faibussowitsch if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n"));
25663a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", (double)error));
257c4762a1bSJed Brown }
2589566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
259c4762a1bSJed Brown
2609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u));
2619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r));
2629566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes));
2639566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
2649566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
265b122ec5aSJacob Faibussowitsch return 0;
266c4762a1bSJed Brown }
267c4762a1bSJed Brown
268c4762a1bSJed Brown /*TEST
269c4762a1bSJed Brown
270c4762a1bSJed Brown build:
271c4762a1bSJed Brown requires: !complex
272c4762a1bSJed Brown
273c4762a1bSJed Brown test:
274c4762a1bSJed Brown suffix: 0
275c4762a1bSJed Brown requires: triangle
276c4762a1bSJed Brown args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -conductivity_petscspace_degree 1 -multiplier_petscspace_degree 2
277c4762a1bSJed Brown
278c4762a1bSJed Brown test:
279c4762a1bSJed Brown suffix: 1
280c4762a1bSJed Brown requires: triangle
281c4762a1bSJed 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
282c4762a1bSJed Brown
283c4762a1bSJed Brown TEST*/
284