]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Bugfix: Tensor<rank,dim>: Always store a tensor type internally
authorMatthias Maier <tamiko@43-1.org>
Fri, 4 Sep 2015 04:01:11 +0000 (23:01 -0500)
committerMatthias Maier <tamiko@43-1.org>
Mon, 7 Sep 2015 18:36:22 +0000 (13:36 -0500)
include/deal.II/base/point.h
include/deal.II/base/tensor.h
include/deal.II/base/tensor_base.h
include/deal.II/fe/fe_nedelec.h

index 66c83dbfb85bc431e8145c4eec9a773b6d074bd2..575637d44ddc054344aea2fb03002a47d30af850 100644 (file)
@@ -19,7 +19,7 @@
 
 #include <deal.II/base/config.h>
 #include <deal.II/base/exceptions.h>
-#include <deal.II/base/tensor_base.h>
+#include <deal.II/base/tensor.h>
 #include <cmath>
 
 DEAL_II_NAMESPACE_OPEN
index 53c1314bcac145b92c4d351d0c38187cd3996346..a338636851643032ee44be2bdb1299a933770ef6 100644 (file)
@@ -81,40 +81,6 @@ std::ostream &operator << (std::ostream &out, const Tensor<0,dim,Number> &p)
 }
 
 
-/**
- * Output operator for tensors of rank 1. Print the elements consecutively,
- * with a space in between.
- *
- * @relates Tensor<1,dim,Number>
- */
-template <int dim, typename Number>
-inline
-std::ostream &operator << (std::ostream &out, const Tensor<1,dim,Number> &p)
-{
-  for (unsigned int i=0; i<dim-1; ++i)
-    out << p[i] << ' ';
-  out << p[dim-1];
-
-  return out;
-}
-
-
-/**
- * Output operator for tensors of rank 1 and dimension 1. This is implemented
- * specialized from the general template in order to avoid a compiler warning
- * that the loop is empty.
- *
- * @relates Tensor<1,dim,Number>
- */
-inline
-std::ostream &operator << (std::ostream &out, const Tensor<1,1,double> &p)
-{
-  out << p[0];
-
-  return out;
-}
-
-
 //@}
 /**
  * @name Vector space operations on Tensor objects:
index 9ec99e4123240d4922a8126526600d67da183972..e9a263450e65f8674a72ad25400f7c2acb455b77 100644 (file)
 #define dealii__tensor_base_h
 
 
-// single this file out from tensor.h, since we want to derive
-// Point<dim,Number> from Tensor<1,dim,Number>. However, the point class will
-// not need all the tensor stuff, so we don't want the whole tensor package to
-// be included every time we use a point.
-
-
 #include <deal.II/base/config.h>
 #include <deal.II/base/exceptions.h>
 #include <deal.II/base/table_indices.h>
@@ -115,14 +109,8 @@ public:
   typedef typename numbers::NumberTraits<Number>::real_type real_type;
 
   /**
-   * The tensor type this object represents. In the special case of a
-   * tensor or rank 0 we strip the tensor class and just store the scalar
-   * object.
-   */
-  typedef Number tensor_type;
-
-  /**
-   * Type of stored objects. This is a Number for a rank 0 tensor.
+   * Type of objects encapsulated by this container and returned by
+   * operator[](). This is a scalar number type for a rank 0 tensor.
    */
   typedef Number value_type;
 
@@ -158,19 +146,22 @@ public:
   Tensor (const OtherNumber initializer);
 
   /**
-   * Conversion to Number. Since rank-0 tensors are scalars, this is a natural
-   * operation.
+   * Return a reference to the encapsulated Number object. Since rank-0
+   * tensors are scalars, this is a natural operation.
+   *
+   * This is the non-const conversion operator that returns a writable
+   * reference.
    */
-  operator Number () const;
+  operator Number &();
 
   /**
-   * Conversion to Number. Since rank-0 tensors are scalars, this is a natural
-   * operation.
+   * Return a reference to the encapsulated Number object. Since rank-0
+   * tensors are scalars, this is a natural operation.
    *
-   * This is the non-const conversion operator that returns a writable
+   * This is the const conversion operator that returns a read-only
    * reference.
    */
-  operator Number &();
+  operator const Number &() const;
 
   /**
    * Copy assignment operator.
@@ -279,6 +270,12 @@ public:
                   << "dim must be positive, but was " << arg1);
 
 private:
+  /**
+   * Internal type declaration that is used to specialize the return type
+   * of operator[]() for Tensor<1,dim,Number>
+   */
+  typedef Number tensor_type;
+
   /**
    * The value of this scalar object.
    */
@@ -377,13 +374,9 @@ public:
   typedef typename numbers::NumberTraits<Number>::real_type real_type;
 
   /**
-   * The tensor type this object represents.
-   */
-  typedef Tensor<rank_,dim,Number> tensor_type;
-
-  /**
-   * Type of stored objects (i.e., the object returned by operator[]()).
-   * This is a tensor of lower rank.
+   * Type of objects encapsulated by this container and returned by
+   * operator[](). This is a tensor of lower rank for a general tensor, and
+   * a scalar number type for Tensor<1,dim,Number>.
    */
   typedef typename Tensor<rank_-1,dim,Number>::tensor_type value_type;
 
@@ -605,10 +598,16 @@ public:
                   << "dim must be positive, but was " << arg1);
 
 private:
+  /**
+   * Internal type declaration that is used to specialize the return type
+   * of operator[]() for Tensor<1,dim,Number>
+   */
+  typedef Tensor<rank_, dim, Number> tensor_type;
+
   /**
    * Array of tensors holding the subelements.
    */
-  value_type values[(dim != 0) ? dim : 1];
+  Tensor<rank_-1, dim, Number> values[(dim != 0) ? dim : 1];
 
   /**
    * Help function for unroll.
@@ -678,7 +677,7 @@ Tensor<0,dim,Number>::Tensor (const Tensor<0,dim,OtherNumber> &p)
 
 template <int dim, typename Number>
 inline
-Tensor<0,dim,Number>::operator Number () const
+Tensor<0,dim,Number>::operator Number &()
 {
   return value;
 }
@@ -686,7 +685,7 @@ Tensor<0,dim,Number>::operator Number () const
 
 template <int dim, typename Number>
 inline
-Tensor<0,dim,Number>::operator Number &()
+Tensor<0,dim,Number>::operator const Number &() const
 {
   return value;
 }
@@ -1154,7 +1153,7 @@ Tensor<rank_,dim,Number>::norm_square () const
 {
   real_type s = 0;
   for (unsigned int i=0; i<dim; ++i)
-    s += static_cast<Tensor<rank_-1,dim,Number> >(values[i]).norm_square();
+    s += values[i].norm_square();
 
   return s;
 }
@@ -1181,8 +1180,7 @@ Tensor<rank_, dim, Number>::unroll_recursion (Vector<OtherNumber> &result,
                                               unsigned int        &index) const
 {
   for (unsigned int i=0; i<dim; ++i)
-    static_cast<Tensor<rank_ - 1, dim, Number> >(values[i]).
-    unroll_recursion(result, index);
+    values[i].unroll_recursion(result, index);
 }
 
 
index 3b57774b4c9e4c5de585e610c46f5bdc90342c7a..06f1c56458e66f7ab65eb218572f52deb133c79c 100644 (file)
@@ -19,7 +19,7 @@
 #include <deal.II/base/config.h>
 #include <deal.II/base/table.h>
 #include <deal.II/base/tensor.h>
-#include <deal.II/base/tensor_base.h>
+#include <deal.II/base/tensor.h>
 #include <deal.II/base/polynomials_nedelec.h>
 #include <deal.II/base/polynomial.h>
 #include <deal.II/base/tensor_product_polynomials.h>

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.