const Number * values,
const bool elide_zero_values = false);
- void
- set(const size_type row,
- const size_type n_cols,
- const size_type * col_indices,
- const TrilinosScalar *values,
- const bool elide_zero_values = false);
-
/**
* Add @p value to the element (<i>i,j</i>).
*
* distributed. Otherwise, an exception will be thrown.
*/
template <typename VectorType>
- void
+ typename std::enable_if<std::is_same<typename VectorType::value_type,
+ TrilinosScalar>::value>::type
+ vmult(VectorType &dst, const VectorType &src) const;
+
+ template <typename VectorType>
+ typename std::enable_if<!std::is_same<typename VectorType::value_type,
+ TrilinosScalar>::value>::type
vmult(VectorType &dst, const VectorType &src) const;
/**
* see the discussion about @p VectorType in vmult().
*/
template <typename VectorType>
- void
+ typename std::enable_if<std::is_same<typename VectorType::value_type,
+ TrilinosScalar>::value>::type
+ Tvmult(VectorType &dst, const VectorType &src) const;
+
+ template <typename VectorType>
+ typename std::enable_if<!std::is_same<typename VectorType::value_type,
+ TrilinosScalar>::value>::type
Tvmult(VectorType &dst, const VectorType &src) const;
/**
// Inline the set() and add() functions, since they will be called
// frequently, and the compiler can optimize away some unnecessary loops
// when the sizes are given at compile time.
+ template <>
+ void
+ SparseMatrix::set<TrilinosScalar>(const size_type row,
+ const size_type n_cols,
+ const size_type * col_indices,
+ const TrilinosScalar *values,
+ const bool elide_zero_values);
+
+
+
+ template <typename Number>
+ void
+ SparseMatrix::set(const size_type row,
+ const size_type n_cols,
+ const size_type *col_indices,
+ const Number * values,
+ const bool elide_zero_values)
+ {
+ std::vector<TrilinosScalar> trilinos_values(n_cols);
+ std::copy(values, values + n_cols, trilinos_values.begin());
+ this->set(
+ row, n_cols, col_indices, trilinos_values.data(), elide_zero_values);
+ }
+
+
+
inline void
SparseMatrix::set(const size_type i,
const size_type j,
} // namespace LinearOperatorImplementation
} // namespace internal
+ template <>
+ void
+ SparseMatrix::set<TrilinosScalar>(const size_type row,
+ const size_type n_cols,
+ const size_type * col_indices,
+ const TrilinosScalar *values,
+ const bool elide_zero_values);
# endif // DOXYGEN
} /* namespace TrilinosWrappers */
- template <typename Number>
+ template <>
void
- SparseMatrix::set(const size_type row,
- const size_type n_cols,
- const size_type *col_indices,
- const Number * values,
- const bool elide_zero_values)
- {
- std::vector<TrilinosScalar> trilinos_values(n_cols);
- std::copy(values, values + n_cols, trilinos_values.begin());
- this->set(
- row, n_cols, col_indices, trilinos_values.data(), elide_zero_values);
- }
-
-
-
- void
- SparseMatrix::set(const size_type row,
- const size_type n_cols,
- const size_type * col_indices,
- const TrilinosScalar *values,
- const bool elide_zero_values)
+ SparseMatrix::set<TrilinosScalar>(const size_type row,
+ const size_type n_cols,
+ const size_type * col_indices,
+ const TrilinosScalar *values,
+ const bool elide_zero_values)
{
AssertIndexRange(row, this->m());
template <typename VectorType>
- void
+ typename std::enable_if<
+ std::is_same<typename VectorType::value_type, TrilinosScalar>::value>::type
SparseMatrix::vmult(VectorType &dst, const VectorType &src) const
{
Assert(&src != &dst, ExcSourceEqualsDestination());
template <typename VectorType>
- void
+ typename std::enable_if<
+ !std::is_same<typename VectorType::value_type, TrilinosScalar>::value>::type
+ SparseMatrix::vmult(VectorType & /*dst*/, const VectorType & /*src*/) const
+ {
+ AssertThrow(false, ExcNotImplemented());
+ }
+
+
+
+ template <typename VectorType>
+ typename std::enable_if<
+ std::is_same<typename VectorType::value_type, TrilinosScalar>::value>::type
SparseMatrix::Tvmult(VectorType &dst, const VectorType &src) const
{
Assert(&src != &dst, ExcSourceEqualsDestination());
+ template <typename VectorType>
+ typename std::enable_if<
+ !std::is_same<typename VectorType::value_type, TrilinosScalar>::value>::type
+ SparseMatrix::Tvmult(VectorType & /*dst*/, const VectorType & /*src*/) const
+ {
+ AssertThrow(false, ExcNotImplemented());
+ }
+
+
+
template <typename VectorType>
void
SparseMatrix::vmult_add(VectorType &dst, const VectorType &src) const
dealii::LinearAlgebra::TpetraWrappers::Vector<double> &,
const dealii::LinearAlgebra::TpetraWrappers::Vector<double> &) const;
- /*template void
- SparseMatrix::vmult(
- dealii::LinearAlgebra::TpetraWrappers::Vector<float> &,
- const dealii::LinearAlgebra::TpetraWrappers::Vector<float> &) const;*/
+ template void
+ SparseMatrix::vmult(
+ dealii::LinearAlgebra::TpetraWrappers::Vector<float> &,
+ const dealii::LinearAlgebra::TpetraWrappers::Vector<float> &) const;
template void
SparseMatrix::vmult(
dealii::LinearAlgebra::TpetraWrappers::Vector<double> &,
const dealii::LinearAlgebra::TpetraWrappers::Vector<double> &) const;
+ template void
+ SparseMatrix::Tvmult(
+ dealii::LinearAlgebra::TpetraWrappers::Vector<float> &,
+ const dealii::LinearAlgebra::TpetraWrappers::Vector<float> &) const;
+
template void
SparseMatrix::Tvmult(
dealii::LinearAlgebra::EpetraWrappers::Vector &,
dealii::LinearAlgebra::TpetraWrappers::Vector<double> &,
const dealii::LinearAlgebra::TpetraWrappers::Vector<double> &) const;
+ template void
+ SparseMatrix::vmult_add(
+ dealii::LinearAlgebra::TpetraWrappers::Vector<float> &,
+ const dealii::LinearAlgebra::TpetraWrappers::Vector<float> &) const;
+
template void
SparseMatrix::vmult_add(
dealii::LinearAlgebra::EpetraWrappers::Vector &,
dealii::LinearAlgebra::TpetraWrappers::Vector<double> &,
const dealii::LinearAlgebra::TpetraWrappers::Vector<double> &) const;
+ template void
+ SparseMatrix::Tvmult_add(
+ dealii::LinearAlgebra::TpetraWrappers::Vector<float> &,
+ const dealii::LinearAlgebra::TpetraWrappers::Vector<float> &) const;
+
template void
SparseMatrix::Tvmult_add(
dealii::LinearAlgebra::EpetraWrappers::Vector &,