From 2a6ff76d2615994544c64294233b3eda811b8d9a Mon Sep 17 00:00:00 2001
From: Wolfgang Bangerth
Date: Fri, 6 Dec 2013 04:10:38 +0000
Subject: [PATCH] In the conversion to 64-bit indices, a lot of classes
acquired a local typedef 'size_type'. In some cases, it was aliased to
types::global_dof_index, but this may have been a mistake. One such example
is FullMatrix, which definitely doesn't need 64-bit indices. Revert this and
go back to 'unsigned int' for the local size_type type.
git-svn-id: https://svn.dealii.org/trunk@31905 0785d39b-7218-0410-832d-ea1e28bc413d
---
deal.II/doc/news/changes.h | 13 ++++
deal.II/include/deal.II/lac/full_matrix.h | 92 +++++++++++------------
deal.II/source/lac/full_matrix.inst.in | 14 ++--
3 files changed, 66 insertions(+), 53 deletions(-)
diff --git a/deal.II/doc/news/changes.h b/deal.II/doc/news/changes.h
index 93c0ecd3a3..de3399598c 100644
--- a/deal.II/doc/news/changes.h
+++ b/deal.II/doc/news/changes.h
@@ -40,6 +40,19 @@ inconvenience this causes.
+ - Changed: During the implementation of the 64-bit features for deal.II
+ 8.0, many linear algebra classes obtained a local typedef
+
size_type
indicating the integer type that is used to index
+ into them. In some instances, this was accidentally set to
+ types::global_dof_index
(which may be a 64-bit data type)
+ even in cases where this is clearly not going to work, for example for
+ FullMatrix::size_type, since we will not be able to store full matrix
+ objects of sizes for which a 32-bit index type is not sufficient. In
+ these cases, the typedef was reverted to just unsigned int
.
+
+ (Wolfgang Bangerth, 2013/12/04)
+
+
- Removed: With the switch of the testsuite to CMake, the old report_features
and build test facilities are removed.
diff --git a/deal.II/include/deal.II/lac/full_matrix.h b/deal.II/include/deal.II/lac/full_matrix.h
index ad8a7fc09a..4829da893e 100644
--- a/deal.II/include/deal.II/lac/full_matrix.h
+++ b/deal.II/include/deal.II/lac/full_matrix.h
@@ -61,17 +61,17 @@ template class LAPACKFullMatrix;
*
* @author Guido Kanschat, Franz-Theo Suttmeier, Wolfgang Bangerth, 1993-2004
*/
-template
+template
class FullMatrix : public Table<2,number>
{
public:
/**
- * Declare type of container size.
+ * 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, unsigned int
is completely sufficient
+ * as an index type.
*/
- // A FullMatrix will never require to use unsigned long long int instead of
- // unsigned but because a ConstraintMatrix can be a SparseMatrix or a
- // FullMatrix, we need have the same interface.
- typedef types::global_dof_index size_type;
+ typedef unsigned int size_type;
/**
* Type of matrix entries. In analogy to
@@ -324,7 +324,7 @@ public:
/**
* Variable assignment operator.
*/
- template
+ template
FullMatrix &
operator = (const FullMatrix &);
@@ -448,8 +448,8 @@ public:
*/
template
void extract_submatrix_from (const MatrixType &matrix,
- const std::vector &row_index_set,
- const std::vector &column_index_set);
+ const std::vector &row_index_set,
+ const std::vector &column_index_set);
/**
* Copy the elements of the current matrix object into a specified
@@ -487,7 +487,7 @@ public:
* determined either by the size
* of this or src.
*/
- template
+ template
void fill (const FullMatrix &src,
const size_type dst_offset_i = 0,
const size_type dst_offset_j = 0,
@@ -499,7 +499,7 @@ public:
* Make function of base class
* available.
*/
- template
+ template
void fill (const number2 *);
/**
@@ -520,7 +520,7 @@ public:
* possible to duplicate rows or
* columns by this method.
*/
- template
+ template
void fill_permutation (const FullMatrix &src,
const std::vector &p_rows,
const std::vector &p_cols);
@@ -612,7 +612,7 @@ public:
* complex-valued, but not mixed, for
* this function to make sense.
*/
- template
+ template
number2 matrix_norm_square (const Vector &v) const;
/**
@@ -630,7 +630,7 @@ public:
* complex-valued, but not mixed, for
* this function to make sense.
*/
- template
+ template
number2 matrix_scalar_product (const Vector &u,
const Vector &v) const;
@@ -825,7 +825,7 @@ public:
* convertible to the data type
* of this matrix.
*/
- template
+ template
void add (const number a,
const FullMatrix &A);
@@ -842,7 +842,7 @@ public:
* convertible to the data type
* of this matrix.
*/
- template
+ template
void add (const number a,
const FullMatrix &A,
const number b,
@@ -861,7 +861,7 @@ public:
* is convertible to the data
* type of this matrix.
*/
- template
+ template
void add (const number a,
const FullMatrix &A,
const number b,
@@ -886,7 +886,7 @@ public:
* size of this or src
* and the given offsets.
*/
- template
+ template
void add (const FullMatrix &src,
const number factor,
const size_type dst_offset_i = 0,
@@ -901,7 +901,7 @@ public:
*
* A += s BT
*/
- template
+ template
void Tadd (const number s,
const FullMatrix &B);
@@ -927,7 +927,7 @@ public:
* of this or
* src.
*/
- template
+ template
void Tadd (const FullMatrix &src,
const number factor,
const size_type dst_offset_i = 0,
@@ -959,10 +959,10 @@ public:
* and are indeed ignored in the
* implementation.
*/
- template
+ template
void add (const size_type row,
- const size_type n_cols,
- const size_type *col_indices,
+ const unsigned int n_cols,
+ const index_type *col_indices,
const number2 *values,
const bool elide_zero_values = true,
const bool col_indices_are_sorted = false);
@@ -1030,7 +1030,7 @@ public:
* Assignment *this =
* a*A.
*/
- template
+ template
void equ (const number a,
const FullMatrix &A);
@@ -1038,7 +1038,7 @@ public:
* Assignment *this = a*A +
* b*B.
*/
- template
+ template
void equ (const number a,
const FullMatrix &A,
const number b,
@@ -1048,7 +1048,7 @@ public:
* Assignment *this = a*A +
* b*B + c*C.
*/
- template
+ template
void equ (const number a,
const FullMatrix &A,
const number b,
@@ -1176,7 +1176,7 @@ public:
* usually results in considerable
* performance gains.
*/
- template
+ template
void mmult (FullMatrix &C,
const FullMatrix &B,
const bool adding=false) const;
@@ -1208,7 +1208,7 @@ public:
* BLAS usually results in considerable
* performance gains.
*/
- template
+ template
void Tmmult (FullMatrix &C,
const FullMatrix &B,
const bool adding=false) const;
@@ -1240,7 +1240,7 @@ public:
* usually results in considerable
* performance gains.
*/
- template
+ template
void mTmult (FullMatrix &C,
const FullMatrix &B,
const bool adding=false) const;
@@ -1273,7 +1273,7 @@ public:
* BLAS usually results in considerable
* performance gains.
*/
- template
+ template
void TmTmult (FullMatrix &C,
const FullMatrix &B,
const bool adding=false) const;
@@ -1322,7 +1322,7 @@ public:
* Source and destination must
* not be the same vector.
*/
- template
+ template
void vmult (Vector &w,
const Vector &v,
const bool adding=false) const;
@@ -1334,7 +1334,7 @@ public:
* Source and destination must
* not be the same vector.
*/
- template
+ template
void vmult_add (Vector &w,
const Vector &v) const;
@@ -1357,7 +1357,7 @@ public:
* Source and destination must
* not be the same vector.
*/
- template
+ template
void Tvmult (Vector &w,
const Vector &v,
const bool adding=false) const;
@@ -1370,7 +1370,7 @@ public:
* Source and destination must
* not be the same vector.
*/
- template
+ template
void Tvmult_add (Vector &w,
const Vector &v) const;
@@ -1398,7 +1398,7 @@ public:
* dst must not be the same
* vector.
*/
- template
+ template
number residual (Vector &dst,
const Vector &x,
const Vector &b) const;
@@ -1420,7 +1420,7 @@ public:
* same object for @p dst and @p
* src.
*/
- template
+ template
void forward (Vector &dst,
const Vector &src) const;
@@ -1434,7 +1434,7 @@ public:
* same object for @p dst and @p
* src.
*/
- template
+ template
void backward (Vector &dst,
const Vector &src) const;
@@ -1576,8 +1576,8 @@ template
inline
void
FullMatrix::extract_submatrix_from (const MatrixType &matrix,
- const std::vector &row_index_set,
- const std::vector &column_index_set)
+ const std::vector &row_index_set,
+ const std::vector &column_index_set)
{
AssertDimension(row_index_set.size(), this->n_rows());
AssertDimension(column_index_set.size(), this->n_cols());
@@ -1627,7 +1627,7 @@ FullMatrix::set (const size_type i,
template
-template
+template
void
FullMatrix::vmult_add (Vector &w,
const Vector &v) const
@@ -1637,7 +1637,7 @@ FullMatrix::vmult_add (Vector &w,
template
-template
+template
void
FullMatrix::Tvmult_add (Vector &w,
const Vector &v) const
@@ -1833,13 +1833,13 @@ FullMatrix::add (const size_type r, const size_type c, const number v)
template
-template
+template
inline
void
-FullMatrix::add (const size_type row,
- const size_type n_cols,
- const size_type *col_indices,
- const number2 *values,
+FullMatrix::add (const size_type row,
+ const unsigned int n_cols,
+ const index_type *col_indices,
+ const number2 *values,
const bool,
const bool)
{
diff --git a/deal.II/source/lac/full_matrix.inst.in b/deal.II/source/lac/full_matrix.inst.in
index f902a72697..7484bbd483 100644
--- a/deal.II/source/lac/full_matrix.inst.in
+++ b/deal.II/source/lac/full_matrix.inst.in
@@ -51,7 +51,7 @@ for (S1, S2 : REAL_SCALARS)
template
void FullMatrix::fill (
- const FullMatrix&, types::global_dof_index, types::global_dof_index, types::global_dof_index, types::global_dof_index);
+ const FullMatrix&, size_type, size_type, size_type, size_type);
template
void FullMatrix::add (const S1, const FullMatrix&);
template
@@ -63,12 +63,12 @@ for (S1, S2 : REAL_SCALARS)
const S1, const FullMatrix&);
template
void FullMatrix::add (
- const FullMatrix&, S1, types::global_dof_index, types::global_dof_index, types::global_dof_index, types::global_dof_index);
+ const FullMatrix&, S1, size_type, size_type, size_type, size_type);
template
void FullMatrix::Tadd (const S1, const FullMatrix&);
template
void FullMatrix::Tadd (
- const FullMatrix&, S1, types::global_dof_index, types::global_dof_index, types::global_dof_index, types::global_dof_index);
+ const FullMatrix&, S1, size_type, size_type, size_type, size_type);
template
void FullMatrix::equ (const S1, const FullMatrix&);
template
@@ -159,7 +159,7 @@ for (S1, S2 : COMPLEX_SCALARS)
{
template
void FullMatrix::fill (
- const FullMatrix&, types::global_dof_index, types::global_dof_index, types::global_dof_index, types::global_dof_index);
+ const FullMatrix&, size_type, size_type, size_type, size_type);
template
void FullMatrix::add (const S1, const FullMatrix&);
template
@@ -171,13 +171,13 @@ for (S1, S2 : COMPLEX_SCALARS)
const S1, const FullMatrix&);
template
void FullMatrix::add (
- const FullMatrix&, S1, types::global_dof_index, types::global_dof_index, types::global_dof_index, types::global_dof_index);
+ const FullMatrix&, S1, size_type, size_type, size_type, size_type);
template
void FullMatrix::Tadd (const S1, const FullMatrix&);
template
void FullMatrix::Tadd (
- const FullMatrix&, S1, types::global_dof_index, types::global_dof_index,
- types::global_dof_index, types::global_dof_index);
+ const FullMatrix&, S1, size_type, size_type,
+ size_type, size_type);
template
void FullMatrix::equ (const S1, const FullMatrix&);
template
--
2.39.5