Actual source code: ex2.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 $y$ and $u$, given by
7: \begin{align}
8: 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)
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 null vector for the reference control $a_r$. Right now $\beta = 1$.
15: The PDE will be the Laplace equation with homogeneous boundary conditions
16: \begin{align}
17: -Delta u = a
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: PetscBool useDualPenalty; /* Penalize deviation from both goals */
33: } AppCtx;
35: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
36: {
37: const char *runTypes[2] = {"full", "test"};
38: PetscInt run;
40: PetscFunctionBeginUser;
41: options->runType = RUN_FULL;
42: options->useDualPenalty = PETSC_FALSE;
43: PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");
44: run = options->runType;
45: PetscCall(PetscOptionsEList("-run_type", "The run type", "ex2.c", runTypes, 2, runTypes[options->runType], &run, NULL));
46: options->runType = (RunType)run;
47: PetscCall(PetscOptionsBool("-use_dual_penalty", "Penalize deviation from both goals", "ex2.c", options->useDualPenalty, &options->useDualPenalty, NULL));
48: PetscOptionsEnd();
49: PetscFunctionReturn(PETSC_SUCCESS);
50: }
52: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
53: {
54: PetscFunctionBeginUser;
55: PetscCall(DMCreate(comm, dm));
56: PetscCall(DMSetType(*dm, DMPLEX));
57: PetscCall(DMSetFromOptions(*dm));
58: PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
59: PetscFunctionReturn(PETSC_SUCCESS);
60: }
62: 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[])
63: {
64: f0[0] = (u[0] - (x[0] * x[0] + x[1] * x[1]));
65: }
66: 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[])
67: {
68: 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])));
69: }
70: 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[])
71: {
72: PetscInt d;
73: for (d = 0; d < dim; ++d) f1[d] = u_x[dim * 2 + d];
74: }
75: 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[])
76: {
77: g0[0] = 1.0;
78: }
79: 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[])
80: {
81: 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]));
82: }
83: 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[])
84: {
85: PetscInt d;
86: for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
87: }
89: 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[])
90: {
91: f0[0] = u[1] - 4.0 /* 0.0 */ + u[2];
92: }
93: 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[])
94: {
95: g0[0] = 1.0;
96: }
97: 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[])
98: {
99: g0[0] = 1.0;
100: }
102: 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[])
103: {
104: f0[0] = u[1];
105: }
106: 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[])
107: {
108: PetscInt d;
109: for (d = 0; d < dim; ++d) f1[d] = u_x[d];
110: }
111: 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[])
112: {
113: g0[0] = 1.0;
114: }
115: 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[])
116: {
117: PetscInt d;
118: for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
119: }
121: /*
122: In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:
124: u = x^2 + y^2
125: a = 4
126: d_A = 4
127: d_B = sin(2*pi*x[0]) * sin(2*pi*x[1])
129: so that
131: -\Delta u + a = -4 + 4 = 0
132: */
133: PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
134: {
135: *u = x[0] * x[0] + x[1] * x[1];
136: return PETSC_SUCCESS;
137: }
138: PetscErrorCode constant_a_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx)
139: {
140: *a = 4;
141: return PETSC_SUCCESS;
142: }
143: PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx)
144: {
145: *l = 0.0;
146: return PETSC_SUCCESS;
147: }
149: PetscErrorCode SetupProblem(DM dm, AppCtx *user)
150: {
151: PetscDS ds;
152: DMLabel label;
153: const PetscInt id = 1;
155: PetscFunctionBeginUser;
156: PetscCall(DMGetDS(dm, &ds));
157: PetscCall(PetscDSSetResidual(ds, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u));
158: PetscCall(PetscDSSetResidual(ds, 1, f0_a, NULL));
159: PetscCall(PetscDSSetResidual(ds, 2, f0_l, f1_l));
160: PetscCall(PetscDSSetJacobian(ds, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL));
161: PetscCall(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul));
162: PetscCall(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL));
163: PetscCall(PetscDSSetJacobian(ds, 1, 2, g0_al, NULL, NULL, NULL));
164: PetscCall(PetscDSSetJacobian(ds, 2, 1, g0_la, NULL, NULL, NULL));
165: PetscCall(PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu));
167: PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL));
168: PetscCall(PetscDSSetExactSolution(ds, 1, constant_a_2d, NULL));
169: PetscCall(PetscDSSetExactSolution(ds, 2, zero, NULL));
170: PetscCall(DMGetLabel(dm, "marker", &label));
171: PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)())quadratic_u_2d, NULL, user, NULL));
172: PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)())constant_a_2d, NULL, user, NULL));
173: PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)())zero, NULL, user, NULL));
174: PetscFunctionReturn(PETSC_SUCCESS);
175: }
177: PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
178: {
179: DM cdm = dm;
180: const PetscInt dim = 2;
181: PetscFE fe[3];
182: PetscInt f;
183: MPI_Comm comm;
185: PetscFunctionBeginUser;
186: /* Create finite element */
187: PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
188: PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]));
189: PetscCall(PetscObjectSetName((PetscObject)fe[0], "potential"));
190: PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "charge_", -1, &fe[1]));
191: PetscCall(PetscObjectSetName((PetscObject)fe[1], "charge"));
192: PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
193: PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]));
194: PetscCall(PetscObjectSetName((PetscObject)fe[2], "multiplier"));
195: PetscCall(PetscFECopyQuadrature(fe[0], fe[2]));
196: /* Set discretization and boundary conditions for each mesh */
197: for (f = 0; f < 3; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe[f]));
198: PetscCall(DMCreateDS(cdm));
199: PetscCall(SetupProblem(dm, user));
200: while (cdm) {
201: PetscCall(DMCopyDisc(dm, cdm));
202: PetscCall(DMGetCoarseDM(cdm, &cdm));
203: }
204: for (f = 0; f < 3; ++f) PetscCall(PetscFEDestroy(&fe[f]));
205: PetscFunctionReturn(PETSC_SUCCESS);
206: }
208: int main(int argc, char **argv)
209: {
210: DM dm;
211: SNES snes;
212: Vec u, r;
213: AppCtx user;
215: PetscFunctionBeginUser;
216: PetscCall(PetscInitialize(&argc, &argv, NULL, help));
217: PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
218: PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
219: PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
220: PetscCall(SNESSetDM(snes, dm));
221: PetscCall(SetupDiscretization(dm, &user));
223: PetscCall(DMCreateGlobalVector(dm, &u));
224: PetscCall(PetscObjectSetName((PetscObject)u, "solution"));
225: PetscCall(VecDuplicate(u, &r));
226: PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
227: PetscCall(SNESSetFromOptions(snes));
229: PetscCall(DMSNESCheckFromOptions(snes, u));
230: if (user.runType == RUN_FULL) {
231: PetscDS ds;
232: PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
233: PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
234: PetscReal error;
236: PetscCall(DMGetDS(dm, &ds));
237: PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL));
238: PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL));
239: PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL));
240: initialGuess[0] = zero;
241: initialGuess[1] = zero;
242: initialGuess[2] = zero;
243: PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u));
244: PetscCall(VecViewFromOptions(u, NULL, "-initial_vec_view"));
245: PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error));
246: if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n"));
247: else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", (double)error));
248: PetscCall(SNESSolve(snes, NULL, u));
249: PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error));
250: if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n"));
251: else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", (double)error));
252: }
253: PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
255: PetscCall(VecDestroy(&u));
256: PetscCall(VecDestroy(&r));
257: PetscCall(SNESDestroy(&snes));
258: PetscCall(DMDestroy(&dm));
259: PetscCall(PetscFinalize());
260: return 0;
261: }
263: /*TEST
265: build:
266: requires: !complex triangle
268: test:
269: suffix: 0
270: args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1
272: test:
273: suffix: 1
274: 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
276: test:
277: suffix: 2
278: 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
280: TEST*/