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";
We are solving the parameter estimation problem for the Laplacian. We will ask to minimize a Lagrangian
function over $a$ and $u$, given by
\begin{align}
L(u, a, \lambda) = \frac{1}{2} || Qu - d ||^2 + \frac{1}{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 exact control for the reference $a_r$.
The PDE will be the Laplace equation with homogeneous boundary conditions
\begin{align}
-nabla \cdot a \nabla u = f
\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: } AppCtx;
31: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
32: {
33: const char *runTypes[2] = {"full", "test"};
34: PetscInt run;
38: options->runType = RUN_FULL;
40: PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");
41: run = options->runType;
42: PetscOptionsEList("-run_type", "The run type", "ex1.c", runTypes, 2, runTypes[options->runType], &run, NULL);
43: options->runType = (RunType) run;
44: PetscOptionsEnd();
45: return(0);
46: }
48: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
49: {
50: DM distributedMesh = NULL;
54: DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);
55: PetscObjectSetName((PetscObject) *dm, "Mesh");
56: DMPlexDistribute(*dm, 0, NULL, &distributedMesh);
57: if (distributedMesh) {
58: DMDestroy(dm);
59: *dm = distributedMesh;
60: }
61: DMSetFromOptions(*dm);
62: DMViewFromOptions(*dm, NULL, "-dm_view");
63: return(0);
64: }
66: /* u - (x^2 + y^2) */
67: void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
68: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
69: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
70: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
71: {
72: f0[0] = u[0] - (x[0]*x[0] + x[1]*x[1]);
73: }
74: /* a \nabla\lambda */
75: void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
76: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
77: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
78: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
79: {
80: PetscInt d;
81: for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[dim*2+d];
82: }
83: /* I */
84: void g0_uu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
88: {
89: g0[0] = 1.0;
90: }
91: /* \nabla */
92: void g2_ua(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 g2[])
96: {
97: PetscInt d;
98: for (d = 0; d < dim; ++d) g2[d] = u_x[dim*2+d];
99: }
100: /* a */
101: void g3_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux,
102: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
103: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
104: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
105: {
106: PetscInt d;
107: for (d = 0; d < dim; ++d) g3[d*dim+d] = u[1];
108: }
109: /* a - (x + y) */
110: void f0_a(PetscInt dim, PetscInt Nf, PetscInt NfAux,
111: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
112: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
113: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
114: {
115: f0[0] = u[1] - (x[0] + x[1]);
116: }
117: /* \lambda \nabla u */
118: void f1_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 f1[])
122: {
123: PetscInt d;
124: for (d = 0; d < dim; ++d) f1[d] = u[2]*u_x[d];
125: }
126: /* I */
127: void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux,
128: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
129: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
130: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
131: {
132: g0[0] = 1.0;
133: }
134: /* 6 (x + y) */
135: void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
136: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
137: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
138: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
139: {
140: f0[0] = 6.0*(x[0] + x[1]);
141: }
142: /* a \nabla u */
143: void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
144: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
145: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
146: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
147: {
148: PetscInt d;
149: for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[d];
150: }
151: /* \nabla u */
152: void g2_la(PetscInt dim, PetscInt Nf, PetscInt NfAux,
153: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
154: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
155: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
156: {
157: PetscInt d;
158: for (d = 0; d < dim; ++d) g2[d] = u_x[d];
159: }
160: /* a */
161: void g3_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
162: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
163: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
164: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
165: {
166: PetscInt d;
167: for (d = 0; d < dim; ++d) g3[d*dim+d] = u[1];
168: }
170: /*
171: In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:
173: u = x^2 + y^2
174: f = 6 (x + y)
175: kappa(a) = a = (x + y)
177: so that
179: -\div \kappa(a) \grad u + f = -6 (x + y) + 6 (x + y) = 0
180: */
181: PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
182: {
183: *u = x[0]*x[0] + x[1]*x[1];
184: return 0;
185: }
186: PetscErrorCode linear_a_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx)
187: {
188: *a = x[0] + x[1];
189: return 0;
190: }
191: PetscErrorCode zero(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx)
192: {
193: *l = 0.0;
194: return 0;
195: }
197: PetscErrorCode SetupProblem(DM dm, AppCtx *user)
198: {
199: PetscDS prob;
200: const PetscInt id = 1;
204: DMGetDS(dm, &prob);
205: PetscDSSetResidual(prob, 0, f0_u, f1_u);
206: PetscDSSetResidual(prob, 1, f0_a, f1_a);
207: PetscDSSetResidual(prob, 2, f0_l, f1_l);
208: PetscDSSetJacobian(prob, 0, 0, g0_uu, NULL, NULL, NULL);
209: PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_ua, NULL);
210: PetscDSSetJacobian(prob, 0, 2, NULL, NULL, NULL, g3_ul);
211: PetscDSSetJacobian(prob, 1, 1, g0_aa, NULL, NULL, NULL);
212: PetscDSSetJacobian(prob, 2, 1, NULL, NULL, g2_la, NULL);
213: PetscDSSetJacobian(prob, 2, 0, NULL, NULL, NULL, g3_lu);
215: PetscDSSetExactSolution(prob, 0, quadratic_u_2d, NULL);
216: PetscDSSetExactSolution(prob, 1, linear_a_2d, NULL);
217: PetscDSSetExactSolution(prob, 2, zero, NULL);
218: DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) quadratic_u_2d, NULL, 1, &id, user);
219: DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 1, 0, NULL, (void (*)(void)) linear_a_2d, NULL, 1, &id, user);
220: DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 2, 0, NULL, (void (*)(void)) zero, NULL, 1, &id, user);
221: return(0);
222: }
224: PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
225: {
226: DM cdm = dm;
227: const PetscInt dim = 2;
228: PetscFE fe[3];
229: PetscInt f;
230: MPI_Comm comm;
231: PetscErrorCode ierr;
234: /* Create finite element */
235: PetscObjectGetComm((PetscObject) dm, &comm);
236: PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]);
237: PetscObjectSetName((PetscObject) fe[0], "potential");
238: PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "conductivity_", -1, &fe[1]);
239: PetscObjectSetName((PetscObject) fe[1], "conductivity");
240: PetscFECopyQuadrature(fe[0], fe[1]);
241: PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]);
242: PetscObjectSetName((PetscObject) fe[2], "multiplier");
243: PetscFECopyQuadrature(fe[0], fe[2]);
244: /* Set discretization and boundary conditions for each mesh */
245: for (f = 0; f < 3; ++f) {DMSetField(dm, f, NULL, (PetscObject) fe[f]);}
246: DMCreateDS(dm);
247: SetupProblem(dm, user);
248: while (cdm) {
249: DMCopyDisc(dm, cdm);
250: DMGetCoarseDM(cdm, &cdm);
251: }
252: for (f = 0; f < 3; ++f) {PetscFEDestroy(&fe[f]);}
253: return(0);
254: }
256: int main(int argc, char **argv)
257: {
258: DM dm;
259: SNES snes;
260: Vec u, r;
261: AppCtx user;
264: PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
265: ProcessOptions(PETSC_COMM_WORLD, &user);
266: SNESCreate(PETSC_COMM_WORLD, &snes);
267: CreateMesh(PETSC_COMM_WORLD, &user, &dm);
268: SNESSetDM(snes, dm);
269: SetupDiscretization(dm, &user);
271: DMCreateGlobalVector(dm, &u);
272: PetscObjectSetName((PetscObject) u, "solution");
273: VecDuplicate(u, &r);
274: DMPlexSetSNESLocalFEM(dm,&user,&user,&user);
275: SNESSetFromOptions(snes);
277: DMSNESCheckFromOptions(snes, u);
278: if (user.runType == RUN_FULL) {
279: PetscDS ds;
280: PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
281: PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
282: PetscReal error;
284: DMGetDS(dm, &ds);
285: PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL);
286: PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL);
287: PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL);
288: initialGuess[0] = zero;
289: initialGuess[1] = zero;
290: initialGuess[2] = zero;
291: DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);
292: VecViewFromOptions(u, NULL, "-initial_vec_view");
293: DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);
294: if (error < 1.0e-11) {PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n");}
295: else {PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", error);}
296: SNESSolve(snes, NULL, u);
297: DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);
298: if (error < 1.0e-11) {PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n");}
299: else {PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", error);}
300: }
301: VecViewFromOptions(u, NULL, "-sol_vec_view");
303: VecDestroy(&u);
304: VecDestroy(&r);
305: SNESDestroy(&snes);
306: DMDestroy(&dm);
307: PetscFinalize();
308: return ierr;
309: }
311: /*TEST
313: build:
314: requires: !complex
316: test:
317: suffix: 0
318: requires: triangle
319: args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -conductivity_petscspace_degree 1 -multiplier_petscspace_degree 2
321: test:
322: suffix: 1
323: requires: triangle
324: 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
326: TEST*/