<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>
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
* 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
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;
};
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());
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());
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());
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
// 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());
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,
// 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;
+ }
}
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());
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());
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());
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());