From: Pasquale Africa Date: Wed, 3 Feb 2021 12:29:07 +0000 (+0000) Subject: Rename PETScWrappers::PreconditionerBase to PETScWrappers::PreconditionBase for consi... X-Git-Tag: v9.3.0-rc1~478^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F11670%2Fhead;p=dealii.git Rename PETScWrappers::PreconditionerBase to PETScWrappers::PreconditionBase for consistency --- diff --git a/doc/news/changes/minor/20210203PasqualeAfrica b/doc/news/changes/minor/20210203PasqualeAfrica new file mode 100644 index 0000000000..42d63052e9 --- /dev/null +++ b/doc/news/changes/minor/20210203PasqualeAfrica @@ -0,0 +1,4 @@ +Update: class PETScWrappers::PreconditionerBase renamed to +PETScWrappers::PreconditionBase for consistency. +
+(Pasquale Claudio Africa, 2021/02/03) diff --git a/include/deal.II/lac/petsc_precondition.h b/include/deal.II/lac/petsc_precondition.h index 7e46a63658..c5a27bbc62 100644 --- a/include/deal.II/lac/petsc_precondition.h +++ b/include/deal.II/lac/petsc_precondition.h @@ -53,18 +53,18 @@ namespace PETScWrappers * * @ingroup PETScWrappers */ - class PreconditionerBase + class PreconditionBase { public: /** * Constructor. */ - PreconditionerBase(); + PreconditionBase(); /** * Destructor. */ - virtual ~PreconditionerBase(); + virtual ~PreconditionBase(); /** * Destroys the preconditioner, leaving an object like just after having @@ -129,12 +129,12 @@ namespace PETScWrappers * preconditioner. * * See the comment in the base class - * @ref PreconditionerBase + * @ref PreconditionBase * for when this preconditioner may or may not work. * * @ingroup PETScWrappers */ - class PreconditionJacobi : public PreconditionerBase + class PreconditionJacobi : public PreconditionBase { public: /** @@ -212,12 +212,12 @@ namespace PETScWrappers * the relevant section of the PETSc manual, but is not implemented here. * * See the comment in the base class - * @ref PreconditionerBase + * @ref PreconditionBase * for when this preconditioner may or may not work. * * @ingroup PETScWrappers */ - class PreconditionBlockJacobi : public PreconditionerBase + class PreconditionBlockJacobi : public PreconditionBase { public: /** @@ -285,7 +285,7 @@ namespace PETScWrappers * * @ingroup PETScWrappers */ - class PreconditionSOR : public PreconditionerBase + class PreconditionSOR : public PreconditionBase { public: /** @@ -345,7 +345,7 @@ namespace PETScWrappers * * @ingroup PETScWrappers */ - class PreconditionSSOR : public PreconditionerBase + class PreconditionSSOR : public PreconditionBase { public: /** @@ -403,12 +403,12 @@ namespace PETScWrappers * each processor. * * See the comment in the base class - * @ref PreconditionerBase + * @ref PreconditionBase * for when this preconditioner may or may not work. * * @ingroup PETScWrappers */ - class PreconditionEisenstat : public PreconditionerBase + class PreconditionEisenstat : public PreconditionBase { public: /** @@ -469,7 +469,7 @@ namespace PETScWrappers * * @ingroup PETScWrappers */ - class PreconditionICC : public PreconditionerBase + class PreconditionICC : public PreconditionBase { public: /** @@ -529,7 +529,7 @@ namespace PETScWrappers * * @ingroup PETScWrappers */ - class PreconditionILU : public PreconditionerBase + class PreconditionILU : public PreconditionBase { public: /** @@ -587,14 +587,14 @@ namespace PETScWrappers * (depending on the settings) performs an exact factorization of the * matrix, so it is not necessary to wrap it in an iterative solver. This * class is typically used with SolverPreOnly to get a direct - * solver. Alternatively, you can use PreconditionerBase::vmult() directly. + * solver. Alternatively, you can use PreconditionBase::vmult() directly. * * @note This is not a parallel preconditioner so it only works in serial * computations with a single processor. * * @ingroup PETScWrappers */ - class PreconditionLU : public PreconditionerBase + class PreconditionLU : public PreconditionBase { public: /** @@ -673,7 +673,7 @@ namespace PETScWrappers * * @ingroup PETScWrappers */ - class PreconditionBoomerAMG : public PreconditionerBase + class PreconditionBoomerAMG : public PreconditionBase { public: /** @@ -803,7 +803,7 @@ namespace PETScWrappers * * @ingroup PETScWrappers */ - class PreconditionParaSails : public PreconditionerBase + class PreconditionParaSails : public PreconditionBase { public: /** @@ -913,7 +913,7 @@ namespace PETScWrappers * * @ingroup PETScWrappers */ - class PreconditionNone : public PreconditionerBase + class PreconditionNone : public PreconditionBase { public: /** @@ -954,6 +954,12 @@ namespace PETScWrappers */ AdditionalData additional_data; }; + + /** + * Alias for backwards-compatibility. + * @deprecated Use PETScWrappers::PreconditionBase instead. + */ + using PreconditionerBase DEAL_II_DEPRECATED_EARLY = PreconditionBase; } // namespace PETScWrappers diff --git a/include/deal.II/lac/petsc_solver.h b/include/deal.II/lac/petsc_solver.h index 8bd7aa7675..46e7a7eea7 100644 --- a/include/deal.II/lac/petsc_solver.h +++ b/include/deal.II/lac/petsc_solver.h @@ -51,7 +51,7 @@ namespace PETScWrappers # ifndef DOXYGEN class MatrixBase; class VectorBase; - class PreconditionerBase; + class PreconditionBase; # endif @@ -136,10 +136,10 @@ namespace PETScWrappers * performance reasons. See class Documentation. */ void - solve(const MatrixBase & A, - VectorBase & x, - const VectorBase & b, - const PreconditionerBase &preconditioner); + solve(const MatrixBase & A, + VectorBase & x, + const VectorBase & b, + const PreconditionBase &preconditioner); /** @@ -169,7 +169,7 @@ namespace PETScWrappers * intended for use with SLEPc spectral transformation class. */ void - initialize(const PreconditionerBase &preconditioner); + initialize(const PreconditionBase &preconditioner); protected: /** diff --git a/source/lac/petsc_precondition.cc b/source/lac/petsc_precondition.cc index de93af0f36..13aa14efb9 100644 --- a/source/lac/petsc_precondition.cc +++ b/source/lac/petsc_precondition.cc @@ -33,12 +33,12 @@ DEAL_II_NAMESPACE_OPEN namespace PETScWrappers { - PreconditionerBase::PreconditionerBase() + PreconditionBase::PreconditionBase() : pc(nullptr) , matrix(nullptr) {} - PreconditionerBase::~PreconditionerBase() + PreconditionBase::~PreconditionBase() { try { @@ -49,7 +49,7 @@ namespace PETScWrappers } void - PreconditionerBase::clear() + PreconditionBase::clear() { matrix = nullptr; @@ -63,7 +63,7 @@ namespace PETScWrappers void - PreconditionerBase::vmult(VectorBase &dst, const VectorBase &src) const + PreconditionBase::vmult(VectorBase &dst, const VectorBase &src) const { AssertThrow(pc != nullptr, StandardExceptions::ExcInvalidState()); @@ -73,7 +73,7 @@ namespace PETScWrappers void - PreconditionerBase::Tvmult(VectorBase &dst, const VectorBase &src) const + PreconditionBase::Tvmult(VectorBase &dst, const VectorBase &src) const { AssertThrow(pc != nullptr, StandardExceptions::ExcInvalidState()); @@ -83,7 +83,7 @@ namespace PETScWrappers void - PreconditionerBase::create_pc() + PreconditionBase::create_pc() { // only allow the creation of the // preconditioner once @@ -110,13 +110,13 @@ namespace PETScWrappers const PC & - PreconditionerBase::get_pc() const + PreconditionBase::get_pc() const { return pc; } - PreconditionerBase::operator Mat() const + PreconditionBase::operator Mat() const { return matrix; } diff --git a/source/lac/petsc_solver.cc b/source/lac/petsc_solver.cc index 5462c3cb51..ad9d181573 100644 --- a/source/lac/petsc_solver.cc +++ b/source/lac/petsc_solver.cc @@ -50,10 +50,10 @@ namespace PETScWrappers void - SolverBase::solve(const MatrixBase & A, - VectorBase & x, - const VectorBase & b, - const PreconditionerBase &preconditioner) + SolverBase::solve(const MatrixBase & A, + VectorBase & x, + const VectorBase & b, + const PreconditionBase &preconditioner) { /* TODO: PETSc duplicates communicators, so this does not work (you put @@ -210,7 +210,7 @@ namespace PETScWrappers } void - SolverBase::initialize(const PreconditionerBase &preconditioner) + SolverBase::initialize(const PreconditionBase &preconditioner) { PetscErrorCode ierr; diff --git a/tests/slepc/step-36_parallel.cc b/tests/slepc/step-36_parallel.cc index 22f10a2fc2..43745d78a9 100644 --- a/tests/slepc/step-36_parallel.cc +++ b/tests/slepc/step-36_parallel.cc @@ -217,7 +217,7 @@ test(std::string solver_name, std::string preconditioner_name) // test SLEPc by { - PETScWrappers::PreconditionerBase *preconditioner; + PETScWrappers::PreconditionBase *preconditioner; dealii::deallog << preconditioner_name << std::endl; if (preconditioner_name == "Jacobi") diff --git a/tests/slepc/step-36_parallel_02.cc b/tests/slepc/step-36_parallel_02.cc index 503bcc5e60..bf0010d0b1 100644 --- a/tests/slepc/step-36_parallel_02.cc +++ b/tests/slepc/step-36_parallel_02.cc @@ -217,7 +217,7 @@ test(std::string solver_name, std::string preconditioner_name) // test SLEPc by { - PETScWrappers::PreconditionerBase *preconditioner; + PETScWrappers::PreconditionBase *preconditioner; dealii::deallog << preconditioner_name << std::endl; if (preconditioner_name == "Jacobi") diff --git a/tests/slepc/step-36_parallel_03.cc b/tests/slepc/step-36_parallel_03.cc index 1a86c9284b..d2b63d9435 100644 --- a/tests/slepc/step-36_parallel_03.cc +++ b/tests/slepc/step-36_parallel_03.cc @@ -217,7 +217,7 @@ test(std::string solver_name, std::string preconditioner_name) // test SLEPc by { - PETScWrappers::PreconditionerBase *preconditioner; + PETScWrappers::PreconditionBase *preconditioner; dealii::deallog << preconditioner_name << std::endl; if (preconditioner_name == "Jacobi")