Skip to content

Commit

Permalink
Merge pull request #15781 from bangerth/no-mpi
Browse files Browse the repository at this point in the history
Allow compilation with PETSc but without MPI.
  • Loading branch information
tamiko committed Jul 24, 2023
2 parents 9b9f5e1 + efb96ff commit db3b59f
Showing 1 changed file with 76 additions and 4 deletions.
80 changes: 76 additions & 4 deletions source/lac/petsc_communication_pattern.cc
Original file line number Diff line number Diff line change
Expand Up @@ -37,11 +37,15 @@ namespace PETScWrappers
: sf(nullptr)
{}



CommunicationPattern::~CommunicationPattern()
{
clear();
}



void
CommunicationPattern::reinit(const types::global_dof_index local_size,
const IndexSet & ghost_indices,
Expand Down Expand Up @@ -78,6 +82,8 @@ namespace PETScWrappers
AssertPETSc(PetscLayoutDestroy(&layout));
}



void
CommunicationPattern::reinit(const IndexSet &locally_owned_indices,
const IndexSet &ghost_indices,
Expand All @@ -96,6 +102,8 @@ namespace PETScWrappers
this->do_reinit(in_petsc, dummy, out_petsc, dummy, communicator);
}



void
CommunicationPattern::reinit(
const std::vector<types::global_dof_index> &indices_has,
Expand Down Expand Up @@ -149,6 +157,8 @@ namespace PETScWrappers
communicator);
}



void
CommunicationPattern::do_reinit(const std::vector<PetscInt> &inidx,
const std::vector<PetscInt> &inloc,
Expand Down Expand Up @@ -243,50 +253,80 @@ namespace PETScWrappers
AssertPETSc(PetscSFDestroy(&sf2));
}



void
CommunicationPattern::clear()
{
AssertPETSc(PetscSFDestroy(&sf));
}



MPI_Comm
CommunicationPattern::get_mpi_communicator() const
{
return PetscObjectComm(reinterpret_cast<PetscObject>(sf));
}



template <typename Number>
void
CommunicationPattern::export_to_ghosted_array_start(
const ArrayView<const Number> &src,
const ArrayView<Number> & dst) const
{
# ifdef DEAL_II_WITH_MPI
auto datatype = Utilities::MPI::mpi_type_id_for_type<Number>;

# if DEAL_II_PETSC_VERSION_LT(3, 15, 0)
# if DEAL_II_PETSC_VERSION_LT(3, 15, 0)
AssertPETSc(PetscSFBcastBegin(sf, datatype, src.data(), dst.data()));
# else
# else
AssertPETSc(
PetscSFBcastBegin(sf, datatype, src.data(), dst.data(), MPI_REPLACE));
# endif

# else

(void)src;
(void)dst;
Assert(false,
ExcMessage("This program is running without MPI. There should "
"not be anything to import or export!"));
# endif
}



template <typename Number>
void
CommunicationPattern::export_to_ghosted_array_finish(
const ArrayView<const Number> &src,
const ArrayView<Number> & dst) const
{
# ifdef DEAL_II_WITH_MPI
auto datatype = Utilities::MPI::mpi_type_id_for_type<Number>;

# if DEAL_II_PETSC_VERSION_LT(3, 15, 0)
# if DEAL_II_PETSC_VERSION_LT(3, 15, 0)
AssertPETSc(PetscSFBcastEnd(sf, datatype, src.data(), dst.data()));
# else
# else
AssertPETSc(
PetscSFBcastEnd(sf, datatype, src.data(), dst.data(), MPI_REPLACE));
# endif

# else

(void)src;
(void)dst;
Assert(false,
ExcMessage("This program is running without MPI. There should "
"not be anything to import or export!"));
# endif
}



template <typename Number>
void
CommunicationPattern::export_to_ghosted_array(
Expand All @@ -297,33 +337,61 @@ namespace PETScWrappers
export_to_ghosted_array_finish(src, dst);
}



template <typename Number>
void
CommunicationPattern::import_from_ghosted_array_start(
const VectorOperation::values op,
const ArrayView<const Number> &src,
const ArrayView<Number> & dst) const
{
# ifdef DEAL_II_WITH_MPI
MPI_Op mpiop = (op == VectorOperation::insert) ? MPI_REPLACE : MPI_SUM;
auto datatype = Utilities::MPI::mpi_type_id_for_type<Number>;

AssertPETSc(
PetscSFReduceBegin(sf, datatype, src.data(), dst.data(), mpiop));

# else

(void)op;
(void)src;
(void)dst;
Assert(false,
ExcMessage("This program is running without MPI. There should "
"not be anything to import or export!"));
# endif
}



template <typename Number>
void
CommunicationPattern::import_from_ghosted_array_finish(
const VectorOperation::values op,
const ArrayView<const Number> &src,
const ArrayView<Number> & dst) const
{
# ifdef DEAL_II_WITH_MPI
MPI_Op mpiop = (op == VectorOperation::insert) ? MPI_REPLACE : MPI_SUM;
auto datatype = Utilities::MPI::mpi_type_id_for_type<Number>;

AssertPETSc(PetscSFReduceEnd(sf, datatype, src.data(), dst.data(), mpiop));

# else

(void)op;
(void)src;
(void)dst;
Assert(false,
ExcMessage("This program is running without MPI. There should "
"not be anything to import or export!"));
# endif
}



template <typename Number>
void
CommunicationPattern::import_from_ghosted_array(
Expand All @@ -335,6 +403,8 @@ namespace PETScWrappers
import_from_ghosted_array_finish(op, src, dst);
}



// Partitioner

Partitioner::Partitioner()
Expand All @@ -345,6 +415,8 @@ namespace PETScWrappers
, n_ghost_indices_larger(numbers::invalid_dof_index)
{}



void
Partitioner::reinit(const IndexSet &locally_owned_indices,
const IndexSet &ghost_indices,
Expand Down

0 comments on commit db3b59f

Please sign in to comment.