xref: /petsc/src/tao/tutorials/ex3.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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.:
2809fb3ddcSToby Isaac   ./ex3 -laplace_ksp_type cg -laplace_pc_type hypre -tao_blmvm_mat_lmvm_J0_ksp_type cg -tao_blmvm_mat_lmvm_J0_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 
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)57d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
58d71ae5a4SJacob Faibussowitsch {
59c4762a1bSJed Brown   PetscBool flg;
6030602db0SMatthew G. Knepley   char      filename[2048];
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   PetscFunctionBeginUser;
63c4762a1bSJed Brown   filename[0]     = '\0';
64c4762a1bSJed Brown   user->use_riesz = PETSC_TRUE;
65c4762a1bSJed Brown 
66d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Poisson mother problem options", "DMPLEX");
679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-use_riesz", "Use the Riesz map to achieve mesh independence", "ex3.c", user->use_riesz, &user->use_riesz, NULL));
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-f", "filename to read", "ex3.c", filename, filename, sizeof(filename), &flg));
69d0609cedSBarry Smith   PetscOptionsEnd();
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   if (!flg) {
729566063dSJacob Faibussowitsch     PetscCall(DMCreate(comm, dm));
739566063dSJacob Faibussowitsch     PetscCall(DMSetType(*dm, DMPLEX));
74c4762a1bSJed Brown   } else {
7530602db0SMatthew G. Knepley     /* TODO Eliminate this in favor of DMLoad() in new code */
76c4762a1bSJed Brown #if defined(PETSC_HAVE_HDF5)
77c4762a1bSJed Brown     const PetscInt vertices_per_cell = 3;
78c4762a1bSJed Brown     PetscViewer    viewer;
79c4762a1bSJed Brown     Vec            coordinates;
80c4762a1bSJed Brown     Vec            topology;
8130602db0SMatthew G. Knepley     PetscInt       dim = 2, numCells;
82c4762a1bSJed Brown     PetscInt       numVertices;
83c4762a1bSJed Brown     PetscScalar   *coords;
84c4762a1bSJed Brown     PetscScalar   *topo_f;
85c4762a1bSJed Brown     PetscInt      *cells;
86c4762a1bSJed Brown     PetscInt       j;
87c4762a1bSJed Brown     DMLabel        label;
88c4762a1bSJed Brown 
89c4762a1bSJed Brown     /* Read in FEniCS HDF5 output */
909566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Open(comm, filename, FILE_MODE_READ, &viewer));
91c4762a1bSJed Brown 
92c4762a1bSJed Brown     /* create Vecs to read in the data from H5 */
939566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm, &coordinates));
949566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
959566063dSJacob Faibussowitsch     PetscCall(VecSetBlockSize(coordinates, dim));
969566063dSJacob Faibussowitsch     PetscCall(VecCreate(comm, &topology));
979566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)topology, "topology"));
989566063dSJacob Faibussowitsch     PetscCall(VecSetBlockSize(topology, vertices_per_cell));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown     /* navigate to the right group */
1019566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5PushGroup(viewer, "/Mesh/mesh"));
102c4762a1bSJed Brown 
103c4762a1bSJed Brown     /* Read the Vecs */
1049566063dSJacob Faibussowitsch     PetscCall(VecLoad(coordinates, viewer));
1059566063dSJacob Faibussowitsch     PetscCall(VecLoad(topology, viewer));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown     /* do some ugly calculations */
1089566063dSJacob Faibussowitsch     PetscCall(VecGetSize(topology, &numCells));
109c4762a1bSJed Brown     numCells = numCells / vertices_per_cell;
1109566063dSJacob Faibussowitsch     PetscCall(VecGetSize(coordinates, &numVertices));
111c4762a1bSJed Brown     numVertices = numVertices / dim;
112c4762a1bSJed Brown 
1139566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &coords));
1149566063dSJacob Faibussowitsch     PetscCall(VecGetArray(topology, &topo_f));
115c4762a1bSJed Brown     /* and now we have to convert the double representation to integers to pass over, argh */
1169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numCells * vertices_per_cell, &cells));
117ad540459SPierre Jolivet     for (j = 0; j < numCells * vertices_per_cell; j++) cells[j] = (PetscInt)topo_f[j];
118c4762a1bSJed Brown 
119c4762a1bSJed Brown     /* Now create the DM */
1209566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, vertices_per_cell, PETSC_TRUE, cells, dim, coords, dm));
121c4762a1bSJed Brown     /* Check for flipped first cell */
122c4762a1bSJed Brown     {
123c4762a1bSJed Brown       PetscReal v0[3], J[9], invJ[9], detJ;
124c4762a1bSJed Brown 
1259566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFEM(*dm, 0, NULL, v0, J, invJ, &detJ));
126c4762a1bSJed Brown       if (detJ < 0) {
1279566063dSJacob Faibussowitsch         PetscCall(DMPlexOrientPoint(*dm, 0, -1));
1289566063dSJacob Faibussowitsch         PetscCall(DMPlexComputeCellGeometryFEM(*dm, 0, NULL, v0, J, invJ, &detJ));
1293c859ba3SBarry Smith         PetscCheck(detJ >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Something is wrong");
130c4762a1bSJed Brown       }
131c4762a1bSJed Brown     }
1329566063dSJacob Faibussowitsch     PetscCall(DMPlexOrient(*dm));
1339566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(*dm, "marker"));
1349566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "marker", &label));
1359566063dSJacob Faibussowitsch     PetscCall(DMPlexMarkBoundaryFaces(*dm, 1, label));
1369566063dSJacob Faibussowitsch     PetscCall(DMPlexLabelComplete(*dm, label));
137c4762a1bSJed Brown 
1389566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
1399566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &coords));
1409566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(topology, &topo_f));
1419566063dSJacob Faibussowitsch     PetscCall(PetscFree(cells));
1429566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&coordinates));
1439566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&topology));
144c4762a1bSJed Brown #else
145c4762a1bSJed Brown     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Reconfigure PETSc with --download-hdf5");
146c4762a1bSJed Brown #endif
147c4762a1bSJed Brown   }
1489566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
1499566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
1503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
mass_kernel(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[])153d71ae5a4SJacob Faibussowitsch void mass_kernel(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[])
154d71ae5a4SJacob Faibussowitsch {
155c4762a1bSJed Brown   g0[0] = 1.0;
156c4762a1bSJed Brown }
157c4762a1bSJed Brown 
laplace_kernel(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[])158d71ae5a4SJacob Faibussowitsch void laplace_kernel(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[])
159d71ae5a4SJacob Faibussowitsch {
160c4762a1bSJed Brown   PetscInt d;
161c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
162c4762a1bSJed Brown }
163c4762a1bSJed Brown 
164c4762a1bSJed Brown /* data we seek to match */
data_kernel(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf,PetscScalar * y,PetscCtx ctx)165*2a8381b2SBarry Smith PetscErrorCode data_kernel(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *y, PetscCtx ctx)
166d71ae5a4SJacob Faibussowitsch {
167c4762a1bSJed Brown   *y = 1.0 / (2 * PETSC_PI * PETSC_PI) * PetscSinReal(PETSC_PI * x[0]) * PetscSinReal(PETSC_PI * x[1]);
168c4762a1bSJed Brown   /* the associated control is sin(pi*x[0])*sin(pi*x[1]) */
1693ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
170c4762a1bSJed Brown }
zero(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)171*2a8381b2SBarry Smith PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
172d71ae5a4SJacob Faibussowitsch {
173c4762a1bSJed Brown   *u = 0.0;
1743ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
175c4762a1bSJed Brown }
176c4762a1bSJed Brown 
CreateCtx(DM dm,AppCtx * user)177d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateCtx(DM dm, AppCtx *user)
178d71ae5a4SJacob Faibussowitsch {
179c4762a1bSJed Brown   DM              dm_mass;
180c4762a1bSJed Brown   DM              dm_laplace;
181c4762a1bSJed Brown   PetscDS         prob_mass;
182c4762a1bSJed Brown   PetscDS         prob_laplace;
183c4762a1bSJed Brown   PetscFE         fe;
184c4762a1bSJed Brown   DMLabel         label;
185c4762a1bSJed Brown   PetscSection    section;
186c4762a1bSJed Brown   PetscInt        n, k, p, d;
187c4762a1bSJed Brown   PetscInt        dof, off;
188c4762a1bSJed Brown   IS              is;
189c4762a1bSJed Brown   const PetscInt *points;
190c4762a1bSJed Brown   const PetscInt  dim = 2;
191*2a8381b2SBarry Smith   PetscErrorCode (**wtf)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx);
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   PetscFunctionBeginUser;
194c4762a1bSJed Brown   /* make the data we seek to match */
1959566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, PETSC_TRUE, NULL, 4, &fe));
196c4762a1bSJed Brown 
1979566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
1989566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
1999566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &user->data));
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /* ugh, this is hideous */
202c4762a1bSJed Brown   /* y_d = interpolate(Expression("sin(x[0]) + .."), V) */
2039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc(1 * sizeof(void (*)(const PetscReal[], PetscScalar *, void *)), &wtf));
204c4762a1bSJed Brown   wtf[0] = data_kernel;
2059566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, wtf, NULL, INSERT_VALUES, user->data));
2069566063dSJacob Faibussowitsch   PetscCall(PetscFree(wtf));
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   /* assemble(inner(u, v)*dx), almost */
2099566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dm_mass));
2109566063dSJacob Faibussowitsch   PetscCall(DMCopyDisc(dm, dm_mass));
2119566063dSJacob Faibussowitsch   PetscCall(DMSetNumFields(dm_mass, 1));
2129566063dSJacob Faibussowitsch   PetscCall(DMPlexCopyCoordinates(dm, dm_mass)); /* why do I have to do this separately? */
2139566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm_mass, &prob_mass));
2149566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob_mass, 0, 0, mass_kernel, NULL, NULL, NULL));
2159566063dSJacob Faibussowitsch   PetscCall(PetscDSSetDiscretization(prob_mass, 0, (PetscObject)fe));
2169566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm_mass, &user->mass));
2179566063dSJacob Faibussowitsch   PetscCall(DMPlexSNESComputeJacobianFEM(dm_mass, user->data, user->mass, user->mass, NULL));
2189566063dSJacob Faibussowitsch   PetscCall(MatSetOption(user->mass, MAT_SYMMETRIC, PETSC_TRUE));
2199566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm_mass));
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   /* inner(grad(u), grad(v))*dx with homogeneous Dirichlet boundary conditions */
2229566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dm_laplace));
2239566063dSJacob Faibussowitsch   PetscCall(DMCopyDisc(dm, dm_laplace));
2249566063dSJacob Faibussowitsch   PetscCall(DMSetNumFields(dm_laplace, 1));
2259566063dSJacob Faibussowitsch   PetscCall(DMPlexCopyCoordinates(dm, dm_laplace));
2269566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm_laplace, &prob_laplace));
2279566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob_laplace, 0, 0, NULL, NULL, NULL, laplace_kernel));
2289566063dSJacob Faibussowitsch   PetscCall(PetscDSSetDiscretization(prob_laplace, 0, (PetscObject)fe));
2299566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm_laplace, &user->laplace));
2309566063dSJacob Faibussowitsch   PetscCall(DMPlexSNESComputeJacobianFEM(dm_laplace, user->data, user->laplace, user->laplace, NULL));
231c4762a1bSJed Brown 
2329566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm_laplace, "marker", &label));
2339566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm_laplace, &section));
2349566063dSJacob Faibussowitsch   PetscCall(DMLabelGetStratumSize(label, 1, &n));
2359566063dSJacob Faibussowitsch   PetscCall(DMLabelGetStratumIS(label, 1, &is));
2369566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &points));
237c4762a1bSJed Brown   user->num_bc_dofs = 0;
238c4762a1bSJed Brown   for (p = 0; p < n; ++p) {
2399566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, points[p], &dof));
240c4762a1bSJed Brown     user->num_bc_dofs += dof;
241c4762a1bSJed Brown   }
2429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->num_bc_dofs, &user->bc_indices));
243c4762a1bSJed Brown   for (p = 0, k = 0; p < n; ++p) {
2449566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, points[p], &dof));
2459566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, points[p], &off));
246c4762a1bSJed Brown     for (d = 0; d < dof; ++d) user->bc_indices[k++] = off + d;
247c4762a1bSJed Brown   }
2489566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &points));
2499566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
2509566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm_laplace));
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /* This is how I handle boundary conditions. I can't figure out how to get
253c4762a1bSJed Brown      plex to play with the way I want to impose the BCs. This loses symmetry,
254c4762a1bSJed Brown      but not in a disastrous way. If someone can improve it, please do! */
2559566063dSJacob Faibussowitsch   PetscCall(MatZeroRows(user->laplace, user->num_bc_dofs, user->bc_indices, 1.0, NULL, NULL));
2569566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(user->num_bc_dofs, &user->bc_values));
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   /* also create the KSP for solving the Laplace system */
2599566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD, &user->ksp_laplace));
2609566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(user->ksp_laplace, user->laplace, user->laplace));
2619566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(user->ksp_laplace, "laplace_"));
2629566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(user->ksp_laplace));
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   /* A bit of setting up the user context */
265c4762a1bSJed Brown   user->dm = dm;
2669566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->data, &user->state));
2679566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->data, &user->adjoint));
2689566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->data, &user->tmp1));
2699566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->data, &user->tmp2));
270c4762a1bSJed Brown 
2719566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
2723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
273c4762a1bSJed Brown }
274c4762a1bSJed Brown 
DestroyCtx(AppCtx * user)275d71ae5a4SJacob Faibussowitsch PetscErrorCode DestroyCtx(AppCtx *user)
276d71ae5a4SJacob Faibussowitsch {
277c4762a1bSJed Brown   PetscFunctionBeginUser;
2789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->mass));
2799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->laplace));
2809566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&user->ksp_laplace));
2819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->data));
2829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->state));
2839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->adjoint));
2849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->tmp1));
2859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->tmp2));
2869566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->bc_indices));
2879566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->bc_values));
2883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
289c4762a1bSJed Brown }
290c4762a1bSJed Brown 
ReducedFunctionGradient(Tao tao,Vec u,PetscReal * func,Vec g,void * userv)291d71ae5a4SJacob Faibussowitsch PetscErrorCode ReducedFunctionGradient(Tao tao, Vec u, PetscReal *func, Vec g, void *userv)
292d71ae5a4SJacob Faibussowitsch {
293c4762a1bSJed Brown   AppCtx         *user  = (AppCtx *)userv;
294c4762a1bSJed Brown   const PetscReal alpha = 1.0e-6; /* regularisation parameter */
295c4762a1bSJed Brown   PetscReal       inner;
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   PetscFunctionBeginUser;
2989566063dSJacob Faibussowitsch   PetscCall(MatMult(user->mass, u, user->tmp1));
2999566063dSJacob Faibussowitsch   PetscCall(VecDot(u, user->tmp1, &inner)); /* regularisation contribution to */
300c4762a1bSJed Brown   *func = alpha * 0.5 * inner;              /* the functional                 */
301c4762a1bSJed Brown 
3029566063dSJacob Faibussowitsch   PetscCall(VecSet(g, 0.0));
3039566063dSJacob Faibussowitsch   PetscCall(VecAXPY(g, alpha, user->tmp1)); /* regularisation contribution to the gradient */
304c4762a1bSJed Brown 
305c4762a1bSJed Brown   /* Now compute the forward state. */
3069566063dSJacob Faibussowitsch   PetscCall(VecSetValues(user->tmp1, user->num_bc_dofs, user->bc_indices, user->bc_values, INSERT_VALUES));
3079566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(user->tmp1));
3089566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(user->tmp1));
3099566063dSJacob Faibussowitsch   PetscCall(KSPSolve(user->ksp_laplace, user->tmp1, user->state)); /* forward solve */
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   /* Now compute the adjoint state also. */
3129566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->state, user->tmp1));
3139566063dSJacob Faibussowitsch   PetscCall(VecAXPY(user->tmp1, -1.0, user->data));
3149566063dSJacob Faibussowitsch   PetscCall(MatMult(user->mass, user->tmp1, user->tmp2));
3159566063dSJacob Faibussowitsch   PetscCall(VecDot(user->tmp1, user->tmp2, &inner)); /* misfit contribution to */
316c4762a1bSJed Brown   *func += 0.5 * inner;                              /* the functional         */
317c4762a1bSJed Brown 
3189566063dSJacob Faibussowitsch   PetscCall(VecSetValues(user->tmp2, user->num_bc_dofs, user->bc_indices, user->bc_values, INSERT_VALUES));
3199566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(user->tmp2));
3209566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(user->tmp2));
3219566063dSJacob Faibussowitsch   PetscCall(KSPSolve(user->ksp_laplace, user->tmp2, user->adjoint)); /* adjoint solve */
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   /* And bring it home with the gradient. */
3249566063dSJacob Faibussowitsch   PetscCall(MatMult(user->mass, user->adjoint, user->tmp1));
3259566063dSJacob Faibussowitsch   PetscCall(VecAXPY(g, 1.0, user->tmp1)); /* adjoint contribution to the gradient */
3263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
327c4762a1bSJed Brown }
328c4762a1bSJed Brown 
main(int argc,char ** argv)329d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
330d71ae5a4SJacob Faibussowitsch {
331c4762a1bSJed Brown   DM     dm;
332c4762a1bSJed Brown   Tao    tao;
333c4762a1bSJed Brown   Vec    u, lb, ub;
334c4762a1bSJed Brown   AppCtx user;
335c4762a1bSJed Brown 
336327415f7SBarry Smith   PetscFunctionBeginUser;
3379566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3389566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
3399566063dSJacob Faibussowitsch   PetscCall(CreateCtx(dm, &user));
340c4762a1bSJed Brown 
3419566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
3429566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 0.0));
3439566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &lb));
3449566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &ub));
3459566063dSJacob Faibussowitsch   PetscCall(VecSet(lb, 0.0)); /* satisfied at the minimum anyway */
3469566063dSJacob Faibussowitsch   PetscCall(VecSet(ub, 0.8)); /* a nontrivial upper bound */
347c4762a1bSJed Brown 
3489566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
3499566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, u));
3509566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, ReducedFunctionGradient, &user));
3519566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao, lb, ub));
3529566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOBLMVM));
3539566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
354c4762a1bSJed Brown 
355c4762a1bSJed Brown   if (user.use_riesz) {
3569566063dSJacob Faibussowitsch     PetscCall(TaoLMVMSetH0(tao, user.mass)); /* crucial for mesh independence */
3579566063dSJacob Faibussowitsch     PetscCall(TaoSetGradientNorm(tao, user.mass));
358c4762a1bSJed Brown   }
359c4762a1bSJed Brown 
3609566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
361c4762a1bSJed Brown 
3629566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
3639566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-sol_view"));
364c4762a1bSJed Brown 
3659566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
3669566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
3689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lb));
3699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ub));
3709566063dSJacob Faibussowitsch   PetscCall(DestroyCtx(&user));
3719566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
372b122ec5aSJacob Faibussowitsch   return 0;
373c4762a1bSJed Brown }
374c4762a1bSJed Brown 
375c4762a1bSJed Brown /*TEST
376c4762a1bSJed Brown 
377c4762a1bSJed Brown     build:
378c4762a1bSJed Brown       requires: !complex !single
379c4762a1bSJed Brown 
380c4762a1bSJed Brown     test:
381342ad626SSatish Balay       requires: hdf5 double datafilespath !defined(PETSC_USE_64BIT_INDICES) hypre !cuda
38209fb3ddcSToby Isaac       args: -laplace_ksp_type cg -laplace_pc_type hypre -laplace_ksp_monitor_true_residual -laplace_ksp_max_it 5 -tao_blmvm_mat_lmvm_J0_ksp_type cg -tao_blmvm_mat_lmvm_J0_ksp_rtol 1e-5 -tao_blmvm_mat_lmvm_J0_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
383c4762a1bSJed Brown       filter: sed -e "s/-nan/nan/g"
384c4762a1bSJed Brown 
385c4762a1bSJed Brown     test:
386c4762a1bSJed Brown       suffix: guess_pod
387c4762a1bSJed Brown       requires: double triangle
38809fb3ddcSToby Isaac       args: -laplace_ksp_type cg -laplace_pc_type gamg -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 -tao_blmvm_mat_lmvm_J0_ksp_type cg -tao_blmvm_mat_lmvm_J0_ksp_rtol 1e-5 -tao_blmvm_mat_lmvm_J0_pc_type gamg -tao_blmvm_mat_lmvm_J0_pc_gamg_esteig_ksp_type cg -tao_blmvm_mat_lmvm_J0_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 -laplace_pc_gamg_aggressive_coarsening 0
38973f7197eSJed Brown       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"
390c4762a1bSJed Brown 
391c4762a1bSJed Brown TEST*/
392