From: Daniel Arndt Date: Tue, 5 Feb 2019 20:56:33 +0000 (+0100) Subject: Make Tpetra optional X-Git-Tag: v9.1.0-rc1~309^2~16 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=c1356b3c4ad5a1e7503166467a1d5e568eb0711d;p=dealii.git Make Tpetra optional --- diff --git a/cmake/configure/configure_2_trilinos.cmake b/cmake/configure/configure_2_trilinos.cmake index c159fe2552..eca8862a0a 100644 --- a/cmake/configure/configure_2_trilinos.cmake +++ b/cmake/configure/configure_2_trilinos.cmake @@ -42,7 +42,7 @@ MACRO(FEATURE_TRILINOS_FIND_EXTERNAL var) ) FOREACH(_module - Amesos Epetra Ifpack AztecOO Teuchos Tpetra ML MueLu + Amesos Epetra Ifpack AztecOO Teuchos ML MueLu ) ITEM_MATCHES(_module_found ${_module} ${Trilinos_PACKAGE_LIST}) IF(_module_found) @@ -146,7 +146,7 @@ MACRO(FEATURE_TRILINOS_FIND_EXTERNAL var) CHECK_MPI_INTERFACE(TRILINOS ${var}) IF (${var}) - FOREACH(_optional_module EpetraExt ROL Sacado Zoltan) + FOREACH(_optional_module EpetraExt ROL Sacado Tpetra Zoltan) ITEM_MATCHES(_module_found ${_optional_module} ${Trilinos_PACKAGE_LIST}) IF(_module_found) MESSAGE(STATUS "Found ${_optional_module}") @@ -227,8 +227,10 @@ MACRO(FEATURE_TRILINOS_CONFIGURE_EXTERNAL) SET(DEAL_II_EXPAND_TRILINOS_MPI_VECTOR "TrilinosWrappers::MPI::Vector") IF (TRILINOS_WITH_MPI) SET(DEAL_II_EXPAND_EPETRA_VECTOR "LinearAlgebra::EpetraWrappers::Vector") - SET(DEAL_II_EXPAND_TPETRA_VECTOR_DOUBLE "LinearAlgebra::TpetraWrappers::Vector") - SET(DEAL_II_EXPAND_TPETRA_VECTOR_FLOAT "LinearAlgebra::TpetraWrappers::Vector") + IF (${DEAL_II_TRILINOS_WITH_TPETRA}) + SET(DEAL_II_EXPAND_TPETRA_VECTOR_DOUBLE "LinearAlgebra::TpetraWrappers::Vector") + SET(DEAL_II_EXPAND_TPETRA_VECTOR_FLOAT "LinearAlgebra::TpetraWrappers::Vector") + ENDIF() ENDIF() IF(${DEAL_II_TRILINOS_WITH_SACADO}) # Note: Only CMake 3.0 and greater support line continuation with the "\" character diff --git a/include/deal.II/base/config.h.in b/include/deal.II/base/config.h.in index 2369050759..6cf6971a90 100644 --- a/include/deal.II/base/config.h.in +++ b/include/deal.II/base/config.h.in @@ -186,6 +186,7 @@ #cmakedefine DEAL_II_TRILINOS_WITH_EPETRAEXT #cmakedefine DEAL_II_TRILINOS_WITH_ROL #cmakedefine DEAL_II_TRILINOS_WITH_SACADO +#cmakedefine DEAL_II_TRILINOS_WITH_TPETRA #cmakedefine DEAL_II_TRILINOS_WITH_ZOLTAN diff --git a/include/deal.II/base/index_set.h b/include/deal.II/base/index_set.h index afb90fe2d7..4800633ada 100644 --- a/include/deal.II/base/index_set.h +++ b/include/deal.II/base/index_set.h @@ -31,7 +31,9 @@ #ifdef DEAL_II_WITH_TRILINOS # include -# include +# ifdef DEAL_II_TRILINOS_WITH_TPETRA +# include +# endif #endif #if defined(DEAL_II_WITH_MPI) || defined(DEAL_II_WITH_PETSC) @@ -448,9 +450,11 @@ public: make_trilinos_map(const MPI_Comm &communicator = MPI_COMM_WORLD, const bool overlapping = false) const; +# ifdef DEAL_II_TRILINOS_WITH_TPETRA Tpetra::Map make_tpetra_map(const MPI_Comm &communicator = MPI_COMM_WORLD, const bool overlapping = false) const; +# endif #endif diff --git a/include/deal.II/dofs/dof_accessor.templates.h b/include/deal.II/dofs/dof_accessor.templates.h index 9b89e6e31f..d7a097f669 100644 --- a/include/deal.II/dofs/dof_accessor.templates.h +++ b/include/deal.II/dofs/dof_accessor.templates.h @@ -1534,6 +1534,7 @@ namespace internal +# ifdef DEAL_II_TRILINOS_WITH_TPETRA template static void extract_subvector_to( @@ -1559,6 +1560,7 @@ namespace internal for (unsigned int i = 0; i < cache_size; ++i, ++local_values_begin) *local_values_begin = read_write_vector[sorted_indices_pos[i]]; } +# endif diff --git a/include/deal.II/fe/fe_tools_extrapolate.templates.h b/include/deal.II/fe/fe_tools_extrapolate.templates.h index ea76059356..eb210626d7 100644 --- a/include/deal.II/fe/fe_tools_extrapolate.templates.h +++ b/include/deal.II/fe/fe_tools_extrapolate.templates.h @@ -1581,6 +1581,7 @@ namespace FETools # ifdef DEAL_II_WITH_MPI +# ifdef DEAL_II_TRILINOS_WITH_TPETRA template void reinit_distributed(const DoFHandler & dh, @@ -1595,6 +1596,7 @@ namespace FETools const IndexSet &locally_owned_dofs = dh.locally_owned_dofs(); vector.reinit(locally_owned_dofs, parallel_tria->get_communicator()); } +# endif template void diff --git a/include/deal.II/fe/fe_tools_interpolate.templates.h b/include/deal.II/fe/fe_tools_interpolate.templates.h index 9699c0a9b4..a453ae0572 100644 --- a/include/deal.II/fe/fe_tools_interpolate.templates.h +++ b/include/deal.II/fe/fe_tools_interpolate.templates.h @@ -513,6 +513,7 @@ namespace FETools AssertThrow(false, ExcNotImplemented()); } +# ifdef DEAL_II_TRILINOS_WITH_TPETRA template void back_interpolate( @@ -527,6 +528,7 @@ namespace FETools { AssertThrow(false, ExcNotImplemented()); } +# endif #endif diff --git a/include/deal.II/lac/read_write_vector.h b/include/deal.II/lac/read_write_vector.h index b72a95eae2..f752ff815a 100644 --- a/include/deal.II/lac/read_write_vector.h +++ b/include/deal.II/lac/read_write_vector.h @@ -335,6 +335,7 @@ namespace LinearAlgebra std::shared_ptr()); # ifdef DEAL_II_WITH_MPI +# ifdef DEAL_II_TRILINOS_WITH_TPETRA /** * Imports all the elements present in the vector's IndexSet from the input * vector @p tpetra_vec. VectorOperation::values @p operation is used to @@ -349,6 +350,7 @@ namespace LinearAlgebra const std::shared_ptr &communication_pattern = std::shared_ptr()); +# endif /** * Imports all the elements present in the vector's IndexSet from the input @@ -609,6 +611,7 @@ namespace LinearAlgebra protected: #ifdef DEAL_II_WITH_TRILINOS +# ifdef DEAL_II_TRILINOS_WITH_TPETRA /** * Import all the elements present in the vector's IndexSet from the input * vector @p tpetra_vector. This is an helper function and it should not be @@ -622,6 +625,7 @@ namespace LinearAlgebra const MPI_Comm & mpi_comm, const std::shared_ptr &communication_pattern); +# endif /** * Import all the elements present in the vector's IndexSet from the input @@ -656,6 +660,7 @@ namespace LinearAlgebra resize_val(const size_type new_allocated_size); #if defined(DEAL_II_WITH_TRILINOS) && defined(DEAL_II_WITH_MPI) +# ifdef DEAL_II_TRILINOS_WITH_TPETRA /** * Return a TpetraWrappers::CommunicationPattern and store it for future * use. @@ -663,6 +668,7 @@ namespace LinearAlgebra TpetraWrappers::CommunicationPattern create_tpetra_comm_pattern(const IndexSet &source_index_set, const MPI_Comm &mpi_comm); +# endif /** * Return a EpetraWrappers::CommunicationPattern and store it for future diff --git a/include/deal.II/lac/read_write_vector.templates.h b/include/deal.II/lac/read_write_vector.templates.h index 572c28f1de..1444b00856 100644 --- a/include/deal.II/lac/read_write_vector.templates.h +++ b/include/deal.II/lac/read_write_vector.templates.h @@ -506,6 +506,7 @@ namespace LinearAlgebra #if defined(DEAL_II_WITH_TRILINOS) && defined(DEAL_II_WITH_MPI) +# ifdef DEAL_II_TRILINOS_WITH_TPETRA template void ReadWriteVector::import( @@ -595,6 +596,7 @@ namespace LinearAlgebra AssertThrow(false, ExcNotImplemented()); } } +# endif @@ -743,6 +745,7 @@ namespace LinearAlgebra +# ifdef DEAL_II_TRILINOS_WITH_TPETRA template void ReadWriteVector::import( @@ -757,6 +760,7 @@ namespace LinearAlgebra trilinos_vec.get_mpi_communicator(), communication_pattern); } +# endif @@ -901,6 +905,7 @@ namespace LinearAlgebra #if defined(DEAL_II_WITH_TRILINOS) && defined(DEAL_II_WITH_MPI) +# ifdef DEAL_II_TRILINOS_WITH_TPETRA template TpetraWrappers::CommunicationPattern ReadWriteVector::create_tpetra_comm_pattern( @@ -915,6 +920,8 @@ namespace LinearAlgebra return epetra_comm_pattern; } +# endif + template diff --git a/include/deal.II/lac/trilinos_tpetra_communication_pattern.h b/include/deal.II/lac/trilinos_tpetra_communication_pattern.h index fd3e074b2a..45c6f1d34d 100644 --- a/include/deal.II/lac/trilinos_tpetra_communication_pattern.h +++ b/include/deal.II/lac/trilinos_tpetra_communication_pattern.h @@ -19,16 +19,14 @@ #include -#ifdef DEAL_II_WITH_TRILINOS +#if defined(DEAL_II_TRILINOS_WITH_TPETRA) && defined(DEAL_II_WITH_MPI) -# ifdef DEAL_II_WITH_MPI +# include -# include +# include +# include -# include -# include - -# include +# include DEAL_II_NAMESPACE_OPEN @@ -102,8 +100,6 @@ namespace LinearAlgebra DEAL_II_NAMESPACE_CLOSE -# endif - #endif #endif diff --git a/include/deal.II/lac/trilinos_tpetra_vector.h b/include/deal.II/lac/trilinos_tpetra_vector.h index 911863c9d9..a91d4e4b55 100644 --- a/include/deal.II/lac/trilinos_tpetra_vector.h +++ b/include/deal.II/lac/trilinos_tpetra_vector.h @@ -19,7 +19,7 @@ #include -#if defined(DEAL_II_WITH_TRILINOS) && defined(DEAL_II_WITH_MPI) +#if defined(DEAL_II_TRILINOS_WITH_TPETRA) && defined(DEAL_II_WITH_MPI) # include # include diff --git a/include/deal.II/lac/vector_element_access.h b/include/deal.II/lac/vector_element_access.h index 4a81d97fdb..236203e198 100644 --- a/include/deal.II/lac/vector_element_access.h +++ b/include/deal.II/lac/vector_element_access.h @@ -125,6 +125,7 @@ namespace internal +# ifdef DEAL_II_TRILINOS_WITH_TPETRA template <> inline void ElementAccess>::add( @@ -270,6 +271,7 @@ namespace internal // We're going to modify the data on host. return vector_1d(trilinos_i); } +# endif #endif } // namespace internal diff --git a/include/deal.II/multigrid/mg_transfer.h b/include/deal.II/multigrid/mg_transfer.h index 961b0223df..82a82150e8 100644 --- a/include/deal.II/multigrid/mg_transfer.h +++ b/include/deal.II/multigrid/mg_transfer.h @@ -132,6 +132,7 @@ namespace internal }; # ifdef DEAL_II_WITH_MPI +# ifdef DEAL_II_TRILINOS_WITH_TPETRA template struct MatrixSelector> { @@ -163,6 +164,7 @@ namespace internal true); } }; +# endif template <> struct MatrixSelector diff --git a/source/base/index_set.cc b/source/base/index_set.cc index d925745485..9e5cc43c27 100644 --- a/source/base/index_set.cc +++ b/source/base/index_set.cc @@ -519,6 +519,7 @@ IndexSet::fill_index_vector(std::vector &indices) const #ifdef DEAL_II_WITH_TRILINOS +# ifdef DEAL_II_TRILINOS_WITH_TPETRA Tpetra::Map IndexSet::make_tpetra_map(const MPI_Comm &communicator, @@ -527,7 +528,7 @@ IndexSet::make_tpetra_map(const MPI_Comm &communicator, compress(); (void)communicator; -# ifdef DEBUG +# ifdef DEBUG if (!overlapping) { const size_type n_global_elements = @@ -549,7 +550,7 @@ IndexSet::make_tpetra_map(const MPI_Comm &communicator, "by any processor, or there are indices that are " "claimed by multiple processors.")); } -# endif +# endif // Find out if the IndexSet is ascending and 1:1. This corresponds to a // linear Tpetra::Map. Overlapping IndexSets are never 1:1. @@ -560,11 +561,11 @@ IndexSet::make_tpetra_map(const MPI_Comm &communicator, size(), n_elements(), 0, -# ifdef DEAL_II_WITH_MPI +# ifdef DEAL_II_WITH_MPI Teuchos::rcp(new Teuchos::MpiComm(communicator)) -# else +# else Teuchos::rcp(new Teuchos::Comm()) -# endif +# endif ); else { @@ -577,14 +578,15 @@ IndexSet::make_tpetra_map(const MPI_Comm &communicator, size(), arr_view, 0, -# ifdef DEAL_II_WITH_MPI +# ifdef DEAL_II_WITH_MPI Teuchos::rcp(new Teuchos::MpiComm(communicator)) -# else +# else Teuchos::rcp(new Teuchos::Comm()) -# endif +# endif ); } } +# endif @@ -657,7 +659,6 @@ IndexSet::make_trilinos_map(const MPI_Comm &communicator, #endif - bool IndexSet::is_ascending_and_one_to_one(const MPI_Comm &communicator) const { diff --git a/source/lac/trilinos_sparse_matrix.cc b/source/lac/trilinos_sparse_matrix.cc index 250a6d0519..d87c7e88c1 100644 --- a/source/lac/trilinos_sparse_matrix.cc +++ b/source/lac/trilinos_sparse_matrix.cc @@ -104,6 +104,7 @@ namespace TrilinosWrappers return V.trilinos_vector()[0] + V.trilinos_vector().MyLength(); } +# ifdef DEAL_II_TRILINOS_WITH_TPETRA template Number * begin(LinearAlgebra::TpetraWrappers::Vector &V) @@ -133,6 +134,7 @@ namespace TrilinosWrappers return V.trilinos_vector().getData().get() + V.trilinos_vector().getLocalLength(); } +# endif # endif } // namespace internal @@ -3510,6 +3512,7 @@ namespace TrilinosWrappers const dealii::LinearAlgebra::distributed::Vector &) const; # ifdef DEAL_II_WITH_MPI +# ifdef DEAL_II_TRILINOS_WITH_TPETRA template void SparseMatrix::vmult( dealii::LinearAlgebra::TpetraWrappers::Vector &, @@ -3519,6 +3522,7 @@ namespace TrilinosWrappers SparseMatrix::vmult( dealii::LinearAlgebra::TpetraWrappers::Vector &, const dealii::LinearAlgebra::TpetraWrappers::Vector &) const; +# endif template void SparseMatrix::vmult( @@ -3539,6 +3543,7 @@ namespace TrilinosWrappers const dealii::LinearAlgebra::distributed::Vector &) const; # ifdef DEAL_II_WITH_MPI +# ifdef DEAL_II_TRILINOS_WITH_TPETRA template void SparseMatrix::Tvmult( dealii::LinearAlgebra::TpetraWrappers::Vector &, @@ -3548,6 +3553,7 @@ namespace TrilinosWrappers SparseMatrix::Tvmult( dealii::LinearAlgebra::TpetraWrappers::Vector &, const dealii::LinearAlgebra::TpetraWrappers::Vector &) const; +# endif template void SparseMatrix::Tvmult( @@ -3568,6 +3574,7 @@ namespace TrilinosWrappers const dealii::LinearAlgebra::distributed::Vector &) const; # ifdef DEAL_II_WITH_MPI +# ifdef DEAL_II_TRILINOS_WITH_TPETRA template void SparseMatrix::vmult_add( dealii::LinearAlgebra::TpetraWrappers::Vector &, @@ -3577,6 +3584,7 @@ namespace TrilinosWrappers SparseMatrix::vmult_add( dealii::LinearAlgebra::TpetraWrappers::Vector &, const dealii::LinearAlgebra::TpetraWrappers::Vector &) const; +# endif template void SparseMatrix::vmult_add( @@ -3597,6 +3605,7 @@ namespace TrilinosWrappers const dealii::LinearAlgebra::distributed::Vector &) const; # ifdef DEAL_II_WITH_MPI +# ifdef DEAL_II_TRILINOS_WITH_TPETRA template void SparseMatrix::Tvmult_add( dealii::LinearAlgebra::TpetraWrappers::Vector &, @@ -3606,6 +3615,7 @@ namespace TrilinosWrappers SparseMatrix::Tvmult_add( dealii::LinearAlgebra::TpetraWrappers::Vector &, const dealii::LinearAlgebra::TpetraWrappers::Vector &) const; +# endif template void SparseMatrix::Tvmult_add( diff --git a/source/lac/trilinos_tpetra_communication_pattern.cc b/source/lac/trilinos_tpetra_communication_pattern.cc index 626622d49d..99d3aa5f57 100644 --- a/source/lac/trilinos_tpetra_communication_pattern.cc +++ b/source/lac/trilinos_tpetra_communication_pattern.cc @@ -17,7 +17,7 @@ #include -#ifdef DEAL_II_WITH_TRILINOS +#ifdef DEAL_II_TRILINOS_WITH_TPETRA # ifdef DEAL_II_WITH_MPI diff --git a/source/lac/trilinos_tpetra_vector.cc b/source/lac/trilinos_tpetra_vector.cc index 6eb2b6e1f7..428182e916 100644 --- a/source/lac/trilinos_tpetra_vector.cc +++ b/source/lac/trilinos_tpetra_vector.cc @@ -17,7 +17,7 @@ #include -#ifdef DEAL_II_WITH_TRILINOS +#ifdef DEAL_II_TRILINOS_WITH_TPETRA # ifdef DEAL_II_WITH_MPI diff --git a/tests/trilinos/tpetra_vector_01.mpirun=2.output b/tests/trilinos/tpetra_vector_01.mpirun=2.with_trilinos_with_tpetra=on.output similarity index 100% rename from tests/trilinos/tpetra_vector_01.mpirun=2.output rename to tests/trilinos/tpetra_vector_01.mpirun=2.with_trilinos_with_tpetra=on.output diff --git a/tests/trilinos/tpetra_vector_02.mpirun=2.output b/tests/trilinos/tpetra_vector_02.mpirun=2.with_trilinos_with_tpetra=on.output similarity index 100% rename from tests/trilinos/tpetra_vector_02.mpirun=2.output rename to tests/trilinos/tpetra_vector_02.mpirun=2.with_trilinos_with_tpetra=on.output diff --git a/tests/trilinos/tpetra_vector_03.with_complex_values=off.mpirun=2.output b/tests/trilinos/tpetra_vector_03.with_trilinos_with_tpetra=on.with_complex_values=off.mpirun=2.output similarity index 100% rename from tests/trilinos/tpetra_vector_03.with_complex_values=off.mpirun=2.output rename to tests/trilinos/tpetra_vector_03.with_trilinos_with_tpetra=on.with_complex_values=off.mpirun=2.output diff --git a/tests/trilinos/tpetra_vector_03.with_complex_values=off.mpirun=4.output b/tests/trilinos/tpetra_vector_03.with_trilinos_with_tpetra=on.with_complex_values=off.mpirun=4.output similarity index 100% rename from tests/trilinos/tpetra_vector_03.with_complex_values=off.mpirun=4.output rename to tests/trilinos/tpetra_vector_03.with_trilinos_with_tpetra=on.with_complex_values=off.mpirun=4.output diff --git a/tests/trilinos/tpetra_vector_03.with_complex_values=on.mpirun=2.output b/tests/trilinos/tpetra_vector_03.with_trilinos_with_tpetra=on.with_complex_values=on.mpirun=2.output similarity index 100% rename from tests/trilinos/tpetra_vector_03.with_complex_values=on.mpirun=2.output rename to tests/trilinos/tpetra_vector_03.with_trilinos_with_tpetra=on.with_complex_values=on.mpirun=2.output diff --git a/tests/trilinos/tpetra_vector_03.with_complex_values=on.mpirun=4.output b/tests/trilinos/tpetra_vector_03.with_trilinos_with_tpetra=on.with_complex_values=on.mpirun=4.output similarity index 100% rename from tests/trilinos/tpetra_vector_03.with_complex_values=on.mpirun=4.output rename to tests/trilinos/tpetra_vector_03.with_trilinos_with_tpetra=on.with_complex_values=on.mpirun=4.output