]> https://gitweb.dealii.org/ - dealii.git/commitdiff
FE_Nothing: pass reference-cell type to constructor 11403/head
authorPeter Munch <peterrmuench@gmail.com>
Sun, 20 Dec 2020 07:34:57 +0000 (08:34 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Tue, 22 Dec 2020 07:21:53 +0000 (08:21 +0100)
include/deal.II/fe/fe_nothing.h
source/fe/fe_nothing.cc

index 3b219b0f0ade9306f07a221f051ed55edd85e583..b9185e372610662a0866fd218c114a8305a75e4e 100644 (file)
@@ -79,16 +79,27 @@ class FE_Nothing : public FiniteElement<dim, spacedim>
 {
 public:
   /**
-   * Constructor. First argument denotes the number of components to give this
+   * Constructor.
+   *
+   * The first argument specifies the reference-cell type.
+   *
+   * The second 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
+   * The third argument decides whether FE_Nothing will dominate any other FE in
    * compare_for_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 ReferenceCell::Type &type,
+             const unsigned int         n_components = 1,
+             const bool                 dominate     = false);
+
+  /**
+   * Same as above but for a hypercube reference-cell type.
+   */
   FE_Nothing(const unsigned int n_components = 1, const bool dominate = false);
 
   virtual std::unique_ptr<FiniteElement<dim, spacedim>>
index ce77a56f79482a7581cd694949e574953858f330..35a1436786fd86f1c209dec03fa2b80b5d7447ab 100644 (file)
@@ -22,10 +22,12 @@ DEAL_II_NAMESPACE_OPEN
 
 
 template <int dim, int spacedim>
-FE_Nothing<dim, spacedim>::FE_Nothing(const unsigned int n_components,
-                                      const bool         dominate)
+FE_Nothing<dim, spacedim>::FE_Nothing(const ReferenceCell::Type &type,
+                                      const unsigned int         n_components,
+                                      const bool                 dominate)
   : FiniteElement<dim, spacedim>(
       FiniteElementData<dim>(std::vector<unsigned>(dim + 1, 0),
+                             type,
                              n_components,
                              0,
                              FiniteElementData<dim>::unknown),
@@ -41,6 +43,17 @@ FE_Nothing<dim, spacedim>::FE_Nothing(const unsigned int n_components,
 }
 
 
+
+template <int dim, int spacedim>
+FE_Nothing<dim, spacedim>::FE_Nothing(const unsigned int n_components,
+                                      const bool         dominate)
+  : FE_Nothing<dim, spacedim>(ReferenceCell::get_hypercube(dim),
+                              n_components,
+                              dominate)
+{}
+
+
+
 template <int dim, int spacedim>
 std::unique_ptr<FiniteElement<dim, spacedim>>
 FE_Nothing<dim, spacedim>::clone() const
@@ -194,10 +207,11 @@ operator==(const FiniteElement<dim, spacedim> &f) const
   // Then make sure the other object is really of type FE_Nothing,
   // and compare the data that has been passed to both objects'
   // constructors.
-  if (const FE_Nothing<dim, spacedim> *f_nothing =
+  if (const FE_Nothing<dim, spacedim> *fe_nothing =
         dynamic_cast<const FE_Nothing<dim, spacedim> *>(&f))
-    return ((dominate == f_nothing->dominate) &&
-            (this->components == f_nothing->components));
+    return ((dominate == fe_nothing->dominate) &&
+            (this->components == fe_nothing->components) &&
+            (this->reference_cell_type() == fe_nothing->reference_cell_type()));
   else
     return false;
 }

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.