Skip to content
Snippets Groups Projects
Select Git revision
  • af1fa0f8303e09f0a0946698c5cae201f75c3c69
  • develop default protected
  • save_clemence
  • feature/gmsh-reader
  • origin/stage/bouguettaia
  • feature/kinetic-schemes
  • feature/reconstruction
  • feature/local-dt-fsi
  • feature/composite-scheme-sources
  • feature/composite-scheme-other-fluxes
  • feature/serraille
  • feature/variational-hydro
  • feature/composite-scheme
  • hyperplastic
  • feature/polynomials
  • feature/gks
  • feature/implicit-solver-o2
  • feature/coupling_module
  • feature/implicit-solver
  • feature/merge-local-dt-fsi
  • master protected
  • v0.5.0 protected
  • v0.4.1 protected
  • v0.4.0 protected
  • v0.3.0 protected
  • v0.2.0 protected
  • v0.1.0 protected
  • Kidder
  • v0.0.4 protected
  • v0.0.3 protected
  • v0.0.2 protected
  • v0 protected
  • v0.0.1 protected
33 results

ConnectivityDispatcher.hpp

Blame
  • ConnectivityDispatcher.cpp 27.52 KiB
    #include <ConnectivityDispatcher.hpp>
    #include <Partitioner.hpp>
    
    #include <ItemOfItemType.hpp>
    
    #include <unordered_map>
    
    template <int Dimension>
    template <ItemType item_type>
    void
    ConnectivityDispatcher<Dimension>::_buildNewOwner()
    {
      if constexpr (item_type == ItemType::cell) {
        CSRGraph connectivity_graph = m_connectivity.cellToCellGraph();
        Partitioner P;
    
        CellValue<int> cell_new_owner(m_connectivity);
        cell_new_owner = P.partition(connectivity_graph);
    
        this->_dispatchedInfo<ItemType::cell>().m_new_owner = cell_new_owner;
      } else {
        const auto& item_to_cell_matrix
            = m_connectivity.template getItemToItemMatrix<item_type,ItemType::cell>();
    
        const auto& cell_number = m_connectivity.cellNumber();
    
        const auto& cell_new_owner = this->_dispatchedInfo<ItemType::cell>().m_new_owner;
    
        using ItemId = ItemIdT<item_type>;
        ItemValue<int, item_type> item_new_owner(m_connectivity);
        parallel_for(item_new_owner.size(), PASTIS_LAMBDA(const ItemId& l) {
            const auto& item_to_cell = item_to_cell_matrix[l];
            CellId Jmin = item_to_cell[0];
    
            for (size_t j=1; j<item_to_cell.size(); ++j) {
              const CellId J = item_to_cell[j];
              if (cell_number[J] < cell_number[Jmin]) {
                Jmin=J;
              }
            }
            item_new_owner[l] = cell_new_owner[Jmin];
          });
    
    #warning Add missing synchronize (fix it!)
        // synchronize(face_new_owner);
        this->_dispatchedInfo<item_type>().m_new_owner = item_new_owner;
      }
    }
    
    template <int Dimension>
    template <ItemType item_type>
    void
    ConnectivityDispatcher<Dimension>::_buildItemListToSend()
    {
      if constexpr (item_type == ItemType::cell) {
        const auto& node_to_cell_matrix
            = m_connectivity.nodeToCellMatrix();
        const auto& cell_to_node_matrix
            = m_connectivity.cellToNodeMatrix();
    
        const auto& cell_new_owner = this->_dispatchedInfo<ItemType::cell>().m_new_owner;
    
        std::vector<std::vector<CellId>> cell_vector_to_send_by_proc(parallel::size());
        Array<bool> send_to_rank(parallel::size());
        for (CellId j=0; j<m_connectivity.numberOfCells(); ++j) {
          send_to_rank.fill(false);
          const auto& cell_to_node = cell_to_node_matrix[j];
    
          for (size_t R=0; R<cell_to_node.size(); ++R) {
            const NodeId& r = cell_to_node[R];
            const auto& node_to_cell = node_to_cell_matrix[r];
            for (size_t K=0; K<node_to_cell.size(); ++K) {
              const CellId& k = node_to_cell[K];
              send_to_rank[cell_new_owner[k]] = true;
            }
          }
    
          for (size_t k=0; k<send_to_rank.size(); ++k) {
            if (send_to_rank[k]) {
              cell_vector_to_send_by_proc[k].push_back(j);
            }
          }
        }
    
        auto& cell_list_to_send_by_proc = this->_dispatchedInfo<ItemType::cell>().m_list_to_send_by_proc;
        cell_list_to_send_by_proc.resize(parallel::size());
        for (size_t i=0; i<parallel::size(); ++i) {
          cell_list_to_send_by_proc[i] = convert_to_array(cell_vector_to_send_by_proc[i]);
        }
      } else {
        const auto& cell_list_to_send_by_proc = this->_dispatchedInfo<ItemType::cell>().m_list_to_send_by_proc;
    
        using ItemId = ItemIdT<item_type>;
        const auto& cell_to_sub_item_matrix = m_connectivity.template getItemToItemMatrix<ItemType::cell,item_type>();
    
        auto& item_list_to_send_by_proc = this->_dispatchedInfo<item_type>().m_list_to_send_by_proc;
        item_list_to_send_by_proc.resize(parallel::size());
    
        for (size_t i_rank=0; i_rank<parallel::size(); ++i_rank) {
          Array<bool> tag(m_connectivity.template numberOf<item_type>());
          tag.fill(false);
          std::vector<ItemId> item_id_vector;
          for (size_t j=0; j<cell_list_to_send_by_proc[i_rank].size(); ++j) {
            const CellId& cell_id = cell_list_to_send_by_proc[i_rank][j];
            const auto& cell_sub_item_list = cell_to_sub_item_matrix[cell_id];
            for (size_t r=0; r<cell_sub_item_list.size(); ++r) {
              const ItemId& item_id = cell_sub_item_list[r];
              if (not tag[item_id]) {
                item_id_vector.push_back(item_id);
                tag[item_id] = true;
              }
            }
          }
          item_list_to_send_by_proc[i_rank] = convert_to_array(item_id_vector);
        }
      }
    }
    
    template <int Dimension>
    Array<const unsigned int>
    ConnectivityDispatcher<Dimension>::_buildNbCellToSend()
    {
      const auto& item_list_to_send_by_proc = this->_dispatchedInfo<ItemType::cell>().m_list_to_send_by_proc;
      Array<unsigned int> nb_cell_to_send_by_proc(parallel::size());
      for (size_t i=0; i<parallel::size(); ++i) {
        nb_cell_to_send_by_proc[i] = item_list_to_send_by_proc[i].size();
      }
      return nb_cell_to_send_by_proc;
    }
    
    template <int Dimension>
    template<typename DataType, ItemType item_type, typename ConnectivityPtr>
    void
    ConnectivityDispatcher<Dimension>::
    _gatherFrom(const ItemValue<DataType, item_type, ConnectivityPtr>& data_to_gather,
                std::vector<std::remove_const_t<DataType>>& gathered_vector)
    {
      std::vector<Array<const DataType>> recv_item_data_by_proc = this->exchange(data_to_gather);
    
      const auto& recv_id_correspondance_by_proc = this->_dispatchedInfo<item_type>().m_recv_id_correspondance_by_proc;
      Assert(recv_id_correspondance_by_proc.size()==parallel::size());
    
      gathered_vector.resize(this->_dispatchedInfo<item_type>().m_number_to_id_map.size());
      for (size_t i_rank=0; i_rank<parallel::size(); ++i_rank) {
        Assert(recv_id_correspondance_by_proc[i_rank].size()==recv_item_data_by_proc[i_rank].size());
        for (size_t r=0; r<recv_id_correspondance_by_proc[i_rank].size(); ++r) {
          const auto& item_id = recv_id_correspondance_by_proc[i_rank][r];
          gathered_vector[item_id] = recv_item_data_by_proc[i_rank][r];
        }
      }
    }
    
    template <int Dimension>
    void
    ConnectivityDispatcher<Dimension>::
    _buildCellNumberIdMap()
    {
      const auto recv_cell_number_by_proc = this->exchange(m_connectivity.template number<ItemType::cell>());
      auto& cell_number_id_map = this->_dispatchedInfo<ItemType::cell>().m_number_to_id_map;
      for (size_t i_rank=0; i_rank<parallel::size(); ++i_rank) {
        CellId cell_id=0;
        for (size_t i=0; i<recv_cell_number_by_proc[i_rank].size(); ++i) {
          const int cell_number = recv_cell_number_by_proc[i_rank][i];
          auto [iterator, inserted] = cell_number_id_map.insert(std::make_pair(cell_number, cell_id));
          if (inserted) ++cell_id;
        }
      }
    }
    
    template <int Dimension>
    template <typename ItemOfItemT>
    void
    ConnectivityDispatcher<Dimension>::
    _buildSubItemNumberToIdMap(const std::vector<Array<const int>>& recv_cell_sub_item_number_by_proc)
    {
      static_assert(ItemOfItemT::item_type == ItemType::cell, "Dispatcher requires to be build using cell as master entities");
    
      auto& sub_item_number_id_map = this->_dispatchedInfo<ItemOfItemT::sub_item_type>().m_number_to_id_map;
      for (size_t i_rank=0; i_rank<parallel::size(); ++i_rank) {
        int sub_item_id=0;
        for (size_t i=0; i<recv_cell_sub_item_number_by_proc[i_rank].size(); ++i) {
          int sub_item_number = recv_cell_sub_item_number_by_proc[i_rank][i];
          auto [iterator, inserted] = sub_item_number_id_map.insert(std::make_pair(sub_item_number, sub_item_id));
          if (inserted) sub_item_id++;
        }
      }
    }
    
    template <int Dimension>
    template <typename SubItemOfItemT>
    std::vector<Array<const int>>
    ConnectivityDispatcher<Dimension>::_getRecvNumberOfSubItemPerItemByProc()
    {
      const auto& item_to_sub_item_matrix
          = m_connectivity.template getItemToItemMatrix<SubItemOfItemT::item_type,
                                                        SubItemOfItemT::sub_item_type>();
    
      ItemValue<int, SubItemOfItemT::item_type> number_of_sub_item_per_item(m_connectivity);
    
      using ItemId = ItemIdT<SubItemOfItemT::item_type>;
      parallel_for(number_of_sub_item_per_item.size(), PASTIS_LAMBDA(const ItemId& j){
          number_of_sub_item_per_item[j] = item_to_sub_item_matrix[j].size();
        });
      return this->exchange(number_of_sub_item_per_item);
    }
    
    template <int Dimension>
    template <typename SubItemOfItemT>
    std::vector<Array<const int>>
    ConnectivityDispatcher<Dimension>::
    _getRecvItemSubItemNumberingByProc(const std::vector<Array<const int>>& recv_number_of_sub_item_per_item_by_proc)
    {
      std::vector<Array<const int>> item_sub_item_numbering_to_send_by_proc =
          [&] () {
            const auto& item_to_sub_item_matrix
                = m_connectivity.template getItemToItemMatrix<SubItemOfItemT::item_type,
                                                              SubItemOfItemT::sub_item_type>();
    
            const ItemValue<const int, SubItemOfItemT::sub_item_type>& sub_item_number =
                m_connectivity.template number<SubItemOfItemT::sub_item_type>();
    
            using ItemId = ItemIdT<SubItemOfItemT::item_type>;
            using SubItemId = ItemIdT<SubItemOfItemT::sub_item_type>;
    
            std::vector<Array<const int>> item_sub_item_numbering_to_send_by_proc(parallel::size());
            for (size_t i_rank=0; i_rank < parallel::size(); ++i_rank) {
              const auto& item_list_to_send_by_proc
                  = this->_dispatchedInfo<SubItemOfItemT::item_type>().m_list_to_send_by_proc;
    
              std::vector<int> sub_item_numbering_by_item_vector;
              for (size_t j=0; j<item_list_to_send_by_proc[i_rank].size(); ++j) {
                const ItemId& item_id = item_list_to_send_by_proc[i_rank][j];
                const auto& sub_item_list = item_to_sub_item_matrix[item_id];
                for (size_t r=0; r<sub_item_list.size(); ++r) {
                  const SubItemId& sub_item_id = sub_item_list[r];
                  sub_item_numbering_by_item_vector.push_back(sub_item_number[sub_item_id]);
                }
              }
              item_sub_item_numbering_to_send_by_proc[i_rank] = convert_to_array(sub_item_numbering_by_item_vector);
            }
            return item_sub_item_numbering_to_send_by_proc;
          } ();
    
      std::vector<Array<int>> recv_item_sub_item_numbering_by_proc(parallel::size());
      for (size_t i_rank=0; i_rank < parallel::size(); ++i_rank) {
        recv_item_sub_item_numbering_by_proc[i_rank]
            = Array<int>(sum(recv_number_of_sub_item_per_item_by_proc[i_rank]));
      }
      parallel::exchange(item_sub_item_numbering_to_send_by_proc, recv_item_sub_item_numbering_by_proc);
    
      std::vector<Array<const int>> const_recv_item_sub_item_numbering_by_proc(parallel::size());
      for (size_t i_rank=0; i_rank < parallel::size(); ++i_rank) {
        const_recv_item_sub_item_numbering_by_proc[i_rank] = recv_item_sub_item_numbering_by_proc[i_rank];
      }
      return const_recv_item_sub_item_numbering_by_proc;
    }
    
    template <int Dimension>
    template <ItemType item_type>
    void
    ConnectivityDispatcher<Dimension>::
    _buildRecvItemIdCorrespondanceByProc()
    {
      const auto& item_list_to_send_by_proc = this->_dispatchedInfo<item_type>().m_list_to_send_by_proc;
      using ItemId = ItemIdT<item_type>;
    
      std::vector<Array<const ItemId>> recv_item_id_correspondance_by_proc(parallel::size());
      const ItemValue<const int,item_type>& item_number =
          m_connectivity.template number<item_type>();
    
      std::vector<Array<const int>> send_item_number_by_proc(parallel::size());
      for (size_t i_rank=0; i_rank<parallel::size(); ++i_rank) {
        Array<int> send_item_number(item_list_to_send_by_proc[i_rank].size());
        const Array<const ItemId> send_item_id = item_list_to_send_by_proc[i_rank];
        parallel_for(send_item_number.size(), PASTIS_LAMBDA(const size_t& j){
            send_item_number[j] = item_number[send_item_id[j]];
          });
        send_item_number_by_proc[i_rank] = send_item_number;
      }
    
      const auto& item_list_to_recv_size_by_proc = this->_dispatchedInfo<item_type>().m_list_to_recv_size_by_proc;
      std::vector<Array<int>> recv_item_number_by_proc(parallel::size());
      for (size_t i_rank=0; i_rank < parallel::size(); ++i_rank) {
        recv_item_number_by_proc[i_rank] = Array<int>(item_list_to_recv_size_by_proc[i_rank]);
      }
      parallel::exchange(send_item_number_by_proc, recv_item_number_by_proc);
    
      const auto& item_number_to_id_map = this->_dispatchedInfo<item_type>().m_number_to_id_map;
      for (size_t i_rank=0; i_rank<item_list_to_recv_size_by_proc.size(); ++i_rank) {
        Array<ItemId> item_id_correspondance(item_list_to_recv_size_by_proc[i_rank]);
        for (size_t l=0; l<item_list_to_recv_size_by_proc[i_rank]; ++l) {
          const int& item_number = recv_item_number_by_proc[i_rank][l];
          const auto& searched_item_id = item_number_to_id_map.find(item_number);
          Assert(searched_item_id != item_number_to_id_map.end());
          item_id_correspondance[l] = searched_item_id->second;
        }
        recv_item_id_correspondance_by_proc[i_rank] = item_id_correspondance;
      }
      this->_dispatchedInfo<item_type>().m_recv_id_correspondance_by_proc = recv_item_id_correspondance_by_proc;
    }
    
    template <int Dimension>
    void
    ConnectivityDispatcher<Dimension>::_dispatchFaces()
    {
      if constexpr (Dimension>1) {
        std::vector<Array<const int>> recv_number_of_face_per_cell_by_proc =
            _getRecvNumberOfSubItemPerItemByProc<FaceOfCell>();
    
        std::vector<Array<const int>> recv_cell_face_numbering_by_proc
          = this->_getRecvItemSubItemNumberingByProc<FaceOfCell>(recv_number_of_face_per_cell_by_proc);
    
        this->_buildSubItemNumberToIdMap<FaceOfCell>(recv_cell_face_numbering_by_proc);
    
        const std::unordered_map<int, int>& face_number_id_map =
            this->_dispatchedInfo<ItemType::face>().m_number_to_id_map;
    
        this->_buildItemListToSend<ItemType::face>();
    
    #warning this patterns repeats for each item type. Should be factorized
        Array<unsigned int> nb_face_to_send_by_proc(parallel::size());
        for (size_t i_rank=0; i_rank<parallel::size(); ++i_rank) {
          nb_face_to_send_by_proc[i_rank] = m_dispatched_face_info.m_list_to_send_by_proc[i_rank].size();
        }
        this->_dispatchedInfo<ItemType::face>().m_list_to_recv_size_by_proc
            = parallel::allToAll(nb_face_to_send_by_proc);
    
        this->_buildRecvItemIdCorrespondanceByProc<ItemType::face>();
        this->_gatherFrom(m_connectivity.template number<ItemType::face>(), m_new_descriptor.face_number_vector);
    
        {
          const auto& cell_list_to_recv_size_by_proc =
              this->_dispatchedInfo<ItemType::cell>().m_list_to_recv_size_by_proc;
    
          for (size_t i_rank=0; i_rank < parallel::size(); ++i_rank) {
            int l=0;
            for (size_t i=0; i<cell_list_to_recv_size_by_proc[i_rank]; ++i) {
              std::vector<unsigned int> face_vector;
              for (int k=0; k<recv_number_of_face_per_cell_by_proc[i_rank][i]; ++k) {
                const auto& searched_face_id = face_number_id_map.find(recv_cell_face_numbering_by_proc[i_rank][l++]);
                Assert(searched_face_id != face_number_id_map.end());
                face_vector.push_back(searched_face_id->second);
              }
              m_new_descriptor.cell_to_face_vector.emplace_back(face_vector);
            }
          }
        }
    
        {
          std::vector<Array<bool>> cell_face_is_reversed_to_send_by_proc(parallel::size());
          {
            const auto& cell_list_to_send_by_proc = this->_dispatchedInfo<ItemType::cell>().m_list_to_send_by_proc;
    
            const auto& cell_face_is_reversed = m_connectivity.cellFaceIsReversed();
            for (size_t i_rank=0; i_rank < parallel::size(); ++i_rank) {
              std::vector<bool> face_is_reversed_by_cell_vector;
              for (size_t j=0; j<cell_list_to_send_by_proc[i_rank].size(); ++j) {
                const CellId& cell_id = cell_list_to_send_by_proc[i_rank][j];
                const auto& face_is_reversed = cell_face_is_reversed.itemValues(cell_id);
                for (size_t L=0; L<face_is_reversed.size(); ++L) {
                  face_is_reversed_by_cell_vector.push_back(face_is_reversed[L]);
                }
              }
              cell_face_is_reversed_to_send_by_proc[i_rank] = convert_to_array(face_is_reversed_by_cell_vector);
            }
          }
    
          std::vector<Array<bool>> recv_cell_face_is_reversed_by_proc(parallel::size());
          for (size_t i_rank=0; i_rank < parallel::size(); ++i_rank) {
            recv_cell_face_is_reversed_by_proc[i_rank]
                = Array<bool>(sum(recv_number_of_face_per_cell_by_proc[i_rank]));
          }
    
          parallel::exchange(cell_face_is_reversed_to_send_by_proc, recv_cell_face_is_reversed_by_proc);
    
          const auto& cell_list_to_recv_size_by_proc =
              this->_dispatchedInfo<ItemType::cell>().m_list_to_recv_size_by_proc;
    
          for (size_t i_rank=0; i_rank < parallel::size(); ++i_rank) {
            int l=0;
            for (size_t i=0; i<cell_list_to_recv_size_by_proc[i_rank]; ++i) {
              std::vector<bool> face_is_reversed_vector;
              for (int k=0; k<recv_number_of_face_per_cell_by_proc[i_rank][i]; ++k) {
                face_is_reversed_vector.push_back(recv_cell_face_is_reversed_by_proc[i_rank][l++]);
              }
              m_new_descriptor.cell_face_is_reversed_vector.emplace_back(face_is_reversed_vector);
            }
          }
        }
    
        this->_gatherFrom(this->_dispatchedInfo<ItemType::face>().m_new_owner, m_new_descriptor.face_owner_vector);
    
        std::vector<Array<const int>> recv_number_of_node_per_face_by_proc =
            _getRecvNumberOfSubItemPerItemByProc<NodeOfFace>();
    
        std::vector<Array<const int>> recv_face_node_numbering_by_proc
            = this->_getRecvItemSubItemNumberingByProc<NodeOfFace>(recv_number_of_node_per_face_by_proc);
    
        {
          const auto& node_number_id_map = this->_dispatchedInfo<ItemType::node>().m_number_to_id_map;
          for (size_t i_rank=0; i_rank < parallel::size(); ++i_rank) {
            int l=0;
            for (size_t i=0; i<recv_number_of_node_per_face_by_proc[i_rank].size(); ++i) {
              std::vector<unsigned int> node_vector;
              for (int k=0; k<recv_number_of_node_per_face_by_proc[i_rank][i]; ++k) {
                const auto& searched_node_id = node_number_id_map.find(recv_face_node_numbering_by_proc[i_rank][l++]);
                Assert(searched_node_id != node_number_id_map.end());
                node_vector.push_back(searched_node_id->second);
              }
              m_new_descriptor.face_to_node_vector.emplace_back(node_vector);
            }
          }
        }
    
        // Getting references
        Array<const size_t> number_of_ref_face_list_per_proc
            = parallel::allGather(m_connectivity.numberOfRefFaceList());
    
        const size_t number_of_face_list_sender
            = [&] () {
                size_t number_of_face_list_sender=0;
                for (size_t i_rank = 0; i_rank < parallel::size(); ++i_rank) {
                  number_of_face_list_sender
                      += (number_of_ref_face_list_per_proc[i_rank] > 0);
                }
                return number_of_face_list_sender;
              }();
    
        if (number_of_face_list_sender > 0) {
          if (number_of_face_list_sender > 1) {
            perr() << __FILE__ << ':' << __LINE__ << ": "
                   << rang::fgB::red
                   <<"need to check that knowing procs know the same ref_face_lists!"
                   << rang::fg::reset << '\n';
          }
    
          if (number_of_face_list_sender < parallel::size()) {
            const size_t sender_rank
                = [&] () {
                    size_t i_rank = 0;
                    for (; i_rank < parallel::size(); ++i_rank) {
                      if (number_of_ref_face_list_per_proc[i_rank] > 0) {
                        break;
                      }
                    }
                    return i_rank;
                  }();
    
            Assert(number_of_face_list_sender < parallel::size());
    
            // sending references tags
            Array<RefId::TagNumberType> ref_tag_list{number_of_ref_face_list_per_proc[sender_rank]};
            if (parallel::rank() == sender_rank){
              for (size_t i_ref_face_list=0; i_ref_face_list<m_connectivity.numberOfRefFaceList();
                   ++i_ref_face_list) {
                auto ref_face_list = m_connectivity.refFaceList(i_ref_face_list);
                ref_tag_list[i_ref_face_list] = ref_face_list.refId().tagNumber();
              }
            }
            parallel::broadcast(ref_tag_list, sender_rank);
    
            // sending references name size
            Array<size_t> ref_name_size_list{number_of_ref_face_list_per_proc[sender_rank]};
            if (parallel::rank() == sender_rank){
              for (size_t i_ref_face_list=0; i_ref_face_list<m_connectivity.numberOfRefFaceList();
                   ++i_ref_face_list) {
                auto ref_face_list = m_connectivity.refFaceList(i_ref_face_list);
                ref_name_size_list[i_ref_face_list] = ref_face_list.refId().tagName().size();
              }
            }
            parallel::broadcast(ref_name_size_list, sender_rank);
    
            // sending references name size
            Array<RefId::TagNameType::value_type> ref_name_cat{sum(ref_name_size_list)};
            if (parallel::rank() == sender_rank){
              size_t i_char=0;
              for (size_t i_ref_face_list=0; i_ref_face_list<m_connectivity.numberOfRefFaceList();
                   ++i_ref_face_list) {
                auto ref_face_list = m_connectivity.refFaceList(i_ref_face_list);
                for (auto c : ref_face_list.refId().tagName()) {
                  ref_name_cat[i_char++] = c;
                }
              }
            }
            parallel::broadcast(ref_name_cat, sender_rank);
    
            std::vector<RefId> ref_id_list
                = [&] () {
                    std::vector<RefId> ref_id_list;
                    ref_id_list.reserve(ref_name_size_list.size());
                    size_t begining=0;
                    for (size_t i_ref=0; i_ref < ref_name_size_list.size(); ++i_ref) {
                      const size_t size = ref_name_size_list[i_ref];
                      ref_id_list.emplace_back(ref_tag_list[i_ref],
                                               std::string{&(ref_name_cat[begining]), size});
                      begining += size;
                    }
                    return ref_id_list;
                  } ();
    
            using block_type = int32_t;
            constexpr size_t block_size = sizeof(block_type);
            const size_t nb_block = ref_id_list.size()/block_size + (ref_id_list.size()%block_size != 0);
            for (size_t i_block=0; i_block<nb_block; ++i_block) {
              FaceValue<block_type> face_references(m_connectivity);
              face_references.fill(0);
    
              if (m_connectivity.numberOfRefFaceList() > 0) {
                const size_t max_i_ref = std::min(ref_id_list.size(), block_size*(i_block+1));
                for (size_t i_ref=block_size*i_block, i=0; i_ref<max_i_ref; ++i_ref, ++i) {
                  block_type ref_bit{1<<i};
                  auto ref_face_list = m_connectivity.refFaceList(i_ref);
    
                  const auto& face_list = ref_face_list.faceList();
                  for (size_t i_face=0; i_face<face_list.size(); ++i_face) {
                    const FaceId& face_id = face_list[i_face];
                    face_references[face_id] |= ref_bit;
                  }
                }
              }
    
              const auto& send_face_id_by_proc = m_dispatched_face_info.m_list_to_send_by_proc;
              std::vector<Array<const block_type>> send_face_refs_by_proc(parallel::size());
              for (size_t i_rank=0; i_rank<parallel::size(); ++i_rank) {
                Array<block_type> send_face_refs(nb_face_to_send_by_proc[i_rank]);
                const Array<const FaceId> send_face_id = send_face_id_by_proc[i_rank];
                parallel_for(send_face_id.size(), PASTIS_LAMBDA(const size_t& l) {
                    const FaceId& face_id = send_face_id[l];
                    send_face_refs[l] = face_references[face_id];
                  });
                send_face_refs_by_proc[i_rank] = send_face_refs;
              }
    
              std::vector<Array<block_type>> recv_face_refs_by_proc(parallel::size());
              for (size_t i_rank=0; i_rank<parallel::size(); ++i_rank) {
                recv_face_refs_by_proc[i_rank] = Array<block_type>(m_dispatched_face_info.m_list_to_recv_size_by_proc[i_rank]);
              }
              parallel::exchange(send_face_refs_by_proc, recv_face_refs_by_proc);
    
              const auto& recv_face_id_correspondance_by_proc = m_dispatched_face_info.m_recv_id_correspondance_by_proc;
              std::vector<block_type> face_refs(m_new_descriptor.face_number_vector.size());
              for (size_t i_rank=0; i_rank<parallel::size(); ++i_rank) {
                for (size_t r=0; r<recv_face_refs_by_proc[i_rank].size(); ++r) {
                  const FaceId& face_id = recv_face_id_correspondance_by_proc[i_rank][r];
                  face_refs[face_id] = recv_face_refs_by_proc[i_rank][r];
                }
              }
    
              const size_t max_i_ref = std::min(ref_id_list.size(), block_size*(i_block+1));
              for (size_t i_ref=block_size*i_block, i=0; i_ref<max_i_ref; ++i_ref, ++i) {
                block_type ref_bit{1<<i};
    
                std::vector<FaceId> face_id_vector;
    
                for (uint32_t i_face=0; i_face<face_refs.size(); ++i_face) {
                  const FaceId face_id{i_face};
                  if (face_refs[face_id] & ref_bit) {
                    face_id_vector.push_back(face_id);
                  }
                }
    
                Array<const FaceId> face_id_array = convert_to_array(face_id_vector);
    
                m_new_descriptor.addRefFaceList(RefFaceList(ref_id_list[i_ref], face_id_array));
              }
    
              pout() << __FILE__ << ':' << __LINE__ << ": remains to build lists\n";
            }
          }
        }
      }
    }
    
    
    template <int Dimension>
    ConnectivityDispatcher<Dimension>::ConnectivityDispatcher(const ConnectivityType& connectivity)
        : m_connectivity(connectivity)
    {
      this->_buildNewOwner<ItemType::cell>();
      this->_buildNewOwner<ItemType::face>();
      // this->_buildNewOwner<ItemType::edge>();
      this->_buildNewOwner<ItemType::node>();
    
      this->_buildItemListToSend<ItemType::cell>();
      this->_dispatchedInfo<ItemType::cell>().m_list_to_recv_size_by_proc
          = parallel::allToAll(this->_buildNbCellToSend());
    
      this->_buildCellNumberIdMap();
    
      const std::vector<Array<const int>> recv_number_of_node_per_cell_by_proc
          = this->_getRecvNumberOfSubItemPerItemByProc<NodeOfCell>();
    
      const std::vector<Array<const int>> recv_cell_node_numbering_by_proc
          = this->_getRecvItemSubItemNumberingByProc<NodeOfCell>(recv_number_of_node_per_cell_by_proc);
    
      this->_buildRecvItemIdCorrespondanceByProc<ItemType::cell>();
      this->_gatherFrom(m_connectivity.template number<ItemType::cell>(), m_new_descriptor.cell_number_vector);
    
      this->_buildSubItemNumberToIdMap<NodeOfCell>(recv_cell_node_numbering_by_proc);
      this->_buildItemListToSend<ItemType::node>();
    
      Array<unsigned int> nb_node_to_send_by_proc(parallel::size());
      for (size_t i_rank=0; i_rank<parallel::size(); ++i_rank) {
        nb_node_to_send_by_proc[i_rank] = m_dispatched_node_info.m_list_to_send_by_proc[i_rank].size();
      }
      this->_dispatchedInfo<ItemType::node>().m_list_to_recv_size_by_proc
          = parallel::allToAll(nb_node_to_send_by_proc);
    
      this->_buildRecvItemIdCorrespondanceByProc<ItemType::node>();
    
      // Fill new descriptor
      this->_gatherFrom(m_connectivity.cellType(), m_new_descriptor.cell_type_vector);
      this->_gatherFrom(this->_dispatchedInfo<ItemType::cell>().m_new_owner, m_new_descriptor.cell_owner_vector);
    
      this->_gatherFrom(m_connectivity.template number<ItemType::node>(), m_new_descriptor.node_number_vector);
      this->_gatherFrom(this->_dispatchedInfo<ItemType::node>().m_new_owner, m_new_descriptor.node_owner_vector);
    
      { // build cells connectivity
        const auto& cell_list_to_recv_size_by_proc =
            this->_dispatchedInfo<ItemType::cell>().m_list_to_recv_size_by_proc;
    
        const auto& node_number_id_map  = this->_dispatchedInfo<ItemType::node>().m_number_to_id_map;
        for (size_t i_rank=0; i_rank < parallel::size(); ++i_rank) {
          int l=0;
          for (size_t i=0; i<cell_list_to_recv_size_by_proc[i_rank]; ++i) {
            std::vector<unsigned int> node_vector;
            for (int k=0; k<recv_number_of_node_per_cell_by_proc[i_rank][i]; ++k) {
              const auto& searched_node_id = node_number_id_map.find(recv_cell_node_numbering_by_proc[i_rank][l++]);
              Assert(searched_node_id != node_number_id_map.end());
              node_vector.push_back(searched_node_id->second);
            }
            m_new_descriptor.cell_by_node_vector.emplace_back(node_vector);
          }
        }
      }
      this->_dispatchFaces();
    
      m_dispatched_connectivity = ConnectivityType::build(m_new_descriptor);
    }
    
    template ConnectivityDispatcher<1>::ConnectivityDispatcher(const ConnectivityType&);
    template ConnectivityDispatcher<2>::ConnectivityDispatcher(const ConnectivityType&);
    template ConnectivityDispatcher<3>::ConnectivityDispatcher(const ConnectivityType&);