From e5e41d50175b29c02b07132a5c6b1608100898a3 Mon Sep 17 00:00:00 2001 From: Stephane Del Pino <stephane.delpino44@gmail.com> Date: Thu, 23 Jul 2020 00:13:15 +0200 Subject: [PATCH] Improve integration of item id mapping when building diamond mesh 1d remains to be treated --- src/mesh/DiamondDualConnectivityBuilder.cpp | 65 +++++++++------------ 1 file changed, 28 insertions(+), 37 deletions(-) diff --git a/src/mesh/DiamondDualConnectivityBuilder.cpp b/src/mesh/DiamondDualConnectivityBuilder.cpp index 666eda055..76adb60de 100644 --- a/src/mesh/DiamondDualConnectivityBuilder.cpp +++ b/src/mesh/DiamondDualConnectivityBuilder.cpp @@ -26,17 +26,18 @@ DiamondDualConnectivityBuilder::_buildDiamondConnectivityDescriptor(const Connec const auto& primal_cell_number = primal_connectivity.cellNumber(); const size_t diamond_number_of_nodes = primal_number_of_cells + primal_number_of_nodes; + const size_t diamond_number_of_cells = primal_number_of_faces; { m_primal_node_to_dual_node_map = NodeIdToNodeIdMap{primal_number_of_nodes}; NodeId diamond_node_id = 0; - for (NodeId primal_node_id = 0; primal_node_id < primal_connectivity.numberOfNodes(); ++primal_node_id) { + for (NodeId primal_node_id = 0; primal_node_id < primal_number_of_nodes; ++primal_node_id) { m_primal_node_to_dual_node_map[primal_node_id] = std::make_pair(primal_node_id, diamond_node_id++); } m_primal_cell_to_dual_node_map = CellIdToNodeIdMap{primal_number_of_cells}; - for (CellId primal_cell_id = 0; primal_cell_id < primal_connectivity.numberOfCells(); ++primal_cell_id) { + for (CellId primal_cell_id = 0; primal_cell_id < primal_number_of_cells; ++primal_cell_id) { m_primal_cell_to_dual_node_map[primal_cell_id] = std::make_pair(primal_cell_id, diamond_node_id++); } } @@ -50,19 +51,26 @@ DiamondDualConnectivityBuilder::_buildDiamondConnectivityDescriptor(const Connec } const size_t cell_number_shift = max(primal_node_number) + 1; - for (size_t i = 0; i < m_primal_cell_to_dual_node_map.size(); ++i) { + for (size_t i = 0; i < primal_number_of_cells; ++i) { const auto [primal_cell_id, diamond_dual_node_id] = m_primal_cell_to_dual_node_map[i]; diamond_descriptor.node_number_vector[diamond_dual_node_id] = primal_cell_number[primal_cell_id] + cell_number_shift; } - const size_t diamond_number_of_cells = primal_number_of_faces; - diamond_descriptor.cell_number_vector.resize(diamond_number_of_cells); + { + m_primal_face_to_dual_cell_map = FaceIdToCellIdMap{primal_number_of_faces}; + CellId diamond_cell_id = 0; + for (FaceId primal_face_id = 0; primal_face_id < primal_number_of_faces; ++primal_face_id) { + m_primal_face_to_dual_cell_map[primal_face_id] = std::make_pair(primal_face_id, diamond_cell_id++); + } + } + diamond_descriptor.cell_number_vector.resize(diamond_number_of_cells); const auto& primal_face_number = primal_connectivity.faceNumber(); - for (FaceId i_primal_face = 0; i_primal_face < primal_number_of_faces; ++i_primal_face) { - diamond_descriptor.cell_number_vector[i_primal_face] = primal_face_number[i_primal_face]; + for (size_t i = 0; i < diamond_number_of_cells; ++i) { + const auto [primal_face_id, dual_cell_id] = m_primal_face_to_dual_cell_map[i]; + diamond_descriptor.cell_number_vector[dual_cell_id] = primal_face_number[primal_face_id]; } if constexpr (Dimension == 3) { @@ -80,7 +88,7 @@ DiamondDualConnectivityBuilder::_buildDiamondConnectivityDescriptor(const Connec const auto& primal_face_to_cell_matrix = primal_connectivity.faceToCellMatrix(); - for (FaceId i_face = 0; i_face < primal_connectivity.numberOfFaces(); ++i_face) { + for (FaceId i_face = 0; i_face < primal_number_of_faces; ++i_face) { const size_t i_cell = i_face; const auto& primal_face_cell_list = primal_face_to_cell_matrix[i_face]; @@ -106,7 +114,7 @@ DiamondDualConnectivityBuilder::_buildDiamondConnectivityDescriptor(const Connec const auto& primal_face_to_node_matrix = primal_connectivity.faceToNodeMatrix(); const auto& primal_face_local_number_in_their_cells = primal_connectivity.faceLocalNumbersInTheirCells(); const auto& cell_face_is_reversed = primal_connectivity.cellFaceIsReversed(); - for (FaceId i_face = 0; i_face < primal_connectivity.numberOfFaces(); ++i_face) { + for (FaceId i_face = 0; i_face < primal_number_of_faces; ++i_face) { const size_t& i_diamond_cell = i_face; const auto& primal_face_cell_list = primal_face_to_cell_matrix[i_face]; const auto& primal_face_node_list = primal_face_to_node_matrix[i_face]; @@ -516,44 +524,24 @@ DiamondDualConnectivityBuilder::_buildDiamondConnectivityFrom(const IConnectivit }(); m_primal_cell_to_dual_node_map = [&]() { - CellIdToNodeIdMap primal_cell_to_dual_node_map{primal_connectivity.numberOfCells()}; + CellIdToNodeIdMap primal_cell_to_dual_node_map{primal_number_of_cells}; for (CellId primal_cell_id = 0; primal_cell_id < primal_cell_to_dual_node_map.size(); ++primal_cell_id) { primal_cell_to_dual_node_map[primal_cell_id] = std::make_pair(primal_cell_id, dual_node_id++); } return primal_cell_to_dual_node_map; }(); - } else { - m_primal_node_to_dual_node_map = [&]() { - NodeIdToNodeIdMap primal_node_to_dual_node_map{primal_connectivity.numberOfNodes()}; - for (NodeId primal_node_id = 0; primal_node_id < primal_node_to_dual_node_map.size(); ++primal_node_id) { - const NodeId dual_node_id = primal_node_id; + m_primal_face_to_dual_cell_map = [&]() { + FaceIdToCellIdMap primal_face_to_dual_cell_map{primal_connectivity.numberOfFaces()}; + for (size_t id = 0; id < primal_face_to_dual_cell_map.size(); ++id) { + const CellId dual_cell_id = id; + const FaceId primal_face_id = id; - primal_node_to_dual_node_map[primal_node_id] = std::make_pair(primal_node_id, dual_node_id); + primal_face_to_dual_cell_map[id] = std::make_pair(primal_face_id, dual_cell_id); } - return primal_node_to_dual_node_map; - }(); - - m_primal_cell_to_dual_node_map = [&]() { - CellIdToNodeIdMap primal_cell_to_dual_node_map{primal_connectivity.numberOfCells()}; - NodeId dual_node_id = m_primal_node_to_dual_node_map.size(); - for (CellId primal_cell_id = 0; primal_cell_id < primal_cell_to_dual_node_map.size(); ++primal_cell_id) { - primal_cell_to_dual_node_map[primal_cell_id] = std::make_pair(primal_cell_id, dual_node_id++); - } - return primal_cell_to_dual_node_map; + return primal_face_to_dual_cell_map; }(); } - - m_primal_face_to_dual_cell_map = [&]() { - FaceIdToCellIdMap primal_face_to_dual_cell_map{primal_connectivity.numberOfFaces()}; - for (size_t id = 0; id < primal_face_to_dual_cell_map.size(); ++id) { - const CellId dual_cell_id = id; - const FaceId primal_face_id = id; - - primal_face_to_dual_cell_map[id] = std::make_pair(primal_face_id, dual_cell_id); - } - return primal_face_to_dual_cell_map; - }(); } m_mapper = @@ -567,6 +555,9 @@ DiamondDualConnectivityBuilder::_buildDiamondConnectivityFrom(const IConnectivit DiamondDualConnectivityBuilder::DiamondDualConnectivityBuilder(const IConnectivity& connectivity) { + if (parallel::size() > 1) { + throw NotImplementedError("Construction of diamond dual mesh is not implemented in parallel"); + } switch (connectivity.dimension()) { case 1: { this->_buildDiamondConnectivityFrom<1>(connectivity); -- GitLab