From: Timo Heister Date: Fri, 21 Mar 2025 18:25:10 +0000 (-0400) Subject: Trilinos: remove <13.2 code X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=b6f89559797e8f016711d21829ec5f84cd070e49;p=dealii.git Trilinos: remove <13.2 code --- diff --git a/doc/news/changes/incompatibilities/20250321tjhei b/doc/news/changes/incompatibilities/20250321tjhei new file mode 100644 index 0000000000..041e2dd9bb --- /dev/null +++ b/doc/news/changes/incompatibilities/20250321tjhei @@ -0,0 +1,3 @@ +Changed: deal.II now requires Trilinos version 13.2. +
+(Timo Heister, 2025/03/21) diff --git a/include/deal.II/lac/trilinos_tpetra_sparse_matrix.h b/include/deal.II/lac/trilinos_tpetra_sparse_matrix.h index a8587ec021..c510db9add 100644 --- a/include/deal.II/lac/trilinos_tpetra_sparse_matrix.h +++ b/include/deal.II/lac/trilinos_tpetra_sparse_matrix.h @@ -2045,17 +2045,12 @@ namespace LinearAlgebra colnum_cache->resize(colnums); } -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) typename TpetraTypes::MatrixType:: nonconst_global_inds_host_view_type col_indices(colnum_cache->data(), colnums); typename TpetraTypes::MatrixType:: nonconst_values_host_view_type values(value_cache->data(), colnums); -# else - Teuchos::ArrayView col_indices( - *colnum_cache); - Teuchos::ArrayView values(*value_cache); -# endif + matrix->trilinos_matrix().getGlobalRowCopy(this->a_row, col_indices, values, diff --git a/include/deal.II/lac/trilinos_tpetra_sparse_matrix.templates.h b/include/deal.II/lac/trilinos_tpetra_sparse_matrix.templates.h index d382c19942..a25d8195ef 100644 --- a/include/deal.II/lac/trilinos_tpetra_sparse_matrix.templates.h +++ b/include/deal.II/lac/trilinos_tpetra_sparse_matrix.templates.h @@ -222,24 +222,17 @@ namespace LinearAlgebra for (size_type row = first_row; row < last_row; ++row) n_entries_per_row[row - first_row] = sparsity_pattern.row_length(row); - // The deal.II notation of a Sparsity pattern corresponds to the - // Tpetra concept of a Graph. Hence, we generate a graph by copying - // the sparsity pattern into it, and then build up the matrix from the - // graph. This is considerable faster than directly filling elements - // into the matrix. Moreover, it consumes less memory, since the - // internal reordering is done on ints only, and we can leave the - // doubles aside. -# if DEAL_II_TRILINOS_VERSION_GTE(12, 16, 0) + // The deal.II notation of a Sparsity pattern corresponds to the + // Tpetra concept of a Graph. Hence, we generate a graph by copying + // the sparsity pattern into it, and then build up the matrix from the + // graph. This is considerable faster than directly filling elements + // into the matrix. Moreover, it consumes less memory, since the + // internal reordering is done on ints only, and we can leave the + // doubles aside. Teuchos::RCP> graph = Utilities::Trilinos::internal::make_rcp< TpetraTypes::GraphType>(row_space_map, n_entries_per_row); -# else - Teuchos::RCP> graph = - Utilities::Trilinos::internal::make_rcp< - TpetraTypes::GraphType>( - row_space_map, Teuchos::arcpFromArray(n_entries_per_row)); -# endif // This functions assumes that the sparsity pattern sits on all // processors (completely). The parallel version uses a Tpetra graph @@ -366,17 +359,10 @@ namespace LinearAlgebra // into the matrix. Moreover, it consumes less memory, since the // internal reordering is done on ints only, and we can leave the // doubles aside. -# if DEAL_II_TRILINOS_VERSION_GTE(12, 16, 0) Teuchos::RCP> graph = Utilities::Trilinos::internal::make_rcp< TpetraTypes::GraphType>(row_space_map, n_entries_per_row); -# else - Teuchos::RCP> graph = - Utilities::Trilinos::internal::make_rcp< - TpetraTypes::GraphType>( - row_space_map, Teuchos::arcpFromArray(n_entries_per_row)); -# endif // This functions assumes that the sparsity pattern sits on all // processors (completely). The parallel version uses a Tpetra graph @@ -509,15 +495,7 @@ namespace LinearAlgebra TpetraTypes::MapType>( m, 0, Utilities::Trilinos::tpetra_comm_self()), column_space_map, -# if DEAL_II_TRILINOS_VERSION_GTE(13, 0, 0) - Teuchos::ArrayView{entries_per_row_size_type} -# else - Teuchos::ArrayRCP(entries_per_row_size_type.data(), - 0, - entries_per_row_size_type.size(), - false) -# endif - ); + Teuchos::ArrayView{entries_per_row_size_type}); } @@ -622,15 +600,9 @@ namespace LinearAlgebra { Teuchos::Array n_entries_per_row_array(n_entries_per_row.begin(), n_entries_per_row.end()); -# if DEAL_II_TRILINOS_VERSION_GTE(12, 16, 0) matrix = Utilities::Trilinos::internal::make_rcp< TpetraTypes::MatrixType>(column_space_map, n_entries_per_row_array); -# else - matrix = Utilities::Trilinos::internal::make_rcp< - TpetraTypes::MatrixType>( - column_space_map, Teuchos::arcpFromArray(n_entries_per_row_array)); -# endif } @@ -667,21 +639,13 @@ namespace LinearAlgebra { Teuchos::Array n_entries_per_row_array(n_entries_per_row.begin(), n_entries_per_row.end()); -# if DEAL_II_TRILINOS_VERSION_GTE(12, 16, 0) + matrix = Utilities::Trilinos::internal::make_rcp< TpetraTypes::MatrixType>( row_parallel_partitioning .template make_tpetra_map_rcp>( communicator, false), n_entries_per_row_array); -# else - matrix = Utilities::Trilinos::internal::make_rcp< - TpetraTypes::MatrixType>( - row_parallel_partitioning - .template make_tpetra_map_rcp>( - communicator, false), - Teuchos::arcpFromArray(n_entries_per_row_array)); -# endif } @@ -1287,17 +1251,14 @@ namespace LinearAlgebra size_t nnz = matrix->getNumEntriesInLocalRow(local_row); col_indices_vector.resize(nnz); values_vector.resize(nnz); -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) + typename TpetraTypes::MatrixType:: nonconst_local_inds_host_view_type col_indices( col_indices_vector.data(), nnz); typename TpetraTypes::MatrixType:: nonconst_values_host_view_type values(values_vector.data(), nnz); -# else - Teuchos::ArrayView col_indices(col_indices_vector); - Teuchos::ArrayView values(values_vector); -# endif + matrix->getLocalRowCopy(local_row, col_indices, values, nnz); const size_t diag_index = std::find(col_indices_vector.begin(), @@ -1339,14 +1300,10 @@ namespace LinearAlgebra // not need to perform a deep copy. // Perform a deep copy -# if DEAL_II_TRILINOS_VERSION_GTE(12, 18, 1) matrix = Utilities::Trilinos::internal::make_rcp< TpetraTypes::MatrixType>(*source.matrix, Teuchos::Copy); -# else - matrix = source.matrix->clone(Utilities::Trilinos::internal::make_rcp< - TpetraTypes::NodeType>()); -# endif + column_space_map = Teuchos::rcp_const_cast>( matrix->getColMap()); @@ -1453,16 +1410,11 @@ namespace LinearAlgebra } else { -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) typename TpetraTypes::MatrixType::values_host_view_type values; typename TpetraTypes::MatrixType:: local_inds_host_view_type indices; -# else - Teuchos::ArrayView values; - Teuchos::ArrayView indices; -# endif # if DEAL_II_TRILINOS_VERSION_GTE(14, 0, 0) for (size_t i = 0; i < matrix->getLocalNumRows(); ++i) @@ -1536,18 +1488,12 @@ namespace LinearAlgebra // Prepare pointers for extraction of a view of the row. size_t nnz_present = matrix->getNumEntriesInLocalRow(trilinos_i); -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) + typename TpetraTypes::MatrixType:: nonconst_local_inds_host_view_type col_indices("indices", nnz_present); typename TpetraTypes::MatrixType:: nonconst_values_host_view_type values("values", nnz_present); -# else - std::vector col_indices_vector(nnz_present); - Teuchos::ArrayView col_indices(col_indices_vector); - std::vector values_vector(nnz_present); - Teuchos::ArrayView values(values_vector); -# endif matrix->getLocalRowCopy(trilinos_i, col_indices, values, nnz_present); diff --git a/include/deal.II/lac/trilinos_tpetra_vector.h b/include/deal.II/lac/trilinos_tpetra_vector.h index f44819042c..1dee27c90f 100644 --- a/include/deal.II/lac/trilinos_tpetra_vector.h +++ b/include/deal.II/lac/trilinos_tpetra_vector.h @@ -1129,14 +1129,8 @@ namespace LinearAlgebra // writing to this vector at all. Assert(!has_ghost_elements(), ExcGhostsPresent()); -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) auto vector_2d_local = vector->template getLocalView( Tpetra::Access::ReadWrite); -# else - vector->template sync(); - - auto vector_2d_local = vector->template getLocalView(); -# endif // Having extracted a view into the multivectors above, now also // extract a view into the one vector we actually store. We can diff --git a/include/deal.II/lac/trilinos_tpetra_vector.templates.h b/include/deal.II/lac/trilinos_tpetra_vector.templates.h index eb55f93aaa..3caea75354 100644 --- a/include/deal.II/lac/trilinos_tpetra_vector.templates.h +++ b/include/deal.II/lac/trilinos_tpetra_vector.templates.h @@ -294,22 +294,8 @@ namespace LinearAlgebra { AssertDimension(indices.size(), elements.size()); -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) auto vector_2d = vector->template getLocalView( Tpetra::Access::ReadOnly); -# else - /* - * For Trilinos older than 13.2 we would normally have to call - * vector.template sync() at this place in order - * to sync between memory spaces. This is necessary for GPU support. - * Unfortunately, we are in a const context here and cannot call to - * sync() (which is a non-const member function). - * - * Let us choose to simply ignore this problem for such an old - * Trilinos version. - */ - auto vector_2d = vector->template getLocalView(); -# endif auto vector_1d = Kokkos::subview(vector_2d, Kokkos::ALL(), 0); for (unsigned int i = 0; i < indices.size(); ++i) @@ -355,42 +341,22 @@ namespace LinearAlgebra if (vector->getMap()->isSameAs(*V.vector->getMap())) { // Create a read-only Kokkos view from the source vector -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) auto source_vector_2d = V.vector->template getLocalView( Tpetra::Access::ReadOnly); -# else - auto source_vector_2d = - V.vector->template getLocalView(); -# endif + auto source_vector_1d = Kokkos::subview(source_vector_2d, Kokkos::ALL(), 0); // Create a read/write Kokkos view from the target vector -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) auto target_vector_2d = vector->template getLocalView( Tpetra::Access::ReadWrite); -# else - vector->template sync(); - auto target_vector_2d = - vector->template getLocalView(); -# endif auto target_vector_1d = Kokkos::subview(target_vector_2d, Kokkos::ALL(), 0); -# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) - vector->template modify(); -# endif // Copy the data Kokkos::deep_copy(target_vector_1d, source_vector_1d); - -# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) - vector->template sync::device_type::memory_space>(); -# endif } else if (size() == V.size()) { @@ -540,24 +506,14 @@ namespace LinearAlgebra tpetra_export->getSourceMap()); { -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) auto x_2d = source_vector.template getLocalView( Tpetra::Access::ReadWrite); -# else - source_vector.template sync(); - auto x_2d = source_vector.template getLocalView(); -# endif auto x_1d = Kokkos::subview(x_2d, Kokkos::ALL(), 0); -# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) - source_vector.template modify(); -# endif + const size_t localLength = source_vector.getLocalLength(); auto values_it = V.begin(); for (size_t k = 0; k < localLength; ++k) x_1d(k) = *values_it++; -# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) - source_vector.template sync(); -# endif } Tpetra::CombineMode tpetra_operation = Tpetra::ZERO; if (operation == VectorOperation::insert) @@ -735,25 +691,15 @@ namespace LinearAlgebra // writing to this vector at all. Assert(!has_ghost_elements(), ExcGhostsPresent()); -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) auto vector_2d = vector->template getLocalView( Tpetra::Access::ReadWrite); -# else - vector->template sync(); - auto vector_2d = vector->template getLocalView(); -# endif auto vector_1d = Kokkos::subview(vector_2d, Kokkos::ALL(), 0); -# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) - vector->template modify(); -# endif + const size_t localLength = vector->getLocalLength(); for (size_t k = 0; k < localLength; ++k) { vector_1d(k) += a; } -# if !DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) - vector->template sync(); -# endif } @@ -1010,19 +956,11 @@ namespace LinearAlgebra if (this_local_length != other_local_length) return false; -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) auto this_vector_2d = vector->template getLocalView( Tpetra::Access::ReadOnly); auto other_vector_2d = v.vector->template getLocalView( Tpetra::Access::ReadOnly); -# else - vector->template sync(); - v.vector->template sync(); - auto this_vector_2d = vector->template getLocalView(); - auto other_vector_2d = - v.vector->template getLocalView(); -# endif auto this_vector_1d = Kokkos::subview(this_vector_2d, Kokkos::ALL(), 0); auto other_vector_1d = Kokkos::subview(other_vector_2d, Kokkos::ALL(), 0); @@ -1207,13 +1145,9 @@ namespace LinearAlgebra else out.setf(std::ios::fixed, std::ios::floatfield); -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) auto vector_2d = vector->template getLocalView( Tpetra::Access::ReadOnly); -# else - vector->template sync(); - auto vector_2d = vector->template getLocalView(); -# endif + auto vector_1d = Kokkos::subview(vector_2d, Kokkos::ALL(), 0); const size_t local_length = vector->getLocalLength(); diff --git a/source/lac/trilinos_tpetra_sparsity_pattern.cc b/source/lac/trilinos_tpetra_sparsity_pattern.cc index cc96178063..f3e34922ef 100644 --- a/source/lac/trilinos_tpetra_sparsity_pattern.cc +++ b/source/lac/trilinos_tpetra_sparsity_pattern.cc @@ -59,14 +59,10 @@ namespace LinearAlgebra { // get a representation of the present row std::size_t ncols; -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) typename TpetraTypes::GraphType< MemorySpace>::nonconst_global_inds_host_view_type column_indices_view(colnum_cache->data(), colnum_cache->size()); -# else - Teuchos::ArrayView - column_indices_view(colnum_cache->data(), colnum_cache->size()); -# endif + sparsity_pattern->graph->getGlobalRowCopy(this->a_row, column_indices_view, ncols); @@ -409,7 +405,6 @@ namespace LinearAlgebra "that. Either use more MPI processes or recompile Trilinos " "with 'local ordinate = long long' ")); -# if DEAL_II_TRILINOS_VERSION_GTE(12, 16, 0) if (row_map->getComm()->getSize() > 1) graph = Utilities::Trilinos::internal::make_rcp< TpetraTypes::GraphType>(row_map, n_entries_per_row); @@ -418,17 +413,6 @@ namespace LinearAlgebra TpetraTypes::GraphType>(row_map, col_map, n_entries_per_row); -# else - if (row_map->getComm()->getSize() > 1) - graph = Utilities::Trilinos::internal::make_rcp< - TpetraTypes::GraphType>( - row_map, Teuchos::arcpFromArray(n_entries_per_row)); - else - graph = Utilities::Trilinos::internal::make_rcp< - TpetraTypes::GraphType>( - row_map, col_map, Teuchos::arcpFromArray(n_entries_per_row)); - -# endif AssertDimension(sp.n_rows(), graph->getGlobalNumRows()); AssertDimension(sp.n_cols(), graph->getGlobalNumEntries()); @@ -829,7 +813,6 @@ namespace LinearAlgebra SparsityPattern::exists(const size_type i, const size_type j) const { -# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0) if (!row_is_stored_locally(i)) return false; @@ -851,12 +834,6 @@ namespace LinearAlgebra col_indices.data(); return static_cast(local_col_index) != col_indices.size(); -# else - (void)i; - (void)j; - DEAL_II_NOT_IMPLEMENTED(); - return false; -# endif } diff --git a/tests/trilinos_tpetra/precondition.cc b/tests/trilinos_tpetra/precondition.cc index 0d7cf4a032..2e64edd394 100644 --- a/tests/trilinos_tpetra/precondition.cc +++ b/tests/trilinos_tpetra/precondition.cc @@ -537,7 +537,6 @@ Step4::solve(int cycle) // deallog.pop(); // } -#if DEAL_II_TRILINOS_VERSION_GTE(13, 0, 0) { // Make sure the general Ifpack preconditioner throws an Exception for an // unsupported solver. @@ -558,7 +557,6 @@ Step4::solve(int cycle) } deallog.pop(); } -#endif deallog.pop(); }