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";

  4: /*F
  5:   We are solving the parameter estimation problem for the Laplacian. We will ask to minimize a Lagrangian
  6: function over $y$ and $u$, given by
  7: \begin{align}
  8:   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)
  9: \end{align}
 10: where $Q$ is a sampling operator, $L$ is a regularization operator, $F$ defines the PDE.

 12: Currently, we have perfect information, meaning $Q = I$, and then we need no regularization, $L = I$. We
 13: also give the null vector for the reference control $a_r$. Right now $\beta = 1$.

 15: The PDE will be the Laplace equation with homogeneous boundary conditions
 16: \begin{align}
 17:   -Delta u = a
 18: \end{align}

 20: F*/

 22: #include <petsc.h>
 23: #include <petscfe.h>

 25: typedef enum {
 26:   RUN_FULL,
 27:   RUN_TEST
 28: } RunType;

 30: typedef struct {
 31:   RunType   runType;        /* Whether to run tests, or solve the full problem */
 32:   PetscBool useDualPenalty; /* Penalize deviation from both goals */
 33: } AppCtx;

 35: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
 36: {
 37:   const char *runTypes[2] = {"full", "test"};
 38:   PetscInt    run;

 40:   PetscFunctionBeginUser;
 41:   options->runType        = RUN_FULL;
 42:   options->useDualPenalty = PETSC_FALSE;
 43:   PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");
 44:   run = options->runType;
 45:   PetscCall(PetscOptionsEList("-run_type", "The run type", "ex2.c", runTypes, 2, runTypes[options->runType], &run, NULL));
 46:   options->runType = (RunType)run;
 47:   PetscCall(PetscOptionsBool("-use_dual_penalty", "Penalize deviation from both goals", "ex2.c", options->useDualPenalty, &options->useDualPenalty, NULL));
 48:   PetscOptionsEnd();
 49:   PetscFunctionReturn(PETSC_SUCCESS);
 50: }

 52: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
 53: {
 54:   PetscFunctionBeginUser;
 55:   PetscCall(DMCreate(comm, dm));
 56:   PetscCall(DMSetType(*dm, DMPLEX));
 57:   PetscCall(DMSetFromOptions(*dm));
 58:   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
 59:   PetscFunctionReturn(PETSC_SUCCESS);
 60: }

 62: void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
 63: {
 64:   f0[0] = (u[0] - (x[0] * x[0] + x[1] * x[1]));
 65: }
 66: void f0_u_full(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 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])) * PetscSqr(u[0] - (sin(2.0 * PETSC_PI * x[0]) * sin(2.0 * PETSC_PI * x[1]))) + 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])));
 69: }
 70: void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
 71: {
 72:   PetscInt d;
 73:   for (d = 0; d < dim; ++d) f1[d] = u_x[dim * 2 + d];
 74: }
 75: void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
 76: {
 77:   g0[0] = 1.0;
 78: }
 79: void g0_uu_full(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
 80: {
 81:   g0[0] = PetscSqr(u[0] - sin(2.0 * PETSC_PI * x[0]) * sin(2.0 * PETSC_PI * x[1])) + PetscSqr(u[0] - (x[0] * x[0] + x[1] * x[1])) - 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] + 4.0 * (x[0] * x[0] + x[1] * x[1]) * (sin(2.0 * PETSC_PI * x[0]) * sin(2.0 * PETSC_PI * x[1]));
 82: }
 83: void g3_ul(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
 84: {
 85:   for (PetscInt d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
 86: }

 88: void f0_a(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
 89: {
 90:   f0[0] = u[1] - 4.0 /* 0.0 */ + u[2];
 91: }
 92: void g0_aa(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
 93: {
 94:   g0[0] = 1.0;
 95: }
 96: void g0_al(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
 97: {
 98:   g0[0] = 1.0;
 99: }

101: void f0_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
102: {
103:   f0[0] = u[1];
104: }
105: void f1_l(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
106: {
107:   for (PetscInt d = 0; d < dim; ++d) f1[d] = u_x[d];
108: }
109: void g0_la(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
110: {
111:   g0[0] = 1.0;
112: }
113: void g3_lu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
114: {
115:   for (PetscInt d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
116: }

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

121:     u   = x^2 + y^2
122:     a   = 4
123:     d_A = 4
124:     d_B = sin(2*pi*x[0]) * sin(2*pi*x[1])

126:   so that

128:     -\Delta u + a = -4 + 4 = 0
129: */
130: PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
131: {
132:   *u = x[0] * x[0] + x[1] * x[1];
133:   return PETSC_SUCCESS;
134: }
135: PetscErrorCode constant_a_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *a, PetscCtx ctx)
136: {
137:   *a = 4;
138:   return PETSC_SUCCESS;
139: }
140: PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *l, PetscCtx ctx)
141: {
142:   *l = 0.0;
143:   return PETSC_SUCCESS;
144: }

146: PetscErrorCode SetupProblem(DM dm, AppCtx *user)
147: {
148:   PetscDS        ds;
149:   DMLabel        label;
150:   const PetscInt id = 1;

152:   PetscFunctionBeginUser;
153:   PetscCall(DMGetDS(dm, &ds));
154:   PetscCall(PetscDSSetResidual(ds, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u));
155:   PetscCall(PetscDSSetResidual(ds, 1, f0_a, NULL));
156:   PetscCall(PetscDSSetResidual(ds, 2, f0_l, f1_l));
157:   PetscCall(PetscDSSetJacobian(ds, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL));
158:   PetscCall(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul));
159:   PetscCall(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL));
160:   PetscCall(PetscDSSetJacobian(ds, 1, 2, g0_al, NULL, NULL, NULL));
161:   PetscCall(PetscDSSetJacobian(ds, 2, 1, g0_la, NULL, NULL, NULL));
162:   PetscCall(PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu));

164:   PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL));
165:   PetscCall(PetscDSSetExactSolution(ds, 1, constant_a_2d, NULL));
166:   PetscCall(PetscDSSetExactSolution(ds, 2, zero, NULL));
167:   PetscCall(DMGetLabel(dm, "marker", &label));
168:   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)quadratic_u_2d, NULL, user, NULL));
169:   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (PetscVoidFn *)constant_a_2d, NULL, user, NULL));
170:   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (PetscVoidFn *)zero, NULL, user, NULL));
171:   PetscFunctionReturn(PETSC_SUCCESS);
172: }

174: PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
175: {
176:   DM             cdm = dm;
177:   const PetscInt dim = 2;
178:   PetscFE        fe[3];
179:   MPI_Comm       comm;

181:   PetscFunctionBeginUser;
182:   /* Create finite element */
183:   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
184:   PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]));
185:   PetscCall(PetscObjectSetName((PetscObject)fe[0], "potential"));
186:   PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "charge_", -1, &fe[1]));
187:   PetscCall(PetscObjectSetName((PetscObject)fe[1], "charge"));
188:   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
189:   PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]));
190:   PetscCall(PetscObjectSetName((PetscObject)fe[2], "multiplier"));
191:   PetscCall(PetscFECopyQuadrature(fe[0], fe[2]));
192:   /* Set discretization and boundary conditions for each mesh */
193:   for (PetscInt f = 0; f < 3; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe[f]));
194:   PetscCall(DMCreateDS(cdm));
195:   PetscCall(SetupProblem(dm, user));
196:   while (cdm) {
197:     PetscCall(DMCopyDisc(dm, cdm));
198:     PetscCall(DMGetCoarseDM(cdm, &cdm));
199:   }
200:   for (PetscInt f = 0; f < 3; ++f) PetscCall(PetscFEDestroy(&fe[f]));
201:   PetscFunctionReturn(PETSC_SUCCESS);
202: }

204: int main(int argc, char **argv)
205: {
206:   DM     dm;
207:   SNES   snes;
208:   Vec    u, r;
209:   AppCtx user;

211:   PetscFunctionBeginUser;
212:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
213:   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
214:   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
215:   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
216:   PetscCall(SNESSetDM(snes, dm));
217:   PetscCall(SetupDiscretization(dm, &user));

219:   PetscCall(DMCreateGlobalVector(dm, &u));
220:   PetscCall(PetscObjectSetName((PetscObject)u, "solution"));
221:   PetscCall(VecDuplicate(u, &r));
222:   PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
223:   PetscCall(SNESSetFromOptions(snes));

225:   PetscCall(DMSNESCheckFromOptions(snes, u));
226:   if (user.runType == RUN_FULL) {
227:     PetscDS ds;
228:     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
229:     PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], PetscCtx ctx);
230:     PetscReal error;

232:     PetscCall(DMGetDS(dm, &ds));
233:     PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL));
234:     PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL));
235:     PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL));
236:     initialGuess[0] = zero;
237:     initialGuess[1] = zero;
238:     initialGuess[2] = zero;
239:     PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u));
240:     PetscCall(VecViewFromOptions(u, NULL, "-initial_vec_view"));
241:     PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error));
242:     if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n"));
243:     else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", (double)error));
244:     PetscCall(SNESSolve(snes, NULL, u));
245:     PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error));
246:     if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n"));
247:     else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", (double)error));
248:   }
249:   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));

251:   PetscCall(VecDestroy(&u));
252:   PetscCall(VecDestroy(&r));
253:   PetscCall(SNESDestroy(&snes));
254:   PetscCall(DMDestroy(&dm));
255:   PetscCall(PetscFinalize());
256:   return 0;
257: }

259: /*TEST

261:   build:
262:     requires: !complex triangle

264:   test:
265:     suffix: 0
266:     args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1

268:   test:
269:     suffix: 1
270:     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

272:   test:
273:     suffix: 2
274:     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

276: TEST*/