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";
We are solving the parameter estimation problem for the Laplacian. We will ask to minimize a Lagrangian
function over $y$ and $u$, given by
\begin{align}
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)
\end{align}
where $Q$ is a sampling operator, $L$ is a regularization operator, $F$ defines the PDE.
Currently, we have perfect information, meaning $Q = I$, and then we need no regularization, $L = I$. We
also give the null vector for the reference control $a_r$. Right now $\beta = 1$.
The PDE will be the Laplace equation with homogeneous boundary conditions
\begin{align}
-Delta u = a
\end{align}
22: #include <petsc.h>
23: #include <petscfe.h>
25: typedef enum {RUN_FULL, RUN_TEST} RunType;
27: typedef struct {
28: RunType runType; /* Whether to run tests, or solve the full problem */
29: PetscBool useDualPenalty; /* Penalize deviation from both goals */
30: } AppCtx;
32: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
33: {
34: const char *runTypes[2] = {"full", "test"};
35: PetscInt run;
39: options->runType = RUN_FULL;
40: options->useDualPenalty = PETSC_FALSE;
42: PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");
43: run = options->runType;
44: PetscOptionsEList("-run_type", "The run type", "ex2.c", runTypes, 2, runTypes[options->runType], &run, NULL);
45: options->runType = (RunType) run;
46: PetscOptionsBool("-use_dual_penalty", "Penalize deviation from both goals", "ex2.c", options->useDualPenalty, &options->useDualPenalty, NULL);
47: PetscOptionsEnd();
48: return(0);
49: }
51: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
52: {
53: DM distributedMesh = NULL;
57: DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);
58: PetscObjectSetName((PetscObject) *dm, "Mesh");
59: DMPlexDistribute(*dm, 0, NULL, &distributedMesh);
60: if (distributedMesh) {
61: DMDestroy(dm);
62: *dm = distributedMesh;
63: }
64: DMSetFromOptions(*dm);
65: DMViewFromOptions(*dm, NULL, "-dm_view");
66: return(0);
67: }
69: void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
70: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
71: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
72: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
73: {
74: f0[0] = (u[0] - (x[0]*x[0] + x[1]*x[1]));
75: }
76: void f0_u_full(PetscInt dim, PetscInt Nf, PetscInt NfAux,
77: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
78: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
79: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
80: {
81: 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]))) +
82: 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])));
83: }
84: void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
85: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
86: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
87: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
88: {
89: PetscInt d;
90: for (d = 0; d < dim; ++d) f1[d] = u_x[dim*2+d];
91: }
92: void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
93: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
94: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
95: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
96: {
97: g0[0] = 1.0;
98: }
99: void g0_uu_full(PetscInt dim, PetscInt Nf, PetscInt NfAux,
100: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
101: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
102: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
103: {
104: g0[0] = PetscSqr(u[0] - sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1]))
105: + PetscSqr(u[0] - (x[0]*x[0] + x[1]*x[1]))
106: - 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]
107: + 4.0*(x[0]*x[0] + x[1]*x[1])*(sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1]));
108: }
109: void g3_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux,
110: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
111: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
112: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
113: {
114: PetscInt d;
115: for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
116: }
118: void f0_a(PetscInt dim, PetscInt Nf, PetscInt NfAux,
119: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
120: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
121: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
122: {
123: f0[0] = u[1] - 4.0 /* 0.0 */ + u[2];
124: }
125: void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux,
126: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
127: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
128: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
129: {
130: g0[0] = 1.0;
131: }
132: void g0_al(PetscInt dim, PetscInt Nf, PetscInt NfAux,
133: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
134: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
135: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
136: {
137: g0[0] = 1.0;
138: }
140: void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
141: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
142: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
143: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
144: {
145: f0[0] = u[1];
146: }
147: void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
148: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
149: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
150: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
151: {
152: PetscInt d;
153: for (d = 0; d < dim; ++d) f1[d] = u_x[d];
154: }
155: void g0_la(PetscInt dim, PetscInt Nf, PetscInt NfAux,
156: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
157: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
158: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
159: {
160: g0[0] = 1.0;
161: }
162: void g3_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
163: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
164: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
165: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
166: {
167: PetscInt d;
168: for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
169: }
171: /*
172: In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:
174: u = x^2 + y^2
175: a = 4
176: d_A = 4
177: d_B = sin(2*pi*x[0]) * sin(2*pi*x[1])
179: so that
181: -\Delta u + a = -4 + 4 = 0
182: */
183: PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
184: {
185: *u = x[0]*x[0] + x[1]*x[1];
186: return 0;
187: }
188: PetscErrorCode constant_a_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx)
189: {
190: *a = 4;
191: return 0;
192: }
193: PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx)
194: {
195: *l = 0.0;
196: return 0;
197: }
199: PetscErrorCode SetupProblem(DM dm, AppCtx *user)
200: {
201: PetscDS prob;
202: const PetscInt id = 1;
206: DMGetDS(dm, &prob);
207: PetscDSSetResidual(prob, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u);
208: PetscDSSetResidual(prob, 1, f0_a, NULL);
209: PetscDSSetResidual(prob, 2, f0_l, f1_l);
210: PetscDSSetJacobian(prob, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL);
211: PetscDSSetJacobian(prob, 0, 2, NULL, NULL, NULL, g3_ul);
212: PetscDSSetJacobian(prob, 1, 1, g0_aa, NULL, NULL, NULL);
213: PetscDSSetJacobian(prob, 1, 2, g0_al, NULL, NULL, NULL);
214: PetscDSSetJacobian(prob, 2, 1, g0_la, NULL, NULL, NULL);
215: PetscDSSetJacobian(prob, 2, 0, NULL, NULL, NULL, g3_lu);
217: PetscDSSetExactSolution(prob, 0, quadratic_u_2d, NULL);
218: PetscDSSetExactSolution(prob, 1, constant_a_2d, NULL);
219: PetscDSSetExactSolution(prob, 2, zero, NULL);
220: DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)()) quadratic_u_2d, NULL, 1, &id, user);
221: DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 1, 0, NULL, (void (*)()) constant_a_2d, NULL, 1, &id, user);
222: DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 2, 0, NULL, (void (*)()) zero, NULL, 1, &id, user);
223: return(0);
224: }
226: PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
227: {
228: DM cdm = dm;
229: const PetscInt dim = 2;
230: PetscFE fe[3];
231: PetscInt f;
232: MPI_Comm comm;
233: PetscErrorCode ierr;
236: /* Create finite element */
237: PetscObjectGetComm((PetscObject) dm, &comm);
238: PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]);
239: PetscObjectSetName((PetscObject) fe[0], "potential");
240: PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "charge_", -1, &fe[1]);
241: PetscObjectSetName((PetscObject) fe[1], "charge");
242: PetscFECopyQuadrature(fe[0], fe[1]);
243: PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]);
244: PetscObjectSetName((PetscObject) fe[2], "multiplier");
245: PetscFECopyQuadrature(fe[0], fe[2]);
246: /* Set discretization and boundary conditions for each mesh */
247: for (f = 0; f < 3; ++f) {DMSetField(dm, f, NULL, (PetscObject) fe[f]);}
248: DMCreateDS(cdm);
249: SetupProblem(dm, user);
250: while (cdm) {
251: DMCopyDisc(dm, cdm);
252: DMGetCoarseDM(cdm, &cdm);
253: }
254: for (f = 0; f < 3; ++f) {PetscFEDestroy(&fe[f]);}
255: return(0);
256: }
258: int main(int argc, char **argv)
259: {
260: DM dm;
261: SNES snes;
262: Vec u, r;
263: AppCtx user;
266: PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
267: ProcessOptions(PETSC_COMM_WORLD, &user);
268: SNESCreate(PETSC_COMM_WORLD, &snes);
269: CreateMesh(PETSC_COMM_WORLD, &user, &dm);
270: SNESSetDM(snes, dm);
271: SetupDiscretization(dm, &user);
273: DMCreateGlobalVector(dm, &u);
274: PetscObjectSetName((PetscObject) u, "solution");
275: VecDuplicate(u, &r);
276: DMPlexSetSNESLocalFEM(dm,&user,&user,&user);
277: SNESSetFromOptions(snes);
279: DMSNESCheckFromOptions(snes, u);
280: if (user.runType == RUN_FULL) {
281: PetscDS ds;
282: PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
283: PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
284: PetscReal error;
286: DMGetDS(dm, &ds);
287: PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL);
288: PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL);
289: PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL);
290: initialGuess[0] = zero;
291: initialGuess[1] = zero;
292: initialGuess[2] = zero;
293: DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);
294: VecViewFromOptions(u, NULL, "-initial_vec_view");
295: DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);
296: if (error < 1.0e-11) {PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n");}
297: else {PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", error);}
298: SNESSolve(snes, NULL, u);
299: DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);
300: if (error < 1.0e-11) {PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n");}
301: else {PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", error);}
302: }
303: VecViewFromOptions(u, NULL, "-sol_vec_view");
305: VecDestroy(&u);
306: VecDestroy(&r);
307: SNESDestroy(&snes);
308: DMDestroy(&dm);
309: PetscFinalize();
310: return ierr;
311: }
313: /*TEST
315: build:
316: requires: !complex triangle
318: test:
319: suffix: 0
320: args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1
322: test:
323: suffix: 1
324: 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
326: test:
327: suffix: 2
328: 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
330: TEST*/