Skip to content

Commit 0f07931

Browse files
committed
Fixes to get libmesh compiling with PETSc master as of late July 2019.
* Fix for SNES_CONVERGED_TR_DELTA. * Fix for dm->ops->getinjection. * Fix for dm->ops->getaggregates. All of these changes are currently in PETSc master but not yet in any PETSc release. For example, the SNES_CONVERGED_TR_DELTA was merged into PETSc master 16 Jul 2019 in: c1c6be5a2e, and the latest PETSc release is v3.11.3 (26 Jun 2019). These changes are therefore likely to be in 3.12 but not in any 3.11.x release. Also, this change will break anyone who is using libmesh with an older version of PETSc master, but there is basically nothing we can do about that, they just have to update.
1 parent 16f585a commit 0f07931

File tree

2 files changed

+14
-0
lines changed

2 files changed

+14
-0
lines changed

src/solvers/petsc_diff_solver.C

+5
Original file line numberDiff line numberDiff line change
@@ -263,7 +263,12 @@ DiffSolver::SolveResult convert_solve_result(SNESConvergedReason r)
263263
#endif
264264
return DiffSolver::CONVERGED_RELATIVE_STEP;
265265
case SNES_CONVERGED_ITS:
266+
// SNES_CONVERGED_TR_DELTA was changed to a diverged condition,
267+
// SNES_DIVERGED_TR_DELTA, in PETSc 1c6b2ff8df. This change will
268+
// likely be in 3.12 and later releases.
269+
#if PETSC_RELEASE_LESS_THAN(3,12,0)
266270
case SNES_CONVERGED_TR_DELTA:
271+
#endif
267272
return DiffSolver::CONVERGED_NO_REASON;
268273
case SNES_DIVERGED_FUNCTION_DOMAIN:
269274
case SNES_DIVERGED_FUNCTION_COUNT:

src/solvers/petscdmlibmeshimpl.C

+9
Original file line numberDiff line numberDiff line change
@@ -1192,8 +1192,17 @@ PetscErrorCode DMCreate_libMesh(DM dm)
11921192

11931193
dm->ops->refine = 0; // DMRefine_libMesh;
11941194
dm->ops->coarsen = 0; // DMCoarsen_libMesh;
1195+
1196+
// * dm->ops->getinjection was renamed to dm->ops->createinjection in PETSc 5a84ad338 (5 Jul 2019)
1197+
// * dm->ops-getaggregates was removed in PETSc 97779f9a (5 Jul 2019)
1198+
// * Both changes were merged into PETSc master in 94aad3ce (7 Jul 2019).
1199+
#if PETSC_RELEASE_LESS_THAN(3,12,0)
11951200
dm->ops->getinjection = 0; // DMGetInjection_libMesh;
11961201
dm->ops->getaggregates = 0; // DMGetAggregates_libMesh;
1202+
#else
1203+
dm->ops->createinjection = 0;
1204+
#endif
1205+
11971206

11981207
#if PETSC_RELEASE_LESS_THAN(3,3,1)
11991208
dm->ops->createfielddecompositiondm = DMCreateFieldDecompositionDM_libMesh;

0 commit comments

Comments
 (0)