/*
* Rather than using ifdefs everywhere, try to wrap older versions of PETSc
* functions in one place.
+ *
+ * Functions that are not inlined are:
+ * - Functions returning PetscErrorCode that are supposed to be called within
+ * PETSc callbacks.
+ * - Functions that need access to internal PETSc headers that we don't want
+ * to expose to deal.II users
*/
#ifndef dealii_petsc_compatibility_h
#define dealii_petsc_compatibility_h
set_matrix_option(matrix, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE);
}
+
+
/**
* Tell PETSc that the status of the vector has changed.
*/
void
petsc_increment_state_counter(Vec v);
+
+
/**
* Tell PETSc that the status of the matrix has changed.
*/
void
petsc_increment_state_counter(Mat A);
+
+ /**
+ * Set the failed reason for the preconditioner.
+ */
+ PetscErrorCode
+ pc_set_failed_reason(PC pc, PCFailedReason reason);
+
+
+
/**
* Resets internal domain error flags in the SNES object.
*/
void
snes_reset_domain_flags(SNES snes);
+
+
+ /**
+ * Resets internal domain error flags in the SNES object.
+ */
+ void
+ snes_set_jacobian_domain_error(SNES snes);
+
+
+
/**
* Tell PETSc nonlinear solver to use matrix free finite differencing (MFFD).
*
void
set_use_matrix_free(SNES snes, const bool mf_operator, const bool mf);
+
+
/**
* Tell PETSc ODE solver to use matrix free finite differencing (MFFD).
*
void
set_use_matrix_free(TS ts, const bool mf_operator, const bool mf);
+
+
/**
* Set final time for ODE integration.
*/
void
ts_set_max_time(TS ts, const PetscReal maxtime);
+
+
/**
* Set maximum number of steps for ODE integration.
*/
void
ts_set_max_steps(TS ts, const PetscInt maxsteps);
+
+
/**
* Return current step number.
*/
unsigned int
ts_get_step_number(TS ts);
+
+
/**
* Return true if the TS has a SNES object.
*/
AssertPETSc(SNESReset(snes));
// By default PETSc does not check for Jacobian errors in optimized
// mode. Here we do it unconditionally.
+# if DEAL_II_PETSC_VERSION_GTE(3, 11, 0)
AssertPETSc(SNESSetCheckJacobianDomainError(snes, PETSC_TRUE));
+# endif
}
const int err = call_and_possibly_capture_snes_exception(
user->jacobian,
user->pending_exception,
- [snes]() -> void { AssertPETSc(SNESSetJacobianDomainError(snes)); },
+ [snes]() -> void { snes_set_jacobian_domain_error(snes); },
xdealii,
Adealii,
Pdealii);
const int err = call_and_possibly_capture_snes_exception(
user->setup_jacobian,
user->pending_exception,
- [snes]() -> void { AssertPETSc(SNESSetJacobianDomainError(snes)); },
+ [snes]() -> void { snes_set_jacobian_domain_error(snes); },
xdealii);
if (err)
return PetscError(
[ts]() -> void {
SNES snes;
AssertPETSc(TSGetSNES(ts, &snes));
- AssertPETSc(SNESSetJacobianDomainError(snes));
+ snes_set_jacobian_domain_error(snes);
},
t,
xdealii,
[ts]() -> void {
SNES snes;
AssertPETSc(TSGetSNES(ts, &snes));
- AssertPETSc(SNESSetJacobianDomainError(snes));
+ snes_set_jacobian_domain_error(snes);
},
t,
xdealii,
[ts]() -> void {
SNES snes;
AssertPETSc(TSGetSNES(ts, &snes));
- AssertPETSc(SNESSetJacobianDomainError(snes));
+ snes_set_jacobian_domain_error(snes);
},
t,
xdealii,
// By default PETSc does not check for Jacobian errors in optimized
// mode. Here we do it unconditionally.
+# if DEAL_II_PETSC_VERSION_GTE(3, 11, 0)
if (ts_has_snes(ts))
{
SNES snes;
AssertPETSc(TSGetSNES(ts, &snes));
AssertPETSc(SNESSetCheckJacobianDomainError(snes, PETSC_TRUE));
}
-
+# endif
// Having set everything up, now do the actual work
// and let PETSc do the time stepping. If there is
// a pending exception, then one of the user callbacks
#ifdef DEAL_II_WITH_PETSC
+# include <petsc/private/pcimpl.h>
# include <petsc/private/petscimpl.h>
# include <petsc/private/snesimpl.h>
# include <petsc/private/tsimpl.h>
AssertPETSc(PetscObjectStateIncrease(reinterpret_cast<PetscObject>(A)));
}
+ PetscErrorCode
+ pc_set_failed_reason(PC pc, PCFailedReason reason)
+ {
+# if DEAL_II_PETSC_VERSION_GTE(3, 14, 0)
+ return PCSetFailedReason(pc, reason);
+# else
+ pc->failedreason = reason;
+ return 0;
+# endif
+ }
+
void
snes_reset_domain_flags(SNES snes)
{
+# if DEAL_II_PETSC_VERSION_GTE(3, 11, 0)
snes->jacobiandomainerror = PETSC_FALSE;
- snes->domainerror = PETSC_FALSE;
+# endif
+ snes->domainerror = PETSC_FALSE;
+ }
+
+ void
+ snes_set_jacobian_domain_error(SNES snes)
+ {
+# if DEAL_II_PETSC_VERSION_GTE(3, 11, 0)
+ snes->jacobiandomainerror = PETSC_TRUE;
+# else
+ // There is no equivalent, and since this used to stop
+ // computations, we opt to set the converged reason
+ snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
+# endif
}
void
// Failed reason is not reset uniformly within the
// interface code of PCSetUp in PETSc.
// We handle it here.
- PetscCall(PCSetFailedReason(ppc, PC_NOERROR));
+ PetscCall(pc_set_failed_reason(ppc, PC_NOERROR));
PetscFunctionReturn(0);
}
}
catch (const RecoverableUserCallbackError &)
{
- PetscCall(PCSetFailedReason(ppc, PC_SUBPC_ERROR));
+ PetscCall(pc_set_failed_reason(ppc, PC_SUBPC_ERROR));
}
catch (...)
{
}
catch (const RecoverableUserCallbackError &)
{
- PetscCall(PCSetFailedReason(ppc, PC_SUBPC_ERROR));
+ PetscCall(pc_set_failed_reason(ppc, PC_SUBPC_ERROR));
}
catch (...)
{
test(v, w);
}
-#if DEAL_II_PETSC_VERSION_LT(3, 2, 0)
- ierr = VecDestroy(vpetsc);
-#else
ierr = VecDestroy(&vpetsc);
-#endif
AssertThrow(ierr == 0, ExcPETScError(ierr));
}
{
const PetscScalar s(1. * k, 2. * k);
v(k) = s;
-#if DEAL_II_PETSC_VERSION_LT(3, 7, 0)
- norm += std::fabs(1. * k) + std::fabs(2. * k);
-#else
norm += std::abs(s);
-#endif
}
v.compress(VectorOperation::insert);
test(v, w);
}
-#if DEAL_II_PETSC_VERSION_LT(3, 2, 0)
- ierr = VecDestroy(vpetsc);
-#else
ierr = VecDestroy(&vpetsc);
-#endif
AssertThrow(ierr == 0, ExcPETScError(ierr));
}
void
check_petsc_allocations()
{
-# if DEAL_II_PETSC_VERSION_GTE(3, 2, 0)
PetscStageLog stageLog;
PetscLogGetStageLog(&stageLog);
if (errors)
throw dealii::ExcMessage("PETSc memory leak");
-# endif
}
} // namespace
#endif