Skip to content

Commit

Permalink
Require C++17 compiler (gcc 7 and later) (#641)
Browse files Browse the repository at this point in the history
* Require C++17 in CMake file.

* Compress some code.

* More cleaning up.

* More C++17.

* Fix tuple return.

* Undo some changes.

* More std::tie updates.

* Fix for unused variable.

* Fix cpp version for unit tests.

* Revert a file.

* Fix some confusing mesh topology code.

* Simplify some mesh topology.
  • Loading branch information
garth-wells committed Nov 15, 2019
1 parent 2b67be3 commit c051360
Show file tree
Hide file tree
Showing 16 changed files with 49 additions and 105 deletions.
6 changes: 3 additions & 3 deletions cpp/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -26,10 +26,10 @@ if (POLICY CMP0075)
endif()

#------------------------------------------------------------------------------
# Use C++14
set(CMAKE_CXX_STANDARD 14)
# Use C++17
set(CMAKE_CXX_STANDARD 17)

# Require C++14
# Require C++17
set(CMAKE_CXX_STANDARD_REQUIRED ON)

# Do not enable compler-specific extensions
Expand Down
11 changes: 3 additions & 8 deletions cpp/dolfin/fem/DofMapBuilder.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -761,14 +761,9 @@ DofMapBuilder::build(const mesh::Mesh& mesh,
// Compute node ownership:
// (a) Number of owned nodes;
// (b) owned and shared nodes (and owned and un-owned):
// -1: unowned, 0: owned and shared, 1: owned and not shared;
// (c) map from shared node to sharing processes; and
std::int32_t num_owned_nodes;
std::vector<ownership> node_ownership0;
std::unordered_map<std::int32_t, std::vector<std::int32_t>>
shared_node_to_processes0;
std::set<std::int32_t> neighbouring_procs;
std::tie(num_owned_nodes, node_ownership0, shared_node_to_processes0)
// -1: unowned, 0: owned and shared, 1: owned and not shared; and
// (c) map from shared node to sharing processes.
auto [num_owned_nodes, node_ownership0, shared_node_to_processes0]
= compute_ownership(node_graph0, shared_nodes, mesh, global_dimension);

// Build re-ordering map for data locality. Owned dofs are re-ordred
Expand Down
4 changes: 1 addition & 3 deletions cpp/dolfin/function/Function.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -100,9 +100,7 @@ Function Function::sub(int i) const
Function Function::collapse() const
{
// Create new collapsed FunctionSpace
std::shared_ptr<const FunctionSpace> function_space_new;
std::vector<PetscInt> collapsed_map;
std::tie(function_space_new, collapsed_map) = _function_space->collapse();
auto [function_space_new, collapsed_map] = _function_space->collapse();

// Create new vector
assert(function_space_new);
Expand Down
10 changes: 2 additions & 8 deletions cpp/dolfin/graph/GraphBuilder.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -429,19 +429,13 @@ graph::GraphBuilder::compute_dual_graph(
{
LOG(INFO) << "Build mesh dual graph";

std::vector<std::vector<std::size_t>> local_graph;
std::int32_t num_ghost_nodes;

// Compute local part of dual graph
graph::GraphBuilder::FacetCellMap facet_cell_map;
std::int32_t num_local_edges;
std::tie(local_graph, facet_cell_map, num_local_edges)
auto [local_graph, facet_cell_map, num_local_edges]
= graph::GraphBuilder::compute_local_dual_graph(mpi_comm, cell_vertices,
cell_type);

// Compute nonlocal part
std::int32_t num_nonlocal_edges;
std::tie(num_ghost_nodes, num_nonlocal_edges) = compute_nonlocal_dual_graph(
auto [num_ghost_nodes, num_nonlocal_edges] = compute_nonlocal_dual_graph(
mpi_comm, cell_vertices, cell_type, facet_cell_map, local_graph);

// Shrink to fit
Expand Down
9 changes: 1 addition & 8 deletions cpp/dolfin/io/HDF5File.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -1366,15 +1366,8 @@ mesh::Mesh HDF5File::read_mesh(const std::string data_path,
bool use_partition_from_file,
const mesh::GhostMode ghost_mode) const
{
mesh::CellType cell_type;
Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> points;
Eigen::Array<std::int64_t, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>
cells;
std::vector<std::int64_t> global_cell_indices;
std::vector<std::int64_t> cell_distribution;

// Read local mesh data
std::tie(cell_type, points, cells, global_cell_indices, cell_distribution)
auto [cell_type, points, cells, global_cell_indices, cell_distribution]
= read_mesh_data(data_path);

if (use_partition_from_file)
Expand Down
4 changes: 1 addition & 3 deletions cpp/dolfin/io/HDF5Utility.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -289,9 +289,7 @@ void HDF5Utility::set_local_vector_values(

// Calculate one (global cell, local_dof_index) to associate with
// each item in the vector on this process
std::vector<std::size_t> global_cells;
std::vector<std::size_t> remote_local_dofi;
std::tie(global_cells, remote_local_dofi) = HDF5Utility::map_gdof_to_cell(
auto [global_cells, remote_local_dofi] = HDF5Utility::map_gdof_to_cell(
mpi_comm, cells, cell_dofs, x_cell_dofs, input_vector_range);

// At this point, each process has a set of data, and for each
Expand Down
9 changes: 1 addition & 8 deletions cpp/dolfin/io/XDMFFile.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -1458,15 +1458,8 @@ XDMFFile::read_mesh_data(MPI_Comm comm) const
//----------------------------------------------------------------------------
mesh::Mesh XDMFFile::read_mesh(const mesh::GhostMode ghost_mode) const
{

mesh::CellType cell_type;
Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> points;
Eigen::Array<std::int64_t, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>
cells;
std::vector<std::int64_t> global_cell_indices;

// Read local mesh data
std::tie(cell_type, points, cells, global_cell_indices)
auto [cell_type, points, cells, global_cell_indices]
= read_mesh_data(_mpi_comm.comm());

// Permute cells to DOLFIN ordering
Expand Down
13 changes: 5 additions & 8 deletions cpp/dolfin/mesh/DistributedMeshTools.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -402,16 +402,13 @@ void DistributedMeshTools::number_entities(const Mesh& mesh, int d)
return;
}

// Get shared entities map
std::map<std::int32_t, std::set<std::int32_t>>& shared_entities
= _mesh.topology().shared_entities(d);

// Number entities
std::vector<std::int64_t> global_entity_indices;
std::size_t num_global_entities;
std::tie(global_entity_indices, shared_entities, num_global_entities)
auto [global_entity_indices, shared_entities, num_global_entities]
= compute_entity_numbering(mesh, d);

// Set shared entities
_mesh.topology().set_shared_entities(d, shared_entities);

// Set global entity numbers in mesh
_mesh.topology().set_num_entities_global(d, num_global_entities);
_mesh.topology().set_global_indices(d, global_entity_indices);
Expand Down Expand Up @@ -439,7 +436,7 @@ void DistributedMeshTools::init_facet_cell_connections(Mesh& mesh)
Eigen::Array<std::int32_t, Eigen::Dynamic, 1> num_global_neighbors(
mesh.num_entities(D - 1));

std::map<std::int32_t, std::set<std::int32_t>>& shared_facets
const std::map<std::int32_t, std::set<std::int32_t>>& shared_facets
= mesh.topology().shared_entities(D - 1);

// Check if no ghost cells
Expand Down
25 changes: 6 additions & 19 deletions cpp/dolfin/mesh/Mesh.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -170,16 +170,11 @@ compute_point_distribution(
// Distribute points to processes that need them, and calculate
// shared points. Points are returned in same order as in global_index_set.
// Sharing information is (global_index -> [remote sharing processes]).
std::map<std::int64_t, std::set<int>> shared_points_global;
Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>
recv_points;
std::tie(shared_points_global, recv_points)
auto [shared_points_global, recv_points]
= Partitioning::distribute_points(mpi_comm, points, global_index_set);

// Get local to global mapping for points
std::vector<std::int64_t> local_to_global;
std::array<int, 4> num_vertices_local;
std::tie(local_to_global, num_vertices_local)
auto [local_to_global, num_vertices_local]
= compute_local_to_global_point_map(mpi_comm, num_vertices_per_cell,
shared_points_global, cell_nodes,
type);
Expand Down Expand Up @@ -259,16 +254,8 @@ Mesh::Mesh(

// Compute node local-to-global map from global indices, and compute
// cell topology using new local indices
std::array<int, 4> num_vertices_local;
std::vector<std::int64_t> node_indices_global;
Eigen::Array<std::int32_t, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>
coordinate_nodes;
std::map<std::int32_t, std::set<std::int32_t>> nodes_shared;
Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>
points_received;

std::tie(node_indices_global, nodes_shared, coordinate_nodes, points_received,
num_vertices_local)
auto [node_indices_global, nodes_shared, coordinate_nodes, points_received,
num_vertices_local]
= compute_point_distribution(comm, num_vertices_per_cell, cells, points,
type);

Expand Down Expand Up @@ -309,14 +296,14 @@ Mesh::Mesh(
_topology = std::make_unique<Topology>(tdim, num_vertices_local[2],
num_vertices_global);
_topology->set_global_indices(0, vertex_indices_global);
_topology->shared_entities(0) = shared_vertices;
_topology->set_shared_entities(0, shared_vertices);
_topology->init_ghost(0, num_vertices_local[1]);

// Set vertex ownership
std::vector<int> vertex_owner;
for (int i = num_vertices_local[1]; i < num_vertices_local[2]; ++i)
vertex_owner.push_back(*(shared_vertices[i].begin()));
_topology->entity_owner(0) = vertex_owner;
_topology->set_entity_owner(0, vertex_owner);

// Initialise cell topology
_topology->set_num_entities_global(tdim, num_cells_global);
Expand Down
17 changes: 6 additions & 11 deletions cpp/dolfin/mesh/Partitioning.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -552,9 +552,7 @@ PartitionData Partitioning::partition_cells(
if (mpi_comm != MPI_COMM_NULL)
{
// Compute dual graph (for this partition)
std::vector<std::vector<std::size_t>> local_graph;
std::tuple<std::int32_t, std::int32_t, std::int32_t> graph_info;
std::tie(local_graph, graph_info) = graph::GraphBuilder::compute_dual_graph(
auto [local_graph, graph_info] = graph::GraphBuilder::compute_dual_graph(
mpi_comm, cell_vertices, cell_type);

const std::size_t global_graph_size
Expand Down Expand Up @@ -685,14 +683,12 @@ mesh::Mesh Partitioning::build_from_partition(
return mesh;

// Copy cell ownership (only needed for ghost cells)
std::vector<std::int32_t>& cell_owner = mesh.topology().entity_owner(tdim);
cell_owner.clear();
cell_owner.insert(cell_owner.begin(),
new_cell_partition.begin() + num_regular_cells,
new_cell_partition.end());
std::vector<std::int32_t> cell_owner(
new_cell_partition.begin() + num_regular_cells, new_cell_partition.end());

// Assign map of shared cells (only needed for ghost cells)
mesh.topology().shared_entities(tdim) = shared_cells;
mesh.topology().set_entity_owner(tdim, cell_owner);
mesh.topology().set_shared_entities(tdim, shared_cells);
DistributedMeshTools::init_facet_cell_connections(mesh);

return mesh;
Expand Down Expand Up @@ -811,8 +807,7 @@ std::map<std::int64_t, std::vector<int>> Partitioning::compute_halo_cells(
{
// Compute dual graph (for this partition)
std::vector<std::vector<std::size_t>> local_graph;
std::tuple<std::int32_t, std::int32_t, std::int32_t> graph_info;
std::tie(local_graph, graph_info) = graph::GraphBuilder::compute_dual_graph(
std::tie(local_graph, std::ignore) = graph::GraphBuilder::compute_dual_graph(
mpi_comm, cell_vertices, cell_type);

graph::CSRGraph<std::int64_t> csr_graph(mpi_comm, local_graph);
Expand Down
11 changes: 6 additions & 5 deletions cpp/dolfin/mesh/Topology.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -103,11 +103,11 @@ const std::vector<std::int64_t>& Topology::global_indices(std::size_t d) const
return _global_indices[d];
}
//-----------------------------------------------------------------------------
std::map<std::int32_t, std::set<std::int32_t>>&
Topology::shared_entities(int dim)
void Topology::set_shared_entities(
int dim, const std::map<std::int32_t, std::set<std::int32_t>>& entities)
{
assert(dim <= this->dim());
return _shared_entities[dim];
_shared_entities[dim] = entities;
}
//-----------------------------------------------------------------------------
const std::map<std::int32_t, std::set<std::int32_t>>&
Expand All @@ -117,9 +117,10 @@ Topology::shared_entities(int dim) const
return _shared_entities[dim];
}
//-----------------------------------------------------------------------------
std::vector<std::int32_t>& Topology::entity_owner(int dim)
void Topology::set_entity_owner(int dim, const std::vector<std::int32_t>& owners)
{
return _entity_owner[dim];
assert(dim <= this->dim());
_entity_owner[dim] = owners;
}
//-----------------------------------------------------------------------------
const std::vector<std::int32_t>& Topology::entity_owner(int dim) const
Expand Down
13 changes: 7 additions & 6 deletions cpp/dolfin/mesh/Topology.h
Original file line number Diff line number Diff line change
Expand Up @@ -80,19 +80,20 @@ class Topology
/// dimension d
const std::vector<std::int64_t>& global_indices(std::size_t d) const;

/// Return map from shared entities (local index) to processes
/// that share the entity
std::map<std::int32_t, std::set<std::int32_t>>& shared_entities(int dim);
/// Set the map from shared entities (local index) to processes that
/// share the entity
void set_shared_entities(
int dim, const std::map<std::int32_t, std::set<std::int32_t>>& entities);

/// Return map from shared entities (local index) to process that
/// share the entity (const version)
const std::map<std::int32_t, std::set<std::int32_t>>&
shared_entities(int dim) const;

/// Return mapping from local ghost cell index to owning process.
/// Since ghost cells are at the end of the range, this is just a vector
/// Set map from local ghost cell index to owning process. Since
/// ghost cells are at the end of the range, this is just a vector
/// over those cells
std::vector<std::int32_t>& entity_owner(int dim);
void set_entity_owner(int dim, const std::vector<std::int32_t>& owners);

/// Return mapping from local ghost cell index to owning process
/// (const version). Since ghost cells are at the end of the range,
Expand Down
8 changes: 2 additions & 6 deletions cpp/dolfin/refinement/PlazaRefinementND.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -379,9 +379,7 @@ mesh::Mesh PlazaRefinementND::refine(const mesh::Mesh& mesh, bool redistribute)
}

common::Timer t0("PLAZA: refine");
std::vector<std::int32_t> long_edge;
std::vector<bool> edge_ratio_ok;
std::tie(long_edge, edge_ratio_ok) = face_long_edge(mesh);
auto [long_edge, edge_ratio_ok] = face_long_edge(mesh);

ParallelRefinement p_ref(mesh);
p_ref.mark_all();
Expand All @@ -402,9 +400,7 @@ PlazaRefinementND::refine(const mesh::Mesh& mesh,
}

common::Timer t0("PLAZA: refine");
std::vector<std::int32_t> long_edge;
std::vector<bool> edge_ratio_ok;
std::tie(long_edge, edge_ratio_ok) = face_long_edge(mesh);
auto [long_edge, edge_ratio_ok] = face_long_edge(mesh);

ParallelRefinement p_ref(mesh);
p_ref.mark(refinement_marker);
Expand Down
3 changes: 3 additions & 0 deletions cpp/test/unit/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -4,6 +4,9 @@ project(dolfin-tests)
find_package(DOLFIN REQUIRED)
include(${DOLFIN_USE_FILE})

set(CMAKE_CXX_STANDARD 17)
set(CMAKE_CXX_STANDARD_REQUIRED ON)

# Prepare "Catch" library for other executables
set(CATCH_INCLUDE_DIR ${CMAKE_CURRENT_SOURCE_DIR}/catch)

Expand Down
8 changes: 1 addition & 7 deletions cpp/test/unit/mesh/DistributedMesh.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -39,15 +39,9 @@ void test_distributed_mesh()
io::XDMFFile file(MPI_COMM_WORLD, "mesh.xdmf");
file.write(*mesh);

mesh::CellType cell_type;
Eigen::Array<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor> points;
Eigen::Array<std::int64_t, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>
cells;
std::vector<std::int64_t> global_cell_indices;

// Read in mesh in mesh data from XDMF file
io::XDMFFile infile(MPI_COMM_WORLD, "mesh.xdmf");
std::tie(cell_type, points, cells, global_cell_indices)
auto [cell_type, points, cells, global_cell_indices]
= infile.read_mesh_data(subset_comm);

// Partition mesh into nparts using local mesh data and subset of
Expand Down
3 changes: 1 addition & 2 deletions python/src/mesh.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -149,8 +149,7 @@ void mesh(py::module& m)
py::none());
},
py::return_value_policy::reference_internal)
.def("shared_entities",
py::overload_cast<int>(&dolfin::mesh::Topology::shared_entities))
.def("shared_entities", &dolfin::mesh::Topology::shared_entities)
.def("str", &dolfin::mesh::Topology::str);

// dolfin::mesh::Mesh
Expand Down

0 comments on commit c051360

Please sign in to comment.