*/
template <int rank, int dim>
class TensorFunction : public FunctionTime,
- public Subscriptor
+ public Subscriptor
{
public:
/**
- * Define typedefs for the return
- * types of the <tt>value</tt>
- * functions.
+ * Define typedefs for the return types of the <tt>value</tt> functions.
*/
typedef Tensor<rank,dim> value_type;
typedef Tensor<rank+1,dim> gradient_type;
/**
- * Constructor. May take an
- * initial value for the time
- * variable, which defaults to
- * zero.
+ * Constructor. May take an initial value for the time variable, which
+ * defaults to zero.
*/
TensorFunction (const double initial_time = 0.0);
/**
- * Virtual destructor; absolutely
- * necessary in this case, as
- * classes are usually not used
- * by their true type, but rather
- * through pointers to this base
- * class.
+ * Virtual destructor; absolutely necessary in this case, as classes are
+ * usually not used by their true type, but rather through pointers to
+ * this base class.
*/
virtual ~TensorFunction ();
/**
- * Return the value of the function
- * at the given point.
+ * Return the value of the function at the given point.
*/
virtual value_type value (const Point<dim> &p) const;
/**
- * Set <tt>values</tt> to the point
- * values of the function at the
- * <tt>points</tt>. It is assumed
- * that <tt>values</tt> already has
- * the right size, i.e. the same
- * size as the <tt>points</tt> array.
+ * Set <tt>values</tt> to the point values of the function at the
+ * <tt>points</tt>. It is assumed that <tt>values</tt> already has the
+ * right size, i.e. the same size as the <tt>points</tt> array.
*/
virtual void value_list (const std::vector<Point<dim> > &points,
std::vector<value_type> &values) const;
/**
- * Return the gradient of the
- * function at the given point.
+ * Return the gradient of the function at the given point.
*/
virtual gradient_type gradient (const Point<dim> &p) const;
/**
- * Set <tt>gradients</tt> to the
- * gradients of the function at
- * the <tt>points</tt>. It is assumed
- * that <tt>values</tt> already has
- * the right size, i.e. the same
- * size as the <tt>points</tt> array.
+ * Set <tt>gradients</tt> to the gradients of the function at the
+ * <tt>points</tt>. It is assumed that <tt>values</tt> already has the
+ * right size, i.e. the same size as the <tt>points</tt> array.
*/
virtual void gradient_list (const std::vector<Point<dim> > &points,
std::vector<gradient_type> &gradients) const;
};
+
+
+/*
+ * Provide a tensor valued function which always returns a constant tensor
+ * value. Obviously, all derivates of this function are zero.
+ *
+ * @ingroup functions
+ * @author Matthias Maier, 2013
+ */
+template <int rank, int dim>
+class ConstantTensorFunction : public TensorFunction<rank, dim>
+{
+public:
+ /**
+ * Constructor; takes the constant tensor value as an argument.
+ * The reference value is copied internally.
+ *
+ * An initial value for the time variable may be specified, otherwise it
+ * defaults to zero.
+ */
+ ConstantTensorFunction (const dealii::Tensor<rank, dim> &value,
+ const double initial_time = 0.0);
+
+ /**
+ * Virtual destructor; absolutely necessary in this case, as classes are
+ * usually not used by their true type, but rather through pointers to
+ * this base class.
+ */
+ virtual ~ConstantTensorFunction ();
+
+ /**
+ * Return the value of the function at the given point.
+ */
+ virtual typename dealii::TensorFunction<rank, dim>::value_type value (const Point<dim> &p) const;
+
+ /**
+ * Set <tt>values</tt> to the point values of the function at the
+ * <tt>points</tt>. It is assumed that <tt>values</tt> already has the
+ * right size, i.e. the same size as the <tt>points</tt> array.
+ */
+ virtual void value_list (const std::vector<Point<dim> > &points,
+ std::vector<typename dealii::TensorFunction<rank, dim>::value_type> &values) const;
+
+ /**
+ * Return the gradient of the function at the given point.
+ */
+ virtual typename dealii::TensorFunction<rank, dim>::gradient_type gradient (const Point<dim> &p) const;
+
+ /**
+ * Set <tt>gradients</tt> to the gradients of the function at the
+ * <tt>points</tt>. It is assumed that <tt>values</tt> already has the
+ * right size, i.e. the same size as the <tt>points</tt> array.
+ */
+ virtual void gradient_list (const std::vector<Point<dim> > &points,
+ std::vector<typename dealii::TensorFunction<rank, dim>::gradient_type> &gradients) const;
+
+private:
+ const dealii::Tensor<rank, dim> _value;
+};
+
+
+
+/*
+ * Provide a tensor valued function which always returns zero.
+ * Obviously, all derivates of this function are zero.
+ *
+ * @ingroup functions
+ * @author Matthias Maier, 2013
+ */
+template <int rank, int dim>
+class ZeroTensorFunction : public ConstantTensorFunction<rank, dim>
+{
+public:
+ /**
+ * Constructor.
+ *
+ * An initial value for the time variable may be specified, otherwise it
+ * defaults to zero.
+ */
+ ZeroTensorFunction (const double initial_time = 0.0);
+};
+
+
DEAL_II_NAMESPACE_CLOSE
#endif
DEAL_II_NAMESPACE_OPEN
+
//////////////////////////////////////////////////////////////////////
// TensorFunction
//////////////////////////////////////////////////////////////////////
}
+
+//////////////////////////////////////////////////////////////////////
+// ConstantTensorFunction
+//////////////////////////////////////////////////////////////////////
+
+
+template <int rank, int dim>
+ConstantTensorFunction<rank, dim>::ConstantTensorFunction (const Tensor<rank, dim> &value,
+ const double initial_time)
+ :
+ TensorFunction<rank, dim> (initial_time),
+ _value(value)
+{}
+
+
+template <int rank, int dim>
+ConstantTensorFunction<rank, dim>::~ConstantTensorFunction ()
+{}
+
+
+template <int rank, int dim>
+typename TensorFunction<rank, dim>::value_type
+ConstantTensorFunction<rank, dim>::value (const Point<dim> &/*point*/) const
+{
+ return _value;
+}
+
+
+template <int rank, int dim>
+void
+ConstantTensorFunction<rank, dim>::value_list (const std::vector<Point<dim> > &points,
+ std::vector<typename TensorFunction<rank, dim>::value_type> &values) const
+{
+ Assert (values.size() == points.size(),
+ ExcDimensionMismatch(values.size(), points.size()));
+
+ for (unsigned int i=0; i<values.size(); ++i)
+ values[i] = _value;
+}
+
+
+template <int rank, int dim>
+typename TensorFunction<rank, dim>::gradient_type
+ConstantTensorFunction<rank, dim>::gradient (const Point<dim> &) const
+{
+ static const Tensor<rank+1,dim> zero;
+
+ return zero;
+}
+
+
+template <int rank, int dim>
+void
+ConstantTensorFunction<rank, dim>::gradient_list (const std::vector<Point<dim> > &points,
+ std::vector<typename TensorFunction<rank, dim>::gradient_type> &gradients) const
+{
+ Assert (gradients.size() == points.size(),
+ ExcDimensionMismatch(gradients.size(), points.size()));
+
+ static const Tensor<rank+1,dim> zero;
+
+ for (unsigned int i=0; i<gradients.size(); ++i)
+ gradients[i] = zero;
+}
+
+
+
+//////////////////////////////////////////////////////////////////////
+// ZeroTensorFunction
+//////////////////////////////////////////////////////////////////////
+
+
+template <int rank, int dim>
+ZeroTensorFunction<rank, dim>::ZeroTensorFunction (const double initial_time)
+ :
+ ConstantTensorFunction<rank, dim> (dealii::Tensor<rank, dim>(), initial_time)
+{}
+
+
+
+//////////////////////////////////////////////////////////////////////
+// Explicit instantiations:
+//////////////////////////////////////////////////////////////////////
+
template class TensorFunction<1,1>;
template class TensorFunction<2,1>;
template class TensorFunction<3,1>;
template class TensorFunction<3,3>;
template class TensorFunction<4,3>;
+template class ConstantTensorFunction<1,1>;
+template class ConstantTensorFunction<2,1>;
+template class ConstantTensorFunction<3,1>;
+template class ConstantTensorFunction<4,1>;
+template class ConstantTensorFunction<1,2>;
+template class ConstantTensorFunction<2,2>;
+template class ConstantTensorFunction<3,2>;
+template class ConstantTensorFunction<4,2>;
+template class ConstantTensorFunction<1,3>;
+template class ConstantTensorFunction<2,3>;
+template class ConstantTensorFunction<3,3>;
+template class ConstantTensorFunction<4,3>;
+
+template class ZeroTensorFunction<1,1>;
+template class ZeroTensorFunction<2,1>;
+template class ZeroTensorFunction<3,1>;
+template class ZeroTensorFunction<4,1>;
+template class ZeroTensorFunction<1,2>;
+template class ZeroTensorFunction<2,2>;
+template class ZeroTensorFunction<3,2>;
+template class ZeroTensorFunction<4,2>;
+template class ZeroTensorFunction<1,3>;
+template class ZeroTensorFunction<2,3>;
+template class ZeroTensorFunction<3,3>;
+template class ZeroTensorFunction<4,3>;
+
+
DEAL_II_NAMESPACE_CLOSE