19#ifdef DEAL_II_WITH_TRILINOS
27# include <Epetra_Export.h>
89 graph = std::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::make_unique<Epetra_Map>(col_map);
268 if (row_map.Comm().NumProc() > 1)
269 graph = std::make_unique<Epetra_FECrsGraph>(
270 Copy, row_map, n_entries_per_row,
false
278 graph = std::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::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::make_unique<Epetra_FECrsGraph>(
315 Copy, row_map, local_entries_per_row.data(),
false
323 graph = std::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::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::make_unique<Epetra_FECrsGraph>(Copy,
366 n_entries_per_row.data(),
369 graph = std::make_unique<Epetra_FECrsGraph>(
370 Copy, row_map, col_map, n_entries_per_row.data(),
false);
374 std::vector<TrilinosWrappers::types::int_type> row_indices;
378 if (exchange_data ==
false)
379 for (
size_type row = first_row; row < last_row; ++row)
386 row_indices.resize(row_length, -1);
388 typename SparsityPatternType::iterator p = sp.begin(row);
391 for (
int col = 0; col < row_length;)
393 row_indices[col++] = p->column();
394 if (col < row_length)
398 graph->Epetra_CrsGraph::InsertGlobalIndices(row,
403 for (
size_type row = 0; row < sp.n_rows(); ++row)
410 row_indices.resize(row_length, -1);
412 typename SparsityPatternType::iterator p = sp.begin(row);
415 for (
int col = 0; col < row_length;)
417 row_indices[col++] = p->column();
418 if (col < row_length)
423 graph->InsertGlobalIndices(1,
430 const auto &range_map =
431 static_cast<const Epetra_Map &
>(graph->RangeMap());
432 int ierr = graph->GlobalAssemble(*column_space_map, range_map,
true);
435 ierr = graph->OptimizeStorage();
458 const std::vector<size_type> &n_entries_per_row)
470 const IndexSet &col_parallel_partitioning,
490 const IndexSet &col_parallel_partitioning,
492 const std::vector<size_type> &n_entries_per_row)
510 const IndexSet &col_parallel_partitioning,
526 IndexSet nonlocal_partitioner = writable_rows;
528 row_parallel_partitioning.
size());
531 IndexSet tmp = writable_rows & row_parallel_partitioning;
532 Assert(tmp == row_parallel_partitioning,
534 "The set of writable rows passed to this method does not "
535 "contain the locally owned rows, which is not allowed."));
538 nonlocal_partitioner.
subtract_set(row_parallel_partitioning);
541 Epetra_Map nonlocal_map =
544 std::make_unique<Epetra_CrsGraph>(Copy, nonlocal_map, 0);
552 template <
typename SparsityPatternType>
555 const IndexSet & row_parallel_partitioning,
556 const IndexSet & col_parallel_partitioning,
557 const SparsityPatternType &nontrilinos_sparsity_pattern,
559 const bool exchange_data)
567 nontrilinos_sparsity_pattern,
576 template <
typename SparsityPatternType>
579 const IndexSet & parallel_partitioning,
580 const SparsityPatternType &nontrilinos_sparsity_pattern,
582 const bool exchange_data)
588 nontrilinos_sparsity_pattern,
610 graph = std::make_unique<Epetra_FECrsGraph>(*sp.
graph);
620 template <
typename SparsityPatternType>
647 graph = std::make_unique<Epetra_FECrsGraph>(
View,
651 graph->FillComplete();
707 const auto &range_map =
708 static_cast<const Epetra_Map &
>(
graph->RangeMap());
713 ierr =
graph->OptimizeStorage();
722 return graph->RowMap().LID(
747 if (
graph->Filled() ==
false)
749 int nnz_present =
graph->NumGlobalIndices(i);
759 int ierr =
graph->ExtractGlobalRowView(trilinos_i,
764 Assert(nnz_present == nnz_extracted,
768 const std::ptrdiff_t local_col_index =
769 std::find(col_indices, col_indices + nnz_present, trilinos_j) -
772 if (local_col_index == nnz_present)
779 int nnz_present =
graph->NumGlobalIndices(i);
787 graph->ExtractMyRowView(trilinos_i, nnz_extracted, col_indices);
791 Assert(nnz_present == nnz_extracted,
795 const std::ptrdiff_t local_col_index =
796 std::find(col_indices, col_indices + nnz_present, trilinos_j) -
799 if (local_col_index == nnz_present)
814 for (
int i = 0; i < static_cast<int>(
local_size()); ++i)
818 graph->ExtractMyRowView(i, num_entries, indices);
819 for (
unsigned int j = 0; j < static_cast<unsigned int>(num_entries);
848 if (
graph->Filled() ==
true)
868 std::pair<SparsityPattern::size_type, SparsityPattern::size_type>
893 int nnz =
graph->MaxNumIndices();
895 return static_cast<unsigned int>(nnz);
913 return graph->NumMyIndices(local_row);
924 const auto &domain_map =
925 static_cast<const Epetra_Map &
>(
graph->DomainMap());
935 const auto &range_map =
936 static_cast<const Epetra_Map &
>(
graph->RangeMap());
945# ifdef DEAL_II_WITH_MPI
947 const Epetra_MpiComm *mpi_comm =
948 dynamic_cast<const Epetra_MpiComm *
>(&
graph->RangeMap().Comm());
950 return mpi_comm->Comm();
953 return MPI_COMM_SELF;
973 const bool write_extended_trilinos_info)
const
975 if (write_extended_trilinos_info)
982 for (
int i = 0; i <
graph->NumMyRows(); ++i)
984 graph->ExtractMyRowView(i, num_entries, indices);
985 for (
int j = 0; j < num_entries; ++j)
989 <<
") " << std::endl;
1006 graph->ExtractMyRowView(row, num_entries, indices);
1017 out <<
static_cast<int>(
1020 << -
static_cast<int>(
1047 const ::SparsityPattern &,
1052 const ::DynamicSparsityPattern &,
1060 const ::SparsityPattern &,
1066 const ::DynamicSparsityPattern &,
size_type n_elements() const
void subtract_set(const IndexSet &other)
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
std::shared_ptr< const std::vector< size_type > > colnum_cache
SparsityPattern * sparsity_pattern
size_type row_length(const size_type row) const
unsigned int max_entries_per_row() const
size_type bandwidth() const
const_iterator end() const
MPI_Comm get_mpi_communicator() const
bool exists(const size_type i, const size_type j) const
const Epetra_Map & domain_partitioner() const
std::pair< size_type, size_type > local_range() const
size_type n_nonzero_elements() const
::types::global_dof_index size_type
bool is_compressed() const
unsigned int local_size() const
void copy_from(const SparsityPattern &input_sparsity_pattern)
const Epetra_Map & range_partitioner() const
const_iterator begin() const
std::size_t memory_consumption() const
SparsityPattern & operator=(const SparsityPattern &input_sparsity_pattern)
bool row_is_stored_locally(const size_type i) const
void reinit(const size_type m, const size_type n, const size_type n_entries_per_row=0)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcIO()
std::unique_ptr< Epetra_FECrsGraph > graph
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
static ::ExceptionBase & ExcNotImplemented()
void print_gnuplot(std::ostream &out) const
#define Assert(cond, exc)
std::unique_ptr< Epetra_CrsGraph > nonlocal_graph
#define AssertDimension(dim1, dim2)
std::unique_ptr< Epetra_Map > column_space_map
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
IndexSet complete_index_set(const IndexSet::size_type N)
typename ::internal::FEValuesViews::ViewType< dim, spacedim, Extractor >::type View
types::global_dof_index size_type
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
TrilinosWrappers::types::int_type n_global_rows(const Epetra_CrsGraph &graph)
TrilinosWrappers::types::int_type min_my_gid(const Epetra_BlockMap &map)
TrilinosWrappers::types::int_type n_global_elements(const Epetra_BlockMap &map)
TrilinosWrappers::types::int_type max_my_gid(const Epetra_BlockMap &map)
TrilinosWrappers::types::int_type n_global_entries(const Epetra_CrsGraph &graph)
TrilinosWrappers::types::int_type n_global_cols(const Epetra_CrsGraph &graph)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
const Epetra_Comm & comm_self()
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index