TRILINOS_TPETRA_IS_FUNCTIONAL
)
- reset_cmake_required()
-
- if(NOT TRILINOS_TPETRA_IS_FUNCTIONAL)
+ if(TRILINOS_TPETRA_IS_FUNCTIONAL)
+ check_cxx_symbol_exists(HAVE_TPETRA_INST_FLOAT "TpetraCore_config.h" DEAL_II_HAVE_TPETRA_INST_FLOAT)
+ check_cxx_symbol_exists(HAVE_TPETRA_INST_DOUBLE "TpetraCore_config.h" DEAL_II_HAVE_TPETRA_INST_DOUBLE)
+ check_cxx_symbol_exists(HAVE_TPETRA_INST_COMPLEX_FLOAT "TpetraCore_config.h" DEAL_II_HAVE_TPETRA_INST_COMPLEX_FLOAT)
+ check_cxx_symbol_exists(HAVE_TPETRA_INST_COMPLEX_DOUBLE "TpetraCore_config.h" DEAL_II_HAVE_TPETRA_INST_COMPLEX_DOUBLE)
+ else()
message(
STATUS
"Tpetra was found but is not usable! Disabling Tpetra support."
)
set(TRILINOS_WITH_TPETRA OFF)
endif()
+ reset_cmake_required()
endif()
if(TRILINOS_WITH_MUELU)
set(DEAL_II_EXPAND_TRILINOS_MPI_BLOCKVECTOR "TrilinosWrappers::MPI::BlockVector")
set(DEAL_II_EXPAND_TRILINOS_MPI_VECTOR "TrilinosWrappers::MPI::Vector")
set(DEAL_II_EXPAND_EPETRA_VECTOR "LinearAlgebra::EpetraWrappers::Vector")
+
if(${DEAL_II_TRILINOS_WITH_TPETRA})
- set(DEAL_II_EXPAND_TPETRA_VECTOR_DOUBLE "LinearAlgebra::TpetraWrappers::Vector<double>")
- set(DEAL_II_EXPAND_TPETRA_VECTOR_FLOAT "LinearAlgebra::TpetraWrappers::Vector<float>")
+ if (DEAL_II_HAVE_TPETRA_INST_DOUBLE)
+ set(DEAL_II_EXPAND_TPETRA_VECTOR_DOUBLE "LinearAlgebra::TpetraWrappers::Vector<double>")
+ endif()
+ if (DEAL_II_HAVE_TPETRA_INST_FLOAT)
+ set(DEAL_II_EXPAND_TPETRA_VECTOR_FLOAT "LinearAlgebra::TpetraWrappers::Vector<float>")
+ endif()
if(${DEAL_II_WITH_COMPLEX_NUMBERS})
- set(DEAL_II_EXPAND_TPETRA_VECTOR_COMPLEX_DOUBLE "LinearAlgebra::TpetraWrappers::Vector<std::complex<double>>")
- set(DEAL_II_EXPAND_TPETRA_VECTOR_COMPLEX_FLOAT "LinearAlgebra::TpetraWrappers::Vector<std::complex<float>>")
+ if(DEAL_II_HAVE_TPETRA_INST_COMPLEX_DOUBLE)
+ set(DEAL_II_EXPAND_TPETRA_VECTOR_COMPLEX_DOUBLE "LinearAlgebra::TpetraWrappers::Vector<std::complex<double>>")
+ endif()
+ if(DEAL_II_HAVE_TPETRA_INST_COMPLEX_FLOAT)
+ set(DEAL_II_EXPAND_TPETRA_VECTOR_COMPLEX_FLOAT "LinearAlgebra::TpetraWrappers::Vector<std::complex<float>>")
+ endif()
endif()
endif()
* communication pattern is used multiple times. This can be used to improve
* performance.
*/
- void
+ template <typename Dummy = Number>
+ std::enable_if_t<std::is_same<Dummy, Number>::value &&
+ dealii::is_tpetra_type<Number>::value>
import(const TpetraWrappers::Vector<Number> &tpetra_vec,
VectorOperation::values operation,
const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
* vector @p tpetra_vector. This is an helper function and it should not be
* used directly.
*/
- void
+ template <typename Dummy = Number>
+ std::enable_if_t<std::is_same<Dummy, Number>::value &&
+ dealii::is_tpetra_type<Number>::value>
import(const Tpetra::Vector<Number, int, types::signed_global_dof_index>
& tpetra_vector,
const IndexSet & locally_owned_elements,
#ifdef DEAL_II_WITH_TRILINOS
# ifdef DEAL_II_TRILINOS_WITH_TPETRA
template <typename Number>
- void
+ template <typename Dummy>
+ std::enable_if_t<std::is_same<Dummy, Number>::value &&
+ dealii::is_tpetra_type<Number>::value>
ReadWriteVector<Number>::import(
const Tpetra::Vector<Number, int, types::signed_global_dof_index> &vector,
const IndexSet & source_elements,
# ifdef DEAL_II_TRILINOS_WITH_TPETRA
template <typename Number>
- void
+ template <typename Dummy>
+ std::enable_if_t<std::is_same<Dummy, Number>::value &&
+ dealii::is_tpetra_type<Number>::value>
ReadWriteVector<Number>::import(
const LinearAlgebra::TpetraWrappers::Vector<Number> &trilinos_vec,
VectorOperation::values operation,
DEAL_II_NAMESPACE_OPEN
+/**
+ * Type trait indicating if a certain number type has been explicitly
+ * instantiated in Tpetra. deal.II only supports those number types in Tpetra
+ * wrapper classes.
+ */
+template <typename Number>
+struct is_tpetra_type : std::false_type
+{};
+
+# ifdef HAVE_TPETRA_INST_FLOAT
+template <>
+struct is_tpetra_type<float> : std::true_type
+{};
+# endif
+
+# ifdef HAVE_TPETRA_INST_DOUBLE
+template <>
+struct is_tpetra_type<double> : std::true_type
+{};
+# endif
+
+# ifdef DEAL_II_WITH_COMPLEX_VALUES
+# ifdef HAVE_TPETRA_INST_COMPLEX_FLOAT
+template <>
+struct is_tpetra_type<std::complex<float>> : std::true_type
+{};
+# endif
+
+# ifdef HAVE_TPETRA_INST_COMPLEX_DOUBLE
+template <>
+struct is_tpetra_type<std::complex<double>> : std::true_type
+{};
+# endif
+# endif
+
namespace LinearAlgebra
{
// Forward declaration
INSTANTIATE_DLTG_MATRIX(TrilinosWrappers::BlockSparseMatrix);
# ifndef DOXYGEN
-# ifdef DEAL_II_TRILINOS_WITH_TPETRA
+# if defined(DEAL_II_TRILINOS_WITH_TPETRA) && defined(HAVE_TPETRA_INST_FLOAT)
// FIXME: This mixed variant is needed for multigrid and matrix free.
template void
dealii::AffineConstraints<double>::distribute<
VectorOperation::values,
const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase> &);
# endif
+
+
+
+# ifdef HAVE_TPETRA_INST_FLOAT
+ template void
+ ReadWriteVector<float>::import(
+ const LinearAlgebra::TpetraWrappers::Vector<float> &,
+ VectorOperation::values,
+ const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase> &);
+# endif
+# ifdef HAVE_TPETRA_INST_DOUBLE
+ template void
+ ReadWriteVector<double>::import(
+ const LinearAlgebra::TpetraWrappers::Vector<double> &,
+ VectorOperation::values,
+ const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase> &);
+# endif
+# ifdef DEAL_II_WITH_COMPLEX_VALUES
+# ifdef HAVE_TPETRA_INST_COMPLEX_FLOAT
+ template void
+ ReadWriteVector<std::complex<float>>::import(
+ const LinearAlgebra::TpetraWrappers::Vector<std::complex<float>> &,
+ VectorOperation::values,
+ const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase> &);
+# endif
+# ifdef HAVE_TPETRA_INST_COMPLEX_DOUBLE
+ template void
+ ReadWriteVector<std::complex<double>>::import(
+ const LinearAlgebra::TpetraWrappers::Vector<std::complex<double>> &,
+ VectorOperation::values,
+ const std::shared_ptr<const Utilities::MPI::CommunicationPatternBase> &);
+# endif
+# endif
+
#endif
} // namespace LinearAlgebra
const dealii::LinearAlgebra::distributed::Vector<double> &) const;
# ifdef DEAL_II_TRILINOS_WITH_TPETRA
+# if defined(HAVE_TPETRA_INST_DOUBLE)
template void
SparseMatrix::vmult(
dealii::LinearAlgebra::TpetraWrappers::Vector<double> &,
const dealii::LinearAlgebra::TpetraWrappers::Vector<double> &) const;
+# endif
+# if defined(HAVE_TPETRA_INST_FLOAT)
template void
SparseMatrix::vmult(
dealii::LinearAlgebra::TpetraWrappers::Vector<float> &,
const dealii::LinearAlgebra::TpetraWrappers::Vector<float> &) const;
+# endif
# endif
template void
const dealii::LinearAlgebra::distributed::Vector<double> &) const;
# ifdef DEAL_II_TRILINOS_WITH_TPETRA
+# if defined(HAVE_TPETRA_INST_DOUBLE)
template void
SparseMatrix::Tvmult(
dealii::LinearAlgebra::TpetraWrappers::Vector<double> &,
const dealii::LinearAlgebra::TpetraWrappers::Vector<double> &) const;
+# endif
+# if defined(HAVE_TPETRA_INST_FLOAT)
template void
SparseMatrix::Tvmult(
dealii::LinearAlgebra::TpetraWrappers::Vector<float> &,
const dealii::LinearAlgebra::TpetraWrappers::Vector<float> &) const;
+# endif
# endif
template void
const dealii::LinearAlgebra::distributed::Vector<double> &) const;
# ifdef DEAL_II_TRILINOS_WITH_TPETRA
+# if defined(HAVE_TPETRA_INST_DOUBLE)
template void
SparseMatrix::vmult_add(
dealii::LinearAlgebra::TpetraWrappers::Vector<double> &,
const dealii::LinearAlgebra::TpetraWrappers::Vector<double> &) const;
+# endif
+# if defined(HAVE_TPETRA_INST_FLOAT)
template void
SparseMatrix::vmult_add(
dealii::LinearAlgebra::TpetraWrappers::Vector<float> &,
const dealii::LinearAlgebra::TpetraWrappers::Vector<float> &) const;
+# endif
# endif
template void
const dealii::LinearAlgebra::distributed::Vector<double> &) const;
# ifdef DEAL_II_TRILINOS_WITH_TPETRA
+# if defined(HAVE_TPETRA_INST_DOUBLE)
template void
SparseMatrix::Tvmult_add(
dealii::LinearAlgebra::TpetraWrappers::Vector<double> &,
const dealii::LinearAlgebra::TpetraWrappers::Vector<double> &) const;
+# endif
+# if defined(HAVE_TPETRA_INST_FLOAT)
template void
SparseMatrix::Tvmult_add(
dealii::LinearAlgebra::TpetraWrappers::Vector<float> &,
const dealii::LinearAlgebra::TpetraWrappers::Vector<float> &) const;
+# endif
# endif
template void
{
namespace TpetraWrappers
{
+# ifdef HAVE_TPETRA_INST_FLOAT
template class Vector<float>;
+# endif
+# ifdef HAVE_TPETRA_INST_DOUBLE
template class Vector<double>;
+# endif
# ifdef DEAL_II_WITH_COMPLEX_VALUES
+# ifdef HAVE_TPETRA_INST_COMPLEX_FLOAT
template class Vector<std::complex<float>>;
+# endif
+# ifdef HAVE_TPETRA_INST_COMPLEX_DOUBLE
template class Vector<std::complex<double>>;
+# endif
# endif
} // namespace TpetraWrappers
} // namespace LinearAlgebra