xref: /petsc/src/tao/tutorials/ex1.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
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 $a$ and $u$, given by
7 \begin{align}
8   L(u, a, \lambda) = \frac{1}{2} || Qu - d ||^2 + \frac{1}{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 exact control for the reference $a_r$.
14 
15 The PDE will be the Laplace equation with homogeneous boundary conditions
16 \begin{align}
17   -nabla \cdot a \nabla u = f
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 } AppCtx;
30 
31 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
32 {
33   const char    *runTypes[2] = {"full", "test"};
34   PetscInt       run;
35   PetscErrorCode ierr;
36 
37   PetscFunctionBeginUser;
38   options->runType = RUN_FULL;
39 
40   ierr = PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");CHKERRQ(ierr);
41   run  = options->runType;
42   ierr = PetscOptionsEList("-run_type", "The run type", "ex1.c", runTypes, 2, runTypes[options->runType], &run, NULL);CHKERRQ(ierr);
43   options->runType = (RunType) run;
44   ierr = PetscOptionsEnd();CHKERRQ(ierr);
45   PetscFunctionReturn(0);
46 }
47 
48 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
49 {
50   PetscErrorCode ierr;
51 
52   PetscFunctionBeginUser;
53   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
54   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
55   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
56   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
57   PetscFunctionReturn(0);
58 }
59 
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 }
163 
164 /*
165   In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:
166 
167     u  = x^2 + y^2
168     f  = 6 (x + y)
169     kappa(a) = a = (x + y)
170 
171   so that
172 
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 }
190 
191 PetscErrorCode SetupProblem(DM dm, AppCtx *user)
192 {
193   PetscDS        ds;
194   DMLabel        label;
195   const PetscInt id = 1;
196   PetscErrorCode ierr;
197 
198   PetscFunctionBeginUser;
199   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
200   ierr = PetscDSSetResidual(ds, 0, f0_u, f1_u);CHKERRQ(ierr);
201   ierr = PetscDSSetResidual(ds, 1, f0_a, f1_a);CHKERRQ(ierr);
202   ierr = PetscDSSetResidual(ds, 2, f0_l, f1_l);CHKERRQ(ierr);
203   ierr = PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, NULL);CHKERRQ(ierr);
204   ierr = PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_ua, NULL);CHKERRQ(ierr);
205   ierr = PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul);CHKERRQ(ierr);
206   ierr = PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL);CHKERRQ(ierr);
207   ierr = PetscDSSetJacobian(ds, 2, 1, NULL, NULL, g2_la, NULL);CHKERRQ(ierr);
208   ierr = PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu);CHKERRQ(ierr);
209 
210   ierr = PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL);CHKERRQ(ierr);
211   ierr = PetscDSSetExactSolution(ds, 1, linear_a_2d, NULL);CHKERRQ(ierr);
212   ierr = PetscDSSetExactSolution(ds, 2, zero, NULL);CHKERRQ(ierr);
213   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
214   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) quadratic_u_2d, NULL, user, NULL);CHKERRQ(ierr);
215   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)(void)) linear_a_2d, NULL, user, NULL);CHKERRQ(ierr);
216   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)(void)) zero, NULL, user, NULL);CHKERRQ(ierr);
217   PetscFunctionReturn(0);
218 }
219 
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;
228 
229   PetscFunctionBeginUser;
230   /* Create finite element */
231   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
232   ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]);CHKERRQ(ierr);
233   ierr = PetscObjectSetName((PetscObject) fe[0], "potential");CHKERRQ(ierr);
234   ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "conductivity_", -1, &fe[1]);CHKERRQ(ierr);
235   ierr = PetscObjectSetName((PetscObject) fe[1], "conductivity");CHKERRQ(ierr);
236   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
237   ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]);CHKERRQ(ierr);
238   ierr = PetscObjectSetName((PetscObject) fe[2], "multiplier");CHKERRQ(ierr);
239   ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr);
240   /* Set discretization and boundary conditions for each mesh */
241   for (f = 0; f < 3; ++f) {ierr = DMSetField(dm, f, NULL, (PetscObject) fe[f]);CHKERRQ(ierr);}
242   ierr = DMCreateDS(dm);CHKERRQ(ierr);
243   ierr = SetupProblem(dm, user);CHKERRQ(ierr);
244   while (cdm) {
245     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
246     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
247   }
248   for (f = 0; f < 3; ++f) {ierr = PetscFEDestroy(&fe[f]);CHKERRQ(ierr);}
249   PetscFunctionReturn(0);
250 }
251 
252 int main(int argc, char **argv)
253 {
254   DM             dm;
255   SNES           snes;
256   Vec            u, r;
257   AppCtx         user;
258   PetscErrorCode ierr;
259 
260   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
261   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
262   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
263   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
264   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
265   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
266 
267   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
268   ierr = PetscObjectSetName((PetscObject) u, "solution");CHKERRQ(ierr);
269   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
270   ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr);
271   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
272 
273   ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr);
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;
279 
280     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
281     ierr = PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL);CHKERRQ(ierr);
282     ierr = PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL);CHKERRQ(ierr);
283     ierr = PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL);CHKERRQ(ierr);
284     initialGuess[0] = zero;
285     initialGuess[1] = zero;
286     initialGuess[2] = zero;
287     ierr = DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr);
288     ierr = VecViewFromOptions(u, NULL, "-initial_vec_view");CHKERRQ(ierr);
289     ierr = DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);CHKERRQ(ierr);
290     if (error < 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
291     else                 {ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", error);CHKERRQ(ierr);}
292     ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
293     ierr = DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error);CHKERRQ(ierr);
294     if (error < 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
295     else                 {ierr = PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", error);CHKERRQ(ierr);}
296   }
297   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
298 
299   ierr = VecDestroy(&u);CHKERRQ(ierr);
300   ierr = VecDestroy(&r);CHKERRQ(ierr);
301   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
302   ierr = DMDestroy(&dm);CHKERRQ(ierr);
303   ierr = PetscFinalize();
304   return ierr;
305 }
306 
307 /*TEST
308 
309   build:
310     requires: !complex
311 
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
316 
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
321 
322 TEST*/
323