Actual source code: ex2.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 $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: PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
31: } AppCtx;
33: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
34: {
35: const char *runTypes[2] = {"full" , "test" };
36: PetscInt run;
40: options->runType = RUN_FULL;
41: options->useDualPenalty = PETSC_FALSE ;
43: PetscOptionsBegin (comm, "" , "Inverse Problem Options" , "DMPLEX " );
44: run = options->runType;
45: PetscOptionsEList ("-run_type" , "The run type" , "ex2.c" , runTypes, 2, runTypes[options->runType], &run, NULL);
46: options->runType = (RunType) run;
47: PetscOptionsBool ("-use_dual_penalty" , "Penalize deviation from both goals" , "ex2.c" , options->useDualPenalty, &options->useDualPenalty, NULL);
48: PetscOptionsEnd ();
49: return (0);
50: }
52: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
53: {
54: DM distributedMesh = NULL;
58: DMPlexCreateBoxMesh (comm, 2, PETSC_TRUE , NULL, NULL, NULL, NULL, PETSC_TRUE , dm);
59: PetscObjectSetName ((PetscObject ) *dm, "Mesh" );
60: DMPlexDistribute (*dm, 0, NULL, &distributedMesh);
61: if (distributedMesh) {
62: DMDestroy (dm);
63: *dm = distributedMesh;
64: }
65: DMSetFromOptions (*dm);
66: DMViewFromOptions (*dm, NULL, "-dm_view" );
67: return (0);
68: }
70: void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
71: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
72: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
73: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
74: {
75: f0[0] = (u[0] - (x[0]*x[0] + x[1]*x[1]));
76: }
77: void f0_u_full(PetscInt dim, PetscInt Nf, PetscInt NfAux,
78: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
79: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
80: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
81: {
82: 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]))) +
83: 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])));
84: }
85: void f1_u(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, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
89: {
90: PetscInt d;
91: for (d = 0; d < dim; ++d) f1[d] = u_x[dim*2+d];
92: }
93: void g0_uu(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 g0[])
97: {
98: g0[0] = 1.0;
99: }
100: void g0_uu_full(PetscInt dim, PetscInt Nf, PetscInt NfAux,
101: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
102: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
103: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
104: {
105: g0[0] = PetscSqr (u[0] - sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1]))
106: + PetscSqr (u[0] - (x[0]*x[0] + x[1]*x[1]))
107: - 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]
108: + 4.0*(x[0]*x[0] + x[1]*x[1])*(sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1]));
109: }
110: void g3_ul(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
114: {
115: PetscInt d;
116: for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
117: }
119: void f0_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 f0[])
123: {
124: f0[0] = u[1] - 4.0 /* 0.0 */ + u[2];
125: }
126: void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux,
127: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
128: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
129: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
130: {
131: g0[0] = 1.0;
132: }
133: void g0_al(PetscInt dim, PetscInt Nf, PetscInt NfAux,
134: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
135: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
136: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
137: {
138: g0[0] = 1.0;
139: }
141: void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
142: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
143: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
144: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
145: {
146: f0[0] = u[1];
147: }
148: void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
149: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
150: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
151: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
152: {
153: PetscInt d;
154: for (d = 0; d < dim; ++d) f1[d] = u_x[d];
155: }
156: void g0_la(PetscInt dim, PetscInt Nf, PetscInt NfAux,
157: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
158: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
159: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
160: {
161: g0[0] = 1.0;
162: }
163: void g3_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
164: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
165: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
166: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
167: {
168: PetscInt d;
169: for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
170: }
172: /*
173: In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:
175: u = x^2 + y^2
176: a = 4
177: d_A = 4
178: d_B = sin(2*pi*x[0]) * sin(2*pi*x[1])
180: so that
182: -\Delta u + a = -4 + 4 = 0
183: */
184: PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
185: {
186: *u = x[0]*x[0] + x[1]*x[1];
187: return 0;
188: }
189: PetscErrorCode constant_a_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx)
190: {
191: *a = 4;
192: return 0;
193: }
194: PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx)
195: {
196: *l = 0.0;
197: return 0;
198: }
200: PetscErrorCode SetupProblem(DM dm, AppCtx *user)
201: {
202: PetscDS prob;
203: const PetscInt id = 1;
207: DMGetDS (dm, &prob);
208: PetscDSSetResidual (prob, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u);
209: PetscDSSetResidual (prob, 1, f0_a, NULL);
210: PetscDSSetResidual (prob, 2, f0_l, f1_l);
211: PetscDSSetJacobian (prob, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL);
212: PetscDSSetJacobian (prob, 0, 2, NULL, NULL, NULL, g3_ul);
213: PetscDSSetJacobian (prob, 1, 1, g0_aa, NULL, NULL, NULL);
214: PetscDSSetJacobian (prob, 1, 2, g0_al, NULL, NULL, NULL);
215: PetscDSSetJacobian (prob, 2, 1, g0_la, NULL, NULL, NULL);
216: PetscDSSetJacobian (prob, 2, 0, NULL, NULL, NULL, g3_lu);
218: user->exactFuncs[0] = quadratic_u_2d;
219: user->exactFuncs[1] = constant_a_2d;
220: user->exactFuncs[2] = zero;
221: PetscDSAddBoundary (prob, DM_BC_ESSENTIAL , "wall" , "marker" , 0, 0, NULL, (void (*)()) user->exactFuncs[0], 1, &id, user);
222: PetscDSAddBoundary (prob, DM_BC_ESSENTIAL , "wall" , "marker" , 1, 0, NULL, (void (*)()) user->exactFuncs[1], 1, &id, user);
223: PetscDSAddBoundary (prob, DM_BC_ESSENTIAL , "wall" , "marker" , 2, 0, NULL, (void (*)()) user->exactFuncs[2], 1, &id, user);
224: return (0);
225: }
227: PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
228: {
229: DM cdm = dm;
230: const PetscInt dim = 2;
231: PetscFE fe[3];
232: PetscInt f;
233: MPI_Comm comm;
234: PetscErrorCode ierr;
237: /* Create finite element */
238: PetscObjectGetComm ((PetscObject ) dm, &comm);
239: PetscFECreateDefault (comm, dim, 1, PETSC_TRUE , "potential_" , -1, &fe[0]);
240: PetscObjectSetName ((PetscObject ) fe[0], "potential" );
241: PetscFECreateDefault (comm, dim, 1, PETSC_TRUE , "charge_" , -1, &fe[1]);
242: PetscObjectSetName ((PetscObject ) fe[1], "charge" );
243: PetscFECopyQuadrature (fe[0], fe[1]);
244: PetscFECreateDefault (comm, dim, 1, PETSC_TRUE , "multiplier_" , -1, &fe[2]);
245: PetscObjectSetName ((PetscObject ) fe[2], "multiplier" );
246: PetscFECopyQuadrature (fe[0], fe[2]);
247: /* Set discretization and boundary conditions for each mesh */
248: for (f = 0; f < 3; ++f) {DMSetField (dm, f, NULL, (PetscObject ) fe[f]);}
249: DMCreateDS (cdm);
250: SetupProblem(dm, user);
251: while (cdm) {
252: DMCopyDisc (dm, cdm);
253: DMGetCoarseDM (cdm, &cdm);
254: }
255: for (f = 0; f < 3; ++f) {PetscFEDestroy (&fe[f]);}
256: return (0);
257: }
259: int main(int argc, char **argv)
260: {
261: DM dm;
262: SNES snes;
263: Vec u, r;
264: AppCtx user;
267: PetscInitialize (&argc, &argv, NULL,help);if (ierr) return ierr;
268: ProcessOptions(PETSC_COMM_WORLD , &user);
269: SNESCreate (PETSC_COMM_WORLD , &snes);
270: CreateMesh(PETSC_COMM_WORLD , &user, &dm);
271: SNESSetDM (snes, dm);
273: PetscMalloc (3 * sizeof (void (*)()), &user.exactFuncs);
274: SetupDiscretization(dm, &user);
276: DMCreateGlobalVector (dm, &u);
277: PetscObjectSetName ((PetscObject ) u, "solution" );
278: VecDuplicate (u, &r);
279: DMPlexSetSNESLocalFEM (dm,&user,&user,&user);
280: SNESSetFromOptions (snes);
282: DMProjectFunction (dm, 0.0, user.exactFuncs, NULL, INSERT_ALL_VALUES , u);
283: DMSNESCheckFromOptions (snes, u, user.exactFuncs, NULL);
284: if (user.runType == RUN_FULL) {
285: PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
286: PetscReal error;
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, user.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, user.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: PetscFree (user.exactFuncs);
308: PetscFinalize ();
309: return ierr;
310: }
312: /*TEST
314: build:
315: requires: !complex triangle
317: test:
318: suffix: 0
319: args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1
321: test:
322: suffix: 1
323: 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
325: test:
326: suffix: 2
327: 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
329: TEST*/