Actual source code: ex1.c
1: static char help[] = "One-Shot Multigrid for Parameter Estimation Problem for the Poisson Equation.\n\
2: Using the Interior Point Method.\n\n\n";
4: /*F
5: We are solving the parameter estimation problem for the Laplacian. We will ask to minimize a Lagrangian
6: function over $a$ and $u$, given by
7: \begin{align}
8: L(u, a, \lambda) = \frac{1}{2} || Qu - d ||^2 + \frac{1}{2} || L (a - a_r) ||^2 + \lambda F(u; a)
9: \end{align}
10: where $Q$ is a sampling operator, $L$ is a regularization operator, $F$ defines the PDE.
12: Currently, we have perfect information, meaning $Q = I$, and then we need no regularization, $L = I$. We
13: also give the exact control for the reference $a_r$.
15: The PDE will be the Laplace equation with homogeneous boundary conditions
16: \begin{align}
17: -nabla \cdot a \nabla u = f
18: \end{align}
20: F*/
22: #include <petsc.h>
23: #include <petscfe.h>
25: typedef enum {
26: RUN_FULL,
27: RUN_TEST
28: } RunType;
30: typedef struct {
31: RunType runType; /* Whether to run tests, or solve the full problem */
32: } AppCtx;
34: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
35: {
36: const char *runTypes[2] = {"full", "test"};
37: PetscInt run;
39: PetscFunctionBeginUser;
40: options->runType = RUN_FULL;
41: PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");
42: run = options->runType;
43: PetscCall(PetscOptionsEList("-run_type", "The run type", "ex1.c", runTypes, 2, runTypes[options->runType], &run, NULL));
44: options->runType = (RunType)run;
45: PetscOptionsEnd();
46: PetscFunctionReturn(PETSC_SUCCESS);
47: }
49: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
50: {
51: PetscFunctionBeginUser;
52: PetscCall(DMCreate(comm, dm));
53: PetscCall(DMSetType(*dm, DMPLEX));
54: PetscCall(DMSetFromOptions(*dm));
55: PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
56: PetscFunctionReturn(PETSC_SUCCESS);
57: }
59: /* u - (x^2 + y^2) */
60: 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[])
61: {
62: f0[0] = u[0] - (x[0] * x[0] + x[1] * x[1]);
63: }
64: /* a \nabla\lambda */
65: 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[])
66: {
67: PetscInt d;
68: for (d = 0; d < dim; ++d) f1[d] = u[1] * u_x[dim * 2 + d];
69: }
70: /* I */
71: 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[])
72: {
73: g0[0] = 1.0;
74: }
75: /* \nabla */
76: 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[])
77: {
78: PetscInt d;
79: for (d = 0; d < dim; ++d) g2[d] = u_x[dim * 2 + d];
80: }
81: /* a */
82: 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[])
83: {
84: PetscInt d;
85: for (d = 0; d < dim; ++d) g3[d * dim + d] = u[1];
86: }
87: /* a - (x + y) */
88: 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[])
89: {
90: f0[0] = u[1] - (x[0] + x[1]);
91: }
92: /* \lambda \nabla u */
93: 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[])
94: {
95: PetscInt d;
96: for (d = 0; d < dim; ++d) f1[d] = u[2] * u_x[d];
97: }
98: /* I */
99: 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[])
100: {
101: g0[0] = 1.0;
102: }
103: /* 6 (x + y) */
104: 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[])
105: {
106: f0[0] = 6.0 * (x[0] + x[1]);
107: }
108: /* a \nabla u */
109: 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[])
110: {
111: PetscInt d;
112: for (d = 0; d < dim; ++d) f1[d] = u[1] * u_x[d];
113: }
114: /* \nabla u */
115: 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[])
116: {
117: PetscInt d;
118: for (d = 0; d < dim; ++d) g2[d] = u_x[d];
119: }
120: /* a */
121: 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[])
122: {
123: PetscInt d;
124: for (d = 0; d < dim; ++d) g3[d * dim + d] = u[1];
125: }
127: /*
128: In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:
130: u = x^2 + y^2
131: f = 6 (x + y)
132: kappa(a) = a = (x + y)
134: so that
136: -\div \kappa(a) \grad u + f = -6 (x + y) + 6 (x + y) = 0
137: */
138: PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
139: {
140: *u = x[0] * x[0] + x[1] * x[1];
141: return PETSC_SUCCESS;
142: }
143: PetscErrorCode linear_a_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx)
144: {
145: *a = x[0] + x[1];
146: return PETSC_SUCCESS;
147: }
148: PetscErrorCode zero(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx)
149: {
150: *l = 0.0;
151: return PETSC_SUCCESS;
152: }
154: PetscErrorCode SetupProblem(DM dm, AppCtx *user)
155: {
156: PetscDS ds;
157: DMLabel label;
158: const PetscInt id = 1;
160: PetscFunctionBeginUser;
161: PetscCall(DMGetDS(dm, &ds));
162: PetscCall(PetscDSSetResidual(ds, 0, f0_u, f1_u));
163: PetscCall(PetscDSSetResidual(ds, 1, f0_a, f1_a));
164: PetscCall(PetscDSSetResidual(ds, 2, f0_l, f1_l));
165: PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, NULL));
166: PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_ua, NULL));
167: PetscCall(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul));
168: PetscCall(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL));
169: PetscCall(PetscDSSetJacobian(ds, 2, 1, NULL, NULL, g2_la, NULL));
170: PetscCall(PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu));
172: PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL));
173: PetscCall(PetscDSSetExactSolution(ds, 1, linear_a_2d, NULL));
174: PetscCall(PetscDSSetExactSolution(ds, 2, zero, NULL));
175: PetscCall(DMGetLabel(dm, "marker", &label));
176: PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))quadratic_u_2d, NULL, user, NULL));
177: PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)(void))linear_a_2d, NULL, user, NULL));
178: PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)(void))zero, NULL, user, NULL));
179: PetscFunctionReturn(PETSC_SUCCESS);
180: }
182: PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
183: {
184: DM cdm = dm;
185: const PetscInt dim = 2;
186: PetscFE fe[3];
187: PetscInt f;
188: MPI_Comm comm;
190: PetscFunctionBeginUser;
191: /* Create finite element */
192: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
193: PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]));
194: PetscCall(PetscObjectSetName((PetscObject)fe[0], "potential"));
195: PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "conductivity_", -1, &fe[1]));
196: PetscCall(PetscObjectSetName((PetscObject)fe[1], "conductivity"));
197: PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
198: PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]));
199: PetscCall(PetscObjectSetName((PetscObject)fe[2], "multiplier"));
200: PetscCall(PetscFECopyQuadrature(fe[0], fe[2]));
201: /* Set discretization and boundary conditions for each mesh */
202: for (f = 0; f < 3; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe[f]));
203: PetscCall(DMCreateDS(dm));
204: PetscCall(SetupProblem(dm, user));
205: while (cdm) {
206: PetscCall(DMCopyDisc(dm, cdm));
207: PetscCall(DMGetCoarseDM(cdm, &cdm));
208: }
209: for (f = 0; f < 3; ++f) PetscCall(PetscFEDestroy(&fe[f]));
210: PetscFunctionReturn(PETSC_SUCCESS);
211: }
213: int main(int argc, char **argv)
214: {
215: DM dm;
216: SNES snes;
217: Vec u, r;
218: AppCtx user;
220: PetscFunctionBeginUser;
221: PetscCall(PetscInitialize(&argc, &argv, NULL, help));
222: PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
223: PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
224: PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
225: PetscCall(SNESSetDM(snes, dm));
226: PetscCall(SetupDiscretization(dm, &user));
228: PetscCall(DMCreateGlobalVector(dm, &u));
229: PetscCall(PetscObjectSetName((PetscObject)u, "solution"));
230: PetscCall(VecDuplicate(u, &r));
231: PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
232: PetscCall(SNESSetFromOptions(snes));
234: PetscCall(DMSNESCheckFromOptions(snes, u));
235: if (user.runType == RUN_FULL) {
236: PetscDS ds;
237: PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
238: PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
239: PetscReal error;
241: PetscCall(DMGetDS(dm, &ds));
242: PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL));
243: PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL));
244: PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL));
245: initialGuess[0] = zero;
246: initialGuess[1] = zero;
247: initialGuess[2] = zero;
248: PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u));
249: PetscCall(VecViewFromOptions(u, NULL, "-initial_vec_view"));
250: PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error));
251: if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n"));
252: else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", (double)error));
253: PetscCall(SNESSolve(snes, NULL, u));
254: PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error));
255: if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n"));
256: else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", (double)error));
257: }
258: PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
260: PetscCall(VecDestroy(&u));
261: PetscCall(VecDestroy(&r));
262: PetscCall(SNESDestroy(&snes));
263: PetscCall(DMDestroy(&dm));
264: PetscCall(PetscFinalize());
265: return 0;
266: }
268: /*TEST
270: build:
271: requires: !complex
273: test:
274: suffix: 0
275: requires: triangle
276: args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -conductivity_petscspace_degree 1 -multiplier_petscspace_degree 2
278: test:
279: suffix: 1
280: requires: triangle
281: 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
283: TEST*/