class Quadrature : public Subscriptor
{
public:
+ /**
+ * Define a typedef for a
+ * quadrature that acts on an
+ * object of one dimension
+ * less. For cells, this would
+ * then be a face quadrature.
+ */
+ typedef Quadrature<dim-1> SubQuadrature;
+
/**
* Number of quadrature points.
*/
* less than the present and a
* formula in one dimension.
*/
- Quadrature (const Quadrature<dim-1> &,
- const Quadrature<1> &);
+ Quadrature (const SubQuadrature &,
+ const Quadrature<1> &);
/**
* Construct a quadrature formula
class QProjector
{
public:
+ /**
+ * Define a typedef for a
+ * quadrature that acts on an
+ * object of one dimension
+ * less. For cells, this would
+ * then be a face quadrature.
+ */
+ typedef Quadrature<dim-1> SubQuadrature;
+
/**
* Compute the quadrature points
* on the cell if the given
* details, see the general doc
* for this class.
*/
- static void project_to_face (const Quadrature<dim-1> &quadrature,
- const unsigned int face_no,
+ static void project_to_face (const SubQuadrature &quadrature,
+ const unsigned int face_no,
typename std::vector<Point<dim> > &q_points);
/**
* further details, see the
* general doc for this class.
*/
- static void project_to_subface (const Quadrature<dim-1> &quadrature,
- const unsigned int face_no,
- const unsigned int subface_no,
+ static void project_to_subface (const SubQuadrature &quadrature,
+ const unsigned int face_no,
+ const unsigned int subface_no,
typename std::vector<Point<dim> > &q_points);
/**
* later.
*/
static Quadrature<dim>
- project_to_all_faces (const Quadrature<dim-1> &quadrature);
+ project_to_all_faces (const SubQuadrature &quadrature);
/**
* This function is alike the
* faces of the unit cell.
*/
static Quadrature<dim>
- project_to_all_subfaces (const Quadrature<dim-1> &quadrature);
+ project_to_all_subfaces (const SubQuadrature &quadrature);
/**
* Project a give quadrature
public Subscriptor
{
public:
+ /**
+ * Define typedefs for the return
+ * types of the @p{value} and the
+ * @p{gradient} functions.
+ */
+ typedef Tensor<rank,dim> value_type;
+ typedef Tensor<rank+1,dim> gradient_type;
+
/**
* Constructor. May take an
* initial value for the time
* Return the value of the function
* at the given point.
*/
- virtual Tensor<rank, dim> value (const Point<dim> &p) const;
+ virtual value_type value (const Point<dim> &p) const;
/**
* Set @p{values} to the point
* size as the @p{points} array.
*/
virtual void value_list (const typename std::vector<Point<dim> > &points,
- typename std::vector<Tensor<rank,dim> > &values) const;
+ typename std::vector<value_type> &values) const;
/**
* Return the gradient of the
* function at the given point.
*/
- virtual Tensor<rank+1,dim> gradient (const Point<dim> &p) const;
+ virtual gradient_type gradient (const Point<dim> &p) const;
/**
* Set @p{gradients} to the
* size as the @p{points} array.
*/
virtual void gradient_list (const typename std::vector<Point<dim> > &points,
- typename std::vector<Tensor<rank+1,dim> > &gradients) const;
+ typename std::vector<gradient_type> &gradients) const;
/**
* Exception
template <int dim>
-Quadrature<dim>::Quadrature (const Quadrature<dim-1> &q1,
- const Quadrature<1> &q2) :
+Quadrature<dim>::Quadrature (const SubQuadrature &q1,
+ const Quadrature<1> &q2) :
n_quadrature_points (q1.n_quadrature_points *
q2.n_quadrature_points),
quadrature_points (n_quadrature_points),
template <int dim>
Quadrature<dim>
-QProjector<dim>::project_to_all_faces (const Quadrature<dim-1> &quadrature)
+QProjector<dim>::project_to_all_faces (const SubQuadrature &quadrature)
{
const unsigned int n_points = quadrature.n_quadrature_points,
n_faces = GeometryInfo<dim>::faces_per_cell;
template <int dim>
Quadrature<dim>
-QProjector<dim>::project_to_all_subfaces (const Quadrature<dim-1> &quadrature)
+QProjector<dim>::project_to_all_subfaces (const SubQuadrature &quadrature)
{
const unsigned int n_points = quadrature.n_quadrature_points,
n_faces = GeometryInfo<dim>::faces_per_cell,
template <int rank, int dim>
-Tensor<rank,dim>
+typename TensorFunction<rank, dim>::value_type
TensorFunction<rank, dim>::value (const Point<dim> &) const
{
Assert (false, ExcPureFunctionCalled());
template <int rank, int dim>
void
TensorFunction<rank, dim>::value_list (const typename std::vector<Point<dim> > &points,
- typename std::vector<Tensor<rank,dim> > &values) const
+ typename std::vector<value_type> &values) const
{
Assert (values.size() == points.size(),
ExcDimensionMismatch(values.size(), points.size()));
template <int rank, int dim>
-Tensor<rank+1,dim>
+typename TensorFunction<rank, dim>::gradient_type
TensorFunction<rank, dim>::gradient (const Point<dim> &) const
{
Assert (false, ExcPureFunctionCalled());
template <int rank, int dim>
void
TensorFunction<rank, dim>::gradient_list (const typename std::vector<Point<dim> > &points,
- typename std::vector<Tensor<rank+1,dim> > &gradients) const
+ typename std::vector<gradient_type> &gradients) const
{
Assert (gradients.size() == points.size(),
ExcDimensionMismatch(gradients.size(), points.size()));
<h3>base</h3>
<ol>
+ <li> <p>
+ Fixed: The class <code class="class">TensorFunction</code>
+ now uses local types <code class="class">value_type</code> and
+ <code class="class">gradient_type</code> as return values of
+ its member functions. This works around a bug in Sun's Forte
+ C++ compilers.
+ <br>
+ (WB 2002/03/20)
+ </p>
+
<li> <p>
Improved: The <code class="member">AssertThrow</code> macro now
uses <code class="member">__builtin_expect</code> if the