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

 53:   DMCreate(comm, dm);
 54:   DMSetType(*dm, DMPLEX);
 55:   DMSetFromOptions(*dm);
 56:   DMViewFromOptions(*dm, NULL, "-dm_view");
 57:   return(0);
 58: }

 60: /* u - (x^2 + y^2) */
 61: void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
 62:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
 63:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
 64:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
 65: {
 66:   f0[0] = u[0] - (x[0]*x[0] + x[1]*x[1]);
 67: }
 68: /* a \nabla\lambda */
 69: void f1_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 f1[])
 73: {
 74:   PetscInt d;
 75:   for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[dim*2+d];
 76: }
 77: /* I */
 78: void g0_uu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
 82: {
 83:   g0[0] = 1.0;
 84: }
 85: /* \nabla */
 86: void g2_ua(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 g2[])
 90: {
 91:   PetscInt d;
 92:   for (d = 0; d < dim; ++d) g2[d] = u_x[dim*2+d];
 93: }
 94: /* a */
 95: void g3_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux,
 96:            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
 97:            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
 98:            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
 99: {
100:   PetscInt d;
101:   for (d = 0; d < dim; ++d) g3[d*dim+d] = u[1];
102: }
103: /* a - (x + y) */
104: void f0_a(PetscInt dim, PetscInt Nf, PetscInt NfAux,
105:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
106:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
107:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
108: {
109:   f0[0] = u[1] - (x[0] + x[1]);
110: }
111: /* \lambda \nabla u */
112: void f1_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 f1[])
116: {
117:   PetscInt d;
118:   for (d = 0; d < dim; ++d) f1[d] = u[2]*u_x[d];
119: }
120: /* I */
121: void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux,
122:            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
123:            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
124:            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
125: {
126:   g0[0] = 1.0;
127: }
128: /* 6 (x + y) */
129: void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
130:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
131:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
132:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
133: {
134:   f0[0] = 6.0*(x[0] + x[1]);
135: }
136: /* a \nabla u */
137: void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
138:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
139:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
140:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
141: {
142:   PetscInt d;
143:   for (d = 0; d < dim; ++d) f1[d] = u[1]*u_x[d];
144: }
145: /* \nabla u */
146: void g2_la(PetscInt dim, PetscInt Nf, PetscInt NfAux,
147:            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
148:            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
149:            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
150: {
151:   PetscInt d;
152:   for (d = 0; d < dim; ++d) g2[d] = u_x[d];
153: }
154: /* a */
155: void g3_lu(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 g3[])
159: {
160:   PetscInt d;
161:   for (d = 0; d < dim; ++d) g3[d*dim+d] = u[1];
162: }

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

167:     u  = x^2 + y^2
168:     f  = 6 (x + y)
169:     kappa(a) = a = (x + y)

171:   so that

173:     -\div \kappa(a) \grad u + f = -6 (x + y) + 6 (x + y) = 0
174: */
175: PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
176: {
177:   *u = x[0]*x[0] + x[1]*x[1];
178:   return 0;
179: }
180: PetscErrorCode linear_a_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx)
181: {
182:   *a = x[0] + x[1];
183:   return 0;
184: }
185: PetscErrorCode zero(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx)
186: {
187:   *l = 0.0;
188:   return 0;
189: }

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

199:   DMGetDS(dm, &ds);
200:   PetscDSSetResidual(ds, 0, f0_u, f1_u);
201:   PetscDSSetResidual(ds, 1, f0_a, f1_a);
202:   PetscDSSetResidual(ds, 2, f0_l, f1_l);
203:   PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, NULL);
204:   PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_ua, NULL);
205:   PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul);
206:   PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL);
207:   PetscDSSetJacobian(ds, 2, 1, NULL, NULL, g2_la, NULL);
208:   PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu);

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

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

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

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

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

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

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

280:     DMGetDS(dm, &ds);
281:     PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL);
282:     PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL);
283:     PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL);
284:     initialGuess[0] = zero;
285:     initialGuess[1] = zero;
286:     initialGuess[2] = zero;
287:     DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);
288:     VecViewFromOptions(u, NULL, "-initial_vec_view");
289:     DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);
290:     if (error < 1.0e-11) {PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n");}
291:     else                 {PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", error);}
292:     SNESSolve(snes, NULL, u);
293:     DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);
294:     if (error < 1.0e-11) {PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n");}
295:     else                 {PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", error);}
296:   }
297:   VecViewFromOptions(u, NULL, "-sol_vec_view");

299:   VecDestroy(&u);
300:   VecDestroy(&r);
301:   SNESDestroy(&snes);
302:   DMDestroy(&dm);
303:   PetscFinalize();
304:   return ierr;
305: }

307: /*TEST

309:   build:
310:     requires: !complex

312:   test:
313:     suffix: 0
314:     requires: triangle
315:     args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -conductivity_petscspace_degree 1 -multiplier_petscspace_degree 2

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

322: TEST*/