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