xref: /petsc/src/tao/tutorials/ex2.c (revision 2da392cc7c10228af19ad9843ce5155178acb644)
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";
3 
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.
11 
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$.
14 
15 The PDE will be the Laplace equation with homogeneous boundary conditions
16 \begin{align}
17   -Delta u = a
18 \end{align}
19 
20 F*/
21 
22 #include <petsc.h>
23 #include <petscfe.h>
24 
25 typedef enum {RUN_FULL, RUN_TEST} RunType;
26 
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;
31 
32 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
33 {
34   const char    *runTypes[2] = {"full", "test"};
35   PetscInt       run;
36   PetscErrorCode ierr;
37 
38   PetscFunctionBeginUser;
39   options->runType        = RUN_FULL;
40   options->useDualPenalty = PETSC_FALSE;
41 
42   ierr = PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");CHKERRQ(ierr);
43   run  = options->runType;
44   ierr = PetscOptionsEList("-run_type", "The run type", "ex2.c", runTypes, 2, runTypes[options->runType], &run, NULL);CHKERRQ(ierr);
45   options->runType = (RunType) run;
46   ierr = PetscOptionsBool("-use_dual_penalty", "Penalize deviation from both goals", "ex2.c", options->useDualPenalty, &options->useDualPenalty, NULL);CHKERRQ(ierr);
47   ierr = PetscOptionsEnd();CHKERRQ(ierr);
48   PetscFunctionReturn(0);
49 }
50 
51 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
52 {
53   DM             distributedMesh = NULL;
54   PetscErrorCode ierr;
55 
56   PetscFunctionBeginUser;
57   ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
58   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
59   ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr);
60   if (distributedMesh) {
61     ierr = DMDestroy(dm);CHKERRQ(ierr);
62     *dm  = distributedMesh;
63   }
64   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
65   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
66   PetscFunctionReturn(0);
67 }
68 
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 }
117 
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 }
139 
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 }
170 
171 /*
172   In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:
173 
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])
178 
179   so that
180 
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 }
198 
199 PetscErrorCode SetupProblem(DM dm, AppCtx *user)
200 {
201   PetscDS        prob;
202   const PetscInt id = 1;
203   PetscErrorCode ierr;
204 
205   PetscFunctionBeginUser;
206   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
207   ierr = PetscDSSetResidual(prob, 0, user->useDualPenalty == PETSC_TRUE ? f0_u_full : f0_u, f1_u);CHKERRQ(ierr);
208   ierr = PetscDSSetResidual(prob, 1, f0_a, NULL);CHKERRQ(ierr);
209   ierr = PetscDSSetResidual(prob, 2, f0_l, f1_l);CHKERRQ(ierr);
210   ierr = PetscDSSetJacobian(prob, 0, 0, user->useDualPenalty == PETSC_TRUE ? g0_uu_full : g0_uu, NULL, NULL, NULL);CHKERRQ(ierr);
211   ierr = PetscDSSetJacobian(prob, 0, 2, NULL, NULL, NULL, g3_ul);CHKERRQ(ierr);
212   ierr = PetscDSSetJacobian(prob, 1, 1, g0_aa, NULL, NULL, NULL);CHKERRQ(ierr);
213   ierr = PetscDSSetJacobian(prob, 1, 2, g0_al, NULL, NULL, NULL);CHKERRQ(ierr);
214   ierr = PetscDSSetJacobian(prob, 2, 1, g0_la, NULL, NULL, NULL);CHKERRQ(ierr);
215   ierr = PetscDSSetJacobian(prob, 2, 0, NULL, NULL, NULL, g3_lu);CHKERRQ(ierr);
216 
217   ierr = PetscDSSetExactSolution(prob, 0, quadratic_u_2d, NULL);CHKERRQ(ierr);
218   ierr = PetscDSSetExactSolution(prob, 1, constant_a_2d, NULL);CHKERRQ(ierr);
219   ierr = PetscDSSetExactSolution(prob, 2, zero, NULL);CHKERRQ(ierr);
220   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)()) quadratic_u_2d, NULL, 1, &id, user);CHKERRQ(ierr);
221   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 1, 0, NULL, (void (*)()) constant_a_2d, NULL, 1, &id, user);CHKERRQ(ierr);
222   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 2, 0, NULL, (void (*)()) zero, NULL, 1, &id, user);CHKERRQ(ierr);
223   PetscFunctionReturn(0);
224 }
225 
226 PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
227 {
228   DM              cdm = dm;
229   const PetscInt  dim = 2;
230   PetscFE         fe[3];
231   PetscInt        f;
232   MPI_Comm        comm;
233   PetscErrorCode  ierr;
234 
235   PetscFunctionBeginUser;
236   /* Create finite element */
237   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
238   ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]);CHKERRQ(ierr);
239   ierr = PetscObjectSetName((PetscObject) fe[0], "potential");CHKERRQ(ierr);
240   ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "charge_", -1, &fe[1]);CHKERRQ(ierr);
241   ierr = PetscObjectSetName((PetscObject) fe[1], "charge");CHKERRQ(ierr);
242   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
243   ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]);CHKERRQ(ierr);
244   ierr = PetscObjectSetName((PetscObject) fe[2], "multiplier");CHKERRQ(ierr);
245   ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr);
246   /* Set discretization and boundary conditions for each mesh */
247   for (f = 0; f < 3; ++f) {ierr = DMSetField(dm, f, NULL, (PetscObject) fe[f]);CHKERRQ(ierr);}
248   ierr = DMCreateDS(cdm);CHKERRQ(ierr);
249   ierr = SetupProblem(dm, user);CHKERRQ(ierr);
250   while (cdm) {
251     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
252     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
253   }
254   for (f = 0; f < 3; ++f) {ierr = PetscFEDestroy(&fe[f]);CHKERRQ(ierr);}
255   PetscFunctionReturn(0);
256 }
257 
258 int main(int argc, char **argv)
259 {
260   DM             dm;
261   SNES           snes;
262   Vec            u, r;
263   AppCtx         user;
264   PetscErrorCode ierr;
265 
266   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
267   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
268   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
269   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
270   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
271   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
272 
273   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
274   ierr = PetscObjectSetName((PetscObject) u, "solution");CHKERRQ(ierr);
275   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
276   ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr);
277   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
278 
279   ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr);
280   if (user.runType == RUN_FULL) {
281     PetscDS          ds;
282     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
283     PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
284     PetscReal        error;
285 
286     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
287     ierr = PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL);CHKERRQ(ierr);
288     ierr = PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL);CHKERRQ(ierr);
289     ierr = PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL);CHKERRQ(ierr);
290     initialGuess[0] = zero;
291     initialGuess[1] = zero;
292     initialGuess[2] = zero;
293     ierr = DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr);
294     ierr = VecViewFromOptions(u, NULL, "-initial_vec_view");CHKERRQ(ierr);
295     ierr = DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);CHKERRQ(ierr);
296     if (error < 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
297     else                 {ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", error);CHKERRQ(ierr);}
298     ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
299     ierr = DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);CHKERRQ(ierr);
300     if (error < 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
301     else                 {ierr = PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", error);CHKERRQ(ierr);}
302   }
303   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
304 
305   ierr = VecDestroy(&u);CHKERRQ(ierr);
306   ierr = VecDestroy(&r);CHKERRQ(ierr);
307   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
308   ierr = DMDestroy(&dm);CHKERRQ(ierr);
309   ierr = PetscFinalize();
310   return ierr;
311 }
312 
313 /*TEST
314 
315   build:
316     requires: !complex triangle
317 
318   test:
319     suffix: 0
320     args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -charge_petscspace_degree 1 -multiplier_petscspace_degree 1
321 
322   test:
323     suffix: 1
324     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 
326   test:
327     suffix: 2
328     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 
330 TEST*/
331