19#ifdef DEAL_II_WITH_TRILINOS
32# include <boost/container/small_vector.hpp>
35# ifdef DEAL_II_TRILINOS_WITH_EPETRAEXT
36# include <EpetraExt_MatrixMatrix.h>
38# include <Epetra_Export.h>
39# include <Teuchos_RCP.hpp>
40# include <ml_epetra_utils.h>
41# include <ml_struct.h>
51 template <
typename VectorType>
52 typename VectorType::value_type *
58 template <
typename VectorType>
59 const typename VectorType::value_type *
65 template <
typename VectorType>
66 typename VectorType::value_type *
72 template <
typename VectorType>
73 const typename VectorType::value_type *
74 end(
const VectorType &V)
107# ifdef DEAL_II_TRILINOS_WITH_TPETRA
108 template <
typename Number>
115 template <
typename Number>
122 template <
typename Number>
130 template <
typename Number>
167 std::make_shared<std::vector<TrilinosScalar>>(
matrix->
n());
209 : column_space_map(new Epetra_Map(0, 0,
Utilities::Trilinos::comm_self()))
211 new Epetra_FECrsMatrix(View, *column_space_map, *column_space_map, 0))
222 const unsigned int n_max_entries_per_row)
237 matrix(new Epetra_FECrsMatrix(
243 n_max_entries_per_row,
253 const std::vector<unsigned int> &n_entries_per_row)
258 , matrix(new Epetra_FECrsMatrix(
264 reinterpret_cast<
int *>(
265 const_cast<unsigned
int *>(n_entries_per_row.data())),
275 const unsigned int n_max_entries_per_row)
276 : column_space_map(new Epetra_Map(
277 parallel_partitioning.make_trilinos_map(communicator, false)))
278 , matrix(new Epetra_FECrsMatrix(Copy,
280 n_max_entries_per_row,
290 const std::vector<unsigned int> &n_entries_per_row)
291 : column_space_map(new Epetra_Map(
292 parallel_partitioning.make_trilinos_map(communicator, false)))
293 , matrix(new Epetra_FECrsMatrix(Copy,
295 reinterpret_cast<
int *>(
296 const_cast<unsigned
int *>(
297 n_entries_per_row.data())),
306 const IndexSet &col_parallel_partitioning,
309 : column_space_map(new Epetra_Map(
310 col_parallel_partitioning.make_trilinos_map(communicator, false)))
311 , matrix(new Epetra_FECrsMatrix(
313 row_parallel_partitioning.make_trilinos_map(communicator, false),
314 n_max_entries_per_row,
323 const IndexSet &col_parallel_partitioning,
325 const std::vector<unsigned int> &n_entries_per_row)
326 : column_space_map(new Epetra_Map(
327 col_parallel_partitioning.make_trilinos_map(communicator, false)))
328 , matrix(new Epetra_FECrsMatrix(
330 row_parallel_partitioning.make_trilinos_map(communicator, false),
331 reinterpret_cast<
int *>(
332 const_cast<unsigned
int *>(n_entries_per_row.data())),
341 : column_space_map(new Epetra_Map(sparsity_pattern.domain_partitioner()))
343 new Epetra_FECrsMatrix(Copy,
344 sparsity_pattern.trilinos_sparsity_pattern(),
351 "The Trilinos sparsity pattern has not been compressed."));
358 : column_space_map(std::move(other.column_space_map))
359 , matrix(std::move(other.matrix))
360 , nonlocal_matrix(std::move(other.nonlocal_matrix))
361 , nonlocal_matrix_exporter(std::move(other.nonlocal_matrix_exporter))
362 , last_action(other.last_action)
363 , compressed(other.compressed)
365 other.last_action = Zero;
366 other.compressed =
false;
383 bool needs_deep_copy =
388 if (!needs_deep_copy)
395 const int row_local =
matrix->RowMap().LID(
399 int n_entries, rhs_n_entries;
401 int * index_ptr, *rhs_index_ptr;
402 int ierr = rhs.
matrix->ExtractMyRowView(row_local,
409 ierr =
matrix->ExtractMyRowView(row_local,
415 if (n_entries != rhs_n_entries ||
416 std::memcmp(
static_cast<void *
>(index_ptr),
417 static_cast<void *
>(rhs_index_ptr),
418 sizeof(
int) * n_entries) != 0)
420 needs_deep_copy =
true;
424 for (
int i = 0; i < n_entries; ++i)
425 value_ptr[i] = rhs_value_ptr[i];
435 matrix = std::make_unique<Epetra_FECrsMatrix>(*rhs.
matrix);
442 std::make_unique<Epetra_CrsMatrix>(Copy, rhs.
nonlocal_matrix->Graph());
451 template <
typename SparsityPatternType>
453 reinit_matrix(
const IndexSet & row_parallel_partitioning,
454 const IndexSet & column_parallel_partitioning,
455 const SparsityPatternType & sparsity_pattern,
456 const bool exchange_data,
458 std::unique_ptr<Epetra_Map> &column_space_map,
459 std::unique_ptr<Epetra_FECrsMatrix> &matrix,
460 std::unique_ptr<Epetra_CrsMatrix> & nonlocal_matrix,
461 std::unique_ptr<Epetra_Export> & nonlocal_matrix_exporter)
465 nonlocal_matrix.reset();
466 nonlocal_matrix_exporter.reset();
468 column_space_map = std::make_unique<Epetra_Map>(
471 if (column_space_map->Comm().MyPID() == 0)
474 row_parallel_partitioning.
size());
476 column_parallel_partitioning.
size());
479 Epetra_Map row_space_map =
489 trilinos_sparsity.
reinit(row_parallel_partitioning,
490 column_parallel_partitioning,
494 matrix = std::make_unique<Epetra_FECrsMatrix>(
495 Copy, trilinos_sparsity.trilinos_sparsity_pattern(),
false);
503 std::vector<int> n_entries_per_row(last_row - first_row);
505 for (size_type row = first_row; row < last_row; ++row)
506 n_entries_per_row[row - first_row] = sparsity_pattern.row_length(row);
522 std::unique_ptr<Epetra_CrsGraph> graph;
523 if (row_space_map.Comm().NumProc() > 1)
524 graph = std::make_unique<Epetra_CrsGraph>(Copy,
526 n_entries_per_row.data(),
529 graph = std::make_unique<Epetra_CrsGraph>(Copy,
532 n_entries_per_row.data(),
540 std::vector<TrilinosWrappers::types::int_type> row_indices;
542 for (size_type row = first_row; row < last_row; ++row)
544 const int row_length = sparsity_pattern.row_length(row);
548 row_indices.resize(row_length, -1);
550 typename SparsityPatternType::iterator p =
551 sparsity_pattern.begin(row);
552 for (size_type col = 0; p != sparsity_pattern.end(row); ++p, ++col)
553 row_indices[col] = p->column();
555 graph->Epetra_CrsGraph::InsertGlobalIndices(row,
564 graph->FillComplete(*column_space_map, row_space_map);
565 graph->OptimizeStorage();
573 matrix = std::make_unique<Epetra_FECrsMatrix>(Copy, *graph,
false);
585 class Epetra_CrsGraphMod :
public Epetra_CrsGraph
588 Epetra_CrsGraphMod(
const Epetra_Map &row_map,
589 const int * n_entries_per_row)
590 : Epetra_CrsGraph(Copy, row_map, n_entries_per_row, true)
594 SetIndicesAreGlobal()
596 this->Epetra_CrsGraph::SetIndicesAreGlobal(
true);
606 reinit_matrix(
const IndexSet & row_parallel_partitioning,
607 const IndexSet & column_parallel_partitioning,
609 const bool exchange_data,
611 std::unique_ptr<Epetra_Map> & column_space_map,
612 std::unique_ptr<Epetra_FECrsMatrix> &matrix,
613 std::unique_ptr<Epetra_CrsMatrix> & nonlocal_matrix,
614 std::unique_ptr<Epetra_Export> & nonlocal_matrix_exporter)
617 nonlocal_matrix.reset();
618 nonlocal_matrix_exporter.reset();
620 column_space_map = std::make_unique<Epetra_Map>(
624 row_parallel_partitioning.
size());
626 column_parallel_partitioning.
size());
628 Epetra_Map row_space_map =
633 if (relevant_rows.size() == 0)
635 relevant_rows.set_size(
637 relevant_rows.add_range(
640 relevant_rows.compress();
641 Assert(relevant_rows.n_elements() >=
642 static_cast<unsigned int>(row_space_map.NumMyElements()),
644 "Locally relevant rows of sparsity pattern must contain "
645 "all locally owned rows"));
650 const bool have_ghost_rows = [&]() {
651 std::vector<::types::global_dof_index> indices;
652 relevant_rows.fill_index_vector(indices);
653 Epetra_Map relevant_map(
661 row_space_map.Comm());
662 return !relevant_map.SameAs(row_space_map);
665 const unsigned int n_rows = relevant_rows.n_elements();
666 std::vector<TrilinosWrappers::types::int_type> ghost_rows;
667 std::vector<int> n_entries_per_row(row_space_map.NumMyElements());
668 std::vector<int> n_entries_per_ghost_row;
669 for (
unsigned int i = 0, own = 0; i < n_rows; ++i)
672 relevant_rows.nth_index_in_set(i);
673 if (row_space_map.MyGID(global_row))
674 n_entries_per_row[own++] = sparsity_pattern.
row_length(global_row);
675 else if (sparsity_pattern.
row_length(global_row) > 0)
677 ghost_rows.push_back(global_row);
678 n_entries_per_ghost_row.push_back(
683 Epetra_Map off_processor_map(-1,
685 (ghost_rows.size() > 0) ?
686 (ghost_rows.data()) :
689 row_space_map.Comm());
691 std::unique_ptr<Epetra_CrsGraph> graph;
692 std::unique_ptr<Epetra_CrsGraphMod> nonlocal_graph;
693 if (row_space_map.Comm().NumProc() > 1)
696 std::make_unique<Epetra_CrsGraph>(Copy,
698 (n_entries_per_row.size() > 0) ?
699 (n_entries_per_row.data()) :
701 exchange_data ?
false :
true);
702 if (have_ghost_rows ==
true)
703 nonlocal_graph = std::make_unique<Epetra_CrsGraphMod>(
704 off_processor_map, n_entries_per_ghost_row.data());
708 std::make_unique<Epetra_CrsGraph>(Copy,
711 (n_entries_per_row.size() > 0) ?
712 (n_entries_per_row.data()) :
717 std::vector<TrilinosWrappers::types::int_type> row_indices;
719 for (
unsigned int i = 0; i < n_rows; ++i)
722 relevant_rows.nth_index_in_set(i);
723 const int row_length = sparsity_pattern.
row_length(global_row);
727 row_indices.resize(row_length, -1);
728 for (
int col = 0; col < row_length; ++col)
729 row_indices[col] = sparsity_pattern.
column_number(global_row, col);
731 if (row_space_map.MyGID(global_row))
732 graph->InsertGlobalIndices(global_row,
738 nonlocal_graph->InsertGlobalIndices(global_row,
745 if (nonlocal_graph.get() !=
nullptr)
751 nonlocal_graph->SetIndicesAreGlobal();
752 Assert(nonlocal_graph->IndicesAreGlobal() ==
true,
754 nonlocal_graph->FillComplete(*column_space_map, row_space_map);
755 nonlocal_graph->OptimizeStorage();
760 Epetra_Export exporter(nonlocal_graph->RowMap(), row_space_map);
761 int ierr = graph->Export(*nonlocal_graph, exporter, Add);
767 std::make_unique<Epetra_CrsMatrix>(Copy, *nonlocal_graph);
770 graph->FillComplete(*column_space_map, row_space_map);
771 graph->OptimizeStorage();
776 matrix = std::make_unique<Epetra_FECrsMatrix>(Copy, *graph,
false);
782 template <
typename SparsityPatternType>
799 template <
typename SparsityPatternType>
800 inline typename std::enable_if<
801 !std::is_same<SparsityPatternType,
804 const IndexSet & col_parallel_partitioning,
805 const SparsityPatternType &sparsity_pattern,
807 const bool exchange_data)
809 reinit_matrix(row_parallel_partitioning,
810 col_parallel_partitioning,
836 matrix = std::make_unique<Epetra_FECrsMatrix>(
841 std::make_unique<Epetra_CrsMatrix>(Copy,
855 if (
this == &sparse_matrix)
859 std::make_unique<Epetra_Map>(sparse_matrix.
trilinos_matrix().DomainMap());
862 matrix = std::make_unique<Epetra_FECrsMatrix>(
877 template <
typename number>
880 const IndexSet & row_parallel_partitioning,
881 const IndexSet & col_parallel_partitioning,
882 const ::SparseMatrix<number> &dealii_sparse_matrix,
884 const double drop_tolerance,
885 const bool copy_values,
886 const ::SparsityPattern * use_this_sparsity)
888 if (copy_values ==
false)
892 if (use_this_sparsity ==
nullptr)
893 reinit(row_parallel_partitioning,
894 col_parallel_partitioning,
895 dealii_sparse_matrix.get_sparsity_pattern(),
899 reinit(row_parallel_partitioning,
900 col_parallel_partitioning,
907 const size_type n_rows = dealii_sparse_matrix.m();
912 const ::SparsityPattern &sparsity_pattern =
913 (use_this_sparsity !=
nullptr) ?
915 dealii_sparse_matrix.get_sparsity_pattern();
917 if (
matrix.get() ==
nullptr ||
m() != n_rows ||
920 reinit(row_parallel_partitioning,
921 col_parallel_partitioning,
933 std::vector<size_type> row_indices(maximum_row_length);
934 std::vector<TrilinosScalar> values(maximum_row_length);
936 for (
size_type row = 0; row < n_rows; ++row)
938 if (row_parallel_partitioning.
is_element(row) ==
true)
941 sparsity_pattern.begin(row);
942 typename ::SparseMatrix<number>::const_iterator it =
943 dealii_sparse_matrix.begin(row);
945 if (sparsity_pattern.n_rows() == sparsity_pattern.n_cols())
949 if (std::fabs(it->value()) > drop_tolerance)
951 values[col] = it->value();
952 row_indices[col++] = it->column();
958 while (it != dealii_sparse_matrix.end(row) &&
959 select_index != sparsity_pattern.end(row))
961 while (select_index->column() < it->column() &&
962 select_index != sparsity_pattern.end(row))
964 while (it->column() < select_index->column() &&
965 it != dealii_sparse_matrix.end(row))
968 if (it == dealii_sparse_matrix.end(row))
970 if (std::fabs(it->value()) > drop_tolerance)
972 values[col] = it->value();
973 row_indices[col++] = it->column();
980 reinterpret_cast<size_type *
>(row_indices.data()),
989 template <
typename number>
992 const ::SparseMatrix<number> &dealii_sparse_matrix,
993 const double drop_tolerance,
994 const bool copy_values,
995 const ::SparsityPattern * use_this_sparsity)
999 dealii_sparse_matrix,
1010 const bool copy_values)
1012 Assert(input_matrix.Filled() ==
true,
1013 ExcMessage(
"Input CrsMatrix has not called FillComplete()!"));
1017 const Epetra_CrsGraph *graph = &input_matrix.Graph();
1022 matrix = std::make_unique<Epetra_FECrsMatrix>(Copy, *graph,
false);
1026 if (copy_values ==
true)
1032 const size_type my_nonzeros = input_matrix.NumMyNonzeros();
1033 std::memcpy(values, in_values, my_nonzeros *
sizeof(
TrilinosScalar));
1057 "compress() can only be called with VectorOperation add, insert, or unknown"));
1065 ExcMessage(
"Operation and argument to compress() do not match"));
1091 ierr =
matrix->OptimizeStorage();
1136 matrix->ExtractMyRowView(local_row, num_entries, values, col_indices);
1141 const std::ptrdiff_t diag_index =
1142 std::find(col_indices, col_indices + num_entries, local_row) -
1146 if (diag_index != j || new_diag_value == 0)
1149 if (diag_index != num_entries)
1150 values[diag_index] = new_diag_value;
1160 for (
const auto row : rows)
1183 if (trilinos_i == -1)
1198 int nnz_present =
matrix->NumMyEntries(trilinos_i);
1207 int ierr =
matrix->ExtractMyRowView(trilinos_i,
1214 Assert(nnz_present == nnz_extracted,
1220 const std::ptrdiff_t local_col_index =
1221 std::find(col_indices, col_indices + nnz_present, trilinos_j) -
1232 if (local_col_index == nnz_present)
1237 value = values[local_col_index];
1263 if ((trilinos_i == -1) || (trilinos_j == -1))
1274 int nnz_present =
matrix->NumMyEntries(trilinos_i);
1282 int ierr =
matrix->ExtractMyRowView(trilinos_i,
1289 Assert(nnz_present == nnz_extracted,
1295 const std::ptrdiff_t local_col_index =
1296 std::find(col_indices, col_indices + nnz_present, trilinos_j) -
1306 if (local_col_index == nnz_present)
1309 value = values[local_col_index];
1355 int ierr =
matrix->NumMyRowEntries(local_row, ncols);
1359 return static_cast<unsigned int>(ncols);
1366 const std::vector<size_type> & col_indices,
1368 const bool elide_zero_values)
1370 Assert(row_indices.size() == values.m(),
1372 Assert(col_indices.size() == values.n(),
1375 for (
size_type i = 0; i < row_indices.size(); ++i)
1387 const std::vector<size_type> & col_indices,
1388 const std::vector<TrilinosScalar> &values,
1389 const bool elide_zero_values)
1391 Assert(col_indices.size() == values.size(),
1405 SparseMatrix::set<TrilinosScalar>(
const size_type row,
1409 const bool elide_zero_values)
1429 boost::container::small_vector<TrilinosScalar, 200> local_value_array(
1431 boost::container::small_vector<TrilinosWrappers::types::int_type, 200>
1432 local_index_array(n_cols);
1437 if (elide_zero_values ==
false)
1442 col_value_ptr = values;
1449 col_index_ptr = local_index_array.data();
1450 col_value_ptr = local_value_array.data();
1455 const double value = values[j];
1459 local_index_array[n_columns] = col_indices[j];
1460 local_value_array[n_columns] = value;
1477 if (
matrix->RowMap().MyGID(
1480 if (
matrix->Filled() ==
false)
1482 ierr =
matrix->Epetra_CrsMatrix::InsertGlobalValues(
1483 row,
static_cast<int>(n_columns), col_value_ptr, col_index_ptr);
1492 ierr =
matrix->Epetra_CrsMatrix::ReplaceGlobalValues(row,
1507 if (
matrix->Filled() ==
false)
1509 ierr =
matrix->InsertGlobalValues(1,
1514 Epetra_FECrsMatrix::ROW_MAJOR);
1519 ierr =
matrix->ReplaceGlobalValues(1,
1524 Epetra_FECrsMatrix::ROW_MAJOR);
1541 const bool elide_zero_values)
1543 Assert(indices.size() == values.m(),
1547 for (
size_type i = 0; i < indices.size(); ++i)
1559 const std::vector<size_type> & col_indices,
1561 const bool elide_zero_values)
1563 Assert(row_indices.size() == values.m(),
1565 Assert(col_indices.size() == values.n(),
1568 for (
size_type i = 0; i < row_indices.size(); ++i)
1580 const std::vector<size_type> & col_indices,
1581 const std::vector<TrilinosScalar> &values,
1582 const bool elide_zero_values)
1584 Assert(col_indices.size() == values.size(),
1601 const bool elide_zero_values,
1623 boost::container::small_vector<TrilinosScalar, 100> local_value_array(
1625 boost::container::small_vector<TrilinosWrappers::types::int_type, 100>
1626 local_index_array(n_cols);
1631 if (elide_zero_values ==
false)
1636 col_value_ptr = values;
1647 col_index_ptr = local_index_array.data();
1648 col_value_ptr = local_value_array.data();
1653 const double value = values[j];
1658 local_index_array[n_columns] = col_indices[j];
1659 local_value_array[n_columns] = value;
1675 if (
matrix->RowMap().MyGID(
1678 ierr =
matrix->Epetra_CrsMatrix::SumIntoGlobalValues(row,
1692 ExcMessage(
"Attempted to write into off-processor matrix row "
1693 "that has not be specified as being writable upon "
1711 ierr =
matrix->SumIntoGlobalValues(1,
1716 Epetra_FECrsMatrix::ROW_MAJOR);
1723 std::cout <<
"------------------------------------------" << std::endl;
1724 std::cout <<
"Got error " << ierr <<
" in row " << row <<
" of proc "
1725 <<
matrix->RowMap().Comm().MyPID()
1726 <<
" when trying to add the columns:" << std::endl;
1728 std::cout << col_index_ptr[i] <<
" ";
1729 std::cout << std::endl << std::endl;
1730 std::cout <<
"Matrix row "
1731 << (
matrix->RowMap().MyGID(
1736 <<
" has the following indices:" << std::endl;
1737 std::vector<TrilinosWrappers::types::int_type> indices;
1738 const Epetra_CrsGraph * graph =
1745 indices.resize(graph->NumGlobalIndices(row));
1747 graph->ExtractGlobalRowCopy(row,
1754 std::cout << indices[i] <<
" ";
1755 std::cout << std::endl << std::endl;
1772 const int ierr =
matrix->PutScalar(d);
1790 ExcMessage(
"Can only add matrices with same distribution of rows"));
1792 ExcMessage(
"Addition of matrices only allowed if matrices are "
1793 "filled, i.e., compress() has been called"));
1795 const bool same_col_map =
matrix->ColMap().SameAs(rhs.
matrix->ColMap());
1799 const int row_local =
matrix->RowMap().LID(
1807 int n_entries, rhs_n_entries;
1809 int * index_ptr, *rhs_index_ptr;
1810 int ierr = rhs.
matrix->ExtractMyRowView(row_local,
1818 matrix->ExtractMyRowView(row_local, n_entries, value_ptr, index_ptr);
1820 bool expensive_checks = (n_entries != rhs_n_entries || !same_col_map);
1821 if (!expensive_checks)
1825 expensive_checks = std::memcmp(
static_cast<void *
>(index_ptr),
1826 static_cast<void *
>(rhs_index_ptr),
1827 sizeof(
int) * n_entries) != 0;
1828 if (!expensive_checks)
1829 for (
int i = 0; i < n_entries; ++i)
1830 value_ptr[i] += rhs_value_ptr[i] * factor;
1836 if (expensive_checks)
1838 for (
int i = 0; i < rhs_n_entries; ++i)
1840 if (rhs_value_ptr[i] == 0.)
1844 int local_col =
matrix->ColMap().LID(rhs_global_col);
1846 index_ptr + n_entries,
1848 Assert(local_index != index_ptr + n_entries &&
1849 *local_index == local_col,
1851 "Adding the entries from the other matrix "
1852 "failed, because the sparsity pattern "
1853 "of that matrix includes more elements than the "
1854 "calling matrix, which is not allowed."));
1855 value_ptr[local_index - index_ptr] += factor * rhs_value_ptr[i];
1873 if (!
matrix->UseTranspose())
1875 ierr =
matrix->SetUseTranspose(
true);
1880 ierr =
matrix->SetUseTranspose(
false);
1890 const int ierr =
matrix->Scale(a);
1906 const int ierr =
matrix->Scale(factor);
1919 return matrix->NormOne();
1928 return matrix->NormInf();
1937 return matrix->NormFrobenius();
1944 namespace SparseMatrixImplementation
1946 template <
typename VectorType>
1959 ExcMessage(
"The column partitioning of a matrix does not match "
1960 "the partitioning of a vector you are trying to "
1961 "multiply it with. Are you multiplying the "
1962 "matrix with a vector that has ghost elements?"));
1964 ExcMessage(
"The row partitioning of a matrix does not match "
1965 "the partitioning of a vector you are trying to "
1966 "put the result of a matrix-vector product in. "
1967 "Are you trying to put the product of the "
1968 "matrix with a vector into a vector that has "
1969 "ghost elements?"));
1978 template <
typename VectorType>
1979 typename std::enable_if<
1980 std::is_same<typename VectorType::value_type, TrilinosScalar>::value>::type
1996 Epetra_MultiVector tril_dst(
1998 Epetra_MultiVector tril_src(View,
2005 const int ierr =
matrix->Multiply(
false, tril_src, tril_dst);
2012 template <
typename VectorType>
2013 typename std::enable_if<
2014 !std::is_same<typename VectorType::value_type, TrilinosScalar>::value>::type
2022 template <
typename VectorType>
2023 typename std::enable_if<
2024 std::is_same<typename VectorType::value_type, TrilinosScalar>::value>::type
2038 Epetra_MultiVector tril_dst(
2040 Epetra_MultiVector tril_src(View,
2046 const int ierr =
matrix->Multiply(
true, tril_src, tril_dst);
2053 template <
typename VectorType>
2054 typename std::enable_if<
2055 !std::is_same<typename VectorType::value_type, TrilinosScalar>::value>::type
2063 template <
typename VectorType>
2071 VectorType tmp_vector;
2072 tmp_vector.reinit(dst,
true);
2073 vmult(tmp_vector, src);
2079 template <
typename VectorType>
2087 VectorType tmp_vector;
2088 tmp_vector.reinit(dst,
true);
2101 temp_vector.
reinit(v,
true);
2103 vmult(temp_vector, v);
2104 return temp_vector * v;
2116 temp_vector.
reinit(v,
true);
2118 vmult(temp_vector, v);
2119 return u * temp_vector;
2147 const bool transpose_left)
2149 const bool use_vector = (V.
size() == inputright.
m() ? true :
false);
2150 if (transpose_left ==
false)
2152 Assert(inputleft.
n() == inputright.
m(),
2156 ExcMessage(
"Parallel partitioning of A and B does not fit."));
2160 Assert(inputleft.
m() == inputright.
m(),
2164 ExcMessage(
"Parallel partitioning of A and B does not fit."));
2175 Teuchos::RCP<Epetra_CrsMatrix> mod_B;
2176 if (use_vector ==
false)
2178 mod_B = Teuchos::rcp(
const_cast<Epetra_CrsMatrix *
>(
2184 mod_B = Teuchos::rcp(
2190 ExcMessage(
"Parallel distribution of matrix B and vector V "
2191 "does not match."));
2194 for (
int i = 0; i < local_N; ++i)
2197 double *new_data, *B_data;
2198 mod_B->ExtractMyRowView(i, N_entries, new_data);
2204 new_data[j] = value * B_data[j];
2215# ifdef DEAL_II_TRILINOS_WITH_EPETRAEXT
2220 const_cast<Epetra_CrsMatrix &
>(
2224 ExcMessage(
"This function requires that the Trilinos "
2225 "installation found while running the deal.II "
2226 "CMake scripts contains the optional Trilinos "
2227 "package 'EpetraExt'. However, this optional "
2228 "part of Trilinos was not found."));
2268 const bool print_detailed_trilinos_information)
const
2270 if (print_detailed_trilinos_information ==
true)
2278 for (
int i = 0; i <
matrix->NumMyRows(); ++i)
2281 matrix->ExtractMyRowView(i, num_entries, values, indices);
2289 <<
") " << values[j] << std::endl;
2302 sizeof(*this) +
sizeof(*matrix) +
sizeof(*
matrix->Graph().DataPtr());
2305 matrix->NumMyNonzeros() +
2314 const Epetra_MpiComm *mpi_comm =
2315 dynamic_cast<const Epetra_MpiComm *
>(&
matrix->RangeMap().Comm());
2317 return mpi_comm->Comm();
2326 namespace LinearOperatorImplementation
2329 : use_transpose(false)
2330 , communicator(MPI_COMM_SELF)
2331 , domain_map(
IndexSet().make_trilinos_map(communicator.Comm()))
2332 , range_map(
IndexSet().make_trilinos_map(communicator.Comm()))
2336 ExcMessage(
"Uninitialized TrilinosPayload::vmult called "
2337 "(Default constructor)"));
2342 ExcMessage(
"Uninitialized TrilinosPayload::Tvmult called "
2343 "(Default constructor)"));
2348 ExcMessage(
"Uninitialized TrilinosPayload::inv_vmult called "
2349 "(Default constructor)"));
2354 ExcMessage(
"Uninitialized TrilinosPayload::inv_Tvmult called "
2355 "(Default constructor)"));
2365 matrix.trilinos_matrix()),
2367 matrix_exemplar.trilinos_matrix().UseTranspose(),
2368 matrix_exemplar.get_mpi_communicator(),
2369 matrix_exemplar.locally_owned_domain_indices(),
2370 matrix_exemplar.locally_owned_range_indices())
2380 matrix.trilinos_matrix()),
2382 payload_exemplar.UseTranspose(),
2383 payload_exemplar.get_mpi_communicator(),
2384 payload_exemplar.locally_owned_domain_indices(),
2385 payload_exemplar.locally_owned_range_indices())
2395 matrix_exemplar.trilinos_matrix().UseTranspose(),
2396 matrix_exemplar.get_mpi_communicator(),
2397 matrix_exemplar.locally_owned_domain_indices(),
2398 matrix_exemplar.locally_owned_range_indices())
2407 preconditioner.trilinos_operator(),
2409 preconditioner_exemplar.trilinos_operator().UseTranspose(),
2410 preconditioner_exemplar.get_mpi_communicator(),
2411 preconditioner_exemplar.locally_owned_domain_indices(),
2412 preconditioner_exemplar.locally_owned_range_indices())
2422 payload_exemplar.UseTranspose(),
2423 payload_exemplar.get_mpi_communicator(),
2424 payload_exemplar.locally_owned_domain_indices(),
2425 payload_exemplar.locally_owned_range_indices())
2431 : vmult(payload.vmult)
2432 , Tvmult(payload.Tvmult)
2433 , inv_vmult(payload.inv_vmult)
2434 , inv_Tvmult(payload.inv_Tvmult)
2435 , use_transpose(payload.use_transpose)
2436 , communicator(payload.communicator)
2437 , domain_map(payload.domain_map)
2438 , range_map(payload.range_map)
2447 : use_transpose(false)
2450 communicator(first_op.communicator)
2451 , domain_map(second_op.domain_map)
2452 , range_map(first_op.range_map)
2463 tril_dst = tril_src;
2467 tril_dst = tril_src;
2471 tril_dst = tril_src;
2475 tril_dst = tril_src;
2489 const int ierr = tril_dst.PutScalar(0.0);
2495 const int ierr = tril_dst.PutScalar(0.0);
2502 ExcMessage(
"Cannot compute inverse of null operator"));
2504 const int ierr = tril_dst.PutScalar(0.0);
2511 ExcMessage(
"Cannot compute inverse of null operator"));
2513 const int ierr = tril_dst.PutScalar(0.0);
2613 return "TrilinosPayload";
2671 "Operators are set to work on incompatible IndexSets."));
2675 "Operators are set to work on incompatible IndexSets."));
2680 return_op.
vmult = [first_op, second_op](Range & tril_dst,
2681 const Domain &tril_src) {
2689 i->reinit(
IndexSet(first_op_init_map),
2694 const size_type i_local_size = i->end() - i->begin();
2698 (void)second_op_init_map;
2699 Intermediate tril_int(View,
2707 second_op.
Apply(tril_src, tril_int);
2708 first_op.
Apply(tril_src, tril_dst);
2709 const int ierr = tril_dst.Update(1.0, tril_int, 1.0);
2713 return_op.
Tvmult = [first_op, second_op](Domain & tril_dst,
2714 const Range &tril_src) {
2729 i->reinit(
IndexSet(first_op_init_map),
2734 const size_type i_local_size = i->end() - i->begin();
2738 (void)second_op_init_map;
2739 Intermediate tril_int(View,
2747 second_op.
Apply(tril_src, tril_int);
2748 first_op.
Apply(tril_src, tril_dst);
2749 const int ierr = tril_dst.Update(1.0, tril_int, 1.0);
2757 return_op.
inv_vmult = [first_op, second_op](Domain & tril_dst,
2758 const Range &tril_src) {
2766 i->reinit(
IndexSet(first_op_init_map),
2771 const size_type i_local_size = i->end() - i->begin();
2775 (void)second_op_init_map;
2776 Intermediate tril_int(View,
2786 const int ierr = tril_dst.Update(1.0, tril_int, 1.0);
2790 return_op.
inv_Tvmult = [first_op, second_op](Range & tril_dst,
2791 const Domain &tril_src) {
2806 i->reinit(
IndexSet(first_op_init_map),
2811 const size_type i_local_size = i->end() - i->begin();
2815 (void)second_op_init_map;
2816 Intermediate tril_int(View,
2826 const int ierr = tril_dst.Update(1.0, tril_int, 1.0);
2851 "Operators are set to work on incompatible IndexSets."));
2856 return_op.
vmult = [first_op, second_op](Range & tril_dst,
2857 const Domain &tril_src) {
2865 i->reinit(
IndexSet(first_op_init_map),
2870 const size_type i_local_size = i->end() - i->begin();
2874 (void)second_op_init_map;
2875 Intermediate tril_int(View,
2883 second_op.
Apply(tril_src, tril_int);
2884 first_op.
Apply(tril_int, tril_dst);
2887 return_op.
Tvmult = [first_op, second_op](Domain & tril_dst,
2888 const Range &tril_src) {
2903 i->reinit(
IndexSet(first_op_init_map),
2908 const size_type i_local_size = i->end() - i->begin();
2912 (void)second_op_init_map;
2913 Intermediate tril_int(View,
2920 first_op.
Apply(tril_src, tril_int);
2921 second_op.
Apply(tril_int, tril_dst);
2928 return_op.
inv_vmult = [first_op, second_op](Domain & tril_dst,
2929 const Range &tril_src) {
2937 i->reinit(
IndexSet(first_op_init_map),
2942 const size_type i_local_size = i->end() - i->begin();
2946 (void)second_op_init_map;
2947 Intermediate tril_int(View,
2959 return_op.
inv_Tvmult = [first_op, second_op](Range & tril_dst,
2960 const Domain &tril_src) {
2975 i->reinit(
IndexSet(first_op_init_map),
2980 const size_type i_local_size = i->end() - i->begin();
2984 (void)second_op_init_map;
2985 Intermediate tril_int(View,
3013# include "trilinos_sparse_matrix.inst"
3028 const ::SparsityPattern &,
3044 const ::Vector<double> &)
const;
3049 const ::LinearAlgebra::distributed::Vector<double> &)
const;
3051# ifdef DEAL_II_TRILINOS_WITH_TPETRA
3055 const ::LinearAlgebra::TpetraWrappers::Vector<double> &)
const;
3060 const ::LinearAlgebra::TpetraWrappers::Vector<float> &)
const;
3066 const ::LinearAlgebra::EpetraWrappers::Vector &)
const;
3073 const ::Vector<double> &)
const;
3078 const ::LinearAlgebra::distributed::Vector<double> &)
const;
3080# ifdef DEAL_II_TRILINOS_WITH_TPETRA
3084 const ::LinearAlgebra::TpetraWrappers::Vector<double> &)
const;
3089 const ::LinearAlgebra::TpetraWrappers::Vector<float> &)
const;
3095 const ::LinearAlgebra::EpetraWrappers::Vector &)
const;
3102 const ::Vector<double> &)
const;
3107 const ::LinearAlgebra::distributed::Vector<double> &)
const;
3109# ifdef DEAL_II_TRILINOS_WITH_TPETRA
3113 const ::LinearAlgebra::TpetraWrappers::Vector<double> &)
const;
3118 const ::LinearAlgebra::TpetraWrappers::Vector<float> &)
const;
3124 const ::LinearAlgebra::EpetraWrappers::Vector &)
const;
3131 const ::Vector<double> &)
const;
3136 const ::LinearAlgebra::distributed::Vector<double> &)
const;
3138# ifdef DEAL_II_TRILINOS_WITH_TPETRA
3142 const ::LinearAlgebra::TpetraWrappers::Vector<double> &)
const;
3147 const ::LinearAlgebra::TpetraWrappers::Vector<float> &)
const;
3153 const ::LinearAlgebra::EpetraWrappers::Vector &)
const;
bool is_element(const size_type index) const
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
const Epetra_FEVector & trilinos_vector() const
const Tpetra::Vector< Number, int, types::global_dof_index > & trilinos_vector() const
virtual void reinit(const size_type m, const size_type n, const ArrayView< const unsigned int > &row_lengths) override
std::shared_ptr< std::vector< TrilinosScalar > > value_cache
std::shared_ptr< std::vector< size_type > > colnum_cache
void set(const size_type i, const size_type j, const TrilinosScalar value)
std::unique_ptr< Epetra_Map > column_space_map
std::enable_if< std::is_same< typenameVectorType::value_type, TrilinosScalar >::value >::type vmult(VectorType &dst, const VectorType &src) const
MPI_Comm get_mpi_communicator() const
SparseMatrix & operator*=(const TrilinosScalar factor)
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
std::unique_ptr< Epetra_Export > nonlocal_matrix_exporter
TrilinosScalar l1_norm() const
std::unique_ptr< Epetra_FECrsMatrix > matrix
::types::global_dof_index size_type
TrilinosScalar linfty_norm() const
const Epetra_CrsMatrix & trilinos_matrix() const
size_type memory_consumption() const
TrilinosScalar matrix_norm_square(const MPI::Vector &v) const
Epetra_CombineMode last_action
IndexSet locally_owned_range_indices() const
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
void clear_row(const size_type row, const TrilinosScalar new_diag_value=0)
void reinit(const SparsityPatternType &sparsity_pattern)
void vmult_add(VectorType &dst, const VectorType &src) const
SparseMatrix & operator=(const SparseMatrix &)=delete
void compress(::VectorOperation::values operation)
SparseMatrix & operator/=(const TrilinosScalar factor)
void Tvmult_add(VectorType &dst, const VectorType &src) const
void clear_rows(const std::vector< size_type > &rows, const TrilinosScalar new_diag_value=0)
TrilinosScalar el(const size_type i, const size_type j) const
IndexSet locally_owned_domain_indices() const
TrilinosScalar residual(MPI::Vector &dst, const MPI::Vector &x, const MPI::Vector &b) const
bool in_local_range(const size_type index) const
void copy_from(const SparseMatrix &source)
unsigned int row_length(const size_type row) const
TrilinosScalar frobenius_norm() const
std::uint64_t n_nonzero_elements() const
const Epetra_CrsGraph & trilinos_sparsity_pattern() const
TrilinosScalar diag_element(const size_type i) const
void add(const size_type i, const size_type j, const TrilinosScalar value)
unsigned int local_size() const
std::enable_if< std::is_same< typenameVectorType::value_type, TrilinosScalar >::value >::type Tvmult(VectorType &dst, const VectorType &src) const
TrilinosScalar operator()(const size_type i, const size_type j) const
TrilinosScalar matrix_scalar_product(const MPI::Vector &u, const MPI::Vector &v) const
std::pair< size_type, size_type > local_range() const
std::unique_ptr< Epetra_CrsMatrix > nonlocal_matrix
const Epetra_Map & domain_partitioner() const
const Epetra_FECrsGraph & trilinos_sparsity_pattern() const
Epetra_MpiComm communicator
std::function< void(VectorType &, const VectorType &)> inv_Tvmult
virtual int SetUseTranspose(bool UseTranspose) override
virtual bool UseTranspose() const override
virtual const Epetra_Map & OperatorDomainMap() const override
IndexSet locally_owned_range_indices() const
TrilinosPayload transpose_payload() const
Epetra_MultiVector VectorType
MPI_Comm get_mpi_communicator() const
virtual int ApplyInverse(const VectorType &Y, VectorType &X) const override
std::function< void(VectorType &, const VectorType &)> Tvmult
virtual int Apply(const VectorType &X, VectorType &Y) const override
std::function< void(VectorType &, const VectorType &)> vmult
virtual const Epetra_Map & OperatorRangeMap() const override
std::function< void(VectorType &, const VectorType &)> inv_vmult
virtual const char * Label() const override
IndexSet locally_owned_domain_indices() const
virtual const Epetra_Comm & Comm() const override
virtual bool HasNormInf() const override
TrilinosPayload identity_payload() const
virtual double NormInf() const override
TrilinosPayload null_payload() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
static ::ExceptionBase & ExcInvalidIndex(size_type arg1, size_type arg2)
#define Assert(cond, exc)
std::unique_ptr< Epetra_CrsGraph > nonlocal_graph
static ::ExceptionBase & ExcAccessToNonPresentElement(size_type arg1, size_type arg2)
#define AssertIsFinite(number)
static ::ExceptionBase & ExcAccessToNonlocalRow(std::size_t arg1)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDivideByZero()
static ::ExceptionBase & ExcSourceEqualsDestination()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcNotQuadratic()
static ::ExceptionBase & ExcMatrixNotCompressed()
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcTrilinosError(int arg1)
#define AssertThrow(cond, exc)
const IndexSet & row_index_set() const
size_type row_length(const size_type row) const
size_type column_number(const size_type row, const size_type index) const
const Epetra_BlockMap & trilinos_partitioner() const
real_type l2_norm() const
const Epetra_MultiVector & trilinos_vector() const
void reinit(const Vector &v, const bool omit_zeroing_entries=false, const bool allow_different_maps=false)
std::pair< size_type, size_type > local_range() const
IndexSet complete_index_set(const IndexSet::size_type N)
@ matrix
Contents is actually a matrix.
types::global_dof_index size_type
TrilinosPayload operator+(const TrilinosPayload &first_op, const TrilinosPayload &second_op)
TrilinosPayload operator*(const TrilinosPayload &first_op, const TrilinosPayload &second_op)
void check_vector_map_equality(const Epetra_CrsMatrix &, const VectorType &, const VectorType &)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
void perform_mmult(const SparseMatrix &inputleft, const SparseMatrix &inputright, SparseMatrix &result, const MPI::Vector &V, const bool transpose_left)
TrilinosWrappers::types::int_type min_my_gid(const Epetra_BlockMap &map)
TrilinosWrappers::types::int64_type n_global_elements(const Epetra_BlockMap &map)
TrilinosWrappers::types::int_type global_column_index(const Epetra_CrsMatrix &matrix, const ::types::global_dof_index i)
TrilinosWrappers::types::int_type max_my_gid(const Epetra_BlockMap &map)
TrilinosWrappers::types::int_type global_row_index(const Epetra_CrsMatrix &matrix, const ::types::global_dof_index i)
TrilinosWrappers::types::int_type n_global_cols(const Epetra_CrsGraph &graph)
const Epetra_Comm & comm_self()
Iterator lower_bound(Iterator first, Iterator last, const T &val)
unsigned int global_dof_index