<h3>Specific improvements</h3>
<ol>
+ <li> New: The new function FiniteElement::get_associated_geometry_primitive() allows to
+ query whether a given degree of freedom is associated with a vertex, line,
+ quad, or hex.
+ <br>
+ (Wolfgang Bangerth, 2014/09/26)
+ </li>
+
<li> Fixed: The vector and array versions of Utilities::MPI::sum() and
Utilities::MPI::max() produced segmentation faults with some MPI implementations
if the input and output arguments were the same. This is now fixed.
// ---------------------------------------------------------------------
//
-// Copyright (C) 1998 - 2013 by the deal.II authors
+// Copyright (C) 1998 - 2014 by the deal.II authors
//
// This file is part of the deal.II library.
//
DEAL_II_NAMESPACE_OPEN
+/**
+ * A class that can represent the kinds of objects a triangulation is
+ * made up of: vertices, lines, quads and hexes.
+ *
+ * The class is rather primitive: it only stores a single integer
+ * that represents the dimensionality of the object represented.
+ * In other words, this class is useful primarily as a way to
+ * pass around an object whose data type explains what it does
+ * (unlike just passing around an integer), and for providing
+ * symbolic names for these objects such as GeometryPrimitive::vertex
+ * instead of an integer zero.
+ *
+ * Since the ability to identify such objects with the integral
+ * dimension of the object represented, this class provides
+ * conversion operators to and from unsigned integers.
+ *
+ * @author Wolfgang Bangerth, 2014
+ */
+class GeometryPrimitive
+{
+public:
+ /**
+ * An enumeration providing symbolic names for the
+ * objects that can be represented by this class.
+ * The numeric values of these symbolic names equal
+ * the geometric dimensionality of the represented
+ * objects to make conversion from and to integer
+ * variables simpler.
+ */
+ enum Object
+ {
+ vertex = 0,
+ line = 1,
+ quad = 2,
+ hex = 3
+ };
+
+ /**
+ * Constructor. Initialize the object with the
+ * given argument representing a vertex, line, etc.
+ */
+ GeometryPrimitive (const Object object);
+
+ /**
+ * Constructor. Initialize the object with an
+ * integer that should represent the dimensionality of
+ * the geometric object in question. This will usually be
+ * a number between zero (a vertex) and three (a hexahedron).
+ */
+ GeometryPrimitive (const unsigned int object_dimension);
+
+ /**
+ * Return the integral dimension of the object
+ * currently represented, i.e. zero for a vertex,
+ * one for a line, etc.
+ */
+ operator unsigned int () const;
+
+private:
+ /**
+ * The object currently represented.
+ */
+ Object object;
+};
+
+
+
/**
* A class that provides possible choices for isotropic and
/* -------------- inline functions ------------- */
+
+inline
+GeometryPrimitive::GeometryPrimitive (const Object object)
+ :
+ object (object)
+{}
+
+
+
+inline
+GeometryPrimitive::GeometryPrimitive (const unsigned int object_dimension)
+ :
+ object (static_cast<Object>(object_dimension))
+{}
+
+
+inline
+GeometryPrimitive::operator unsigned int () const
+{
+ return static_cast<unsigned int>(object);
+}
+
+
+
namespace internal
{
bool
has_generalized_face_support_points () const;
+ /**
+ * For a given degree of freedom, return whether it is logically
+ * associated with a vertex, line, quad or hex.
+ *
+ * For instance, for continuous finite elements this coincides with
+ * the lowest dimensional object the support point of the degree of
+ * freedom lies on. To give an example, for $Q_1$ elements in 3d,
+ * every degree of freedom is defined by a shape function that we
+ * get by interpolating using support points that lie on the
+ * vertices of the cell. The support of these points of course
+ * extends to all edges connected to this vertex, as well
+ * as the adjacent faces and the cell interior, but we say that
+ * logically the degree of freedom is associated with the vertex as
+ * this is the lowest-dimensional object it is associated with.
+ * Likewise, for $Q_2$ elements in 3d, the degrees of freedom
+ * with support points at edge midpoints would yield a value of
+ * GeometryPrimitive::line from this function, whereas those on
+ * the centers of faces in 3d would return GeometryPrimitive::quad.
+ *
+ * To make this more formal, the kind of object returned by
+ * this function represents the object so that the support of the
+ * shape function corresponding to the degree of freedom, (i.e., that
+ * part of the domain where the function "lives") is the union
+ * of all of the cells sharing this object. To return to the example
+ * above, for $Q_2$ in 3d, the shape function with support point at
+ * an edge midpoint has support on all cells that share the edge and
+ * not only the cells that share the adjacent faces, and consequently
+ * the function will return GeometryPrimitive::line.
+ *
+ * On the other hand, for discontinuous elements of type $DGQ_2$,
+ * a degree of freedom associated with an interpolation polynomial
+ * that has its support point physically located at a line bounding
+ * a cell, but is nonzero only on one cell. Consequently, it is
+ * logically associated with the interior of that cell (i.e., with a
+ * GeometryPrimitive::quad in 2d and a GeometryPrimitive::hex in 3d).
+ *
+ * @param[in] cell_dof_index The index of a shape function or
+ * degree of freedom. This index must be in the range
+ * <code>[0,dofs_per_cell)</code>.
+ *
+ * @note The integer value of the object returned by this function
+ * equals the dimensionality of the object it describes, and can
+ * consequently be used in generic programming paradigms. For
+ * example, if a degree of freedom is associated with a vertex,
+ * then this function returns GeometryPrimitive::vertex, which has
+ * a numeric value of zero (the dimensionality of a vertex).
+ */
+ GeometryPrimitive
+ get_associated_geometry_primitive (const unsigned int cell_dof_index) const;
+
/**
* Interpolate a set of scalar values, computed in the generalized support
* points.
+template <int dim, int spacedim>
+inline
+GeometryPrimitive
+FiniteElement<dim,spacedim>::get_associated_geometry_primitive (const unsigned int cell_dof_index) const
+{
+ Assert (cell_dof_index < this->dofs_per_cell,
+ ExcIndexRange (cell_dof_index, 0, this->dofs_per_cell));
+
+ // just go through the usual cases, taking into account how DoFs
+ // are enumerated on the reference cell
+ if (cell_dof_index < this->first_line_index)
+ return GeometryPrimitive::vertex;
+ else if (cell_dof_index < this->first_quad_index)
+ return GeometryPrimitive::line;
+ else if (cell_dof_index < this->first_hex_index)
+ return GeometryPrimitive::quad;
+ else
+ return GeometryPrimitive::hex;
+}
+
+
+
DEAL_II_NAMESPACE_CLOSE
#endif
--- /dev/null
+// ---------------------------------------------------------------------
+// $Id$
+//
+// Copyright (C) 1998 - 2014 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+// test FiniteElement::get_dof_association() with a couple of elements
+
+
+#include "../tests.h"
+#include <iostream>
+#include <fstream>
+
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_nedelec.h>
+
+
+template <int dim>
+void test (const FiniteElement<dim> &fe)
+{
+ deallog << fe.get_name() << std::endl;
+
+ for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
+ {
+ switch (fe.get_associated_geometry_primitive(i))
+ {
+ case GeometryPrimitive::vertex:
+ deallog << 'v';
+ break;
+ case GeometryPrimitive::line:
+ deallog << 'l';
+ break;
+ case GeometryPrimitive::quad:
+ deallog << 'q';
+ break;
+ case GeometryPrimitive::hex:
+ deallog << 'h';
+ break;
+ default:
+ Assert (false, ExcInternalError());
+ }
+ }
+ deallog << std::endl;
+}
+
+int main()
+{
+ std::ofstream logfile("output");
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ test(FE_Q<1>(2));
+ test(FE_Q<2>(2));
+ test(FE_Q<3>(2));
+
+ test(FE_Nedelec<2>(0));
+ test(FE_Nedelec<3>(0));
+
+ test(FE_Nedelec<2>(1));
+ test(FE_Nedelec<3>(1));
+}
--- /dev/null
+
+DEAL::FE_Q<1>(2)
+DEAL::vvl
+DEAL::FE_Q<2>(2)
+DEAL::vvvvllllq
+DEAL::FE_Q<3>(2)
+DEAL::vvvvvvvvllllllllllllqqqqqqh
+DEAL::FE_Nedelec<2>(0)
+DEAL::llll
+DEAL::FE_Nedelec<3>(0)
+DEAL::llllllllllll
+DEAL::FE_Nedelec<2>(1)
+DEAL::llllllllqqqq
+DEAL::FE_Nedelec<3>(1)
+DEAL::llllllllllllllllllllllllqqqqqqqqqqqqqqqqqqqqqqqqhhhhhh