From: Wolfgang Bangerth Date: Sat, 27 Jun 2015 00:10:27 +0000 (-0500) Subject: Instantiate SparseMatrix for complex scalar types. X-Git-Tag: v8.3.0-rc1~88^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=adbef06366acdca80651c96a31ee3b0d2bd050f5;p=dealii.git Instantiate SparseMatrix for complex scalar types. --- diff --git a/doc/news/changes.h b/doc/news/changes.h index 25191ecb07..d7ba82cc23 100644 --- a/doc/news/changes.h +++ b/doc/news/changes.h @@ -532,12 +532,18 @@ inconvenience this causes.

Specific improvements

    +
  1. Improved: The SparseMatrix class can now also use std::complex + scalars for its elements. +
    + (Wolfgang Bangerth, 2015/06/26) +
  2. +
  3. Improved: FE_DGQArbitraryNodes::get_name() now also detects if the quadrature rule was Gauss points.
    (Guido Kanschat, 2015/06/22) -
  4. - + +
  5. Improved: DoFRenumbering::Cuthill_McKee() can now also use starting indices for parallel triangulations.
    diff --git a/include/deal.II/lac/sparse_matrix.h b/include/deal.II/lac/sparse_matrix.h index 4ae6850bdb..2d86e0b03b 100644 --- a/include/deal.II/lac/sparse_matrix.h +++ b/include/deal.II/lac/sparse_matrix.h @@ -1633,7 +1633,7 @@ SparseMatrix::set (const size_type i, if (index == SparsityPattern::invalid_entry) { Assert ((index != SparsityPattern::invalid_entry) || - (value == 0.), + (value == number()), ExcInvalidIndex(i, j)); return; } @@ -1710,7 +1710,7 @@ SparseMatrix::add (const size_type i, { AssertIsFinite(value); - if (value == 0) + if (value == number()) return; const size_type index = cols->operator()(i, j); @@ -1720,7 +1720,7 @@ SparseMatrix::add (const size_type i, if (index == SparsityPattern::invalid_entry) { Assert ((index != SparsityPattern::invalid_entry) || - (value == 0.), + (value == number()), ExcInvalidIndex(i, j)); return; } @@ -1814,9 +1814,9 @@ SparseMatrix::operator /= (const number factor) { Assert (cols != 0, ExcNotInitialized()); Assert (val != 0, ExcNotInitialized()); - Assert (factor !=0, ExcDivideByZero()); + Assert (factor != number(), ExcDivideByZero()); - const number factor_inv = 1. / factor; + const number factor_inv = number(1.) / factor; number *val_ptr = &val[0]; const number *const end_ptr = &val[cols->n_nonzero_elements()]; diff --git a/include/deal.II/lac/sparse_matrix.templates.h b/include/deal.II/lac/sparse_matrix.templates.h index ee62ee33fc..c3dc4260cb 100644 --- a/include/deal.II/lac/sparse_matrix.templates.h +++ b/include/deal.II/lac/sparse_matrix.templates.h @@ -310,7 +310,7 @@ SparseMatrix::symmetrize () // compute the mean of this // and the transpose value const number mean_value = (*val_ptr + - val[(*cols)(*colnum_ptr,row)]) / 2.0; + val[(*cols)(*colnum_ptr,row)]) / number(2.0); // set this value and the // transpose one to the // mean @@ -354,8 +354,8 @@ SparseMatrix::copy_from (const FullMatrix &matrix) // then copy old matrix for (size_type row=0; row::add (const number factor, const number *const end_ptr = &val[cols->n_nonzero_elements()]; while (val_ptr != end_ptr) - *val_ptr++ += factor **matrix_ptr++; + *val_ptr++ += factor * number(*matrix_ptr++); } @@ -461,7 +461,7 @@ namespace internal typename OutVector::value_type s = 0.; const number *const val_end_of_row = &values[rowstart[row+1]]; while (val_ptr != val_end_of_row) - s += *val_ptr++ * src(*colnum_ptr++); + s += typename OutVector::value_type(*val_ptr++) * typename OutVector::value_type(src(*colnum_ptr++)); *dst_ptr++ = s; } else @@ -470,7 +470,7 @@ namespace internal typename OutVector::value_type s = *dst_ptr; const number *const val_end_of_row = &values[rowstart[row+1]]; while (val_ptr != val_end_of_row) - s += *val_ptr++ * src(*colnum_ptr++); + s += typename OutVector::value_type(*val_ptr++) * typename OutVector::value_type(src(*colnum_ptr++)); *dst_ptr++ = s; } } @@ -535,7 +535,9 @@ SparseMatrix::add (const size_type row, while (this_cols[counter]::add (const size_type row, while (this_cols[counter]::add (const size_type row, while (this_cols[counter]::add (const size_type row, for (size_type j=0; j::add (const size_type row, // value we add is zero if (index == SparsityPattern::invalid_entry) { - Assert (value == 0., ExcInvalidIndex(row,col_indices[j])); + Assert (value == number(), ExcInvalidIndex(row,col_indices[j])); continue; } @@ -645,10 +651,10 @@ SparseMatrix::set (const size_type row, { for (size_type j=0; j::Tvmult (OutVector &dst, for (size_type j=cols->rowstart[i]; jrowstart[i+1] ; j++) { const size_type p = cols->colnums[j]; - dst(p) += val[j] * src(i); + dst(p) += typename OutVector::value_type(val[j]) * typename OutVector::value_type(src(i)); } } } @@ -805,7 +811,7 @@ SparseMatrix::Tvmult_add (OutVector &dst, for (size_type j=cols->rowstart[i]; jrowstart[i+1] ; j++) { const size_type p = cols->colnums[j]; - dst(p) += val[j] * src(i); + dst(p) += typename OutVector::value_type(val[j]) * typename OutVector::value_type(src(i)); } } @@ -824,21 +830,22 @@ namespace internal */ template - number matrix_norm_sqr_on_subrange (const size_type begin_row, - const size_type end_row, - const number *values, - const std::size_t *rowstart, - const size_type *colnums, - const InVector &v) + typename InVector::value_type + matrix_norm_sqr_on_subrange (const size_type begin_row, + const size_type end_row, + const number *values, + const std::size_t *rowstart, + const size_type *colnums, + const InVector &v) { - number norm_sqr=0.; + typename InVector::value_type norm_sqr=0.; for (size_type i=begin_row; i::conjugate(s); + s += typename InVector::value_type(values[j]) * v(colnums[j]); + norm_sqr += v(i) * numbers::NumberTraits::conjugate(s); } return norm_sqr; } @@ -858,7 +865,7 @@ SparseMatrix::matrix_norm_square (const Vector &v) const Assert(n() == v.size(), ExcDimensionMismatch(n(),v.size())); return - parallel::accumulate_from_subranges + parallel::accumulate_from_subranges (std_cxx11::bind (&internal::SparseMatrix::matrix_norm_sqr_on_subrange >, std_cxx11::_1, std_cxx11::_2, @@ -884,22 +891,23 @@ namespace internal */ template - number matrix_scalar_product_on_subrange (const size_type begin_row, - const size_type end_row, - const number *values, - const std::size_t *rowstart, - const size_type *colnums, - const InVector &u, - const InVector &v) + typename InVector::value_type + matrix_scalar_product_on_subrange (const size_type begin_row, + const size_type end_row, + const number *values, + const std::size_t *rowstart, + const size_type *colnums, + const InVector &u, + const InVector &v) { - number norm_sqr=0.; + typename InVector::value_type norm_sqr=0.; for (size_type i=begin_row; i::conjugate(s); + s += typename InVector::value_type(values[j]) * v(colnums[j]); + norm_sqr += u(i) * numbers::NumberTraits::conjugate(s); } return norm_sqr; } @@ -920,7 +928,7 @@ SparseMatrix::matrix_scalar_product (const Vector &u, Assert(n() == v.size(), ExcDimensionMismatch(n(),v.size())); return - parallel::accumulate_from_subranges + parallel::accumulate_from_subranges (std_cxx11::bind (&internal::SparseMatrix::matrix_scalar_product_on_subrange >, std_cxx11::_1, std_cxx11::_2, @@ -1027,7 +1035,7 @@ SparseMatrix::mmult (SparseMatrix &C, const size_type *const end_rows = &sp_A.colnums[sp_A.rowstart[i+1]]; for (; rows != end_rows; ++rows) { - const double A_val = val[rows-&sp_A.colnums[sp_A.rowstart[0]]]; + const number A_val = val[rows-&sp_A.colnums[sp_A.rowstart[0]]]; const size_type col = *rows; const size_type *new_cols = (&sp_B.colnums[sp_B.rowstart[col]]); @@ -1035,9 +1043,9 @@ SparseMatrix::mmult (SparseMatrix &C, // special treatment for diagonal if (sp_B.n_rows() == sp_B.n_cols()) { - C.add (i, *new_cols, A_val * - B.val[new_cols-&sp_B.colnums[sp_B.rowstart[0]]] * - (use_vector ? V(col) : 1)); + C.add (i, *new_cols, + numberC(A_val) * numberC(B.val[new_cols-&sp_B.colnums[sp_B.rowstart[0]]]) * + numberC(use_vector ? V(col) : 1)); ++new_cols; } @@ -1049,7 +1057,7 @@ SparseMatrix::mmult (SparseMatrix &C, &B.val[new_cols-&sp_B.colnums[sp_B.rowstart[0]]]; const numberB *const end_cols = &B.val[sp_B.rowstart[col+1]]; for (; B_val_ptr != end_cols; ++B_val_ptr) - *new_ptr++ = A_val **B_val_ptr * (use_vector ? V(col) : 1); + *new_ptr++ = numberC(A_val) * numberC(*B_val_ptr) * numberC(use_vector ? V(col) : 1); C.add (i, new_ptr-&new_entries[0], new_cols, &new_entries[0], false, true); @@ -1164,13 +1172,13 @@ SparseMatrix::Tmmult (SparseMatrix &C, for (; rows != end_rows; ++rows) { const size_type row = *rows; - const double A_val = val[rows-&sp_A.colnums[sp_A.rowstart[0]]]; + const number A_val = val[rows-&sp_A.colnums[sp_A.rowstart[0]]]; // special treatment for diagonal if (sp_B.n_rows () == sp_B.n_cols()) - C.add (row, i, A_val * - B.val[new_cols-1-&sp_B.colnums[sp_B.rowstart[0]]] * - (use_vector ? V(i) : 1)); + C.add (row, i, + numberC(A_val) * numberC(B.val[new_cols-1-&sp_B.colnums[sp_B.rowstart[0]]]) * + numberC(use_vector ? V(i) : 1)); // now the innermost loop that goes over all the elements in row // 'col' of matrix B. Cache the elements, and then write them into C @@ -1179,7 +1187,7 @@ SparseMatrix::Tmmult (SparseMatrix &C, const numberB *B_val_ptr = &B.val[new_cols-&sp_B.colnums[sp_B.rowstart[0]]]; for (; B_val_ptr != end_cols; ++B_val_ptr) - *new_ptr++ = A_val **B_val_ptr * (use_vector ? V(i) : 1); + *new_ptr++ = numberC(A_val) * numberC(*B_val_ptr) * numberC(use_vector ? V(i) : 1); C.add (row, new_ptr-&new_entries[0], new_cols, &new_entries[0], false, true); @@ -1265,24 +1273,25 @@ namespace internal template - number residual_sqr_on_subrange (const size_type begin_row, - const size_type end_row, - const number *values, - const std::size_t *rowstart, - const size_type *colnums, - const InVector &u, - const InVector &b, - OutVector &dst) + typename OutVector::value_type + residual_sqr_on_subrange (const size_type begin_row, + const size_type end_row, + const number *values, + const std::size_t *rowstart, + const size_type *colnums, + const InVector &u, + const InVector &b, + OutVector &dst) { - number norm_sqr=0.; + typename OutVector::value_type norm_sqr=0.; for (size_type i=begin_row; i::conjugate(s); + norm_sqr += s * numbers::NumberTraits::conjugate(s); } return norm_sqr; } @@ -1306,7 +1315,7 @@ SparseMatrix::residual (Vector &dst, Assert (&u != &dst, ExcSourceEqualsDestination()); return - std::sqrt (parallel::accumulate_from_subranges + std::sqrt (parallel::accumulate_from_subranges (std_cxx11::bind (&internal::SparseMatrix::residual_sqr_on_subrange ,Vector >, std_cxx11::_1, std_cxx11::_2, @@ -1328,7 +1337,7 @@ namespace { #ifdef DEBUG for (typename SparseMatrix::size_type row=0; row::precondition_Jacobi (Vector &dst, // in each row, i.e. at index // rowstart[i]. and we do have a // square matrix by above assertion - if (om != 1.) + if (om != number(1.)) for (size_type i=0; i::precondition_SSOR (Vector &dst, ExcInternalError()); number s = 0; for (size_type j=(*rowstart_ptr)+1; jcolnums[j]); + s += val[j] * number(dst(cols->colnums[j])); // divide by diagonal element *dst_ptr -= s * om; @@ -1443,7 +1452,7 @@ SparseMatrix::precondition_SSOR (Vector &dst, rowstart_ptr = &cols->rowstart[0]; dst_ptr = &dst(0); for ( ; rowstart_ptr!=&cols->rowstart[n]; ++rowstart_ptr, ++dst_ptr) - *dst_ptr *= om*(2.-om)*val[*rowstart_ptr]; + *dst_ptr *= somenumber(om*(number(2.)-om)) * somenumber(val[*rowstart_ptr]); // backward sweep rowstart_ptr = &cols->rowstart[n-1]; @@ -1455,7 +1464,7 @@ SparseMatrix::precondition_SSOR (Vector &dst, = pos_right_of_diagonal[row]; number s = 0; for (size_type j=first_right_of_diagonal_index; jcolnums[j]); + s += val[j] * number(dst(cols->colnums[j])); *dst_ptr -= s * om; *dst_ptr /= val[*rowstart_ptr]; @@ -1486,18 +1495,18 @@ SparseMatrix::precondition_SSOR (Vector &dst, number s = 0; for (size_type j=(*rowstart_ptr)+1; jcolnums[j]); + s += val[j] * number(dst(cols->colnums[j])); // divide by diagonal element *dst_ptr -= s * om; - Assert(val[*rowstart_ptr]!= 0., ExcDivideByZero()); + Assert(val[*rowstart_ptr] != number(), ExcDivideByZero()); *dst_ptr /= val[*rowstart_ptr]; }; rowstart_ptr = &cols->rowstart[0]; dst_ptr = &dst(0); for (size_type row=0; rowrowstart[n-1]; @@ -1512,9 +1521,9 @@ SparseMatrix::precondition_SSOR (Vector &dst, &cols->colnums[0]); number s = 0; for (size_type j=first_right_of_diagonal_index; jcolnums[j]); + s += val[j] * number(dst(cols->colnums[j])); *dst_ptr -= s * om; - Assert(val[*rowstart_ptr]!= 0., ExcDivideByZero()); + Assert(val[*rowstart_ptr] != number(), ExcDivideByZero()); *dst_ptr /= val[*rowstart_ptr]; }; } @@ -1570,10 +1579,10 @@ SparseMatrix::SOR (Vector &dst, { const size_type col = cols->colnums[j]; if (col < row) - s -= val[j] * dst(col); + s -= somenumber(val[j]) * dst(col); } - dst(row) = s * om / val[cols->rowstart[row]]; + dst(row) = s * somenumber(om) / somenumber(val[cols->rowstart[row]]); } } @@ -1597,9 +1606,9 @@ SparseMatrix::TSOR (Vector &dst, somenumber s = dst(row); for (size_type j=cols->rowstart[row]; jrowstart[row+1]; ++j) if (cols->colnums[j] > row) - s -= val[j] * dst(cols->colnums[j]); + s -= somenumber(val[j]) * dst(cols->colnums[j]); - dst(row) = s * om / val[cols->rowstart[row]]; + dst(row) = s * somenumber(om) / somenumber(val[cols->rowstart[row]]); if (row == 0) break; @@ -1639,11 +1648,11 @@ SparseMatrix::PSOR (Vector &dst, const size_type col = cols->colnums[j]; if (inverse_permutation[col] < urow) { - s -= val[j] * dst(col); + s -= somenumber(val[j]) * dst(col); } } - dst(row) = s * om / val[cols->rowstart[row]]; + dst(row) = s * somenumber(om) / somenumber(val[cols->rowstart[row]]); } } @@ -1677,10 +1686,10 @@ SparseMatrix::TPSOR (Vector &dst, { const size_type col = cols->colnums[j]; if (inverse_permutation[col] > urow) - s -= val[j] * dst(col); + s -= somenumber(val[j]) * dst(col); } - dst(row) = s * om / val[cols->rowstart[row]]; + dst(row) = s * somenumber(om) / somenumber(val[cols->rowstart[row]]); } } @@ -1737,9 +1746,9 @@ SparseMatrix::SOR_step (Vector &v, somenumber s = b(row); for (size_type j=cols->rowstart[row]; jrowstart[row+1]; ++j) { - s -= val[j] * v(cols->colnums[j]); + s -= somenumber(val[j]) * v(cols->colnums[j]); } - v(row) += s * om / val[cols->rowstart[row]]; + v(row) += s * somenumber(om) / somenumber(val[cols->rowstart[row]]); } } @@ -1765,9 +1774,9 @@ SparseMatrix::TSOR_step (Vector &v, somenumber s = b(row); for (size_type j=cols->rowstart[row]; jrowstart[row+1]; ++j) { - s -= val[j] * v(cols->colnums[j]); + s -= somenumber(val[j]) * v(cols->colnums[j]); } - v(row) += s * om / val[cols->rowstart[row]]; + v(row) += s * somenumber(om) / somenumber(val[cols->rowstart[row]]); } } @@ -1814,11 +1823,12 @@ SparseMatrix::SSOR (Vector &dst, const size_type p = cols->colnums[j]; if (p != SparsityPattern::invalid_entry) { - if (i>j) s += val[j] * dst(p); + if (i>j) + s += somenumber(val[j]) * dst(p); } } - dst(i) -= s * om; - dst(i) /= val[cols->rowstart[i]]; + dst(i) -= s * somenumber(om); + dst(i) /= somenumber(val[cols->rowstart[i]]); } for (int i=n-1; i>=0; i--) // this time, i is signed, but always positive! @@ -1829,10 +1839,11 @@ SparseMatrix::SSOR (Vector &dst, const size_type p = cols->colnums[j]; if (p != SparsityPattern::invalid_entry) { - if (static_cast(i)(i) < j) + s += somenumber(val[j]) * dst(p); } } - dst(i) -= s * om / val[cols->rowstart[i]]; + dst(i) -= s * somenumber(om) / somenumber(val[cols->rowstart[i]]); } } @@ -1882,7 +1893,7 @@ void SparseMatrix::print_formatted (std::ostream &out, for (size_type j=0; joperator()(i,j)] * denominator << ' '; + << val[cols->operator()(i,j)] * number(denominator) << ' '; else out << std::setw(width) << zero_string << ' '; out << std::endl; diff --git a/source/lac/sparse_matrix.inst.in b/source/lac/sparse_matrix.inst.in index 0e8e7ea5f0..950bf0dfba 100644 --- a/source/lac/sparse_matrix.inst.in +++ b/source/lac/sparse_matrix.inst.in @@ -1,6 +1,6 @@ // --------------------------------------------------------------------- // -// Copyright (C) 1998 - 2014 by the deal.II authors +// Copyright (C) 1998 - 2015 by the deal.II authors // // This file is part of the deal.II library. // @@ -151,105 +151,132 @@ for (S1, S2, S3: REAL_SCALARS) // complex instantiations -// for (S : COMPLEX_SCALARS) -// { -// template class SparseMatrix; -// } - - - -// for (S1, S2 : COMPLEX_SCALARS) -// { -// template SparseMatrix & -// SparseMatrix::copy_from (const SparseMatrix &); - -// template -// void SparseMatrix::copy_from (const FullMatrix &); - -// template void SparseMatrix::add (const S1, -// const SparseMatrix &); -// } - - -// for (S1, S2 : COMPLEX_SCALARS) -// { -// template S2 -// SparseMatrix:: -// matrix_norm_square (const Vector &) const; - -// template S2 -// SparseMatrix:: -// matrix_scalar_product (const Vector &, -// const Vector &) const; - -// template S2 SparseMatrix:: -// residual (Vector &, -// const Vector &, -// const Vector &) const; - -// template void SparseMatrix:: -// precondition_SSOR (Vector &, -// const Vector &, -// const S1) const; - -// template void SparseMatrix:: -// precondition_SOR (Vector &, -// const Vector &, -// const S1) const; - -// template void SparseMatrix:: -// precondition_TSOR (Vector &, -// const Vector &, -// const S1) const; - -// template void SparseMatrix:: -// precondition_Jacobi (Vector &, -// const Vector &, -// const S1) const; - -// template void SparseMatrix:: -// SOR (Vector &, -// const S1) const; -// template void SparseMatrix:: -// TSOR (Vector &, -// const S1) const; -// template void SparseMatrix:: -// SSOR (Vector &, -// const S1) const; -// template void SparseMatrix:: -// PSOR (Vector &, -// const std::vector&, -// const std::vector&, -// const S1) const; -// template void SparseMatrix:: -// TPSOR (Vector &, -// const std::vector&, -// const std::vector&, -// const S1) const; -// template void SparseMatrix:: -// SOR_step (Vector &, -// const Vector &, -// const S1) const; -// template void SparseMatrix:: -// TSOR_step (Vector &, -// const Vector &, -// const S1) const; -// template void SparseMatrix:: -// SSOR_step (Vector &, -// const Vector &, -// const S1) const; -// } - - -// for (S1, S2, S3 : COMPLEX_SCALARS; -// V1, V2 : DEAL_II_VEC_TEMPLATES) -// { -// template void SparseMatrix:: -// vmult (V1 &, const V2 &) const; -// template void SparseMatrix:: -// Tvmult (V1 &, const V2 &) const; -// template void SparseMatrix:: -// vmult_add (V1 &, const V2 &) const; -// template void SparseMatrix:: -// Tvmult_add (V1 &, const V2 &) const; -// } +for (S : COMPLEX_SCALARS) + { + template class SparseMatrix; + } + + + +for (S1, S2 : COMPLEX_SCALARS) + { + template SparseMatrix & + SparseMatrix::copy_from (const SparseMatrix &); + + template + void SparseMatrix::copy_from (const FullMatrix &); + + template void SparseMatrix::add (const S1, + const SparseMatrix &); + + template void SparseMatrix::add (const size_type, + const size_type, + const size_type *, + const S2 *, + const bool, + const bool); + + template void SparseMatrix::set (const size_type, + const size_type, + const size_type *, + const S2 *, + const bool); + } + + +for (S1, S2 : COMPLEX_SCALARS) + { + template S2 + SparseMatrix:: + matrix_norm_square (const Vector &) const; + + template S2 + SparseMatrix:: + matrix_scalar_product (const Vector &, + const Vector &) const; + + template S2 SparseMatrix:: + residual (Vector &, + const Vector &, + const Vector &) const; + + template void SparseMatrix:: + precondition_SSOR (Vector &, + const Vector &, + const S1, + const std::vector&) const; + + template void SparseMatrix:: + precondition_SOR (Vector &, + const Vector &, + const S1) const; + + template void SparseMatrix:: + precondition_TSOR (Vector &, + const Vector &, + const S1) const; + + template void SparseMatrix:: + precondition_Jacobi (Vector &, + const Vector &, + const S1) const; + + template void SparseMatrix:: + SOR (Vector &, + const S1) const; + template void SparseMatrix:: + TSOR (Vector &, + const S1) const; + template void SparseMatrix:: + SSOR (Vector &, + const S1) const; + template void SparseMatrix:: + PSOR (Vector &, + const std::vector&, + const std::vector&, + const S1) const; + template void SparseMatrix:: + TPSOR (Vector &, + const std::vector&, + const std::vector&, + const S1) const; + template void SparseMatrix:: + Jacobi_step (Vector &, + const Vector &, + const S1) const; + template void SparseMatrix:: + SOR_step (Vector &, + const Vector &, + const S1) const; + template void SparseMatrix:: + TSOR_step (Vector &, + const Vector &, + const S1) const; + template void SparseMatrix:: + SSOR_step (Vector &, + const Vector &, + const S1) const; + } + +for (S1, S2, S3 : COMPLEX_SCALARS; + V1, V2 : DEAL_II_VEC_TEMPLATES) + { + template void SparseMatrix:: + vmult (V1 &, const V2 &) const; + template void SparseMatrix:: + Tvmult (V1 &, const V2 &) const; + template void SparseMatrix:: + vmult_add (V1 &, const V2 &) const; + template void SparseMatrix:: + Tvmult_add (V1 &, const V2 &) const; + } + +for (S1, S2, S3: COMPLEX_SCALARS) + { + template void SparseMatrix:: + mmult (SparseMatrix &, const SparseMatrix &, const Vector&, + const bool) const; + template void SparseMatrix:: + Tmmult (SparseMatrix &, const SparseMatrix &, const Vector&, + const bool) const; + } diff --git a/source/lac/vector_memory.inst.in b/source/lac/vector_memory.inst.in index 08981d5cea..37fc030807 100644 --- a/source/lac/vector_memory.inst.in +++ b/source/lac/vector_memory.inst.in @@ -1,6 +1,6 @@ // --------------------------------------------------------------------- // -// Copyright (C) 1999 - 2014 by the deal.II authors +// Copyright (C) 1999 - 2015 by the deal.II authors // // This file is part of the deal.II library. // @@ -20,3 +20,12 @@ for (VECTOR : SERIAL_VECTORS) template class VectorMemory; template class GrowingVectorMemory; } + + for (SCALAR : COMPLEX_SCALARS) + { + template class VectorMemory >; + template class GrowingVectorMemory >; + + template class VectorMemory >; + template class GrowingVectorMemory >; + }