xref: /petsc/src/tao/tutorials/ex3.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 static char help[] = "Reduced formulation of the mother problem of PDE-constrained optimisation.\n";
2 
3 /*F
4   We solve the mother problem
5 
6   min 1/2 || y - y_d ||^2_2 + \alpha/2 || u ||^2_2
7 
8   subject to
9 
10           - \laplace y = u          on \Omega
11                      y = 0          on \Gamma
12 
13   where u is in L^2 and y is in H^1_0.
14 
15   We formulate the reduced problem solely in terms of the control
16   by using the state equation to express y in terms of u, and then
17   apply LMVM/BLMVM to the resulting reduced problem.
18 
19   Mesh independence is achieved by configuring the Riesz map for the control
20   space.
21 
22   Example meshes where the Riesz map is crucial can be downloaded from the
23   http://people.maths.ox.ac.uk/~farrellp/meshes.tar.gz
24 
25   Contributed by: Patrick Farrell <patrick.farrell@maths.ox.ac.uk>
26 
27   Run with e.g.:
28   ./ex3 -laplace_ksp_type cg -laplace_pc_type hypre -mat_lmvm_ksp_type cg -mat_lmvm_pc_type gamg -laplace_ksp_monitor_true_residual -tao_monitor -petscspace_degree 1 -tao_converged_reason -tao_gatol 1.0e-9 -dm_view hdf5:solution.h5 -sol_view hdf5:solution.h5::append -use_riesz 1 -f $DATAFILESPATH/meshes/mesh-1.h5
29 
30   and visualise in paraview with ../../../../petsc_gen_xdmf.py solution.h5.
31 
32   Toggle the Riesz map (-use_riesz 0) to see the difference setting the Riesz maps makes.
33 
34   TODO: broken for parallel runs
35 F*/
36 
37 #include <petsc.h>
38 #include <petscfe.h>
39 #include <petscviewerhdf5.h>
40 
41 typedef struct {
42   DM  dm;
43   Mat mass;
44   Vec data;
45   Vec state;
46   Vec tmp1;
47   Vec tmp2;
48   Vec adjoint;
49   Mat laplace;
50   KSP ksp_laplace;
51   PetscInt  num_bc_dofs;
52   PetscInt* bc_indices;
53   PetscScalar* bc_values;
54   PetscBool use_riesz;
55 } AppCtx;
56 
57 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
58 {
59   PetscBool      flg;
60   char           filename[2048];
61   PetscErrorCode ierr;
62 
63   PetscFunctionBeginUser;
64   filename[0] = '\0';
65   user->use_riesz = PETSC_TRUE;
66 
67   ierr = PetscOptionsBegin(comm, "", "Poisson mother problem options", "DMPLEX");CHKERRQ(ierr);
68   CHKERRQ(PetscOptionsBool("-use_riesz", "Use the Riesz map to achieve mesh independence", "ex3.c", user->use_riesz, &user->use_riesz, NULL));
69   CHKERRQ(PetscOptionsString("-f", "filename to read", "ex3.c", filename, filename, sizeof(filename), &flg));
70   ierr = PetscOptionsEnd();CHKERRQ(ierr);
71 
72   if (!flg) {
73     CHKERRQ(DMCreate(comm, dm));
74     CHKERRQ(DMSetType(*dm, DMPLEX));
75   } else {
76     /* TODO Eliminate this in favor of DMLoad() in new code */
77 #if defined(PETSC_HAVE_HDF5)
78     const PetscInt vertices_per_cell = 3;
79     PetscViewer    viewer;
80     Vec            coordinates;
81     Vec            topology;
82     PetscInt       dim = 2, numCells;
83     PetscInt       numVertices;
84     PetscScalar*   coords;
85     PetscScalar*   topo_f;
86     PetscInt*      cells;
87     PetscInt       j;
88     DMLabel        label;
89 
90     /* Read in FEniCS HDF5 output */
91     CHKERRQ(PetscViewerHDF5Open(comm, filename, FILE_MODE_READ, &viewer));
92 
93     /* create Vecs to read in the data from H5 */
94     CHKERRQ(VecCreate(comm, &coordinates));
95     CHKERRQ(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
96     CHKERRQ(VecSetBlockSize(coordinates, dim));
97     CHKERRQ(VecCreate(comm, &topology));
98     CHKERRQ(PetscObjectSetName((PetscObject)topology, "topology"));
99     CHKERRQ(VecSetBlockSize(topology, vertices_per_cell));
100 
101     /* navigate to the right group */
102     CHKERRQ(PetscViewerHDF5PushGroup(viewer, "/Mesh/mesh"));
103 
104     /* Read the Vecs */
105     CHKERRQ(VecLoad(coordinates, viewer));
106     CHKERRQ(VecLoad(topology, viewer));
107 
108     /* do some ugly calculations */
109     CHKERRQ(VecGetSize(topology, &numCells));
110     numCells = numCells / vertices_per_cell;
111     CHKERRQ(VecGetSize(coordinates, &numVertices));
112     numVertices = numVertices / dim;
113 
114     CHKERRQ(VecGetArray(coordinates, &coords));
115     CHKERRQ(VecGetArray(topology, &topo_f));
116     /* and now we have to convert the double representation to integers to pass over, argh */
117     CHKERRQ(PetscMalloc1(numCells*vertices_per_cell, &cells));
118     for (j = 0; j < numCells*vertices_per_cell; j++) {
119       cells[j] = (PetscInt) topo_f[j];
120     }
121 
122     /* Now create the DM */
123     CHKERRQ(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, vertices_per_cell, PETSC_TRUE, cells, dim, coords, dm));
124     /* Check for flipped first cell */
125     {
126       PetscReal v0[3], J[9], invJ[9], detJ;
127 
128       CHKERRQ(DMPlexComputeCellGeometryFEM(*dm, 0, NULL, v0, J, invJ, &detJ));
129       if (detJ < 0) {
130         CHKERRQ(DMPlexOrientPoint(*dm, 0, -1));
131         CHKERRQ(DMPlexComputeCellGeometryFEM(*dm, 0, NULL, v0, J, invJ, &detJ));
132         PetscCheck(detJ >= 0,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Something is wrong");
133       }
134     }
135     CHKERRQ(DMPlexOrient(*dm));
136     CHKERRQ(DMCreateLabel(*dm, "marker"));
137     CHKERRQ(DMGetLabel(*dm, "marker", &label));
138     CHKERRQ(DMPlexMarkBoundaryFaces(*dm, 1, label));
139     CHKERRQ(DMPlexLabelComplete(*dm, label));
140 
141     CHKERRQ(PetscViewerDestroy(&viewer));
142     CHKERRQ(VecRestoreArray(coordinates, &coords));
143     CHKERRQ(VecRestoreArray(topology, &topo_f));
144     CHKERRQ(PetscFree(cells));
145     CHKERRQ(VecDestroy(&coordinates));
146     CHKERRQ(VecDestroy(&topology));
147 #else
148     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Reconfigure PETSc with --download-hdf5");
149 #endif
150   }
151   CHKERRQ(DMSetFromOptions(*dm));
152   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
153   PetscFunctionReturn(0);
154 }
155 
156 void mass_kernel(PetscInt dim, PetscInt Nf, PetscInt NfAux,
157            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
158            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
159            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
160 {
161   g0[0] = 1.0;
162 }
163 
164 void laplace_kernel(PetscInt dim, PetscInt Nf, PetscInt NfAux,
165            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
166            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
167            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
168 {
169   PetscInt d;
170   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
171 }
172 
173 /* data we seek to match */
174 PetscErrorCode data_kernel(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *y, void *ctx)
175 {
176   *y = 1.0/(2*PETSC_PI*PETSC_PI) * PetscSinReal(PETSC_PI*x[0]) * PetscSinReal(PETSC_PI*x[1]);
177   /* the associated control is sin(pi*x[0])*sin(pi*x[1]) */
178   return 0;
179 }
180 PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
181 {
182   *u = 0.0;
183   return 0;
184 }
185 
186 PetscErrorCode CreateCtx(DM dm, AppCtx* user)
187 {
188 
189   DM             dm_mass;
190   DM             dm_laplace;
191   PetscDS        prob_mass;
192   PetscDS        prob_laplace;
193   PetscFE        fe;
194   DMLabel        label;
195   PetscSection   section;
196   PetscInt       n, k, p, d;
197   PetscInt       dof, off;
198   IS             is;
199   const PetscInt* points;
200   const PetscInt dim = 2;
201   const PetscInt id  = 1;
202   PetscErrorCode (**wtf)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
203 
204   PetscFunctionBeginUser;
205 
206   /* make the data we seek to match */
207   CHKERRQ(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, PETSC_TRUE, NULL, 4, &fe));
208 
209   CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe));
210   CHKERRQ(DMCreateDS(dm));
211   CHKERRQ(DMCreateGlobalVector(dm, &user->data));
212 
213   /* ugh, this is hideous */
214   /* y_d = interpolate(Expression("sin(x[0]) + .."), V) */
215   CHKERRQ(PetscMalloc(1 * sizeof(void (*)(const PetscReal[], PetscScalar *, void *)), &wtf));
216   wtf[0] = data_kernel;
217   CHKERRQ(DMProjectFunction(dm, 0.0, wtf, NULL, INSERT_VALUES, user->data));
218   CHKERRQ(PetscFree(wtf));
219 
220   /* assemble(inner(u, v)*dx), almost */
221   CHKERRQ(DMClone(dm, &dm_mass));
222   CHKERRQ(DMCopyDisc(dm, dm_mass));
223   CHKERRQ(DMSetNumFields(dm_mass, 1));
224   CHKERRQ(DMPlexCopyCoordinates(dm, dm_mass)); /* why do I have to do this separately? */
225   CHKERRQ(DMGetDS(dm_mass, &prob_mass));
226   CHKERRQ(PetscDSSetJacobian(prob_mass, 0, 0, mass_kernel, NULL, NULL, NULL));
227   CHKERRQ(PetscDSSetDiscretization(prob_mass, 0, (PetscObject) fe));
228   CHKERRQ(DMCreateMatrix(dm_mass, &user->mass));
229   CHKERRQ(DMPlexSNESComputeJacobianFEM(dm_mass, user->data, user->mass, user->mass, NULL));
230   CHKERRQ(MatSetOption(user->mass, MAT_SYMMETRIC, PETSC_TRUE));
231   CHKERRQ(DMDestroy(&dm_mass));
232 
233   /* inner(grad(u), grad(v))*dx with homogeneous Dirichlet boundary conditions */
234   CHKERRQ(DMClone(dm, &dm_laplace));
235   CHKERRQ(DMCopyDisc(dm, dm_laplace));
236   CHKERRQ(DMSetNumFields(dm_laplace, 1));
237   CHKERRQ(DMPlexCopyCoordinates(dm, dm_laplace));
238   CHKERRQ(DMGetDS(dm_laplace, &prob_laplace));
239   CHKERRQ(PetscDSSetJacobian(prob_laplace, 0, 0, NULL, NULL, NULL, laplace_kernel));
240   CHKERRQ(PetscDSSetDiscretization(prob_laplace, 0, (PetscObject) fe));
241   CHKERRQ(DMCreateMatrix(dm_laplace, &user->laplace));
242   CHKERRQ(DMPlexSNESComputeJacobianFEM(dm_laplace, user->data, user->laplace, user->laplace, NULL));
243 
244   /* Code from Matt to get the indices associated with the boundary dofs */
245   CHKERRQ(DMGetLabel(dm_laplace, "marker", &label));
246   CHKERRQ(DMAddBoundary(dm_laplace, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) zero, NULL, NULL, NULL));
247   CHKERRQ(DMGetLocalSection(dm_laplace, &section));
248   CHKERRQ(DMLabelGetStratumSize(label, 1, &n));
249   CHKERRQ(DMLabelGetStratumIS(label, 1, &is));
250   CHKERRQ(ISGetIndices(is, &points));
251   user->num_bc_dofs = 0;
252   for (p = 0; p < n; ++p) {
253     CHKERRQ(PetscSectionGetDof(section, points[p], &dof));
254     user->num_bc_dofs += dof;
255   }
256   CHKERRQ(PetscMalloc1(user->num_bc_dofs, &user->bc_indices));
257   for (p = 0, k = 0; p < n; ++p) {
258     CHKERRQ(PetscSectionGetDof(section, points[p], &dof));
259     CHKERRQ(PetscSectionGetOffset(section, points[p], &off));
260     for (d = 0; d < dof; ++d) user->bc_indices[k++] = off+d;
261   }
262   CHKERRQ(ISRestoreIndices(is, &points));
263   CHKERRQ(ISDestroy(&is));
264   CHKERRQ(DMDestroy(&dm_laplace));
265 
266   /* This is how I handle boundary conditions. I can't figure out how to get
267      plex to play with the way I want to impose the BCs. This loses symmetry,
268      but not in a disastrous way. If someone can improve it, please do! */
269   CHKERRQ(MatZeroRows(user->laplace, user->num_bc_dofs, user->bc_indices, 1.0, NULL, NULL));
270   CHKERRQ(PetscCalloc1(user->num_bc_dofs, &user->bc_values));
271 
272   /* also create the KSP for solving the Laplace system */
273   CHKERRQ(KSPCreate(PETSC_COMM_WORLD, &user->ksp_laplace));
274   CHKERRQ(KSPSetOperators(user->ksp_laplace, user->laplace, user->laplace));
275   CHKERRQ(KSPSetOptionsPrefix(user->ksp_laplace, "laplace_"));
276   CHKERRQ(KSPSetFromOptions(user->ksp_laplace));
277 
278   /* A bit of setting up the user context */
279   user->dm = dm;
280   CHKERRQ(VecDuplicate(user->data, &user->state));
281   CHKERRQ(VecDuplicate(user->data, &user->adjoint));
282   CHKERRQ(VecDuplicate(user->data, &user->tmp1));
283   CHKERRQ(VecDuplicate(user->data, &user->tmp2));
284 
285   CHKERRQ(PetscFEDestroy(&fe));
286 
287   PetscFunctionReturn(0);
288 }
289 
290 PetscErrorCode DestroyCtx(AppCtx* user)
291 {
292   PetscFunctionBeginUser;
293 
294   CHKERRQ(MatDestroy(&user->mass));
295   CHKERRQ(MatDestroy(&user->laplace));
296   CHKERRQ(KSPDestroy(&user->ksp_laplace));
297   CHKERRQ(VecDestroy(&user->data));
298   CHKERRQ(VecDestroy(&user->state));
299   CHKERRQ(VecDestroy(&user->adjoint));
300   CHKERRQ(VecDestroy(&user->tmp1));
301   CHKERRQ(VecDestroy(&user->tmp2));
302   CHKERRQ(PetscFree(user->bc_indices));
303   CHKERRQ(PetscFree(user->bc_values));
304 
305   PetscFunctionReturn(0);
306 }
307 
308 PetscErrorCode ReducedFunctionGradient(Tao tao, Vec u, PetscReal* func, Vec g, void* userv)
309 {
310   AppCtx* user = (AppCtx*) userv;
311   const PetscReal alpha = 1.0e-6; /* regularisation parameter */
312   PetscReal inner;
313 
314   PetscFunctionBeginUser;
315 
316   CHKERRQ(MatMult(user->mass, u, user->tmp1));
317   CHKERRQ(VecDot(u, user->tmp1, &inner));               /* regularisation contribution to */
318   *func = alpha * 0.5 * inner;                                      /* the functional                 */
319 
320   CHKERRQ(VecSet(g, 0.0));
321   CHKERRQ(VecAXPY(g, alpha, user->tmp1));               /* regularisation contribution to the gradient */
322 
323   /* Now compute the forward state. */
324   CHKERRQ(VecSetValues(user->tmp1, user->num_bc_dofs, user->bc_indices, user->bc_values, INSERT_VALUES));
325   CHKERRQ(VecAssemblyBegin(user->tmp1));
326   CHKERRQ(VecAssemblyEnd(user->tmp1));
327   CHKERRQ(KSPSolve(user->ksp_laplace, user->tmp1, user->state)); /* forward solve */
328 
329   /* Now compute the adjoint state also. */
330   CHKERRQ(VecCopy(user->state, user->tmp1));
331   CHKERRQ(VecAXPY(user->tmp1, -1.0, user->data));
332   CHKERRQ(MatMult(user->mass, user->tmp1, user->tmp2));
333   CHKERRQ(VecDot(user->tmp1, user->tmp2, &inner));      /* misfit contribution to */
334   *func += 0.5 * inner;                                             /* the functional         */
335 
336   CHKERRQ(VecSetValues(user->tmp2, user->num_bc_dofs, user->bc_indices, user->bc_values, INSERT_VALUES));
337   CHKERRQ(VecAssemblyBegin(user->tmp2));
338   CHKERRQ(VecAssemblyEnd(user->tmp2));
339   CHKERRQ(KSPSolve(user->ksp_laplace, user->tmp2, user->adjoint)); /* adjoint solve */
340 
341   /* And bring it home with the gradient. */
342   CHKERRQ(MatMult(user->mass, user->adjoint, user->tmp1));
343   CHKERRQ(VecAXPY(g, 1.0, user->tmp1));                 /* adjoint contribution to the gradient */
344 
345   PetscFunctionReturn(0);
346 }
347 
348 int main(int argc, char **argv)
349 {
350   DM             dm;
351   Tao            tao;
352   Vec            u, lb, ub;
353   AppCtx         user;
354   PetscErrorCode ierr;
355 
356   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
357   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
358   CHKERRQ(CreateCtx(dm, &user));
359 
360   CHKERRQ(DMCreateGlobalVector(dm, &u));
361   CHKERRQ(VecSet(u, 0.0));
362   CHKERRQ(VecDuplicate(u, &lb));
363   CHKERRQ(VecDuplicate(u, &ub));
364   CHKERRQ(VecSet(lb, 0.0)); /* satisfied at the minimum anyway */
365   CHKERRQ(VecSet(ub, 0.8)); /* a nontrivial upper bound */
366 
367   CHKERRQ(TaoCreate(PETSC_COMM_WORLD, &tao));
368   CHKERRQ(TaoSetSolution(tao, u));
369   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL, ReducedFunctionGradient, &user));
370   CHKERRQ(TaoSetVariableBounds(tao, lb, ub));
371   CHKERRQ(TaoSetType(tao, TAOBLMVM));
372   CHKERRQ(TaoSetFromOptions(tao));
373 
374   if (user.use_riesz) {
375     CHKERRQ(TaoLMVMSetH0(tao, user.mass));       /* crucial for mesh independence */
376     CHKERRQ(TaoSetGradientNorm(tao, user.mass));
377   }
378 
379   CHKERRQ(TaoSolve(tao));
380 
381   CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view"));
382   CHKERRQ(VecViewFromOptions(u, NULL, "-sol_view"));
383 
384   CHKERRQ(TaoDestroy(&tao));
385   CHKERRQ(DMDestroy(&dm));
386   CHKERRQ(VecDestroy(&u));
387   CHKERRQ(VecDestroy(&lb));
388   CHKERRQ(VecDestroy(&ub));
389   CHKERRQ(DestroyCtx(&user));
390   ierr = PetscFinalize();
391   return ierr;
392 }
393 
394 /*TEST
395 
396     build:
397       requires: !complex !single
398 
399     test:
400       requires: hdf5 double datafilespath !defined(PETSC_USE_64BIT_INDICES) hypre
401       args: -laplace_ksp_type cg -laplace_pc_type hypre -laplace_ksp_monitor_true_residual -laplace_ksp_max_it 5 -mat_lmvm_ksp_type cg -mat_lmvm_ksp_rtol 1e-5 -mat_lmvm_pc_type gamg -tao_monitor -petscspace_degree 1 -tao_converged_reason -tao_gatol 1.0e-6 -dm_view hdf5:solution.h5 -sol_view hdf5:solution.h5::append -use_riesz 1 -f $DATAFILESPATH/meshes/mesh-1.h5
402       filter: sed -e "s/-nan/nan/g"
403 
404     test:
405       suffix: guess_pod
406       requires: double triangle
407       args: -laplace_ksp_type cg -laplace_pc_type gamg -laplace_ksp_monitor_true_residual -laplace_ksp_max_it 8 -laplace_pc_gamg_esteig_ksp_type cg -laplace_pc_gamg_esteig_ksp_max_it 5 -laplace_mg_levels_ksp_chebyshev_esteig 0,0.25,0,1.0 -laplace_ksp_converged_reason -mat_lmvm_ksp_type cg -mat_lmvm_ksp_rtol 1e-5 -mat_lmvm_pc_type gamg -mat_lmvm_pc_gamg_esteig_ksp_type cg -mat_lmvm_pc_gamg_esteig_ksp_max_it 3 -tao_monitor -petscspace_degree 1 -tao_converged_reason -dm_refine 0 -laplace_ksp_guess_type pod -tao_gatol 1e-6
408       filter: sed -e "s/-nan/nan/g" -e "s/-NaN/nan/g" -e "s/NaN/nan/g" -e "s/CONVERGED_RTOL iterations 9/CONVERGED_RTOL iterations 8/g" -e "s/CONVERGED_RTOL iterations 4/CONVERGED_RTOL iterations 3/g"
409 
410 TEST*/
411