]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Change size_type for FullMatrix
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Fri, 2 Feb 2018 10:15:04 +0000 (11:15 +0100)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Sat, 3 Feb 2018 00:42:51 +0000 (01:42 +0100)
include/deal.II/lac/full_matrix.h
include/deal.II/lac/full_matrix.templates.h
source/lac/full_matrix.inst.in

index 9b3e315a5962a2ce3fc25d440b58daee87ae4c9c..f76bdde8580b486a43ba4a2728eed3b6d6e520a4 100644 (file)
@@ -65,12 +65,9 @@ class FullMatrix : public Table<2,number>
 {
 public:
   /**
-   * A type of used to index into this container. Because we can not expect to
-   * store matrices bigger than what can be indexed by a regular unsigned
-   * integer, <code>unsigned int</code> is completely sufficient as an index
-   * type.
+   * A type of used to index into this container.
    */
-  typedef unsigned int size_type;
+  typedef std::size_t size_type;
 
   /**
    * Type of matrix entries. This typedef is analogous to <tt>value_type</tt>
@@ -317,10 +314,10 @@ public:
   template <int dim>
   void
   copy_from (const Tensor<2,dim> &T,
-             const size_type src_r_i=0,
-             const size_type src_r_j=dim-1,
-             const size_type src_c_i=0,
-             const size_type src_c_j=dim-1,
+             const unsigned int src_r_i=0,
+             const unsigned int src_r_j=dim-1,
+             const unsigned int src_c_i=0,
+             const unsigned int src_c_j=dim-1,
              const size_type dst_r=0,
              const size_type dst_c=0);
 
@@ -338,8 +335,8 @@ public:
           const size_type src_r_j=dim-1,
           const size_type src_c_i=0,
           const size_type src_c_j=dim-1,
-          const size_type dst_r=0,
-          const size_type dst_c=0) const;
+          const unsigned int dst_r=0,
+          const unsigned int dst_c=0) const;
 
   /**
    * Copy a subset of the rows and columns of another matrix into the current
@@ -731,7 +728,7 @@ public:
    */
   template <typename number2, typename index_type>
   void add (const size_type     row,
-            const unsigned int  n_cols,
+            const size_type     n_cols,
             const index_type   *col_indices,
             const number2      *values,
             const bool          elide_zero_values = true,
@@ -1475,7 +1472,7 @@ inline
 typename FullMatrix<number>::const_iterator
 FullMatrix<number>::begin (const size_type r) const
 {
-  Assert (r<m(), ExcIndexRange(r,0,m()));
+  AssertIndexRange(r,m());
   return const_iterator(this, r, 0);
 }
 
@@ -1486,7 +1483,7 @@ inline
 typename FullMatrix<number>::const_iterator
 FullMatrix<number>::end (const size_type r) const
 {
-  Assert (r<m(), ExcIndexRange(r,0,m()));
+  AssertIndexRange(r,m());
   return const_iterator(this, r+1, 0);
 }
 
@@ -1509,10 +1506,10 @@ template <typename number>
 template <typename number2, typename index_type>
 inline
 void
-FullMatrix<number>::add (const size_type    row,
-                         const unsigned int n_cols,
-                         const index_type  *col_indices,
-                         const number2     *values,
+FullMatrix<number>::add (const size_type   row,
+                         const size_type   n_cols,
+                         const index_type *col_indices,
+                         const number2    *values,
                          const bool,
                          const bool)
 {
@@ -1536,8 +1533,8 @@ FullMatrix<number>::print (StreamType         &s,
   Assert (!this->empty(), ExcEmptyMatrix());
 
   // save the state of out stream
-  const unsigned int old_precision = s.precision (p);
-  const unsigned int old_width = s.width (w);
+  const std::streamsize old_precision = s.precision (p);
+  const std::streamsize old_width = s.width (w);
 
   for (size_type i=0; i<this->m(); ++i)
     {
index 06375b3af4638ee5ac5a0b8bfa914289bcf54470..8e5433191000c96f551ecfacdfff3af51ab646ed 100644 (file)
@@ -301,7 +301,7 @@ void FullMatrix<number>::backward (Vector<number2>       &dst,
 
   size_type j;
   size_type nu = (m()<n() ? m() : n());
-  for (int i=nu-1; i>=0; --i)
+  for (std::make_signed<size_type>::type i=nu-1; i>=0; --i)
     {
       typename ProductType<number,number2>::type s = src(i);
       for (j=i+1; j<nu; ++j)
@@ -321,14 +321,10 @@ void FullMatrix<number>::fill (const FullMatrix<number2> &src,
                                const size_type src_offset_i,
                                const size_type src_offset_j)
 {
-  Assert (dst_offset_i < m(),
-          ExcIndexRange (dst_offset_i, 0, m()));
-  Assert (dst_offset_j < n(),
-          ExcIndexRange (dst_offset_j, 0, n()));
-  Assert (src_offset_i < src.m(),
-          ExcIndexRange (src_offset_i, 0, src.m()));
-  Assert (src_offset_j < src.n(),
-          ExcIndexRange (src_offset_j, 0, src.n()));
+  AssertIndexRange (dst_offset_i, m());
+  AssertIndexRange (dst_offset_j, n());
+  AssertIndexRange (src_offset_i, src.m());
+  AssertIndexRange (src_offset_j, src.n());
 
   // Compute maximal size of copied block
   const size_type rows = std::min (m() - dst_offset_i,
@@ -1085,14 +1081,10 @@ void FullMatrix<number>::add (const FullMatrix<number2> &src,
                               const size_type src_offset_i,
                               const size_type src_offset_j)
 {
-  Assert (dst_offset_i < m(),
-          ExcIndexRange (dst_offset_i, 0, m()));
-  Assert (dst_offset_j < n(),
-          ExcIndexRange (dst_offset_j, 0, n()));
-  Assert (src_offset_i < src.m(),
-          ExcIndexRange (src_offset_i, 0, src.m()));
-  Assert (src_offset_j < src.n(),
-          ExcIndexRange (src_offset_j, 0, src.n()));
+  AssertIndexRange (dst_offset_i, m());
+  AssertIndexRange (dst_offset_j, n());
+  AssertIndexRange (src_offset_i, src.m());
+  AssertIndexRange (src_offset_j, src.n());
 
   // Compute maximal size of copied block
   const size_type rows = std::min (m() - dst_offset_i, src.m() - src_offset_i);
@@ -1115,14 +1107,10 @@ void FullMatrix<number>::Tadd (const FullMatrix<number2> &src,
                                const size_type src_offset_i,
                                const size_type src_offset_j)
 {
-  Assert (dst_offset_i < m(),
-          ExcIndexRange (dst_offset_i, 0, m()));
-  Assert (dst_offset_j < n(),
-          ExcIndexRange (dst_offset_j, 0, n()));
-  Assert (src_offset_i < src.n(),
-          ExcIndexRange (src_offset_i, 0, src.n()));
-  Assert (src_offset_j < src.m(),
-          ExcIndexRange (src_offset_j, 0, src.m()));
+  AssertIndexRange (dst_offset_i, m());
+  AssertIndexRange (dst_offset_j, n());
+  AssertIndexRange (src_offset_i, src.n());
+  AssertIndexRange (src_offset_j, src.m());
 
   // Compute maximal size of copied block
   const size_type rows = std::min (m() - dst_offset_i, src.n() - src_offset_j);
@@ -1193,7 +1181,11 @@ namespace internal
 #ifdef DEAL_II_WITH_LAPACK
       static number value (const FullMatrix<number> &A)
       {
-        LAPACKFullMatrix<number> lp_A (A.m(), A.n());
+        using s_type = typename LAPACKFullMatrix<number>::size_type;
+        AssertIndexRange (A.m(), std::numeric_limits<s_type>::max()+1);
+        AssertIndexRange (A.n(), std::numeric_limits<s_type>::max()+1);
+        LAPACKFullMatrix<number> lp_A (static_cast<s_type>(A.m()),
+                                       static_cast<s_type>(A.n()));
         lp_A = A;
         lp_A.compute_lu_factorization();
         return lp_A.determinant();
@@ -1618,27 +1610,28 @@ template <typename number>
 template <int dim>
 void
 FullMatrix<number>::copy_from (const Tensor<2,dim> &T,
-                               const size_type src_r_i,
-                               const size_type src_r_j,
-                               const size_type src_c_i,
-                               const size_type src_c_j,
+                               const unsigned int src_r_i,
+                               const unsigned int src_r_j,
+                               const unsigned int src_c_i,
+                               const unsigned int src_c_j,
                                const size_type dst_r,
                                const size_type dst_c)
 {
-
   Assert (!this->empty(), ExcEmptyMatrix());
-  Assert(this->m()-dst_r>src_r_j-src_r_i,
-         ExcIndexRange(this->m()-dst_r,0,src_r_j-src_r_i));
-  Assert(this->n()-dst_c>src_c_j-src_c_i,
-         ExcIndexRange(this->n()-dst_c,0,src_c_j-src_c_i));
-  Assert(dim>src_r_j, ExcIndexRange(dim,0,src_r_j));
-  Assert(dim>src_c_j, ExcIndexRange(dim,0,src_r_j));
-  Assert(src_r_j>=src_r_i, ExcIndexRange(src_r_j,0,src_r_i));
-  Assert(src_c_j>=src_c_i, ExcIndexRange(src_r_j,0,src_r_i));
+  AssertIndexRange(src_r_j-src_r_i, this->m()-dst_r);
+  AssertIndexRange(src_c_j-src_c_i, this->n()-dst_c);
+  AssertIndexRange(src_r_j,dim);
+  AssertIndexRange(src_c_j,dim);
+  AssertIndexRange(src_r_i,src_r_j+1);
+  AssertIndexRange(src_c_i,src_c_j+1);
 
   for (size_type i=0; i<src_r_j-src_r_i+1; i++)
     for (size_type j=0; j<src_c_j-src_c_i+1; j++)
-      (*this)(i+dst_r,j+dst_c) = number(T[i+src_r_i][j+src_c_i]);
+      {
+        const unsigned int src_r_index = static_cast<unsigned int>(i+src_r_i);
+        const unsigned int src_c_index = static_cast<unsigned int>(j+src_c_i);
+        (*this)(i+dst_r,j+dst_c) = number(T[src_r_index][src_c_index]);
+      }
 
 }
 
@@ -1651,23 +1644,24 @@ FullMatrix<number>::copy_to (Tensor<2,dim> &T,
                              const size_type src_r_j,
                              const size_type src_c_i,
                              const size_type src_c_j,
-                             const size_type dst_r,
-                             const size_type dst_c) const
+                             const unsigned int dst_r,
+                             const unsigned int dst_c) const
 {
   Assert (!this->empty(), ExcEmptyMatrix());
-  Assert(dim-dst_r>src_r_j-src_r_i,
-         ExcIndexRange(dim-dst_r,0,src_r_j-src_r_i));
-  Assert(dim-dst_c>src_c_j-src_c_i,
-         ExcIndexRange(dim-dst_c,0,src_c_j-src_c_i));
-  Assert(this->m()>src_r_j, ExcIndexRange(dim,0,src_r_j));
-  Assert(this->n()>src_c_j, ExcIndexRange(dim,0,src_r_j));
-  Assert(src_r_j>=src_r_i, ExcIndexRange(src_r_j,0,src_r_i));
-  Assert(src_c_j>=src_c_i, ExcIndexRange(src_r_j,0,src_r_i));
-
+  AssertIndexRange(src_r_j-src_r_i,dim-dst_r);
+  AssertIndexRange(src_c_j-src_c_i,dim-dst_c);
+  AssertIndexRange(src_r_j,this->m());
+  AssertIndexRange(src_r_j,this->n());
+  AssertIndexRange(src_r_i,src_r_j+1);
+  AssertIndexRange(src_c_j,src_c_j+1);
 
   for (size_type i=0; i<src_r_j-src_r_i+1; i++)
     for (size_type j=0; j<src_c_j-src_c_i+1; j++)
-      T[i+dst_r][j+dst_c] = double ((*this)(i+src_r_i,j+src_c_i));
+      {
+        const unsigned int dst_r_index = static_cast<unsigned int>(i+dst_r);
+        const unsigned int dst_c_index = static_cast<unsigned int>(j+dst_c);
+        T[dst_r_index][dst_c_index] = double ((*this)(i+src_r_i,j+src_c_i));
+      }
 }
 
 
@@ -1712,7 +1706,7 @@ FullMatrix<number>::print_formatted (
   // set output format, but store old
   // state
   std::ios::fmtflags old_flags = out.flags();
-  unsigned int old_precision = out.precision (precision);
+  std::streamsize old_precision = out.precision (precision);
 
   if (scientific)
     {
index 389b72e8fe91715e2ed7cc9d0f1c60fb053fd818..0bb473eb4334f6a5291f4907d38b9ba803286753 100644 (file)
@@ -25,26 +25,26 @@ for (S : REAL_AND_COMPLEX_SCALARS)
         std::ostream&, const unsigned int, const unsigned int) const;
 
     template void FullMatrix<S>::copy_from<1>(
-        const Tensor<2,1>&, const size_type, const size_type, const size_type, const size_type, const size_type, const size_type);
+        const Tensor<2,1>&, const unsigned int, const unsigned int, unsigned int, const unsigned int, const size_type, const size_type);
 
     template void FullMatrix<S>::copy_from<2>(
-        const Tensor<2,2>&, const size_type, const size_type, const size_type, const size_type, const size_type, const size_type);
+        const Tensor<2,2>&, const unsigned int, const unsigned int, unsigned int, const unsigned int, const size_type, const size_type);
 
     template void FullMatrix<S>::copy_from<3>(
-        const Tensor<2,3>&, const size_type, const size_type, const size_type, const size_type, const size_type, const size_type);
+        const Tensor<2,3>&, const unsigned int, const unsigned int, unsigned int, const unsigned int, const size_type, const size_type);
 }
 
 
 for (S : REAL_SCALARS)
 {
     template void FullMatrix<S>::copy_to<1>(
-        Tensor<2,1>&, const size_type, const size_type, const size_type, const size_type, const size_type, const size_type) const;
+        Tensor<2,1>&, const size_type, const size_type, const size_type, const size_type, const unsigned int, const unsigned int) const;
 
     template void FullMatrix<S>::copy_to<2>(
-        Tensor<2,2>&, const size_type, const size_type, const size_type, const size_type, const size_type, const size_type) const;
+        Tensor<2,2>&, const size_type, const size_type, const size_type, const size_type, const unsigned int, const unsigned int) const;
 
     template void FullMatrix<S>::copy_to<3>(
-        Tensor<2,3>&, const size_type, const size_type, const size_type, const size_type, const size_type, const size_type) const;
+        Tensor<2,3>&, const size_type, const size_type, const size_type, const size_type, const unsigned int, const unsigned int) const;
 }
 
 

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.