]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Adjust FE_DGPNonparametric to not return values on the reference cell.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 20 Nov 2015 22:44:24 +0000 (16:44 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 20 Nov 2015 22:44:24 +0000 (16:44 -0600)
This is because the element does not define its shape functions by mapping
from the reference cell, so evaluating anything on the reference cell
makes no sense. Fix this by throwing an exception, as discussed in
the documentation of the corresponding functions in the base class.

include/deal.II/fe/fe_dgp_nonparametric.h
source/fe/fe_dgp_nonparametric.cc

index b695616e82313079cb8e3a94696c3c70e6023cbc..15c315824421be16dbc74361cb7b24db9b4cafe5 100644 (file)
@@ -289,63 +289,81 @@ public:
   requires_update_flags (const UpdateFlags update_flags) const;
 
   /**
-   * Return the value of the @p ith shape function at the point @p p. See the
-   * FiniteElement base class for more information about the semantics of this
-   * function.
+   * This function is intended to return the value of a shape function
+   * at a point on the reference cell. However, since the current element
+   * does not implement shape functions by mapping from a reference cell,
+   * no shape functions exist on the reference cell.
+   *
+   * Consequently, as discussed in the corresponding function in the base
+   * class, FiniteElement::shape_value(), this function throws an exception
+   * of type FiniteElement::ExcUnitShapeValuesDoNotExist.
    */
   virtual double shape_value (const unsigned int i,
                               const Point<dim> &p) const;
 
   /**
-   * Return the value of the @p componentth vector component of the @p ith
-   * shape function at the point @p p. See the FiniteElement base class for
-   * more information about the semantics of this function.
+   * This function is intended to return the value of a shape function
+   * at a point on the reference cell. However, since the current element
+   * does not implement shape functions by mapping from a reference cell,
+   * no shape functions exist on the reference cell.
    *
-   * Since this element is scalar, the returned value is the same as if the
-   * function without the @p _component suffix were called, provided that the
-   * specified component is zero.
+   * Consequently, as discussed in the corresponding function in the base
+   * class, FiniteElement::shape_value_component(), this function throws an exception
+   * of type FiniteElement::ExcUnitShapeValuesDoNotExist.
    */
   virtual double shape_value_component (const unsigned int i,
                                         const Point<dim> &p,
                                         const unsigned int component) const;
 
   /**
-   * Return the gradient of the @p ith shape function at the point @p p. See
-   * the FiniteElement base class for more information about the semantics of
-   * this function.
+   * This function is intended to return the gradient of a shape function
+   * at a point on the reference cell. However, since the current element
+   * does not implement shape functions by mapping from a reference cell,
+   * no shape functions exist on the reference cell.
+   *
+   * Consequently, as discussed in the corresponding function in the base
+   * class, FiniteElement::shape_grad(), this function throws an exception
+   * of type FiniteElement::ExcUnitShapeValuesDoNotExist.
    */
   virtual Tensor<1,dim> shape_grad (const unsigned int  i,
                                     const Point<dim>   &p) const;
 
   /**
-   * Return the gradient of the @p componentth vector component of the @p ith
-   * shape function at the point @p p. See the FiniteElement base class for
-   * more information about the semantics of this function.
+   * This function is intended to return the gradient of a shape function
+   * at a point on the reference cell. However, since the current element
+   * does not implement shape functions by mapping from a reference cell,
+   * no shape functions exist on the reference cell.
    *
-   * Since this element is scalar, the returned value is the same as if the
-   * function without the @p _component suffix were called, provided that the
-   * specified component is zero.
+   * Consequently, as discussed in the corresponding function in the base
+   * class, FiniteElement::shape_grad_component(), this function throws an exception
+   * of type FiniteElement::ExcUnitShapeValuesDoNotExist.
    */
   virtual Tensor<1,dim> shape_grad_component (const unsigned int i,
                                               const Point<dim> &p,
                                               const unsigned int component) const;
 
   /**
-   * Return the tensor of second derivatives of the @p ith shape function at
-   * point @p p on the unit cell.  See the FiniteElement base class for more
-   * information about the semantics of this function.
+   * This function is intended to return the Hessian of a shape function
+   * at a point on the reference cell. However, since the current element
+   * does not implement shape functions by mapping from a reference cell,
+   * no shape functions exist on the reference cell.
+   *
+   * Consequently, as discussed in the corresponding function in the base
+   * class, FiniteElement::shape_grad_grad(), this function throws an exception
+   * of type FiniteElement::ExcUnitShapeValuesDoNotExist.
    */
   virtual Tensor<2,dim> shape_grad_grad (const unsigned int  i,
                                          const Point<dim> &p) const;
 
   /**
-   * Return the second derivative of the @p componentth vector component of
-   * the @p ith shape function at the point @p p. See the FiniteElement base
-   * class for more information about the semantics of this function.
+   * This function is intended to return the Hessian of a shape function
+   * at a point on the reference cell. However, since the current element
+   * does not implement shape functions by mapping from a reference cell,
+   * no shape functions exist on the reference cell.
    *
-   * Since this element is scalar, the returned value is the same as if the
-   * function without the @p _component suffix were called, provided that the
-   * specified component is zero.
+   * Consequently, as discussed in the corresponding function in the base
+   * class, FiniteElement::shape_grad_grad_component(), this function throws an exception
+   * of type FiniteElement::ExcUnitShapeValuesDoNotExist.
    */
   virtual Tensor<2,dim> shape_grad_grad_component (const unsigned int i,
                                                    const Point<dim> &p,
index b12b7dd6bd5a143c26e821ae2ce6b60e1ffaec85..0d34f424c379c5803e8a67c565245c9e12adbb82 100644 (file)
@@ -128,8 +128,11 @@ double
 FE_DGPNonparametric<dim,spacedim>::shape_value (const unsigned int i,
                                                 const Point<dim> &p) const
 {
+  (void)i;
+  (void)p;
   Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
-  return polynomial_space.compute_value(i, p);
+  AssertThrow (false, (typename FiniteElement<dim>::ExcUnitShapeValuesDoNotExist()));
+  return 0;
 }
 
 
@@ -140,10 +143,13 @@ FE_DGPNonparametric<dim,spacedim>::shape_value_component (const unsigned int i,
                                                           const Point<dim> &p,
                                                           const unsigned int component) const
 {
+  (void)i;
+  (void)p;
   (void)component;
   Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
   Assert (component == 0, ExcIndexRange (component, 0, 1));
-  return polynomial_space.compute_value(i, p);
+  AssertThrow (false, (typename FiniteElement<dim>::ExcUnitShapeValuesDoNotExist()));
+  return 0;
 }
 
 
@@ -153,8 +159,11 @@ Tensor<1,dim>
 FE_DGPNonparametric<dim,spacedim>::shape_grad (const unsigned int i,
                                                const Point<dim> &p) const
 {
+  (void)i;
+  (void)p;
   Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
-  return polynomial_space.compute_grad(i, p);
+  AssertThrow (false, (typename FiniteElement<dim>::ExcUnitShapeValuesDoNotExist()));
+  return Tensor<1,dim>();
 }
 
 
@@ -164,10 +173,13 @@ FE_DGPNonparametric<dim,spacedim>::shape_grad_component (const unsigned int i,
                                                          const Point<dim> &p,
                                                          const unsigned int component) const
 {
+  (void)i;
+  (void)p;
   (void)component;
   Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
   Assert (component == 0, ExcIndexRange (component, 0, 1));
-  return polynomial_space.compute_grad(i, p);
+  AssertThrow (false, (typename FiniteElement<dim>::ExcUnitShapeValuesDoNotExist()));
+  return Tensor<1,dim>();
 }
 
 
@@ -177,8 +189,11 @@ Tensor<2,dim>
 FE_DGPNonparametric<dim,spacedim>::shape_grad_grad (const unsigned int i,
                                                     const Point<dim> &p) const
 {
+  (void)i;
+  (void)p;
   Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
-  return polynomial_space.compute_grad_grad(i, p);
+  AssertThrow (false, (typename FiniteElement<dim>::ExcUnitShapeValuesDoNotExist()));
+  return Tensor<2,dim>();
 }
 
 
@@ -189,10 +204,13 @@ FE_DGPNonparametric<dim,spacedim>::shape_grad_grad_component (const unsigned int
     const Point<dim> &p,
     const unsigned int component) const
 {
+  (void)i;
+  (void)p;
   (void)component;
   Assert (i<this->dofs_per_cell, ExcIndexRange(i, 0, this->dofs_per_cell));
   Assert (component == 0, ExcIndexRange (component, 0, 1));
-  return polynomial_space.compute_grad_grad(i, p);
+  AssertThrow (false, (typename FiniteElement<dim>::ExcUnitShapeValuesDoNotExist()));
+  return Tensor<2,dim>();
 }
 
 

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.