xref: /petsc/src/tao/tutorials/ex1.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
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 {
26   RUN_FULL,
27   RUN_TEST
28 } RunType;
29 
30 typedef struct {
31   RunType runType; /* Whether to run tests, or solve the full problem */
32 } AppCtx;
33 
34 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
35   const char *runTypes[2] = {"full", "test"};
36   PetscInt    run;
37 
38   PetscFunctionBeginUser;
39   options->runType = RUN_FULL;
40   PetscOptionsBegin(comm, "", "Inverse Problem Options", "DMPLEX");
41   run = options->runType;
42   PetscCall(PetscOptionsEList("-run_type", "The run type", "ex1.c", runTypes, 2, runTypes[options->runType], &run, NULL));
43   options->runType = (RunType)run;
44   PetscOptionsEnd();
45   PetscFunctionReturn(0);
46 }
47 
48 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
49   PetscFunctionBeginUser;
50   PetscCall(DMCreate(comm, dm));
51   PetscCall(DMSetType(*dm, DMPLEX));
52   PetscCall(DMSetFromOptions(*dm));
53   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
54   PetscFunctionReturn(0);
55 }
56 
57 /* u - (x^2 + y^2) */
58 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[]) {
59   f0[0] = u[0] - (x[0] * x[0] + x[1] * x[1]);
60 }
61 /* a \nabla\lambda */
62 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[]) {
63   PetscInt d;
64   for (d = 0; d < dim; ++d) f1[d] = u[1] * u_x[dim * 2 + d];
65 }
66 /* I */
67 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[]) {
68   g0[0] = 1.0;
69 }
70 /* \nabla */
71 void g2_ua(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 g2[]) {
72   PetscInt d;
73   for (d = 0; d < dim; ++d) g2[d] = u_x[dim * 2 + d];
74 }
75 /* a */
76 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[]) {
77   PetscInt d;
78   for (d = 0; d < dim; ++d) g3[d * dim + d] = u[1];
79 }
80 /* a - (x + y) */
81 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[]) {
82   f0[0] = u[1] - (x[0] + x[1]);
83 }
84 /* \lambda \nabla u */
85 void f1_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 f1[]) {
86   PetscInt d;
87   for (d = 0; d < dim; ++d) f1[d] = u[2] * u_x[d];
88 }
89 /* I */
90 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[]) {
91   g0[0] = 1.0;
92 }
93 /* 6 (x + y) */
94 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[]) {
95   f0[0] = 6.0 * (x[0] + x[1]);
96 }
97 /* a \nabla u */
98 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[]) {
99   PetscInt d;
100   for (d = 0; d < dim; ++d) f1[d] = u[1] * u_x[d];
101 }
102 /* \nabla u */
103 void g2_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 g2[]) {
104   PetscInt d;
105   for (d = 0; d < dim; ++d) g2[d] = u_x[d];
106 }
107 /* a */
108 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[]) {
109   PetscInt d;
110   for (d = 0; d < dim; ++d) g3[d * dim + d] = u[1];
111 }
112 
113 /*
114   In 2D for Dirichlet conditions with a variable coefficient, we use exact solution:
115 
116     u  = x^2 + y^2
117     f  = 6 (x + y)
118     kappa(a) = a = (x + y)
119 
120   so that
121 
122     -\div \kappa(a) \grad u + f = -6 (x + y) + 6 (x + y) = 0
123 */
124 PetscErrorCode quadratic_u_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) {
125   *u = x[0] * x[0] + x[1] * x[1];
126   return 0;
127 }
128 PetscErrorCode linear_a_2d(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *a, void *ctx) {
129   *a = x[0] + x[1];
130   return 0;
131 }
132 PetscErrorCode zero(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *l, void *ctx) {
133   *l = 0.0;
134   return 0;
135 }
136 
137 PetscErrorCode SetupProblem(DM dm, AppCtx *user) {
138   PetscDS        ds;
139   DMLabel        label;
140   const PetscInt id = 1;
141 
142   PetscFunctionBeginUser;
143   PetscCall(DMGetDS(dm, &ds));
144   PetscCall(PetscDSSetResidual(ds, 0, f0_u, f1_u));
145   PetscCall(PetscDSSetResidual(ds, 1, f0_a, f1_a));
146   PetscCall(PetscDSSetResidual(ds, 2, f0_l, f1_l));
147   PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, NULL, NULL, NULL));
148   PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_ua, NULL));
149   PetscCall(PetscDSSetJacobian(ds, 0, 2, NULL, NULL, NULL, g3_ul));
150   PetscCall(PetscDSSetJacobian(ds, 1, 1, g0_aa, NULL, NULL, NULL));
151   PetscCall(PetscDSSetJacobian(ds, 2, 1, NULL, NULL, g2_la, NULL));
152   PetscCall(PetscDSSetJacobian(ds, 2, 0, NULL, NULL, NULL, g3_lu));
153 
154   PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u_2d, NULL));
155   PetscCall(PetscDSSetExactSolution(ds, 1, linear_a_2d, NULL));
156   PetscCall(PetscDSSetExactSolution(ds, 2, zero, NULL));
157   PetscCall(DMGetLabel(dm, "marker", &label));
158   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))quadratic_u_2d, NULL, user, NULL));
159   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 1, 0, NULL, (void (*)(void))linear_a_2d, NULL, user, NULL));
160   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 2, 0, NULL, (void (*)(void))zero, NULL, user, NULL));
161   PetscFunctionReturn(0);
162 }
163 
164 PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) {
165   DM             cdm = dm;
166   const PetscInt dim = 2;
167   PetscFE        fe[3];
168   PetscInt       f;
169   MPI_Comm       comm;
170 
171   PetscFunctionBeginUser;
172   /* Create finite element */
173   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
174   PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "potential_", -1, &fe[0]));
175   PetscCall(PetscObjectSetName((PetscObject)fe[0], "potential"));
176   PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "conductivity_", -1, &fe[1]));
177   PetscCall(PetscObjectSetName((PetscObject)fe[1], "conductivity"));
178   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
179   PetscCall(PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, "multiplier_", -1, &fe[2]));
180   PetscCall(PetscObjectSetName((PetscObject)fe[2], "multiplier"));
181   PetscCall(PetscFECopyQuadrature(fe[0], fe[2]));
182   /* Set discretization and boundary conditions for each mesh */
183   for (f = 0; f < 3; ++f) PetscCall(DMSetField(dm, f, NULL, (PetscObject)fe[f]));
184   PetscCall(DMCreateDS(dm));
185   PetscCall(SetupProblem(dm, user));
186   while (cdm) {
187     PetscCall(DMCopyDisc(dm, cdm));
188     PetscCall(DMGetCoarseDM(cdm, &cdm));
189   }
190   for (f = 0; f < 3; ++f) PetscCall(PetscFEDestroy(&fe[f]));
191   PetscFunctionReturn(0);
192 }
193 
194 int main(int argc, char **argv) {
195   DM     dm;
196   SNES   snes;
197   Vec    u, r;
198   AppCtx user;
199 
200   PetscFunctionBeginUser;
201   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
202   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
203   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
204   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
205   PetscCall(SNESSetDM(snes, dm));
206   PetscCall(SetupDiscretization(dm, &user));
207 
208   PetscCall(DMCreateGlobalVector(dm, &u));
209   PetscCall(PetscObjectSetName((PetscObject)u, "solution"));
210   PetscCall(VecDuplicate(u, &r));
211   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
212   PetscCall(SNESSetFromOptions(snes));
213 
214   PetscCall(DMSNESCheckFromOptions(snes, u));
215   if (user.runType == RUN_FULL) {
216     PetscDS ds;
217     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
218     PetscErrorCode (*initialGuess[3])(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx);
219     PetscReal error;
220 
221     PetscCall(DMGetDS(dm, &ds));
222     PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], NULL));
223     PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], NULL));
224     PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], NULL));
225     initialGuess[0] = zero;
226     initialGuess[1] = zero;
227     initialGuess[2] = zero;
228     PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u));
229     PetscCall(VecViewFromOptions(u, NULL, "-initial_vec_view"));
230     PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error));
231     if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: < 1.0e-11\n"));
232     else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial L_2 Error: %g\n", (double)error));
233     PetscCall(SNESSolve(snes, NULL, u));
234     PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, NULL, u, &error));
235     if (error < 1.0e-11) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: < 1.0e-11\n"));
236     else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final L_2 Error: %g\n", (double)error));
237   }
238   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
239 
240   PetscCall(VecDestroy(&u));
241   PetscCall(VecDestroy(&r));
242   PetscCall(SNESDestroy(&snes));
243   PetscCall(DMDestroy(&dm));
244   PetscCall(PetscFinalize());
245   return 0;
246 }
247 
248 /*TEST
249 
250   build:
251     requires: !complex
252 
253   test:
254     suffix: 0
255     requires: triangle
256     args: -run_type test -dmsnes_check -potential_petscspace_degree 2 -conductivity_petscspace_degree 1 -multiplier_petscspace_degree 2
257 
258   test:
259     suffix: 1
260     requires: triangle
261     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
262 
263 TEST*/
264