Commits

Garth Wells  committed 7c8cc31 Merge

Merge branch 'garth/fix-const-dof-to-vertex'

  • Participants
  • Parent commits f68f0fb, df16b23

Comments (0)

Files changed (4)

File dolfin/fem/DofMap.cpp

   return x;
 }
 //-----------------------------------------------------------------------------
-std::vector<dolfin::la_index> DofMap::dof_to_vertex_map(Mesh& mesh) const
+std::vector<dolfin::la_index> DofMap::dof_to_vertex_map(const Mesh& mesh) const
 {
   // Check that we only have dofs living on vertices
   assert(_ufc_dofmap);
     const std::vector<dolfin::la_index>& _cell_dofs = cell_dofs(cell.index());
 
     // Tabulate local to local map of dofs on local vertex
-    _ufc_dofmap->tabulate_entity_dofs(local_to_local_map.data(), 0, local_vertex_ind);
+    _ufc_dofmap->tabulate_entity_dofs(local_to_local_map.data(), 0,
+                                      local_vertex_ind);
 
     // Fill local dofs for the vertex
     for (std::size_t local_dof = 0; local_dof < dofs_per_vertex; local_dof++)
   return dof_map;
 }
 //-----------------------------------------------------------------------------
-std::vector<std::size_t> DofMap::vertex_to_dof_map(Mesh& mesh) const
+std::vector<std::size_t> DofMap::vertex_to_dof_map(const Mesh& mesh) const
 {
   // Get dof to vertex map
   const std::vector<dolfin::la_index> dof_map = dof_to_vertex_map(mesh);

File dolfin/fem/DofMap.h

     /// *Returns*
     ///     std::vector<dolfin::la_index>
     ///         The dof to vertex map
-    std::vector<dolfin::la_index> dof_to_vertex_map(Mesh& mesh) const;
+    std::vector<dolfin::la_index> dof_to_vertex_map(const Mesh& mesh) const;
 
     /// Return a map between vertices and dofs
     /// (vert_ind*dofs_per_vertex + local_dof = vertex_to_dof_map[dof_ind],
     /// *Returns*
     ///     std::vector<std::size_t>
     ///         The vertex to dof map
-    std::vector<std::size_t> vertex_to_dof_map(Mesh& mesh) const;
+    std::vector<std::size_t> vertex_to_dof_map(const Mesh& mesh) const;
 
     /// Create a copy of the dof map
     ///

File dolfin/fem/GenericDofMap.h

 
     /// Return a map between vertices and dofs
     virtual std::vector<dolfin::la_index>
-      dof_to_vertex_map(Mesh& mesh) const = 0;
+      dof_to_vertex_map(const Mesh& mesh) const = 0;
 
     /// Return a map between vertices and dofs
-    virtual std::vector<std::size_t> vertex_to_dof_map(Mesh& mesh) const = 0;
+    virtual std::vector<std::size_t>
+      vertex_to_dof_map(const Mesh& mesh) const = 0;
 
     /// Tabulate the coordinates of all dofs on a cell (UFC cell version)
     virtual

File dolfin/geometry/ImplicitSurface.h

     virtual double operator()(const Point& point) const = 0;
 
     /// Signed distance function surface. If f0(p) = 0, the point p is
-    /// possibly on the surface, which case ImplicitSphere::f1 can be
+    /// possibly on the surface, which case ImplicitSurface::f1 can be
     /// called to check.
     virtual double f0(const Point& point) const
     {