]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove implementations for FE_FaceP and update documentation
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Tue, 11 Jul 2017 21:37:37 +0000 (23:37 +0200)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Wed, 12 Jul 2017 09:25:28 +0000 (11:25 +0200)
doc/news/changes/minor/20170711Simfeld-DanielArndt [new file with mode: 0644]
include/deal.II/fe/fe_face.h
source/fe/fe_face.cc

diff --git a/doc/news/changes/minor/20170711Simfeld-DanielArndt b/doc/news/changes/minor/20170711Simfeld-DanielArndt
new file mode 100644 (file)
index 0000000..5060064
--- /dev/null
@@ -0,0 +1,5 @@
+New: Adding the functions FE_FaceQ::hp_vertex_dof_identities(),
+FE_FaceQ::hp_line_dof_identities() and FE_FaceQ::hp_quad_dof_identities()
+allows FE_FaceQ to be used in combination with a hp::DoFHandler.
+<br>
+(Simfeld, Daniel Arndt, 2017/07/11)
index 4c1d42bc87267d4d3213e87182efda27ac7c7d80..0f436a7ff7e56e7dd99d0d33f88f21233d0fcb9b 100644 (file)
@@ -122,8 +122,7 @@ public:
    * the vertex dofs of the present element, whereas the second is the
    * corresponding index of the other finite element.
    *
-   * This being a discontinuous element, the set of such constraints is of
-   * course empty.
+   * The set of such constraints is non-empty only for dim==1.
    */
   virtual
   std::vector<std::pair<unsigned int, unsigned int> >
@@ -133,8 +132,7 @@ public:
    * Same as hp_vertex_dof_indices(), except that the function treats degrees
    * of freedom on lines.
    *
-   * This being a discontinuous element, the set of such constraints is of
-   * course empty.
+   * The set of such constraints is non-empty only for dim==2.
    */
   virtual
   std::vector<std::pair<unsigned int, unsigned int> >
@@ -144,8 +142,7 @@ public:
    * Same as hp_vertex_dof_indices(), except that the function treats degrees
    * of freedom on quads.
    *
-   * This being a discontinuous element, the set of such constraints is of
-   * course empty.
+   * The set of such constraints is non-empty only for dim==3.
    */
   virtual
   std::vector<std::pair<unsigned int, unsigned int> >
@@ -276,8 +273,7 @@ public:
    * the vertex dofs of the present element, whereas the second is the
    * corresponding index of the other finite element.
    *
-   * This being a discontinuous element, the set of such constraints is of
-   * course empty.
+   * The set of such constraints is non-empty only for dim==1.
    */
   virtual
   std::vector<std::pair<unsigned int, unsigned int> >
@@ -287,8 +283,7 @@ public:
    * Same as hp_vertex_dof_indices(), except that the function treats degrees
    * of freedom on lines.
    *
-   * This being a discontinuous element, the set of such constraints is of
-   * course empty.
+   * The set of such constraints is non-empty only for dim==2.
    */
   virtual
   std::vector<std::pair<unsigned int, unsigned int> >
@@ -298,8 +293,7 @@ public:
    * Same as hp_vertex_dof_indices(), except that the function treats degrees
    * of freedom on quads.
    *
-   * This being a discontinuous element, the set of such constraints is of
-   * course empty.
+   * The set of such constraints is non-empty only for dim==3.
    */
   virtual
   std::vector<std::pair<unsigned int, unsigned int> >
@@ -501,50 +495,6 @@ public:
    */
   virtual bool hp_constraints_are_implemented () const;
 
-  /**
-   * If, on a vertex, several finite elements are active, the hp code first
-   * assigns the degrees of freedom of each of these FEs different global
-   * indices. It then calls this function to find out which of them should get
-   * identical values, and consequently can receive the same global DoF index.
-   * This function therefore returns a list of identities between DoFs of the
-   * present finite element object with the DoFs of @p fe_other, which is a
-   * reference to a finite element object representing one of the other finite
-   * elements active on this particular vertex. The function computes which of
-   * the degrees of freedom of the two finite element objects are equivalent,
-   * both numbered between zero and the corresponding value of dofs_per_vertex
-   * of the two finite elements. The first index of each pair denotes one of
-   * the vertex dofs of the present element, whereas the second is the
-   * corresponding index of the other finite element.
-   *
-   * This being a discontinuous element, the set of such constraints is of
-   * course empty.
-   */
-  virtual
-  std::vector<std::pair<unsigned int, unsigned int> >
-  hp_vertex_dof_identities (const FiniteElement<dim, spacedim> &fe_other) const;
-
-  /**
-   * Same as hp_vertex_dof_indices(), except that the function treats degrees
-   * of freedom on lines.
-   *
-   * This being a discontinuous element, the set of such constraints is of
-   * course empty.
-   */
-  virtual
-  std::vector<std::pair<unsigned int, unsigned int> >
-  hp_line_dof_identities (const FiniteElement<dim, spacedim> &fe_other) const;
-
-  /**
-   * Same as hp_vertex_dof_indices(), except that the function treats degrees
-   * of freedom on quads.
-   *
-   * This being a discontinuous element, the set of such constraints is of
-   * course empty.
-   */
-  virtual
-  std::vector<std::pair<unsigned int, unsigned int> >
-  hp_quad_dof_identities (const FiniteElement<dim, spacedim> &fe_other) const;
-
   /**
    * Return whether this element dominates the one given as argument when they
    * meet at a common face, whether it is the other way around, whether
index ed73ac5a99c7008e94b7b3347c7f1797ab591f1c..6106fe38264c7c37b52abc9461299e1271fa5eb7 100644 (file)
@@ -767,7 +767,6 @@ FE_FaceP<dim,spacedim>::get_dpo_vector (const unsigned int deg)
 
 
 
-
 template <int dim, int spacedim>
 bool
 FE_FaceP<dim,spacedim>::hp_constraints_are_implemented () const
@@ -777,51 +776,6 @@ FE_FaceP<dim,spacedim>::hp_constraints_are_implemented () const
 
 
 
-template <int dim, int spacedim>
-std::vector<std::pair<unsigned int, unsigned int> >
-FE_FaceP<dim, spacedim>::
-hp_vertex_dof_identities (const FiniteElement<dim, spacedim> &/*fe_other*/) const
-{
-  // this element is discontinuous, so by definition there can
-  // be no identities between its dofs and those of any neighbor
-  // (of whichever type the neighbor may be -- after all, we have
-  // no face dofs on this side to begin with)
-  return
-    std::vector<std::pair<unsigned int, unsigned int> > ();
-}
-
-
-
-template <int dim, int spacedim>
-std::vector<std::pair<unsigned int, unsigned int> >
-FE_FaceP<dim, spacedim>::
-hp_line_dof_identities (const FiniteElement<dim, spacedim> &/*fe_other*/) const
-{
-  // this element is discontinuous, so by definition there can
-  // be no identities between its dofs and those of any neighbor
-  // (of whichever type the neighbor may be -- after all, we have
-  // no face dofs on this side to begin with)
-  return
-    std::vector<std::pair<unsigned int, unsigned int> > ();
-}
-
-
-
-template <int dim, int spacedim>
-std::vector<std::pair<unsigned int, unsigned int> >
-FE_FaceP<dim, spacedim>::
-hp_quad_dof_identities (const FiniteElement<dim, spacedim> &/*fe_other*/) const
-{
-  // this element is discontinuous, so by definition there can
-  // be no identities between its dofs and those of any neighbor
-  // (of whichever type the neighbor may be -- after all, we have
-  // no face dofs on this side to begin with)
-  return
-    std::vector<std::pair<unsigned int, unsigned int> > ();
-}
-
-
-
 template <int dim, int spacedim>
 FiniteElementDomination::Domination
 FE_FaceP<dim,spacedim>::

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.