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