diff --git a/src/main.cpp b/src/main.cpp
index 762cfb573f01dad21d1688d706591f2fc8941a25..2bc9dd2d294d7105a8ef8050386d58af4f57e03f 100644
--- a/src/main.cpp
+++ b/src/main.cpp
@@ -224,7 +224,8 @@ int main(int argc, char *argv[])
 
           std::vector<BoundaryConditionHandler> bc_list;
           {
-            for (size_t i_ref_face_list=0; i_ref_face_list<mesh.connectivity().numberOfRefFaceList(); ++i_ref_face_list) {
+            for (size_t i_ref_face_list=0; i_ref_face_list<mesh.connectivity().numberOfRefFaceList();
+                 ++i_ref_face_list) {
               const RefFaceList& ref_face_list = mesh.connectivity().refFaceList(i_ref_face_list);
               const RefId& ref = ref_face_list.refId();
 
@@ -238,10 +239,31 @@ int main(int argc, char *argv[])
                 break;
               }
 
-              const Kokkos::View<const unsigned int*> nodes_ids = ref_face_list.nodeList();
-              std::vector<unsigned int> node_boundary_vector(nodes_ids.extent(0));
-              for (size_t r=0; r<nodes_ids.extent(0); ++r) {
-                node_boundary_vector[r] = nodes_ids[r];
+              std::set<unsigned int> node_id_set;
+              const Kokkos::View<const unsigned int**> face_nodes = mesh.connectivity().faceNodes();
+
+              const Kokkos::View<const unsigned int*>& face_id_list = ref_face_list.faceList();
+
+              for (unsigned int l=0; l<face_id_list.extent(0); ++l) {
+                for (unsigned short r=0; r<2; ++r) {
+                  node_id_set.insert(face_nodes(face_id_list[l],r));
+                }
+              }
+
+              Kokkos::View<unsigned int*> node_id_list("node_id_list", node_id_set.size());
+              {
+                int r=0;
+                for (auto node_id : node_id_set) {
+                  node_id_list[r] = node_id;
+                  ++r;
+                }
+              }
+
+              RefNodeList ref_node_list(ref, node_id_list);
+
+              std::vector<unsigned int> node_boundary_vector(node_id_list.extent(0));
+              for (size_t r=0; r<node_id_list.extent(0); ++r) {
+                node_boundary_vector[r] = node_id_list[r];
               }
               SymmetryBoundaryCondition<MeshType::dimension>* sym_bc
                   = new SymmetryBoundaryCondition<MeshType::dimension>(node_boundary_vector, normal);
diff --git a/src/mesh/Connectivity2D.hpp b/src/mesh/Connectivity2D.hpp
index 09a164906ecd547ed4fec02b04a084f817b4dcf4..290e826aa8cd51cd39b68f2a795978ec79a35ca1 100644
--- a/src/mesh/Connectivity2D.hpp
+++ b/src/mesh/Connectivity2D.hpp
@@ -11,6 +11,7 @@
 #include <algorithm>
 
 #include <RefId.hpp>
+#include <RefNodeList.hpp>
 #include <RefFaceList.hpp>
 
 class Connectivity2D
@@ -19,6 +20,7 @@ public:
   static constexpr size_t dimension = 2;
 
   std::vector<RefFaceList> m_ref_face_list;
+  std::vector<RefNodeList> m_ref_node_list;
 
 private:
   const size_t  m_number_of_cells;
@@ -194,6 +196,21 @@ private:
     return m_ref_face_list[i];
   }
 
+  void addRefNodeList(const RefNodeList& ref_node_list)
+  {
+    m_ref_node_list.push_back(ref_node_list);
+  }
+
+  size_t numberOfRefNodeList() const
+  {
+    return m_ref_node_list.size();
+  }
+
+  const RefNodeList& refNodeList(const size_t& i) const
+  {
+    return m_ref_node_list[i];
+  }
+
   const size_t& numberOfNodes() const
   {
     return m_number_of_nodes;
diff --git a/src/mesh/GmshReader.cpp b/src/mesh/GmshReader.cpp
index b658806de4a1f3b86222c9f29edaa7259806aadb..4eb7dde8099d0a17efb79b655e751bc13b78edc0 100644
--- a/src/mesh/GmshReader.cpp
+++ b/src/mesh/GmshReader.cpp
@@ -838,7 +838,7 @@ GmshReader::__proceedData()
         face_list[j]=ref_face_list.second[j];
       }
       const PhysicalRefId& physical_ref_id = m_physical_ref_map.at(ref_face_list.first);
-      connectivity.addRefFaceList(RefFaceList(connectivity, physical_ref_id.refId(), face_list));
+      connectivity.addRefFaceList(RefFaceList(physical_ref_id.refId(), face_list));
     }
 
     const Kokkos::View<const unsigned int**> face_nodes = connectivity.faceNodes();
@@ -856,7 +856,7 @@ GmshReader::__proceedData()
         point_list[j]=ref_point_list.second[j];
       }
       const PhysicalRefId& physical_ref_id = m_physical_ref_map.at(ref_point_list.first);
-      std::cout << rang::fgB::red << "point list (corners?) not handled anymore" << rang::style::reset << '\n';
+      connectivity.addRefNodeList(RefNodeList(physical_ref_id.refId(), point_list));
     }
 
     typedef Mesh<Connectivity2D> MeshType;
diff --git a/src/mesh/RefFaceList.hpp b/src/mesh/RefFaceList.hpp
index 95223cdb28eaee16c4b6e3b27a5ee21d197866cf..e891ea2581bc1dbb19a0e82f4178a2a172112574 100644
--- a/src/mesh/RefFaceList.hpp
+++ b/src/mesh/RefFaceList.hpp
@@ -4,14 +4,11 @@
 #include <Kokkos_Core.hpp>
 #include <RefId.hpp>
 
-#include <set>
-
 class RefFaceList
 {
  private:
-  RefId m_ref_id;
-  Kokkos::View<const unsigned int*> m_face_id_list;
-  Kokkos::View<const unsigned int*> m_node_id_list;
+  const RefId m_ref_id;
+  const Kokkos::View<const unsigned int*> m_face_id_list;
 
  public:
   const RefId& refId() const
@@ -24,37 +21,12 @@ class RefFaceList
     return m_face_id_list;
   }
 
-  [[deprecated("referenced face list should not know about nodes")]]
-  const Kokkos::View<const unsigned int*>&  nodeList() const
-  {
-    return m_node_id_list;
-  }
-
-  template<typename ConnectivityType>
-  RefFaceList(const ConnectivityType& connectivity,
-           const RefId& ref_id,
-           const Kokkos::View<const unsigned int*>& face_id_list)
+  RefFaceList(const RefId& ref_id,
+              const Kokkos::View<const unsigned int*>& face_id_list)
       : m_ref_id(ref_id),
         m_face_id_list(face_id_list)
   {
-    std::set<unsigned int> node_id_set;
-    const Kokkos::View<const unsigned int**> face_nodes = connectivity.faceNodes();
-
-    for (unsigned int l=0; l<m_face_id_list.extent(0); ++l) {
-      for (unsigned short r=0; r<2; ++r) {
-        node_id_set.insert(face_nodes(m_face_id_list[l],r));
-      }
-    }
-
-    Kokkos::View<unsigned int*> node_id_list("node_id_list", node_id_set.size());
-    {
-      int r=0;
-      for (auto node_id : node_id_set) {
-        node_id_list[r] = node_id;
-        ++r;
-      }
-    }
-    m_node_id_list = node_id_list;
+    ;
   }
 
   RefFaceList& operator=(const RefFaceList&) = default;
diff --git a/src/mesh/RefNodeList.hpp b/src/mesh/RefNodeList.hpp
new file mode 100644
index 0000000000000000000000000000000000000000..7bdc36fac84b1dd9cb713be5a5a64e9d992232c0
--- /dev/null
+++ b/src/mesh/RefNodeList.hpp
@@ -0,0 +1,40 @@
+#ifndef REF_NODE_LIST_HPP
+#define REF_NODE_LIST_HPP
+
+#include <Kokkos_Core.hpp>
+#include <RefId.hpp>
+
+class RefNodeList
+{
+ private:
+  RefId m_ref_id;
+  Kokkos::View<const unsigned int*> m_node_id_list;
+
+ public:
+  const RefId& refId() const
+  {
+    return m_ref_id;
+  }
+
+  const Kokkos::View<const unsigned int*>& nodeList() const
+  {
+    return m_node_id_list;
+  }
+
+  RefNodeList(const RefId& ref_id,
+              const Kokkos::View<const unsigned int*>& node_id_list)
+      : m_ref_id(ref_id),
+        m_node_id_list(node_id_list)
+  {
+    ;
+  }
+
+  RefNodeList& operator=(const RefNodeList&) = default;
+  RefNodeList& operator=(RefNodeList&&) = default;
+
+  RefNodeList() = default;
+  RefNodeList(const RefNodeList&) = default;
+  ~RefNodeList() = default;
+};
+
+#endif // REF_NODE_LIST_HPP