16 #include <deal.II/lac/trilinos_sparsity_pattern.h> 18 #ifdef DEAL_II_WITH_TRILINOS 20 # include <deal.II/base/utilities.h> 21 # include <deal.II/lac/sparsity_pattern.h> 22 # include <deal.II/lac/dynamic_sparsity_pattern.h> 24 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
25 # include <Epetra_Export.h> 28 DEAL_II_NAMESPACE_OPEN
38 #ifndef DEAL_II_WITH_64BIT_INDICES 39 int n_global_elements (
const Epetra_BlockMap &map)
41 return map.NumGlobalElements();
44 int min_my_gid(
const Epetra_BlockMap &map)
46 return map.MinMyGID();
49 int max_my_gid(
const Epetra_BlockMap &map)
51 return map.MaxMyGID();
54 int n_global_rows(
const Epetra_CrsGraph &graph)
56 return graph.NumGlobalRows();
59 int n_global_cols(
const Epetra_CrsGraph &graph)
61 return graph.NumGlobalCols();
64 int n_global_entries(
const Epetra_CrsGraph &graph)
66 return graph.NumGlobalEntries();
69 int global_row_index(
const Epetra_CrsGraph &graph,
int i)
74 long long int n_global_elements (
const Epetra_BlockMap &map)
76 return map.NumGlobalElements64();
79 long long int min_my_gid(
const Epetra_BlockMap &map)
81 return map.MinMyGID64();
84 long long int max_my_gid(
const Epetra_BlockMap &map)
86 return map.MaxMyGID64();
89 long long int n_global_rows(
const Epetra_CrsGraph &graph)
91 return graph.NumGlobalRows64();
94 long long int n_global_cols(
const Epetra_CrsGraph &graph)
96 return graph.NumGlobalCols64();
99 long long int n_global_entries(
const Epetra_CrsGraph &graph)
101 return graph.NumGlobalEntries64();
104 long long int global_row_index(
const Epetra_CrsGraph &graph,
int i)
106 return graph.GRID64(i);
141 (TrilinosWrappers::types::int_type *)&(*
colnum_cache)[0]);
150 Assert (ncols != 0, ExcInternalError());
151 colnum_cache.reset (
new std::vector<size_type> (colnums,
170 column_space_map.reset(
new Epetra_Map (TrilinosWrappers::types::int_type(0),
171 TrilinosWrappers::types::int_type(0),
173 graph.reset (
new Epetra_FECrsGraph(View,
177 graph->FillComplete();
184 reinit (input_map, input_map, n_entries_per_row);
190 const std::vector<size_type> &n_entries_per_row)
192 reinit (input_map, input_map, n_entries_per_row);
198 const Epetra_Map &input_col_map,
201 reinit (input_row_map, input_col_map, n_entries_per_row);
207 const Epetra_Map &input_col_map,
208 const std::vector<size_type> &n_entries_per_row)
210 reinit (input_row_map, input_col_map, n_entries_per_row);
219 reinit (m, n, n_entries_per_row);
226 const std::vector<size_type> &n_entries_per_row)
228 reinit (m, n, n_entries_per_row);
240 graph (new Epetra_FECrsGraph(View,
245 (void)input_sparsity;
247 ExcMessage (
"Copy constructor only works for empty sparsity patterns."));
253 const MPI_Comm &communicator,
256 reinit (parallel_partitioning, parallel_partitioning, communicator,
263 const MPI_Comm &communicator,
264 const std::vector<size_type> &n_entries_per_row)
266 reinit (parallel_partitioning, parallel_partitioning, communicator,
273 const IndexSet &col_parallel_partitioning,
274 const MPI_Comm &communicator,
277 reinit (row_parallel_partitioning, col_parallel_partitioning,
278 communicator, n_entries_per_row);
285 const IndexSet &col_parallel_partitioning,
286 const MPI_Comm &communicator,
287 const std::vector<size_type> &n_entries_per_row)
289 reinit (row_parallel_partitioning, col_parallel_partitioning,
290 communicator, n_entries_per_row);
297 const IndexSet &col_parallel_partitioning,
299 const MPI_Comm &communicator,
302 reinit (row_parallel_partitioning, col_parallel_partitioning,
303 writable_rows, communicator, n_max_entries_per_row);
318 reinit (complete_index_set(m), complete_index_set(n), MPI_COMM_SELF,
327 const std::vector<size_type> &n_entries_per_row)
329 reinit (complete_index_set(m), complete_index_set(n), MPI_COMM_SELF,
340 reinit_sp (
const Epetra_Map &row_map,
341 const Epetra_Map &col_map,
342 const size_type n_entries_per_row,
344 std_cxx11::shared_ptr<Epetra_FECrsGraph> &graph,
347 Assert(row_map.IsOneToOne(),
348 ExcMessage(
"Row map must be 1-to-1, i.e., no overlap between " 349 "the maps of different processors."));
350 Assert(col_map.IsOneToOne(),
351 ExcMessage(
"Column map must be 1-to-1, i.e., no overlap between " 352 "the maps of different processors."));
354 nonlocal_graph.reset();
356 column_space_map.reset (
new Epetra_Map (col_map));
366 if (row_map.Comm().NumProc() > 1)
367 graph.reset (
new Epetra_FECrsGraph(Copy, row_map,
368 n_entries_per_row,
false 378 graph.reset (
new Epetra_FECrsGraph(Copy, row_map, col_map,
379 n_entries_per_row,
false));
385 reinit_sp (
const Epetra_Map &row_map,
386 const Epetra_Map &col_map,
387 const std::vector<size_type> &n_entries_per_row,
388 std_cxx11::shared_ptr<Epetra_Map> &column_space_map,
389 std_cxx11::shared_ptr<Epetra_FECrsGraph> &graph,
390 std_cxx11::shared_ptr<Epetra_CrsGraph> &nonlocal_graph)
392 Assert(row_map.IsOneToOne(),
393 ExcMessage(
"Row map must be 1-to-1, i.e., no overlap between " 394 "the maps of different processors."));
395 Assert(col_map.IsOneToOne(),
396 ExcMessage(
"Column map must be 1-to-1, i.e., no overlap between " 397 "the maps of different processors."));
400 nonlocal_graph.reset();
403 static_cast<size_type
>(n_global_elements(row_map)));
405 column_space_map.reset (
new Epetra_Map (col_map));
406 std::vector<int> local_entries_per_row(max_my_gid(row_map)-
407 min_my_gid(row_map));
408 for (
unsigned int i=0; i<local_entries_per_row.size(); ++i)
409 local_entries_per_row[i] = n_entries_per_row[min_my_gid(row_map)+i];
411 if (row_map.Comm().NumProc() > 1)
412 graph.reset(
new Epetra_FECrsGraph(Copy, row_map,
413 &local_entries_per_row[0],
424 graph.reset(
new Epetra_FECrsGraph(Copy, row_map, col_map,
425 &local_entries_per_row[0],
431 template <
typename SparsityPatternType>
433 reinit_sp (
const Epetra_Map &row_map,
434 const Epetra_Map &col_map,
435 const SparsityPatternType &sp,
436 const bool exchange_data,
437 std_cxx11::shared_ptr<Epetra_Map> &column_space_map,
438 std_cxx11::shared_ptr<Epetra_FECrsGraph> &graph,
439 std_cxx11::shared_ptr<Epetra_CrsGraph> &nonlocal_graph)
441 nonlocal_graph.reset ();
445 static_cast<size_type
>(n_global_elements(row_map)));
447 static_cast<size_type
>(n_global_elements(col_map)));
449 column_space_map.reset (
new Epetra_Map (col_map));
451 Assert (row_map.LinearMap() ==
true,
452 ExcMessage (
"This function only works if the row map is contiguous."));
454 const size_type first_row = min_my_gid(row_map),
455 last_row = max_my_gid(row_map)+1;
456 std::vector<int> n_entries_per_row(last_row - first_row);
460 for (size_type row=first_row; row<last_row; ++row)
461 n_entries_per_row[row-first_row] = static_cast<int>(sp.row_length(row));
463 if (row_map.Comm().NumProc() > 1)
464 graph.reset(
new Epetra_FECrsGraph(Copy, row_map,
465 &n_entries_per_row[0],
468 graph.reset (
new Epetra_FECrsGraph(Copy, row_map, col_map,
469 &n_entries_per_row[0],
473 static_cast<size_type
>(n_global_rows(*graph)));
475 std::vector<TrilinosWrappers::types::int_type> row_indices;
479 if (exchange_data==
false)
480 for (size_type row=first_row; row<last_row; ++row)
482 const TrilinosWrappers::types::int_type
row_length = sp.row_length(row);
486 row_indices.resize (row_length, -1);
488 typename SparsityPatternType::iterator p = sp.begin(row);
493 row_indices[col++] = p->column();
494 if (col < row_length)
498 graph->Epetra_CrsGraph::InsertGlobalIndices (row, row_length,
502 for (size_type row=0; row<sp.n_rows(); ++row)
504 const TrilinosWrappers::types::int_type
row_length = sp.row_length(row);
508 row_indices.resize (row_length, -1);
510 typename SparsityPatternType::iterator p = sp.begin(row);
515 row_indices[col++] = p->column();
516 if (col < row_length)
520 graph->InsertGlobalIndices (1,
521 reinterpret_cast<TrilinosWrappers::types::int_type *>(&row),
522 row_length, &row_indices[0]);
526 graph->GlobalAssemble (*column_space_map,
527 static_cast<const Epetra_Map &>(graph->RangeMap()),
531 ierr = graph->OptimizeStorage ();
541 reinit_sp (input_map, input_map, n_entries_per_row,
549 const Epetra_Map &input_col_map,
552 reinit_sp (input_row_map, input_col_map, n_entries_per_row,
560 const std::vector<size_type> &n_entries_per_row)
562 reinit_sp (input_map, input_map, n_entries_per_row,
570 const Epetra_Map &input_col_map,
571 const std::vector<size_type> &n_entries_per_row)
573 reinit_sp (input_row_map, input_col_map, n_entries_per_row,
581 const MPI_Comm &communicator,
586 reinit_sp (map, map, n_entries_per_row,
593 const MPI_Comm &communicator,
594 const std::vector<size_type> &n_entries_per_row)
598 reinit_sp (map, map, n_entries_per_row,
605 const IndexSet &col_parallel_partitioning,
606 const MPI_Comm &communicator,
613 reinit_sp (row_map, col_map, n_entries_per_row,
621 const IndexSet &col_parallel_partitioning,
622 const MPI_Comm &communicator,
623 const std::vector<size_type> &n_entries_per_row)
629 reinit_sp (row_map, col_map, n_entries_per_row,
637 const IndexSet &col_parallel_partitioning,
639 const MPI_Comm &communicator,
646 reinit_sp (row_map, col_map, n_entries_per_row,
649 IndexSet nonlocal_partitioner = writable_rows;
653 IndexSet tmp = writable_rows & row_parallel_partitioning;
654 Assert (tmp == row_parallel_partitioning,
655 ExcMessage(
"The set of writable rows passed to this method does not " 656 "contain the locally owned rows, which is not allowed."));
659 nonlocal_partitioner.
subtract_set(row_parallel_partitioning);
662 Epetra_Map nonlocal_map =
672 template<
typename SparsityPatternType>
675 const IndexSet &col_parallel_partitioning,
676 const SparsityPatternType &nontrilinos_sparsity_pattern,
677 const MPI_Comm &communicator,
678 const bool exchange_data)
684 reinit_sp (row_map, col_map, nontrilinos_sparsity_pattern, exchange_data,
690 template<
typename SparsityPatternType>
693 const SparsityPatternType &nontrilinos_sparsity_pattern,
694 const MPI_Comm &communicator,
695 const bool exchange_data)
699 reinit_sp (map, map, nontrilinos_sparsity_pattern, exchange_data,
705 template <
typename SparsityPatternType>
708 const SparsityPatternType &sp,
709 const bool exchange_data)
711 reinit_sp (input_map, input_map, sp, exchange_data,
717 template <
typename SparsityPatternType>
720 const Epetra_Map &input_col_map,
721 const SparsityPatternType &sp,
722 const bool exchange_data)
724 reinit_sp (input_row_map, input_col_map, sp, exchange_data,
735 Assert (
false, ExcNotImplemented());
741 template <
typename SparsityPatternType>
745 const Epetra_Map rows (TrilinosWrappers::types::int_type(sp.n_rows()), 0,
747 const Epetra_Map columns (TrilinosWrappers::types::int_type(sp.n_cols()), 0,
750 reinit_sp (rows, columns, sp,
false,
762 column_space_map.reset (
new Epetra_Map (TrilinosWrappers::types::int_type(0),
763 TrilinosWrappers::types::int_type(0),
767 graph->FillComplete();
785 TrilinosWrappers::types::int_type row =
nonlocal_graph->RowMap().MyGID(
786 static_cast<TrilinosWrappers::types::int_type> (0));
793 static_cast<const Epetra_Map &>(graph->RangeMap()));
795 Epetra_Export exporter(
nonlocal_graph->RowMap(), graph->RowMap());
800 static_cast<const Epetra_Map &>(graph->RangeMap()));
804 static_cast<const Epetra_Map &>(graph->RangeMap()),
809 ierr = graph->OptimizeStorage ();
821 int trilinos_i = graph->LRID(static_cast<TrilinosWrappers::types::int_type>(i)),
822 trilinos_j = graph->LCID(static_cast<TrilinosWrappers::types::int_type>(j));
830 if (trilinos_i == -1)
839 if (graph->Filled() ==
false)
841 int nnz_present = graph->NumGlobalIndices(i);
843 TrilinosWrappers::types::int_type *col_indices;
851 int ierr = graph->ExtractGlobalRowView(
852 static_cast<TrilinosWrappers::types::int_type>(trilinos_i),
853 nnz_extracted, col_indices);
855 Assert (ierr==0, ExcTrilinosError(ierr));
856 Assert (nnz_present == nnz_extracted,
857 ExcDimensionMismatch(nnz_present, nnz_extracted));
860 TrilinosWrappers::types::int_type *el_find =
861 std::find(col_indices, col_indices + nnz_present, trilinos_j);
863 TrilinosWrappers::types::int_type local_col_index =
864 (TrilinosWrappers::types::int_type)(el_find - col_indices);
866 if (local_col_index == nnz_present)
873 int nnz_present = graph->NumGlobalIndices(
874 static_cast<TrilinosWrappers::types::int_type>(i));
881 int ierr = graph->ExtractMyRowView(trilinos_i,
882 nnz_extracted, col_indices);
884 Assert (ierr==0, ExcTrilinosError(ierr));
886 Assert (nnz_present == nnz_extracted,
887 ExcDimensionMismatch(nnz_present, nnz_extracted));
890 int *el_find = std::find(col_indices, col_indices + nnz_present,
891 static_cast<int>(trilinos_j));
893 int local_col_index = (int)(el_find - col_indices);
895 if (local_col_index == nnz_present)
909 TrilinosWrappers::types::int_type global_b=0;
914 graph->ExtractMyRowView(i, num_entries, indices);
915 for (
unsigned int j=0; j<(
unsigned int)num_entries; ++j)
917 if (static_cast<size_type>(std::abs(static_cast<TrilinosWrappers::types::int_type>(i-indices[j]))) > local_b)
918 local_b = std::abs(static_cast<TrilinosWrappers::types::int_type>(i-indices[j]));
921 graph->Comm().MaxAll((TrilinosWrappers::types::int_type *)&local_b, &global_b, 1);
930 const TrilinosWrappers::types::int_type
n_rows = n_global_rows(*graph);
939 TrilinosWrappers::types::int_type
n_cols;
940 if (graph->Filled() ==
true)
941 n_cols = n_global_cols(*graph);
953 int n_rows = graph -> NumMyRows();
960 std::pair<SparsityPattern::size_type, SparsityPattern::size_type>
964 begin = min_my_gid(graph->RowMap());
965 end = max_my_gid(graph->RowMap())+1;
967 return std::make_pair (begin, end);
975 TrilinosWrappers::types::int_type nnz = n_global_entries(*graph);
985 int nnz = graph->MaxNumIndices();
987 return static_cast<unsigned int>(nnz);
999 TrilinosWrappers::types::int_type ncols = -1;
1000 TrilinosWrappers::types::int_type local_row =
1001 graph->LRID(static_cast<TrilinosWrappers::types::int_type>(row));
1007 ncols = graph->NumMyIndices (local_row);
1017 return static_cast<const Epetra_Map &
>(graph->DomainMap());
1025 return static_cast<const Epetra_Map &
>(graph->RangeMap());
1033 return static_cast<const Epetra_Map &
>(graph->RowMap());
1041 return static_cast<const Epetra_Map &
>(graph->ColMap());
1049 return graph->RangeMap().Comm();
1058 #ifdef DEAL_II_WITH_MPI 1060 const Epetra_MpiComm *mpi_comm
1061 =
dynamic_cast<const Epetra_MpiComm *
>(&graph->RangeMap().Comm());
1062 return mpi_comm->Comm();
1065 return MPI_COMM_SELF;
1076 Assert (
false, ExcNotImplemented());
1086 const bool write_extended_trilinos_info)
const 1088 if (write_extended_trilinos_info)
1095 for (
int i=0; i<graph->NumMyRows(); ++i)
1097 graph->ExtractMyRowView (i, num_entries, indices);
1098 for (
int j=0; j<num_entries; ++j)
1099 out <<
"(" << i <<
"," << indices[global_row_index(*graph,j)] <<
") " 1112 Assert (graph->Filled() ==
true, ExcInternalError());
1113 for (
unsigned int row=0; row<
local_size(); ++row)
1117 graph->ExtractMyRowView (row, num_entries, indices);
1119 for (
unsigned int j=0; j<(
unsigned int)num_entries; ++j)
1125 out << indices[global_row_index(*graph,static_cast<int>(j))]
1126 <<
" " << -static_cast<signed int>(row) << std::endl;
1136 Assert(
false, ExcNotImplemented());
1151 const ::SparsityPattern &,
1155 const ::DynamicSparsityPattern &,
1161 const ::SparsityPattern &,
1166 const ::DynamicSparsityPattern &,
1172 const ::SparsityPattern &,
1177 const ::DynamicSparsityPattern &,
1185 const ::SparsityPattern &,
1191 const ::DynamicSparsityPattern &,
1197 DEAL_II_NAMESPACE_CLOSE
1199 #endif // DEAL_II_WITH_TRILINOS unsigned int local_size() const
#define AssertDimension(dim1, dim2)
::types::global_dof_index size_type
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
size_type row_length(const size_type row) const
::ExceptionBase & ExcMessage(std::string arg1)
std_cxx11::shared_ptr< Epetra_FECrsGraph > graph
size_type n_nonzero_elements() const
std_cxx11::shared_ptr< Epetra_CrsGraph > nonlocal_graph
#define AssertThrow(cond, exc)
const Epetra_Comm & comm_self()
const_iterator begin() const
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
const Epetra_Comm & trilinos_communicator() const DEAL_II_DEPRECATED
bool exists(const size_type i, const size_type j) const
SparsityPattern * sparsity_pattern
std_cxx11::shared_ptr< Epetra_Map > column_space_map
const Epetra_Map & row_partitioner() const DEAL_II_DEPRECATED
std_cxx11::shared_ptr< const std::vector< size_type > > colnum_cache
MPI_Comm get_mpi_communicator() const
void subtract_set(const IndexSet &other)
#define Assert(cond, exc)
void copy_from(const SparsityPattern &input_sparsity_pattern)
size_type bandwidth() const
std::size_t memory_consumption() const
void print_gnuplot(std::ostream &out) const
const_iterator end() const
const Epetra_Map & range_partitioner() const DEAL_II_DEPRECATED
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
unsigned int max_entries_per_row() const
const Epetra_Map & domain_partitioner() const DEAL_II_DEPRECATED
const Epetra_Map & col_partitioner() const DEAL_II_DEPRECATED
SparsityPattern & operator=(const SparsityPattern &input_sparsity_pattern)
void reinit(const size_type m, const size_type n, const size_type n_entries_per_row=0)
size_type n_elements() const
virtual ~SparsityPattern()
std::pair< size_type, size_type > local_range() const