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: {

 56:   DMCreate(comm, dm);
 57:   DMSetType(*dm, DMPLEX);
 58:   DMSetFromOptions(*dm);
 59:   DMViewFromOptions(*dm, NULL, "-dm_view");
 60:   return(0);
 61: }

 63: void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
 64:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
 65:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
 66:           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]));
 69: }
 70: void f0_u_full(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]))*PetscSqr(u[0] - (sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1]))) +
 76:     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])));
 77: }
 78: void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
 79:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
 80:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
 81:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
 82: {
 83:   PetscInt d;
 84:   for (d = 0; d < dim; ++d) f1[d] = u_x[dim*2+d];
 85: }
 86: void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
 87:            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
 88:            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
 89:            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
 90: {
 91:   g0[0] = 1.0;
 92: }
 93: void g0_uu_full(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] = PetscSqr(u[0] - sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1]))
 99:     + PetscSqr(u[0] - (x[0]*x[0] + x[1]*x[1]))
100:     - 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]
101:     + 4.0*(x[0]*x[0] + x[1]*x[1])*(sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1]));
102: }
103: void g3_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux,
104:            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
105:            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
106:            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
107: {
108:   PetscInt d;
109:   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
110: }

112: void f0_a(PetscInt dim, PetscInt Nf, PetscInt NfAux,
113:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
114:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
115:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
116: {
117:   f0[0] = u[1] - 4.0 /* 0.0 */ + u[2];
118: }
119: void g0_aa(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
123: {
124:   g0[0] = 1.0;
125: }
126: void g0_al(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: }

134: void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
135:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
136:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
137:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
138: {
139:   f0[0] = u[1];
140: }
141: void f1_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 f1[])
145: {
146:   PetscInt d;
147:   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
148: }
149: void g0_la(PetscInt dim, PetscInt Nf, PetscInt NfAux,
150:            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
151:            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
152:            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
153: {
154:   g0[0] = 1.0;
155: }
156: void g3_lu(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 g3[])
160: {
161:   PetscInt d;
162:   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
163: }

165: /*
166:   In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:

168:     u   = x^2 + y^2
169:     a   = 4
170:     d_A = 4
171:     d_B = sin(2*pi*x[0]) * sin(2*pi*x[1])

173:   so that

175:     -\Delta u + a = -4 + 4 = 0
176: */
177: PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
178: {
179:   *u = x[0]*x[0] + x[1]*x[1];
180:   return 0;
181: }
182: PetscErrorCode constant_a_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx)
183: {
184:   *a = 4;
185:   return 0;
186: }
187: PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx)
188: {
189:   *l = 0.0;
190:   return 0;
191: }

193: PetscErrorCode SetupProblem(DM dm, AppCtx *user)
194: {
195:   PetscDS        ds;
196:   DMLabel        label;
197:   const PetscInt id = 1;

201:   DMGetDS(dm, &ds);
202:   PetscDSSetResidual(ds, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u);
203:   PetscDSSetResidual(ds, 1, f0_a, NULL);
204:   PetscDSSetResidual(ds, 2, f0_l, f1_l);
205:   PetscDSSetJacobian(ds, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL);
206:   PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul);
207:   PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL);
208:   PetscDSSetJacobian(ds, 1, 2, g0_al, NULL, NULL, NULL);
209:   PetscDSSetJacobian(ds, 2, 1, g0_la, NULL, NULL, NULL);
210:   PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu);

212:   PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL);
213:   PetscDSSetExactSolution(ds, 1, constant_a_2d, NULL);
214:   PetscDSSetExactSolution(ds, 2, zero, NULL);
215:   DMGetLabel(dm, "marker", &label);
216:   DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)()) quadratic_u_2d, NULL, user, NULL);
217:   DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)()) constant_a_2d, NULL, user, NULL);
218:   DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)()) zero, NULL, user, NULL);
219:   return(0);
220: }

222: PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
223: {
224:   DM              cdm = dm;
225:   const PetscInt  dim = 2;
226:   PetscFE         fe[3];
227:   PetscInt        f;
228:   MPI_Comm        comm;
229:   PetscErrorCode  ierr;

232:   /* Create finite element */
233:   PetscObjectGetComm((PetscObject) dm, &comm);
234:   PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]);
235:   PetscObjectSetName((PetscObject) fe[0], "potential");
236:   PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "charge_", -1, &fe[1]);
237:   PetscObjectSetName((PetscObject) fe[1], "charge");
238:   PetscFECopyQuadrature(fe[0], fe[1]);
239:   PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]);
240:   PetscObjectSetName((PetscObject) fe[2], "multiplier");
241:   PetscFECopyQuadrature(fe[0], fe[2]);
242:   /* Set discretization and boundary conditions for each mesh */
243:   for (f = 0; f < 3; ++f) {DMSetField(dm, f, NULL, (PetscObject) fe[f]);}
244:   DMCreateDS(cdm);
245:   SetupProblem(dm, user);
246:   while (cdm) {
247:     DMCopyDisc(dm, cdm);
248:     DMGetCoarseDM(cdm, &cdm);
249:   }
250:   for (f = 0; f < 3; ++f) {PetscFEDestroy(&fe[f]);}
251:   return(0);
252: }

254: int main(int argc, char **argv)
255: {
256:   DM             dm;
257:   SNES           snes;
258:   Vec            u, r;
259:   AppCtx         user;

262:   PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
263:   ProcessOptions(PETSC_COMM_WORLD, &user);
264:   SNESCreate(PETSC_COMM_WORLD, &snes);
265:   CreateMesh(PETSC_COMM_WORLD, &user, &dm);
266:   SNESSetDM(snes, dm);
267:   SetupDiscretization(dm, &user);

269:   DMCreateGlobalVector(dm, &u);
270:   PetscObjectSetName((PetscObject) u, "solution");
271:   VecDuplicate(u, &r);
272:   DMPlexSetSNESLocalFEM(dm,&user,&user,&user);
273:   SNESSetFromOptions(snes);

275:   DMSNESCheckFromOptions(snes, u);
276:   if (user.runType == RUN_FULL) {
277:     PetscDS          ds;
278:     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
279:     PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
280:     PetscReal        error;

282:     DMGetDS(dm, &ds);
283:     PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL);
284:     PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL);
285:     PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL);
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, 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, 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:   PetscFinalize();
306:   return ierr;
307: }

309: /*TEST

311:   build:
312:     requires: !complex triangle

314:   test:
315:     suffix: 0
316:     args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1

318:   test:
319:     suffix: 1
320:     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

322:   test:
323:     suffix: 2
324:     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

326: TEST*/