<h3>Specific improvements</h3>
<ol>
+<li> New: There is a class VectorFunctionFromTensorFunction that converts
+between objects of type TensorFunction and Function.
+<br>
+(Spencer Patty, 2013/4/2)
+</li>
+
<li> Fixed: The ParameterHandler class could not deal with parameters named
<code>"value"</code> (and a few other names). This is now fixed.
#include <vector>
DEAL_II_NAMESPACE_OPEN
+
+
template <typename number> class Vector;
+template <int rank, int dim> class TensorFunction;
/**
* This class is a model for a general function. It serves the purpose
};
+/**
+ * This class is built as a means of translating
+ * the <code>Tensor<1,dim> </code> values produced by objects of type
+ * TensorFunction and returning them as a multiple
+ * component version of the same thing as a Vector for use in, for example,
+ * the VectorTools::interpolate or the many other functions taking Function objects.
+ * It allows the user to place the desired components into an <tt>n_components</tt>
+ * long vector starting at the <tt>selected_component</tt> location in that vector and have all other
+ * components be 0.
+ *
+ * For example: Say you created a class called
+ * @code
+ * class RightHandSide : public TensorFunction<rank,dim>
+ * @endcode
+ * which extends the TensorFunction class and you have an object
+ * @code
+ * RightHandSide<1,dim> rhs;
+ * @endcode
+ * of that class which you want to interpolate onto your mesh using the
+ * VectorTools::interpolate function, but the finite element you use for the DoFHandler
+ * object has 3 copies of a finite element with <tt>dim</tt> components, for a total of 3*dim
+ * components. To interpolate onto that DoFHandler, you need an object of type Function
+ * that has 3*dim vector components. Creating such an object from the existing <code>rhs</code>
+ * object is done using this piece of code:
+ * @code
+ * RighHandSide<1,dim> rhs;
+ * VectorFunctionFromTensorFunction<dim> rhs_vector_function (rhs, 0, 3*dim);
+ * @endcode
+ * where the <code>dim</code> components of the tensor function are placed into
+ * the first <code>dim</code> components of the function object.
+ *
+ * @author Spencer Patty, 2013
+ */
+template <int dim>
+class VectorFunctionFromTensorFunction : public Function<dim>
+{
+public:
+ /** Given a TensorFunction object that takes a <tt>Point</tt> and returns a <tt>Tensor<1,dim></tt>
+ * value, convert this into an object that matches the Function@<dim@>
+ * interface.
+ *
+ * By default, create a Vector object of the same size as <tt>tensor_function</tt>
+ * returns, i.e., with <tt>dim</tt> components.
+ *
+ * @param tensor_function The TensorFunction that will form one component
+ * of the resulting Vector Function object.
+ * @param n_components The total number of vector components of the
+ * resulting TensorFunction object.
+ * @param selected_component The first component that should be
+ * filled by the first argument. This should be such that the entire tensor_function
+ * fits inside the <tt>n_component</tt> length return vector.
+ */
+ VectorFunctionFromTensorFunction (const TensorFunction<1,dim> &tensor_function,
+ const unsigned int selected_component=0,
+ const unsigned int n_components=dim);
+
+ /**
+ * This destructor is defined as
+ * virtual so as to coincide with all
+ * other aspects of class.
+ */
+ virtual ~VectorFunctionFromTensorFunction();
+
+ /**
+ * Return a single component of a
+ * vector-valued function at a
+ * given point.
+ */
+ virtual double value (const Point<dim> &p,
+ const unsigned int component = 0) const;
+
+ /**
+ * Return all components of a
+ * vector-valued function at a
+ * given point.
+ *
+ * <tt>values</tt> shall have the right
+ * size beforehand,
+ * i.e. #n_components.
+ */
+ virtual void vector_value (const Point<dim> &p,
+ Vector<double> &values) const;
+
+ /**
+ * Return all components of a
+ * vector-valued function at a
+ * list of points.
+ *
+ * <tt>value_list</tt> shall be the same
+ * size as <tt>points</tt> and each element
+ * of the vector will be passed to
+ * vector_value() to evaluate
+ * the function
+ */
+ virtual void vector_value_list (const std::vector<Point<dim> > &points,
+ std::vector<Vector<double> > &value_list) const;
+
+private:
+ /**
+ * The TensorFunction object which we call when this class's vector_value()
+ * or vector_value_list() functions are called.
+ **/
+ const TensorFunction<1,dim> & tensor_function;
+
+ /**
+ * The first vector component whose value is to be filled by the
+ * given TensorFunction. The values will be placed in components
+ * selected_component to selected_component+dim-1 for a <tt>TensorFunction<1,dim></tt> object.
+ */
+ const unsigned int selected_component;
+};
+
+
DEAL_II_NAMESPACE_CLOSE
#endif
#include <deal.II/base/function.h>
+#include <deal.II/base/tensor_function.h>
#include <deal.II/base/point.h>
#include <deal.II/lac/vector.h>
#include <vector>
+/**
+ * The constructor for <tt>VectorFunctionFromTensorFunction</tt> which initiates the return vector
+ * to be size <tt>n_components</tt>.
+ */
+template <int dim>
+VectorFunctionFromTensorFunction<dim>::VectorFunctionFromTensorFunction (const TensorFunction<1,dim>& tensor_function,
+ const unsigned int selected_component,
+ const unsigned int n_components)
+:
+Function<dim> (n_components),
+tensor_function (tensor_function),
+selected_component (selected_component)
+{
+
+ // Verify that the Tensor<1,dim> will fit in the given length selected_components
+ // and not hang over the end of the vector.
+ Assert (0 <= selected_component,
+ ExcIndexRange (selected_component,0,0));
+ Assert (selected_component + dim - 1 < this->n_components,
+ ExcIndexRange (selected_component, 0, this->n_components));
+}
+
+
+template <int dim>
+VectorFunctionFromTensorFunction<dim>::~VectorFunctionFromTensorFunction ()
+{}
+
+
+template <int dim>
+inline
+double VectorFunctionFromTensorFunction<dim>::value (const Point<dim> &p,
+ const unsigned int component) const
+{
+ Assert (component<this->n_components,
+ ExcIndexRange(component, 0, this->n_components));
+
+ // if the requested component is out of the range selected, then we
+ // can return early
+ if ((component < selected_component)
+ ||
+ (component >= selected_component+dim))
+ return 0;
+
+ // otherwise retrieve the values from the <tt>tensor_function</tt> to be
+ // placed at the <tt>selected_component</tt> to <tt>selected_component + dim - 1</tt>
+ // elements of the <tt>Vector</tt> values and pick the correct one
+ const Tensor<1,dim> tensor_value = tensor_function.value (p);
+
+ return tensor_value[component-selected_component];
+}
+
+
+template <int dim>
+inline
+void VectorFunctionFromTensorFunction<dim>::vector_value (const Point<dim> &p,
+ Vector<double> &values) const
+{
+ Assert(values.size() == this->n_components,
+ ExcDimensionMismatch(values.size(),this->n_components));
+
+ // Retrieve the values from the <tt>tensor_function</tt> to be
+ // placed at the <tt>selected_component</tt> to <tt>selected_component + dim - 1</tt>
+ // elements of the <tt>Vector</tt> values.
+ const Tensor<1,dim> tensor_value = tensor_function.value (p);
+
+ // First we make all elements of values = 0
+ values = 0;
+
+ // Second we adjust the desired components to take on the values in <tt>tensor_value</tt>.
+ for (unsigned int i=0; i<dim; ++i)
+ values(i+selected_component) = tensor_value[i];
+}
+
+
+/** Member function <tt>vector_value_list </tt> is the interface for giving a list of
+ * points (<code>vector<Point<dim> ></code>) of which to evaluate using the <tt>vector_value</tt>
+ * member function. Again, this function is written so as to not replicate the function definition
+ * but passes each point on to <tt>vector_value</tt> to be evaluated.
+ */
+template <int dim>
+void VectorFunctionFromTensorFunction<dim>::vector_value_list (const std::vector<Point<dim> > &points,
+ std::vector<Vector<double> > &value_list) const
+{
+ Assert (value_list.size() == points.size(),
+ ExcDimensionMismatch (value_list.size(), points.size()));
+
+ const unsigned int n_points = points.size();
+
+ for (unsigned int p=0; p<n_points; ++p)
+ VectorFunctionFromTensorFunction<dim>::vector_value (points[p], value_list[p]);
+}
+
+
+
+
// explicit instantiations
template class Function<1>;
template class ComponentSelectFunction<1>;
template class ScalarFunctionFromFunctionObject<1>;
template class VectorFunctionFromScalarFunctionObject<1>;
+template class VectorFunctionFromTensorFunction<1>;
template class Function<2>;
template class ZeroFunction<2>;
template class ComponentSelectFunction<2>;
template class ScalarFunctionFromFunctionObject<2>;
template class VectorFunctionFromScalarFunctionObject<2>;
+template class VectorFunctionFromTensorFunction<2>;
template class Function<3>;
template class ZeroFunction<3>;
template class ComponentSelectFunction<3>;
template class ScalarFunctionFromFunctionObject<3>;
template class VectorFunctionFromScalarFunctionObject<3>;
-
+template class VectorFunctionFromTensorFunction<3>;
DEAL_II_NAMESPACE_CLOSE
--- /dev/null
+//-----------------------------------------------------------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2007, 2008, 2010, 2011, 2013 by the deal.II authors
+//
+// This file is subject to QPL and may not be distributed
+// without copyright and license information. Please refer
+// to the file deal.II/doc/license.html for the text and
+// further information on this license.
+//
+//-----------------------------------------------------------------------------
+
+// Test VectorFunctionTensorFunction
+
+#include "../tests.h"
+#include <deal.II/base/function.h>
+#include <deal.II/base/tensor_function.h>
+#include <deal.II/lac/vector.h>
+
+
+template <int dim>
+class X : public TensorFunction<1,dim>
+{
+public:
+ virtual Tensor<1,dim> value (const Point<dim> &p) const
+ {
+ return p;
+ }
+};
+
+
+template <int dim>
+void check1 ()
+{
+ X<dim> x;
+ VectorFunctionFromTensorFunction<dim>
+ object (x, 1, dim+2);
+
+ Assert (object.n_components == dim+2, ExcInternalError());
+
+ for (unsigned int i=0; i<10; ++i)
+ {
+ Point<dim> p;
+ for (unsigned int d=0; d<dim; ++d)
+ p[d] = i+d;
+
+ for (unsigned int c=0; c<dim+2; ++c)
+ if (c==0 || c==dim+1)
+ Assert (object.value(p,c) == 0,
+ ExcInternalError())
+ else
+ Assert (object.value(p,c) == p[c-1],
+ ExcInternalError());
+
+ Vector<double> v(dim+2);
+ object.vector_value (p, v);
+ for (unsigned int c=0; c<dim+2; ++c)
+ if (c==0 || c==dim+1)
+ Assert (v(c) == 0,
+ ExcInternalError())
+ else
+ Assert (v(c) == p[c-1],
+ ExcInternalError());
+ }
+
+ deallog << "OK" << std::endl;
+}
+
+
+
+
+int main()
+{
+ std::string logname = JobIdentifier::base_name(__FILE__) + std::string("/output");
+ std::ofstream logfile(logname.c_str());
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ check1<1> ();
+ check1<2> ();
+ check1<3> ();
+}
+
+
--- /dev/null
+
+DEAL::OK
+DEAL::OK
+DEAL::OK