19 #ifdef DEAL_II_WITH_TRILINOS
27 # include <Epetra_Export.h>
89 graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
View,
93 graph->FillComplete();
102 reinit(m, n, n_entries_per_row);
110 const std::vector<size_type> &n_entries_per_row)
112 reinit(m, n, n_entries_per_row);
119 , column_space_map(std::move(other.column_space_map))
120 , graph(std::move(other.graph))
121 , nonlocal_graph(std::move(other.nonlocal_graph))
134 new Epetra_FECrsGraph(
View, *column_space_map, *column_space_map, 0))
136 (void)input_sparsity;
139 "Copy constructor only works for empty sparsity patterns."));
148 reinit(parallel_partitioning,
149 parallel_partitioning,
157 const IndexSet & parallel_partitioning,
159 const std::vector<size_type> &n_entries_per_row)
161 reinit(parallel_partitioning,
162 parallel_partitioning,
170 const IndexSet &col_parallel_partitioning,
174 reinit(row_parallel_partitioning,
175 col_parallel_partitioning,
183 const IndexSet & row_parallel_partitioning,
184 const IndexSet & col_parallel_partitioning,
186 const std::vector<size_type> &n_entries_per_row)
188 reinit(row_parallel_partitioning,
189 col_parallel_partitioning,
197 const IndexSet &col_parallel_partitioning,
202 reinit(row_parallel_partitioning,
203 col_parallel_partitioning,
206 n_max_entries_per_row);
227 const std::vector<size_type> &n_entries_per_row)
242 reinit_sp(
const Epetra_Map & row_map,
243 const Epetra_Map & col_map,
245 std::unique_ptr<Epetra_Map> & column_space_map,
246 std::unique_ptr<Epetra_FECrsGraph> &graph,
247 std::unique_ptr<Epetra_CrsGraph> & nonlocal_graph)
249 Assert(row_map.IsOneToOne(),
250 ExcMessage(
"Row map must be 1-to-1, i.e., no overlap between "
251 "the maps of different processors."));
252 Assert(col_map.IsOneToOne(),
253 ExcMessage(
"Column map must be 1-to-1, i.e., no overlap between "
254 "the maps of different processors."));
256 nonlocal_graph.reset();
258 column_space_map = std_cxx14::make_unique<Epetra_Map>(col_map);
268 if (row_map.Comm().NumProc() > 1)
269 graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
270 Copy, row_map, n_entries_per_row,
false
278 graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
279 Copy, row_map, col_map, n_entries_per_row,
false);
285 reinit_sp(
const Epetra_Map & row_map,
286 const Epetra_Map & col_map,
287 const std::vector<size_type> & n_entries_per_row,
288 std::unique_ptr<Epetra_Map> & column_space_map,
289 std::unique_ptr<Epetra_FECrsGraph> &graph,
290 std::unique_ptr<Epetra_CrsGraph> & nonlocal_graph)
292 Assert(row_map.IsOneToOne(),
293 ExcMessage(
"Row map must be 1-to-1, i.e., no overlap between "
294 "the maps of different processors."));
295 Assert(col_map.IsOneToOne(),
296 ExcMessage(
"Column map must be 1-to-1, i.e., no overlap between "
297 "the maps of different processors."));
300 nonlocal_graph.reset();
305 column_space_map = std_cxx14::make_unique<Epetra_Map>(col_map);
306 std::vector<int> local_entries_per_row(
309 for (
unsigned int i = 0; i < local_entries_per_row.size(); ++i)
310 local_entries_per_row[i] =
313 if (row_map.Comm().NumProc() > 1)
314 graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
315 Copy, row_map, local_entries_per_row.data(),
false
323 graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
324 Copy, row_map, col_map, local_entries_per_row.data(),
false);
329 template <
typename SparsityPatternType>
331 reinit_sp(
const Epetra_Map & row_map,
332 const Epetra_Map & col_map,
333 const SparsityPatternType & sp,
334 const bool exchange_data,
335 std::unique_ptr<Epetra_Map> & column_space_map,
336 std::unique_ptr<Epetra_FECrsGraph> &graph,
337 std::unique_ptr<Epetra_CrsGraph> & nonlocal_graph)
339 nonlocal_graph.reset();
347 column_space_map = std_cxx14::make_unique<Epetra_Map>(col_map);
349 Assert(row_map.LinearMap() ==
true,
351 "This function only works if the row map is contiguous."));
355 std::vector<int> n_entries_per_row(last_row - first_row);
359 for (
size_type row = first_row; row < last_row; ++row)
360 n_entries_per_row[row - first_row] =
361 static_cast<int>(sp.row_length(row));
363 if (row_map.Comm().NumProc() > 1)
364 graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
365 Copy, row_map, n_entries_per_row.data(),
false);
367 graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
368 Copy, row_map, col_map, n_entries_per_row.data(),
false);
372 std::vector<TrilinosWrappers::types::int_type> row_indices;
376 if (exchange_data ==
false)
377 for (
size_type row = first_row; row < last_row; ++row)
384 row_indices.resize(row_length, -1);
386 typename SparsityPatternType::iterator p = sp.begin(row);
389 for (
int col = 0; col < row_length;)
391 row_indices[col++] = p->column();
392 if (col < row_length)
396 graph->Epetra_CrsGraph::InsertGlobalIndices(row,
401 for (
size_type row = 0; row < sp.n_rows(); ++row)
408 row_indices.resize(row_length, -1);
410 typename SparsityPatternType::iterator p = sp.begin(row);
413 for (
int col = 0; col < row_length;)
415 row_indices[col++] = p->column();
416 if (col < row_length)
421 graph->InsertGlobalIndices(1,
428 const auto &range_map =
429 static_cast<const Epetra_Map &
>(graph->RangeMap());
430 int ierr = graph->GlobalAssemble(*column_space_map, range_map,
true);
433 ierr = graph->OptimizeStorage();
456 const std::vector<size_type> &n_entries_per_row)
468 const IndexSet &col_parallel_partitioning,
488 const IndexSet &col_parallel_partitioning,
490 const std::vector<size_type> &n_entries_per_row)
508 const IndexSet &col_parallel_partitioning,
524 IndexSet nonlocal_partitioner = writable_rows;
526 row_parallel_partitioning.
size());
529 IndexSet tmp = writable_rows & row_parallel_partitioning;
530 Assert(tmp == row_parallel_partitioning,
532 "The set of writable rows passed to this method does not "
533 "contain the locally owned rows, which is not allowed."));
536 nonlocal_partitioner.
subtract_set(row_parallel_partitioning);
539 Epetra_Map nonlocal_map =
542 std_cxx14::make_unique<Epetra_CrsGraph>(Copy, nonlocal_map, 0);
550 template <
typename SparsityPatternType>
553 const IndexSet & row_parallel_partitioning,
554 const IndexSet & col_parallel_partitioning,
555 const SparsityPatternType &nontrilinos_sparsity_pattern,
557 const bool exchange_data)
565 nontrilinos_sparsity_pattern,
574 template <
typename SparsityPatternType>
577 const IndexSet & parallel_partitioning,
578 const SparsityPatternType &nontrilinos_sparsity_pattern,
580 const bool exchange_data)
586 nontrilinos_sparsity_pattern,
608 graph = std_cxx14::make_unique<Epetra_FECrsGraph>(*sp.
graph);
619 template <
typename SparsityPatternType>
646 graph = std_cxx14::make_unique<Epetra_FECrsGraph>(
View,
650 graph->FillComplete();
706 const auto &range_map =
707 static_cast<const Epetra_Map &
>(
graph->RangeMap());
712 ierr =
graph->OptimizeStorage();
721 return graph->RowMap().LID(
746 if (
graph->Filled() ==
false)
748 int nnz_present =
graph->NumGlobalIndices(i);
758 int ierr =
graph->ExtractGlobalRowView(trilinos_i,
763 Assert(nnz_present == nnz_extracted,
767 const std::ptrdiff_t local_col_index =
768 std::find(col_indices, col_indices + nnz_present, trilinos_j) -
771 if (local_col_index == nnz_present)
778 int nnz_present =
graph->NumGlobalIndices(i);
786 graph->ExtractMyRowView(trilinos_i, nnz_extracted, col_indices);
790 Assert(nnz_present == nnz_extracted,
794 const std::ptrdiff_t local_col_index =
795 std::find(col_indices, col_indices + nnz_present, trilinos_j) -
798 if (local_col_index == nnz_present)
813 for (
int i = 0; i < static_cast<int>(
local_size()); ++i)
817 graph->ExtractMyRowView(i, num_entries, indices);
818 for (
unsigned int j = 0; j < static_cast<unsigned int>(num_entries);
821 if (
static_cast<size_type>(std::abs(i - indices[j])) > local_b)
822 local_b = std::abs(i - indices[j]);
847 if (
graph->Filled() ==
true)
867 std::pair<SparsityPattern::size_type, SparsityPattern::size_type>
892 int nnz =
graph->MaxNumIndices();
894 return static_cast<unsigned int>(nnz);
912 return graph->NumMyIndices(local_row);
923 const auto &domain_map =
924 static_cast<const Epetra_Map &
>(
graph->DomainMap());
934 const auto &range_map =
935 static_cast<const Epetra_Map &
>(
graph->RangeMap());
944 # ifdef DEAL_II_WITH_MPI
946 const Epetra_MpiComm *mpi_comm =
947 dynamic_cast<const Epetra_MpiComm *
>(&
graph->RangeMap().Comm());
949 return mpi_comm->Comm();
952 return MPI_COMM_SELF;
972 const bool write_extended_trilinos_info)
const
974 if (write_extended_trilinos_info)
981 for (
int i = 0; i <
graph->NumMyRows(); ++i)
983 graph->ExtractMyRowView(i, num_entries, indices);
984 for (
int j = 0; j < num_entries; ++j)
988 <<
") " << std::endl;
1005 graph->ExtractMyRowView(row, num_entries, indices);
1016 out <<
static_cast<int>(
1019 << -
static_cast<int>(
1046 const ::SparsityPattern &,
1051 const ::DynamicSparsityPattern &,
1059 const ::SparsityPattern &,
1065 const ::DynamicSparsityPattern &,
1074 #endif // DEAL_II_WITH_TRILINOS