Actual source code: ex1.c
petsc-3.13.6 2020-09-29
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: PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
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;
41: PetscOptionsBegin (comm, "" , "Inverse Problem Options" , "DMPLEX " );
42: run = options->runType;
43: PetscOptionsEList ("-run_type" , "The run type" , "ex1.c" , runTypes, 2, runTypes[options->runType], &run, NULL);
44: options->runType = (RunType) run;
45: PetscOptionsEnd ();
46: return (0);
47: }
49: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
50: {
51: DM distributedMesh = NULL;
55: DMPlexCreateBoxMesh (comm, 2, PETSC_TRUE , NULL, NULL, NULL, NULL, PETSC_TRUE , dm);
56: PetscObjectSetName ((PetscObject ) *dm, "Mesh" );
57: DMPlexDistribute (*dm, 0, NULL, &distributedMesh);
58: if (distributedMesh) {
59: DMDestroy (dm);
60: *dm = distributedMesh;
61: }
62: DMSetFromOptions (*dm);
63: DMViewFromOptions (*dm, NULL, "-dm_view" );
64: return (0);
65: }
67: /* u - (x^2 + y^2) */
68: void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
69: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
70: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
71: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
72: {
73: f0[0] = u[0] - (x[0]*x[0] + x[1]*x[1]);
74: }
75: /* a \nabla\lambda */
76: void f1_u(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 f1[])
80: {
81: PetscInt d;
82: for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[dim*2+d];
83: }
84: /* I */
85: void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
86: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
87: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
88: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
89: {
90: g0[0] = 1.0;
91: }
92: /* \nabla */
93: void g2_ua(PetscInt dim, PetscInt Nf, PetscInt NfAux,
94: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
95: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
96: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
97: {
98: PetscInt d;
99: for (d = 0; d < dim; ++d) g2[d] = u_x[dim*2+d];
100: }
101: /* a */
102: void g3_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux,
103: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
104: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
105: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
106: {
107: PetscInt d;
108: for (d = 0; d < dim; ++d) g3[d*dim+d] = u[1];
109: }
110: /* a - (x + y) */
111: void f0_a(PetscInt dim, PetscInt Nf, PetscInt NfAux,
112: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
113: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
114: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
115: {
116: f0[0] = u[1] - (x[0] + x[1]);
117: }
118: /* \lambda \nabla u */
119: void f1_a(PetscInt dim, PetscInt Nf, PetscInt NfAux,
120: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
121: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
122: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
123: {
124: PetscInt d;
125: for (d = 0; d < dim; ++d) f1[d] = u[2]*u_x[d];
126: }
127: /* I */
128: void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux,
129: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
130: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
131: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
132: {
133: g0[0] = 1.0;
134: }
135: /* 6 (x + y) */
136: void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
137: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
138: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
139: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
140: {
141: f0[0] = 6.0*(x[0] + x[1]);
142: }
143: /* a \nabla u */
144: void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
145: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
146: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
147: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
148: {
149: PetscInt d;
150: for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[d];
151: }
152: /* \nabla u */
153: void g2_la(PetscInt dim, PetscInt Nf, PetscInt NfAux,
154: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
155: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
156: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
157: {
158: PetscInt d;
159: for (d = 0; d < dim; ++d) g2[d] = u_x[d];
160: }
161: /* a */
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] = u[1];
169: }
171: /*
172: In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:
174: u = x^2 + y^2
175: f = 6 (x + y)
176: kappa(a) = a = (x + y)
178: so that
180: -\div \kappa(a) \grad u + f = -6 (x + y) + 6 (x + y) = 0
181: */
182: PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
183: {
184: *u = x[0]*x[0] + x[1]*x[1];
185: return 0;
186: }
187: PetscErrorCode linear_a_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx)
188: {
189: *a = x[0] + x[1];
190: return 0;
191: }
192: PetscErrorCode zero(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx)
193: {
194: *l = 0.0;
195: return 0;
196: }
198: PetscErrorCode SetupProblem(DM dm, AppCtx *user)
199: {
200: PetscDS prob;
201: const PetscInt id = 1;
205: DMGetDS (dm, &prob);
206: PetscDSSetResidual (prob, 0, f0_u, f1_u);
207: PetscDSSetResidual (prob, 1, f0_a, f1_a);
208: PetscDSSetResidual (prob, 2, f0_l, f1_l);
209: PetscDSSetJacobian (prob, 0, 0, g0_uu, NULL, NULL, NULL);
210: PetscDSSetJacobian (prob, 0, 1, NULL, NULL, g2_ua, NULL);
211: PetscDSSetJacobian (prob, 0, 2, NULL, NULL, NULL, g3_ul);
212: PetscDSSetJacobian (prob, 1, 1, g0_aa, NULL, NULL, NULL);
213: PetscDSSetJacobian (prob, 2, 1, NULL, NULL, g2_la, NULL);
214: PetscDSSetJacobian (prob, 2, 0, NULL, NULL, NULL, g3_lu);
216: user->exactFuncs[0] = quadratic_u_2d;
217: user->exactFuncs[1] = linear_a_2d;
218: user->exactFuncs[2] = zero;
219: PetscDSAddBoundary (prob, DM_BC_ESSENTIAL , "wall" , "marker" , 0, 0, NULL, (void (*)(void)) user->exactFuncs[0], 1, &id, user);
220: PetscDSAddBoundary (prob, DM_BC_ESSENTIAL , "wall" , "marker" , 1, 0, NULL, (void (*)(void)) user->exactFuncs[1], 1, &id, user);
221: PetscDSAddBoundary (prob, DM_BC_ESSENTIAL , "wall" , "marker" , 2, 0, NULL, (void (*)(void)) user->exactFuncs[2], 1, &id, user);
222: return (0);
223: }
225: PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
226: {
227: DM cdm = dm;
228: const PetscInt dim = 2;
229: PetscFE fe[3];
230: PetscInt f;
231: MPI_Comm comm;
232: PetscErrorCode ierr;
235: /* Create finite element */
236: PetscObjectGetComm ((PetscObject ) dm, &comm);
237: PetscFECreateDefault (comm, dim, 1, PETSC_TRUE , "potential_" , -1, &fe[0]);
238: PetscObjectSetName ((PetscObject ) fe[0], "potential" );
239: PetscFECreateDefault (comm, dim, 1, PETSC_TRUE , "conductivity_" , -1, &fe[1]);
240: PetscObjectSetName ((PetscObject ) fe[1], "conductivity" );
241: PetscFECopyQuadrature (fe[0], fe[1]);
242: PetscFECreateDefault (comm, dim, 1, PETSC_TRUE , "multiplier_" , -1, &fe[2]);
243: PetscObjectSetName ((PetscObject ) fe[2], "multiplier" );
244: PetscFECopyQuadrature (fe[0], fe[2]);
245: /* Set discretization and boundary conditions for each mesh */
246: for (f = 0; f < 3; ++f) {DMSetField (dm, f, NULL, (PetscObject ) fe[f]);}
247: DMCreateDS (dm);
248: SetupProblem(dm, user);
249: while (cdm) {
250: DMCopyDisc (dm, cdm);
251: DMGetCoarseDM (cdm, &cdm);
252: }
253: for (f = 0; f < 3; ++f) {PetscFEDestroy (&fe[f]);}
254: return (0);
255: }
257: int main(int argc, char **argv)
258: {
259: DM dm;
260: SNES snes;
261: Vec u, r;
262: AppCtx user;
265: PetscInitialize (&argc, &argv, NULL,help);if (ierr) return ierr;
266: ProcessOptions(PETSC_COMM_WORLD , &user);
267: SNESCreate (PETSC_COMM_WORLD , &snes);
268: CreateMesh(PETSC_COMM_WORLD , &user, &dm);
269: SNESSetDM (snes, dm);
271: PetscMalloc (3 * sizeof (void (*)()), &user.exactFuncs);
272: SetupDiscretization(dm, &user);
274: DMCreateGlobalVector (dm, &u);
275: PetscObjectSetName ((PetscObject ) u, "solution" );
276: VecDuplicate (u, &r);
277: DMPlexSetSNESLocalFEM (dm,&user,&user,&user);
278: SNESSetFromOptions (snes);
280: DMProjectFunction (dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES , u);
281: DMSNESCheckFromOptions (snes, u, user.exactFuncs, NULL);
282: if (user.runType == RUN_FULL) {
283: PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
284: PetscReal error;
286: initialGuess[0] = zero;
287: initialGuess[1] = zero;
288: initialGuess[2] = zero;
289: DMProjectFunction (dm, 0.0, initialGuess, NULL, INSERT_VALUES , u);
290: VecViewFromOptions (u, NULL, "-initial_vec_view" );
291: DMComputeL2Diff (dm, 0.0, user.exactFuncs, NULL, u, &error);
292: if (error < 1.0e-11) {PetscPrintf (PETSC_COMM_WORLD , "Initial L_2 Error: < 1.0e-11\n" );}
293: else {PetscPrintf (PETSC_COMM_WORLD , "Initial L_2 Error: %g\n" , error);}
294: SNESSolve (snes, NULL, u);
295: DMComputeL2Diff (dm, 0.0, user.exactFuncs, NULL, u, &error);
296: if (error < 1.0e-11) {PetscPrintf (PETSC_COMM_WORLD , "Final L_2 Error: < 1.0e-11\n" );}
297: else {PetscPrintf (PETSC_COMM_WORLD , "Final L_2 Error: %g\n" , error);}
298: }
299: VecViewFromOptions (u, NULL, "-sol_vec_view" );
301: VecDestroy (&u);
302: VecDestroy (&r);
303: SNESDestroy (&snes);
304: DMDestroy (&dm);
305: PetscFree (user.exactFuncs);
306: PetscFinalize ();
307: return ierr;
308: }
310: /*TEST
312: build:
313: requires: !complex
315: test:
316: suffix: 0
317: requires: triangle
318: args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -conductivity_petscspace_degree 1 -multiplier_petscspace_degree 2
320: test:
321: suffix: 1
322: requires: triangle
323: 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
325: TEST*/