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 PetscBool flg; 6030602db0SMatthew G. Knepley char filename[2048]; 6130602db0SMatthew G. Knepley PetscErrorCode ierr; 62c4762a1bSJed Brown 63c4762a1bSJed Brown PetscFunctionBeginUser; 64c4762a1bSJed Brown filename[0] = '\0'; 65c4762a1bSJed Brown user->use_riesz = PETSC_TRUE; 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(comm, "", "Poisson mother problem options", "DMPLEX");PetscCall(ierr); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_riesz", "Use the Riesz map to achieve mesh independence", "ex3.c", user->use_riesz, &user->use_riesz, NULL)); 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-f", "filename to read", "ex3.c", filename, filename, sizeof(filename), &flg)); 709566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 71c4762a1bSJed Brown 72c4762a1bSJed Brown if (!flg) { 739566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 749566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 75c4762a1bSJed Brown } else { 7630602db0SMatthew G. Knepley /* TODO Eliminate this in favor of DMLoad() in new code */ 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; 8230602db0SMatthew G. Knepley PetscInt dim = 2, 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 */ 919566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(comm, filename, FILE_MODE_READ, &viewer)); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* create Vecs to read in the data from H5 */ 949566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &coordinates)); 959566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 969566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(coordinates, dim)); 979566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &topology)); 989566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)topology, "topology")); 999566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(topology, vertices_per_cell)); 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* navigate to the right group */ 1029566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5PushGroup(viewer, "/Mesh/mesh")); 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* Read the Vecs */ 1059566063dSJacob Faibussowitsch PetscCall(VecLoad(coordinates, viewer)); 1069566063dSJacob Faibussowitsch PetscCall(VecLoad(topology, viewer)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* do some ugly calculations */ 1099566063dSJacob Faibussowitsch PetscCall(VecGetSize(topology, &numCells)); 110c4762a1bSJed Brown numCells = numCells / vertices_per_cell; 1119566063dSJacob Faibussowitsch PetscCall(VecGetSize(coordinates, &numVertices)); 112c4762a1bSJed Brown numVertices = numVertices / dim; 113c4762a1bSJed Brown 1149566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 1159566063dSJacob Faibussowitsch PetscCall(VecGetArray(topology, &topo_f)); 116c4762a1bSJed Brown /* and now we have to convert the double representation to integers to pass over, argh */ 1179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells*vertices_per_cell, &cells)); 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 */ 1239566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(comm, dim, numCells, numVertices, vertices_per_cell, PETSC_TRUE, cells, dim, coords, dm)); 124c4762a1bSJed Brown /* Check for flipped first cell */ 125c4762a1bSJed Brown { 126c4762a1bSJed Brown PetscReal v0[3], J[9], invJ[9], detJ; 127c4762a1bSJed Brown 1289566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(*dm, 0, NULL, v0, J, invJ, &detJ)); 129c4762a1bSJed Brown if (detJ < 0) { 1309566063dSJacob Faibussowitsch PetscCall(DMPlexOrientPoint(*dm, 0, -1)); 1319566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(*dm, 0, NULL, v0, J, invJ, &detJ)); 1323c859ba3SBarry Smith PetscCheck(detJ >= 0,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Something is wrong"); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown } 1359566063dSJacob Faibussowitsch PetscCall(DMPlexOrient(*dm)); 1369566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dm, "marker")); 1379566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "marker", &label)); 1389566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(*dm, 1, label)); 1399566063dSJacob Faibussowitsch PetscCall(DMPlexLabelComplete(*dm, label)); 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 1439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(topology, &topo_f)); 1449566063dSJacob Faibussowitsch PetscCall(PetscFree(cells)); 1459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 1469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&topology)); 147c4762a1bSJed Brown #else 148c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Reconfigure PETSc with --download-hdf5"); 149c4762a1bSJed Brown #endif 150c4762a1bSJed Brown } 1519566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1529566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 153c4762a1bSJed Brown PetscFunctionReturn(0); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156c4762a1bSJed Brown void mass_kernel(PetscInt dim, PetscInt Nf, PetscInt NfAux, 157c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 158c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 159c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 160c4762a1bSJed Brown { 161c4762a1bSJed Brown g0[0] = 1.0; 162c4762a1bSJed Brown } 163c4762a1bSJed Brown 164c4762a1bSJed Brown void laplace_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 g3[]) 168c4762a1bSJed Brown { 169c4762a1bSJed Brown PetscInt d; 170c4762a1bSJed Brown for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0; 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* data we seek to match */ 174c4762a1bSJed Brown PetscErrorCode data_kernel(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *y, void *ctx) 175c4762a1bSJed Brown { 176c4762a1bSJed Brown *y = 1.0/(2*PETSC_PI*PETSC_PI) * PetscSinReal(PETSC_PI*x[0]) * PetscSinReal(PETSC_PI*x[1]); 177c4762a1bSJed Brown /* the associated control is sin(pi*x[0])*sin(pi*x[1]) */ 178c4762a1bSJed Brown return 0; 179c4762a1bSJed Brown } 180c4762a1bSJed Brown PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 181c4762a1bSJed Brown { 182c4762a1bSJed Brown *u = 0.0; 183c4762a1bSJed Brown return 0; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown 186c4762a1bSJed Brown PetscErrorCode CreateCtx(DM dm, AppCtx* user) 187c4762a1bSJed Brown { 188c4762a1bSJed Brown 189c4762a1bSJed Brown DM dm_mass; 190c4762a1bSJed Brown DM dm_laplace; 191c4762a1bSJed Brown PetscDS prob_mass; 192c4762a1bSJed Brown PetscDS prob_laplace; 193c4762a1bSJed Brown PetscFE fe; 194c4762a1bSJed Brown DMLabel label; 195c4762a1bSJed Brown PetscSection section; 196c4762a1bSJed Brown PetscInt n, k, p, d; 197c4762a1bSJed Brown PetscInt dof, off; 198c4762a1bSJed Brown IS is; 199c4762a1bSJed Brown const PetscInt* points; 200c4762a1bSJed Brown const PetscInt dim = 2; 201c4762a1bSJed Brown const PetscInt id = 1; 202c4762a1bSJed Brown PetscErrorCode (**wtf)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 203c4762a1bSJed Brown 204c4762a1bSJed Brown PetscFunctionBeginUser; 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* make the data we seek to match */ 2079566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, PETSC_TRUE, NULL, 4, &fe)); 208c4762a1bSJed Brown 2099566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 2109566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 2119566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &user->data)); 212c4762a1bSJed Brown 213c4762a1bSJed Brown /* ugh, this is hideous */ 214c4762a1bSJed Brown /* y_d = interpolate(Expression("sin(x[0]) + .."), V) */ 2159566063dSJacob Faibussowitsch PetscCall(PetscMalloc(1 * sizeof(void (*)(const PetscReal[], PetscScalar *, void *)), &wtf)); 216c4762a1bSJed Brown wtf[0] = data_kernel; 2179566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, wtf, NULL, INSERT_VALUES, user->data)); 2189566063dSJacob Faibussowitsch PetscCall(PetscFree(wtf)); 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* assemble(inner(u, v)*dx), almost */ 2219566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dm_mass)); 2229566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, dm_mass)); 2239566063dSJacob Faibussowitsch PetscCall(DMSetNumFields(dm_mass, 1)); 2249566063dSJacob Faibussowitsch PetscCall(DMPlexCopyCoordinates(dm, dm_mass)); /* why do I have to do this separately? */ 2259566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm_mass, &prob_mass)); 2269566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob_mass, 0, 0, mass_kernel, NULL, NULL, NULL)); 2279566063dSJacob Faibussowitsch PetscCall(PetscDSSetDiscretization(prob_mass, 0, (PetscObject) fe)); 2289566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm_mass, &user->mass)); 2299566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm_mass, user->data, user->mass, user->mass, NULL)); 2309566063dSJacob Faibussowitsch PetscCall(MatSetOption(user->mass, MAT_SYMMETRIC, PETSC_TRUE)); 2319566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm_mass)); 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* inner(grad(u), grad(v))*dx with homogeneous Dirichlet boundary conditions */ 2349566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dm_laplace)); 2359566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, dm_laplace)); 2369566063dSJacob Faibussowitsch PetscCall(DMSetNumFields(dm_laplace, 1)); 2379566063dSJacob Faibussowitsch PetscCall(DMPlexCopyCoordinates(dm, dm_laplace)); 2389566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm_laplace, &prob_laplace)); 2399566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob_laplace, 0, 0, NULL, NULL, NULL, laplace_kernel)); 2409566063dSJacob Faibussowitsch PetscCall(PetscDSSetDiscretization(prob_laplace, 0, (PetscObject) fe)); 2419566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm_laplace, &user->laplace)); 2429566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm_laplace, user->data, user->laplace, user->laplace, NULL)); 243c4762a1bSJed Brown 244c4762a1bSJed Brown /* Code from Matt to get the indices associated with the boundary dofs */ 2459566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm_laplace, "marker", &label)); 2469566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm_laplace, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) zero, NULL, NULL, NULL)); 2479566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm_laplace, §ion)); 2489566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumSize(label, 1, &n)); 2499566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, 1, &is)); 2509566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &points)); 251c4762a1bSJed Brown user->num_bc_dofs = 0; 252c4762a1bSJed Brown for (p = 0; p < n; ++p) { 2539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, points[p], &dof)); 254c4762a1bSJed Brown user->num_bc_dofs += dof; 255c4762a1bSJed Brown } 2569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->num_bc_dofs, &user->bc_indices)); 257c4762a1bSJed Brown for (p = 0, k = 0; p < n; ++p) { 2589566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, points[p], &dof)); 2599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, points[p], &off)); 260c4762a1bSJed Brown for (d = 0; d < dof; ++d) user->bc_indices[k++] = off+d; 261c4762a1bSJed Brown } 2629566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &points)); 2639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 2649566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm_laplace)); 265c4762a1bSJed Brown 266c4762a1bSJed Brown /* This is how I handle boundary conditions. I can't figure out how to get 267c4762a1bSJed Brown plex to play with the way I want to impose the BCs. This loses symmetry, 268c4762a1bSJed Brown but not in a disastrous way. If someone can improve it, please do! */ 2699566063dSJacob Faibussowitsch PetscCall(MatZeroRows(user->laplace, user->num_bc_dofs, user->bc_indices, 1.0, NULL, NULL)); 2709566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(user->num_bc_dofs, &user->bc_values)); 271c4762a1bSJed Brown 272c4762a1bSJed Brown /* also create the KSP for solving the Laplace system */ 2739566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &user->ksp_laplace)); 2749566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(user->ksp_laplace, user->laplace, user->laplace)); 2759566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(user->ksp_laplace, "laplace_")); 2769566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(user->ksp_laplace)); 277c4762a1bSJed Brown 278c4762a1bSJed Brown /* A bit of setting up the user context */ 279c4762a1bSJed Brown user->dm = dm; 2809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->data, &user->state)); 2819566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->data, &user->adjoint)); 2829566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->data, &user->tmp1)); 2839566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user->data, &user->tmp2)); 284c4762a1bSJed Brown 2859566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 286c4762a1bSJed Brown 287c4762a1bSJed Brown PetscFunctionReturn(0); 288c4762a1bSJed Brown } 289c4762a1bSJed Brown 290c4762a1bSJed Brown PetscErrorCode DestroyCtx(AppCtx* user) 291c4762a1bSJed Brown { 292c4762a1bSJed Brown PetscFunctionBeginUser; 293c4762a1bSJed Brown 2949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->mass)); 2959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user->laplace)); 2969566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&user->ksp_laplace)); 2979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->data)); 2989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->state)); 2999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->adjoint)); 3009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->tmp1)); 3019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user->tmp2)); 3029566063dSJacob Faibussowitsch PetscCall(PetscFree(user->bc_indices)); 3039566063dSJacob Faibussowitsch PetscCall(PetscFree(user->bc_values)); 304c4762a1bSJed Brown 305c4762a1bSJed Brown PetscFunctionReturn(0); 306c4762a1bSJed Brown } 307c4762a1bSJed Brown 308c4762a1bSJed Brown PetscErrorCode ReducedFunctionGradient(Tao tao, Vec u, PetscReal* func, Vec g, void* userv) 309c4762a1bSJed Brown { 310c4762a1bSJed Brown AppCtx* user = (AppCtx*) userv; 311c4762a1bSJed Brown const PetscReal alpha = 1.0e-6; /* regularisation parameter */ 312c4762a1bSJed Brown PetscReal inner; 313c4762a1bSJed Brown 314c4762a1bSJed Brown PetscFunctionBeginUser; 315c4762a1bSJed Brown 3169566063dSJacob Faibussowitsch PetscCall(MatMult(user->mass, u, user->tmp1)); 3179566063dSJacob Faibussowitsch PetscCall(VecDot(u, user->tmp1, &inner)); /* regularisation contribution to */ 318c4762a1bSJed Brown *func = alpha * 0.5 * inner; /* the functional */ 319c4762a1bSJed Brown 3209566063dSJacob Faibussowitsch PetscCall(VecSet(g, 0.0)); 3219566063dSJacob Faibussowitsch PetscCall(VecAXPY(g, alpha, user->tmp1)); /* regularisation contribution to the gradient */ 322c4762a1bSJed Brown 323c4762a1bSJed Brown /* Now compute the forward state. */ 3249566063dSJacob Faibussowitsch PetscCall(VecSetValues(user->tmp1, user->num_bc_dofs, user->bc_indices, user->bc_values, INSERT_VALUES)); 3259566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(user->tmp1)); 3269566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(user->tmp1)); 3279566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->ksp_laplace, user->tmp1, user->state)); /* forward solve */ 328c4762a1bSJed Brown 329c4762a1bSJed Brown /* Now compute the adjoint state also. */ 3309566063dSJacob Faibussowitsch PetscCall(VecCopy(user->state, user->tmp1)); 3319566063dSJacob Faibussowitsch PetscCall(VecAXPY(user->tmp1, -1.0, user->data)); 3329566063dSJacob Faibussowitsch PetscCall(MatMult(user->mass, user->tmp1, user->tmp2)); 3339566063dSJacob Faibussowitsch PetscCall(VecDot(user->tmp1, user->tmp2, &inner)); /* misfit contribution to */ 334c4762a1bSJed Brown *func += 0.5 * inner; /* the functional */ 335c4762a1bSJed Brown 3369566063dSJacob Faibussowitsch PetscCall(VecSetValues(user->tmp2, user->num_bc_dofs, user->bc_indices, user->bc_values, INSERT_VALUES)); 3379566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(user->tmp2)); 3389566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(user->tmp2)); 3399566063dSJacob Faibussowitsch PetscCall(KSPSolve(user->ksp_laplace, user->tmp2, user->adjoint)); /* adjoint solve */ 340c4762a1bSJed Brown 341c4762a1bSJed Brown /* And bring it home with the gradient. */ 3429566063dSJacob Faibussowitsch PetscCall(MatMult(user->mass, user->adjoint, user->tmp1)); 3439566063dSJacob Faibussowitsch PetscCall(VecAXPY(g, 1.0, user->tmp1)); /* adjoint contribution to the gradient */ 344c4762a1bSJed Brown 345c4762a1bSJed Brown PetscFunctionReturn(0); 346c4762a1bSJed Brown } 347c4762a1bSJed Brown 348c4762a1bSJed Brown int main(int argc, char **argv) 349c4762a1bSJed Brown { 350c4762a1bSJed Brown DM dm; 351c4762a1bSJed Brown Tao tao; 352c4762a1bSJed Brown Vec u, lb, ub; 353c4762a1bSJed Brown AppCtx user; 354c4762a1bSJed Brown 3559566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 3569566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 3579566063dSJacob Faibussowitsch PetscCall(CreateCtx(dm, &user)); 358c4762a1bSJed Brown 3599566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 3609566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 3619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &lb)); 3629566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &ub)); 3639566063dSJacob Faibussowitsch PetscCall(VecSet(lb, 0.0)); /* satisfied at the minimum anyway */ 3649566063dSJacob Faibussowitsch PetscCall(VecSet(ub, 0.8)); /* a nontrivial upper bound */ 365c4762a1bSJed Brown 3669566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 3679566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, u)); 3689566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao,NULL, ReducedFunctionGradient, &user)); 3699566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(tao, lb, ub)); 3709566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBLMVM)); 3719566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 372c4762a1bSJed Brown 373c4762a1bSJed Brown if (user.use_riesz) { 3749566063dSJacob Faibussowitsch PetscCall(TaoLMVMSetH0(tao, user.mass)); /* crucial for mesh independence */ 3759566063dSJacob Faibussowitsch PetscCall(TaoSetGradientNorm(tao, user.mass)); 376c4762a1bSJed Brown } 377c4762a1bSJed Brown 3789566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 379c4762a1bSJed Brown 3809566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 3819566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_view")); 382c4762a1bSJed Brown 3839566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 3849566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 3859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 3869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lb)); 3879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ub)); 3889566063dSJacob Faibussowitsch PetscCall(DestroyCtx(&user)); 3899566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 390b122ec5aSJacob Faibussowitsch return 0; 391c4762a1bSJed Brown } 392c4762a1bSJed Brown 393c4762a1bSJed Brown /*TEST 394c4762a1bSJed Brown 395c4762a1bSJed Brown build: 396c4762a1bSJed Brown requires: !complex !single 397c4762a1bSJed Brown 398c4762a1bSJed Brown test: 399*342ad626SSatish Balay requires: hdf5 double datafilespath !defined(PETSC_USE_64BIT_INDICES) hypre !cuda 40073f7197eSJed Brown 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 401c4762a1bSJed Brown filter: sed -e "s/-nan/nan/g" 402c4762a1bSJed Brown 403c4762a1bSJed Brown test: 404c4762a1bSJed Brown suffix: guess_pod 405c4762a1bSJed Brown requires: double triangle 406074aec55SMark Adams 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 40773f7197eSJed 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" 408c4762a1bSJed Brown 409c4762a1bSJed Brown TEST*/ 410