]> https://gitweb.dealii.org/ - dealii.git/commitdiff
introduce a parameter to allow FE_Nothing dominate other FEs 1484/head
authorDenis Davydov <davydden@gmail.com>
Mon, 24 Aug 2015 15:31:57 +0000 (17:31 +0200)
committerDenis Davydov <davydden@gmail.com>
Mon, 31 Aug 2015 21:59:57 +0000 (23:59 +0200)
made other FEs to comply with the new FE_Nothing;
corrected documentation of FE_Nothing::compare_for_face_domination()

doc/news/changes.h
include/deal.II/fe/fe_nothing.h
source/fe/fe_bernstein.cc
source/fe/fe_face.cc
source/fe/fe_nedelec.cc
source/fe/fe_nothing.cc
source/fe/fe_q_base.cc
source/fe/fe_q_hierarchical.cc
source/fe/fe_q_iso_q1.cc
source/fe/fe_trace.cc

index f78659debbeac43296c6dbce14f2923c10ca44f3..1828acb7af32a162e0cb3bd9ccc5dc9ccd5dd53f 100644 (file)
@@ -235,6 +235,15 @@ inconvenience this causes.
 
 
 <ol>
+  <li> New: Introduce an option for FE_Nothing to dominate any other FE.
+  Therefore at interfaces where, for example, a Q1 meets an FE_Nothing,
+  we will force the traces of the two functions to be the same. Because the
+  FE_Nothing encodes a space that is zero everywhere, this means that the Q1
+  field will be forced to become zero at this interface.
+  <br>
+  (Denis Davydov, 2015/08/31)
+  </li>
+
   <li> New: Jacobian second and third derivatives are now computed by the mapping classes and can be
   accessed through FEValues in much the same way as the Jacobian and Jacobian gradient.
   <br>
index 91b8212547ab4a7729984d16b7ea64b85b6752e9..65bf2297aa0deda5e0c2bdb72a644790b26463a5 100644 (file)
@@ -82,10 +82,18 @@ class FE_Nothing : public FiniteElement<dim,spacedim>
 public:
 
   /**
-   * Constructor. Argument denotes the number of components to give this
+   * Constructor. First argument denotes the number of components to give this
    * finite element (default = 1).
+   *
+   * Second argument decides whether FE_Nothing
+   * will dominate any other FE in compare_for_face_domination() (default = false).
+   * Therefore at interfaces where, for example, a Q1 meets an FE_Nothing,
+   * we will force the traces of the two functions to be the same. Because the
+   * FE_Nothing encodes a space that is zero everywhere, this means that the Q1
+   * field will be forced to become zero at this interface.
    */
-  FE_Nothing (const unsigned int n_components = 1);
+  FE_Nothing (const unsigned int n_components = 1,
+              const bool dominate = false);
 
   /**
    * A sort of virtual copy constructor. Some places in the library, for
@@ -205,9 +213,10 @@ public:
    * particular the
    * @ref hp_paper "hp paper".
    *
-   * In the current case, this element is always assumed to dominate, unless
-   * it is also of type FE_Nothing().  In that situation, either element can
-   * dominate.
+   * In the current case, this element is assumed to dominate if the second
+   * argument in the constructor @p dominate is true. When this argument is false
+   * and @p fe_other is also of type FE_Nothing(), either element can dominate.
+   * Otherwise there are no_requirements.
    */
   virtual
   FiniteElementDomination::Domination
@@ -261,6 +270,16 @@ public:
                                     const unsigned int index,
                                     FullMatrix<double>  &interpolation_matrix) const;
 
+  /**
+   * @return true if the FE dominates any other.
+   */
+  bool is_dominating() const;
+
+private:
+  /**
+   * If true, this element will dominate any other apart from itself in compare_for_face_domination();
+   */
+  const bool dominate;
 
 };
 
index 72bb8a5918e481e72199c15a9cac978abda8f77d..fdba3c32af1628a5ccf85d2a94079532b4d6d3ee 100644 (file)
@@ -236,11 +236,18 @@ FE_Bernstein<dim,spacedim>::compare_for_face_domination (const FiniteElement<dim
       else
         return FiniteElementDomination::other_element_dominates;
     }
-  else if (dynamic_cast<const FE_Nothing<dim>*>(&fe_other) != 0)
+  else if (const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(&fe_other))
     {
-      // the FE_Nothing has no degrees of freedom and it is typically used in
-      // a context where we don't require any continuity along the interface
-      return FiniteElementDomination::no_requirements;
+      if (fe_nothing->is_dominating())
+        {
+          return FiniteElementDomination::other_element_dominates;
+        }
+      else
+        {
+          // the FE_Nothing has no degrees of freedom and it is typically used in
+          // a context where we don't require any continuity along the interface
+          return FiniteElementDomination::no_requirements;
+        }
     }
 
   Assert (false, ExcNotImplemented());
index b3d9e2db48479417a6c10ffdf38f54d915b3438f..9732fab950b1810a41dbbc097d1c2b75fb74d7cc 100644 (file)
@@ -276,11 +276,18 @@ compare_for_face_domination (const FiniteElement<dim,spacedim> &fe_other) const
       else
         return FiniteElementDomination::other_element_dominates;
     }
-  else if (dynamic_cast<const FE_Nothing<dim>*>(&fe_other) != 0)
+  else if (const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(&fe_other))
     {
-      // the FE_Nothing has no degrees of freedom and it is typically used in
-      // a context where we don't require any continuity along the interface
-      return FiniteElementDomination::no_requirements;
+      if (fe_nothing->is_dominating())
+        {
+          return FiniteElementDomination::other_element_dominates;
+        }
+      else
+        {
+          // the FE_Nothing has no degrees of freedom and it is typically used in
+          // a context where we don't require any continuity along the interface
+          return FiniteElementDomination::no_requirements;
+        }
     }
 
   Assert (false, ExcNotImplemented());
@@ -608,11 +615,18 @@ compare_for_face_domination (const FiniteElement<dim,spacedim> &fe_other) const
       else
         return FiniteElementDomination::other_element_dominates;
     }
-  else if (dynamic_cast<const FE_Nothing<dim>*>(&fe_other) != 0)
+  else if (const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(&fe_other))
     {
-      // the FE_Nothing has no degrees of freedom and it is typically used in
-      // a context where we don't require any continuity along the interface
-      return FiniteElementDomination::no_requirements;
+      if (fe_nothing->is_dominating())
+        {
+          return FiniteElementDomination::other_element_dominates;
+        }
+      else
+        {
+          // the FE_Nothing has no degrees of freedom and it is typically used in
+          // a context where we don't require any continuity along the interface
+          return FiniteElementDomination::no_requirements;
+        }
     }
 
   Assert (false, ExcNotImplemented());
index ac21bd9e8f6dd7330958047d6da1f2eb85fb0579..8934aaad2483b945f749cd02e93d9d7239daacfb 100644 (file)
@@ -2296,8 +2296,9 @@ FE_Nedelec<dim>::compare_for_face_domination (const FiniteElement<dim> &fe_other
       else
         return FiniteElementDomination::other_element_dominates;
     }
-  else if (dynamic_cast<const FE_Nothing<dim>*>(&fe_other) != 0)
+  else if (const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(&fe_other))
     {
+      // TODO: ???
       // the FE_Nothing has no
       // degrees of
       // freedom. nevertheless, we
@@ -2308,7 +2309,17 @@ FE_Nedelec<dim>::compare_for_face_domination (const FiniteElement<dim> &fe_other
       // and rather allow the
       // function to be discontinuous
       // along the interface
-      return FiniteElementDomination::other_element_dominates;
+//      return FiniteElementDomination::other_element_dominates;
+      if (fe_nothing->is_dominating())
+        {
+          return FiniteElementDomination::other_element_dominates;
+        }
+      else
+        {
+          // the FE_Nothing has no degrees of freedom and it is typically used in
+          // a context where we don't require any continuity along the interface
+          return FiniteElementDomination::no_requirements;
+        }
     }
 
   Assert (false, ExcNotImplemented());
index fb3ea4497c7294fbe344d431b02ba9170dee6f4d..d122101c4d97c174249c6afcca3cd574fc878d4f 100644 (file)
@@ -28,14 +28,16 @@ namespace
 
 
 template <int dim, int spacedim>
-FE_Nothing<dim,spacedim>::FE_Nothing (const unsigned int n_components)
+FE_Nothing<dim,spacedim>::FE_Nothing (const unsigned int n_components,
+                                      const bool dominate)
   :
   FiniteElement<dim,spacedim>
   (FiniteElementData<dim>(std::vector<unsigned>(dim+1,0),
                           n_components, 0,
                           FiniteElementData<dim>::unknown),
    std::vector<bool>(),
-   std::vector<ComponentMask>() )
+   std::vector<ComponentMask>() ),
+  dominate(dominate)
 {
 // in most other elements we have to set up all sorts of stuff
 // here. there isn't much that we have to do here; in particular,
@@ -162,13 +164,34 @@ fill_fe_subface_values (const Mapping<dim,spacedim> & /*mapping*/,
   // leave data fields empty
 }
 
+template <int dim, int spacedim>
+bool
+FE_Nothing<dim,spacedim>::is_dominating() const
+{
+  return dominate;
+}
+
 
 template <int dim, int spacedim>
 FiniteElementDomination::Domination
 FE_Nothing<dim,spacedim> ::
-compare_for_face_domination (const FiniteElement<dim,spacedim> &) const
+compare_for_face_domination (const FiniteElement<dim,spacedim> &fe) const
 {
-  return FiniteElementDomination::no_requirements;
+  // if FE_Nothing does not dominate, there are no requirements
+  if (!dominate)
+    {
+      return FiniteElementDomination::no_requirements;
+    }
+  // if it does and the other is FE_Nothing, either can dominate
+  else if (dynamic_cast<const FE_Nothing<dim>*>(&fe) != 0)
+    {
+      return FiniteElementDomination::either_element_can_dominate;
+    }
+  // otherwise we dominate whatever fe is provided
+  else
+    {
+      return FiniteElementDomination::this_element_dominates;
+    }
 }
 
 
index f6eb9fdbdc3042695cc3284ac57e24e88d91ee54..3c2b3d39e232b7fa0c471f74e0e0e9dd83b8c0db 100644 (file)
@@ -859,11 +859,18 @@ compare_for_face_domination (const FiniteElement<dim,spacedim> &fe_other) const
       else
         return FiniteElementDomination::other_element_dominates;
     }
-  else if (dynamic_cast<const FE_Nothing<dim>*>(&fe_other) != 0)
+  else if (const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(&fe_other))
     {
-      // the FE_Nothing has no degrees of freedom and it is typically used in
-      // a context where we don't require any continuity along the interface
-      return FiniteElementDomination::no_requirements;
+      if (fe_nothing->is_dominating())
+        {
+          return FiniteElementDomination::other_element_dominates;
+        }
+      else
+        {
+          // the FE_Nothing has no degrees of freedom and it is typically used in
+          // a context where we don't require any continuity along the interface
+          return FiniteElementDomination::no_requirements;
+        }
     }
 
   Assert (false, ExcNotImplemented());
index 98b151c4e4e6dd1b75bb83b8abbab4409c30f8d8..a3a7641eaa16e2fb4b207758fc17a56f7852424b 100644 (file)
@@ -364,15 +364,18 @@ compare_for_face_domination (const FiniteElement<dim> &fe_other) const
       else
         return FiniteElementDomination::other_element_dominates;
     }
-  else if (dynamic_cast<const FE_Nothing<dim>*>(&fe_other) != 0)
+  else if (const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(&fe_other))
     {
-      // the FE_Nothing has no
-      // degrees of freedom and it is
-      // typically used in a context
-      // where we don't require any
-      // continuity along the
-      // interface
-      return FiniteElementDomination::no_requirements;
+      if (fe_nothing->is_dominating())
+        {
+          return FiniteElementDomination::other_element_dominates;
+        }
+      else
+        {
+          // the FE_Nothing has no degrees of freedom and it is typically used in
+          // a context where we don't require any continuity along the interface
+          return FiniteElementDomination::no_requirements;
+        }
     }
 
   Assert (false, ExcNotImplemented());
index 5a154c4ddeb944194cdf195c7a4f86fb5cfb3323..35747b15e6cc12b509e94ea04c6097ca48075737 100644 (file)
@@ -98,11 +98,18 @@ compare_for_face_domination (const FiniteElement<dim,spacedim> &fe_other) const
       else
         return FiniteElementDomination::neither_element_dominates;
     }
-  else if (dynamic_cast<const FE_Nothing<dim>*>(&fe_other) != 0)
+  else if (const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(&fe_other))
     {
-      // the FE_Nothing has no degrees of freedom and it is typically used in
-      // a context where we don't require any continuity along the interface
-      return FiniteElementDomination::no_requirements;
+      if (fe_nothing->is_dominating())
+        {
+          return FiniteElementDomination::other_element_dominates;
+        }
+      else
+        {
+          // the FE_Nothing has no degrees of freedom and it is typically used in
+          // a context where we don't require any continuity along the interface
+          return FiniteElementDomination::no_requirements;
+        }
     }
 
   Assert (false, ExcNotImplemented());
index 3d68c93ed0211f4acc5f3f345939633eaccbb3f1..87ed272a22b12abc0689a7fa8bab5326e56f3f1a 100644 (file)
@@ -156,11 +156,18 @@ compare_for_face_domination (const FiniteElement<dim,spacedim> &fe_other) const
       else
         return FiniteElementDomination::other_element_dominates;
     }
-  else if (dynamic_cast<const FE_Nothing<dim>*>(&fe_other) != 0)
+  else if (const FE_Nothing<dim> *fe_nothing = dynamic_cast<const FE_Nothing<dim>*>(&fe_other))
     {
-      // the FE_Nothing has no degrees of freedom and it is typically used in
-      // a context where we don't require any continuity along the interface
-      return FiniteElementDomination::no_requirements;
+      if (fe_nothing->is_dominating())
+        {
+          return FiniteElementDomination::other_element_dominates;
+        }
+      else
+        {
+          // the FE_Nothing has no degrees of freedom and it is typically used in
+          // a context where we don't require any continuity along the interface
+          return FiniteElementDomination::no_requirements;
+        }
     }
 
   Assert (false, ExcNotImplemented());

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.