Actual source code: ex2.c

petsc-main 2021-04-20
Report Typos and Errors
  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: {
 53:   DM             distributedMesh = NULL;

 57:   DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);
 58:   PetscObjectSetName((PetscObject) *dm, "Mesh");
 59:   DMPlexDistribute(*dm, 0, NULL, &distributedMesh);
 60:   if (distributedMesh) {
 61:     DMDestroy(dm);
 62:     *dm  = distributedMesh;
 63:   }
 64:   DMSetFromOptions(*dm);
 65:   DMViewFromOptions(*dm, NULL, "-dm_view");
 66:   return(0);
 67: }

 69: void f0_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 f0[])
 73: {
 74:   f0[0] = (u[0] - (x[0]*x[0] + x[1]*x[1]));
 75: }
 76: void f0_u_full(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 f0[])
 80: {
 81:   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]))) +
 82:     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])));
 83: }
 84: void f1_u(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, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
 88: {
 89:   PetscInt d;
 90:   for (d = 0; d < dim; ++d) f1[d] = u_x[dim*2+d];
 91: }
 92: void g0_uu(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 g0[])
 96: {
 97:   g0[0] = 1.0;
 98: }
 99: void g0_uu_full(PetscInt dim, PetscInt Nf, PetscInt NfAux,
100:                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
101:                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
102:                 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
103: {
104:   g0[0] = PetscSqr(u[0] - sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1]))
105:     + PetscSqr(u[0] - (x[0]*x[0] + x[1]*x[1]))
106:     - 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]
107:     + 4.0*(x[0]*x[0] + x[1]*x[1])*(sin(2.0*PETSC_PI*x[0]) * sin(2.0*PETSC_PI*x[1]));
108: }
109: void g3_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux,
110:            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
111:            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
112:            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
113: {
114:   PetscInt d;
115:   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
116: }

118: void f0_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 f0[])
122: {
123:   f0[0] = u[1] - 4.0 /* 0.0 */ + u[2];
124: }
125: void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux,
126:            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
127:            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
128:            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
129: {
130:   g0[0] = 1.0;
131: }
132: void g0_al(PetscInt dim, PetscInt Nf, PetscInt NfAux,
133:            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
134:            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
135:            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
136: {
137:   g0[0] = 1.0;
138: }

140: void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
141:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
142:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
143:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
144: {
145:   f0[0] = u[1];
146: }
147: void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux,
148:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
149:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
150:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
151: {
152:   PetscInt d;
153:   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
154: }
155: void g0_la(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 g0[])
159: {
160:   g0[0] = 1.0;
161: }
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] = 1.0;
169: }

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

174:     u   = x^2 + y^2
175:     a   = 4
176:     d_A = 4
177:     d_B = sin(2*pi*x[0]) * sin(2*pi*x[1])

179:   so that

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

199: PetscErrorCode SetupProblem(DM dm, AppCtx *user)
200: {
201:   PetscDS        ds;
202:   DMLabel        label;
203:   const PetscInt id = 1;

207:   DMGetDS(dm, &ds);
208:   PetscDSSetResidual(ds, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u);
209:   PetscDSSetResidual(ds, 1, f0_a, NULL);
210:   PetscDSSetResidual(ds, 2, f0_l, f1_l);
211:   PetscDSSetJacobian(ds, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL);
212:   PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul);
213:   PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL);
214:   PetscDSSetJacobian(ds, 1, 2, g0_al, NULL, NULL, NULL);
215:   PetscDSSetJacobian(ds, 2, 1, g0_la, NULL, NULL, NULL);
216:   PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu);

218:   PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL);
219:   PetscDSSetExactSolution(ds, 1, constant_a_2d, NULL);
220:   PetscDSSetExactSolution(ds, 2, zero, NULL);
221:   DMGetLabel(dm, "marker", &label);
222:   DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)()) quadratic_u_2d, NULL, user, NULL);
223:   DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)()) constant_a_2d, NULL, user, NULL);
224:   DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)()) zero, NULL, user, NULL);
225:   return(0);
226: }

228: PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
229: {
230:   DM              cdm = dm;
231:   const PetscInt  dim = 2;
232:   PetscFE         fe[3];
233:   PetscInt        f;
234:   MPI_Comm        comm;
235:   PetscErrorCode  ierr;

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

260: int main(int argc, char **argv)
261: {
262:   DM             dm;
263:   SNES           snes;
264:   Vec            u, r;
265:   AppCtx         user;

268:   PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
269:   ProcessOptions(PETSC_COMM_WORLD, &user);
270:   SNESCreate(PETSC_COMM_WORLD, &snes);
271:   CreateMesh(PETSC_COMM_WORLD, &user, &dm);
272:   SNESSetDM(snes, dm);
273:   SetupDiscretization(dm, &user);

275:   DMCreateGlobalVector(dm, &u);
276:   PetscObjectSetName((PetscObject) u, "solution");
277:   VecDuplicate(u, &r);
278:   DMPlexSetSNESLocalFEM(dm,&user,&user,&user);
279:   SNESSetFromOptions(snes);

281:   DMSNESCheckFromOptions(snes, u);
282:   if (user.runType == RUN_FULL) {
283:     PetscDS          ds;
284:     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
285:     PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
286:     PetscReal        error;

288:     DMGetDS(dm, &ds);
289:     PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL);
290:     PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL);
291:     PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL);
292:     initialGuess[0] = zero;
293:     initialGuess[1] = zero;
294:     initialGuess[2] = zero;
295:     DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);
296:     VecViewFromOptions(u, NULL, "-initial_vec_view");
297:     DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);
298:     if (error < 1.0e-11) {PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n");}
299:     else                 {PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", error);}
300:     SNESSolve(snes, NULL, u);
301:     DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);
302:     if (error < 1.0e-11) {PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n");}
303:     else                 {PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", error);}
304:   }
305:   VecViewFromOptions(u, NULL, "-sol_vec_view");

307:   VecDestroy(&u);
308:   VecDestroy(&r);
309:   SNESDestroy(&snes);
310:   DMDestroy(&dm);
311:   PetscFinalize();
312:   return ierr;
313: }

315: /*TEST

317:   build:
318:     requires: !complex triangle

320:   test:
321:     suffix: 0
322:     args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1

324:   test:
325:     suffix: 1
326:     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

328:   test:
329:     suffix: 2
330:     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

332: TEST*/