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