16 #include <deal.II/lac/trilinos_sparse_matrix.h> 18 #ifdef DEAL_II_WITH_TRILINOS 20 # include <deal.II/base/utilities.h> 21 # include <deal.II/lac/sparse_matrix.h> 22 # include <deal.II/lac/trilinos_sparsity_pattern.h> 23 # include <deal.II/lac/sparsity_pattern.h> 24 # include <deal.II/lac/dynamic_sparsity_pattern.h> 25 # include <deal.II/lac/sparsity_tools.h> 26 # include <deal.II/lac/parallel_vector.h> 28 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
29 # include <Epetra_Export.h> 30 # include <ml_epetra_utils.h> 31 # include <ml_struct.h> 32 # include <Teuchos_RCP.hpp> 35 DEAL_II_NAMESPACE_OPEN
41 #ifndef DEAL_II_WITH_64BIT_INDICES 46 int n_global_elements (
const Epetra_BlockMap &map)
48 return map.NumGlobalElements();
51 int min_my_gid(
const Epetra_BlockMap &map)
53 return map.MinMyGID();
56 int max_my_gid(
const Epetra_BlockMap &map)
58 return map.MaxMyGID();
61 int n_global_cols(
const Epetra_CrsGraph &graph)
63 return graph.NumGlobalCols();
66 int global_column_index(
const Epetra_CrsMatrix &matrix,
int i)
68 return matrix.GCID(i);
71 int global_row_index(
const Epetra_CrsMatrix &matrix,
int i)
73 return matrix.GRID(i);
80 long long int n_global_elements (
const Epetra_BlockMap &map)
82 return map.NumGlobalElements64();
85 long long int min_my_gid(
const Epetra_BlockMap &map)
87 return map.MinMyGID64();
90 long long int max_my_gid(
const Epetra_BlockMap &map)
92 return map.MaxMyGID64();
95 long long int n_global_cols(
const Epetra_CrsGraph &graph)
97 return graph.NumGlobalCols64();
100 long long int global_column_index(
const Epetra_CrsMatrix &matrix,
int i)
102 return matrix.GCID64(i);
105 long long int global_row_index(
const Epetra_CrsMatrix &matrix,
int i)
107 return matrix.GRID64(i);
124 if ((this->
a_row == matrix->m())
126 (matrix->in_local_range (this->a_row) ==
false))
136 TrilinosWrappers::types::int_type colnums = matrix->n();
139 value_cache.reset (
new std::vector<TrilinosScalar> (matrix->n()));
140 colnum_cache.reset (
new std::vector<size_type> (matrix->n()));
148 int ierr = matrix->trilinos_matrix().
149 ExtractGlobalRowCopy((TrilinosWrappers::types::int_type)this->
a_row,
152 reinterpret_cast<TrilinosWrappers::types::int_type *>(&((*
colnum_cache)[0])));
180 column_space_map (new Epetra_Map (0, 0,
182 matrix (new Epetra_FECrsMatrix(View, *column_space_map,
183 *column_space_map, 0)),
187 matrix->FillComplete();
205 const std::vector<unsigned int> &n_entries_per_row)
208 matrix (new Epetra_FECrsMatrix
210 (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
219 const Epetra_Map &input_col_map,
223 matrix (new Epetra_FECrsMatrix(Copy, input_row_map,
232 const Epetra_Map &input_col_map,
233 const std::vector<unsigned int> &n_entries_per_row)
236 matrix (new Epetra_FECrsMatrix(Copy, input_row_map,
237 (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
247 const unsigned int n_max_entries_per_row)
260 matrix (new Epetra_FECrsMatrix(Copy,
264 n_max_entries_per_row,
274 const std::vector<unsigned int> &n_entries_per_row)
278 matrix (new Epetra_FECrsMatrix(Copy,
282 (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
291 const MPI_Comm &communicator,
292 const unsigned int n_max_entries_per_row)
295 make_trilinos_map(communicator, false))),
296 matrix (new Epetra_FECrsMatrix(Copy,
298 n_max_entries_per_row,
307 const MPI_Comm &communicator,
308 const std::vector<unsigned int> &n_entries_per_row)
311 make_trilinos_map(communicator, false))),
312 matrix (new Epetra_FECrsMatrix(Copy,
314 (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
323 const IndexSet &col_parallel_partitioning,
324 const MPI_Comm &communicator,
328 make_trilinos_map(communicator, false))),
329 matrix (new Epetra_FECrsMatrix(Copy,
330 row_parallel_partitioning.
331 make_trilinos_map(communicator, false),
332 n_max_entries_per_row,
341 const IndexSet &col_parallel_partitioning,
342 const MPI_Comm &communicator,
343 const std::vector<unsigned int> &n_entries_per_row)
346 make_trilinos_map(communicator, false))),
347 matrix (new Epetra_FECrsMatrix(Copy,
348 row_parallel_partitioning.
349 make_trilinos_map(communicator, false),
350 (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
361 matrix (new Epetra_FECrsMatrix(Copy,
368 ExcMessage(
"The Trilinos sparsity pattern has not been compressed."));
391 bool needs_deep_copy =
392 !matrix->RowMap().SameAs(rhs.
matrix->RowMap()) ||
393 !matrix->ColMap().SameAs(rhs.
matrix->ColMap()) ||
394 !matrix->DomainMap().SameAs(rhs.
matrix->DomainMap()) ||
396 if (!needs_deep_copy)
398 const std::pair<size_type, size_type>
405 for (
size_type row=local_range.first; row < local_range.second; ++row)
407 const int row_local =
408 matrix->RowMap().LID(static_cast<TrilinosWrappers::types::int_type>(row));
410 int n_entries, rhs_n_entries;
411 TrilinosScalar *value_ptr, *rhs_value_ptr;
412 int *index_ptr, *rhs_index_ptr;
413 ierr = rhs.
matrix->ExtractMyRowView (row_local, rhs_n_entries,
414 rhs_value_ptr, rhs_index_ptr);
416 Assert (ierr == 0, ExcTrilinosError(ierr));
418 ierr = matrix->ExtractMyRowView (row_local, n_entries, value_ptr,
420 Assert (ierr == 0, ExcTrilinosError(ierr));
422 if (n_entries != rhs_n_entries ||
423 std::memcmp(static_cast<void *>(index_ptr),
424 static_cast<void *>(rhs_index_ptr),
425 sizeof(
int)*n_entries) != 0)
427 needs_deep_copy =
true;
431 for (
int i=0; i<n_entries; ++i)
432 value_ptr[i] = rhs_value_ptr[i];
442 matrix.reset (
new Epetra_FECrsMatrix(*rhs.
matrix));
457 template <
typename SparsityPatternType>
459 reinit_matrix (
const Epetra_Map &input_row_map,
460 const Epetra_Map &input_col_map,
461 const SparsityPatternType &sparsity_pattern,
462 const bool exchange_data,
464 std_cxx11::shared_ptr<Epetra_FECrsMatrix> &matrix,
470 nonlocal_matrix.reset();
471 nonlocal_matrix_exporter.reset();
473 if (input_row_map.Comm().MyPID() == 0)
476 static_cast<size_type
>(n_global_elements(input_row_map)));
478 static_cast<size_type
>(n_global_elements(input_col_map)));
481 column_space_map.reset (
new Epetra_Map (input_col_map));
490 trilinos_sparsity.
reinit (input_row_map, input_col_map,
491 sparsity_pattern, exchange_data);
492 matrix.reset (
new Epetra_FECrsMatrix
498 const size_type first_row = min_my_gid(input_row_map),
499 last_row = max_my_gid(input_row_map)+1;
500 std::vector<int> n_entries_per_row(last_row-first_row);
502 for (size_type row=first_row; row<last_row; ++row)
503 n_entries_per_row[row-first_row] = sparsity_pattern.row_length(row);
519 std_cxx11::shared_ptr<Epetra_CrsGraph> graph;
520 if (input_row_map.Comm().NumProc() > 1)
521 graph.reset (
new Epetra_CrsGraph (Copy, input_row_map,
522 &n_entries_per_row[0],
true));
524 graph.reset (
new Epetra_CrsGraph (Copy, input_row_map, input_col_map,
525 &n_entries_per_row[0],
true));
532 std::vector<TrilinosWrappers::types::int_type> row_indices;
534 for (size_type row=first_row; row<last_row; ++row)
536 const int row_length = sparsity_pattern.row_length(row);
540 row_indices.resize (row_length, -1);
542 typename SparsityPatternType::iterator p = sparsity_pattern.begin(row);
543 for (size_type col=0; p != sparsity_pattern.end(row); ++p, ++col)
544 row_indices[col] = p->column();
546 graph->Epetra_CrsGraph::InsertGlobalIndices (row, row_length,
554 graph->FillComplete(input_col_map, input_row_map);
555 graph->OptimizeStorage();
559 static_cast<size_type
>(n_global_cols(*graph)));
563 matrix.reset (
new Epetra_FECrsMatrix(Copy, *graph,
false));
575 class Epetra_CrsGraphMod :
public Epetra_CrsGraph
578 Epetra_CrsGraphMod (
const Epetra_Map &row_map,
579 const int *n_entries_per_row)
581 Epetra_CrsGraph(Copy, row_map, n_entries_per_row,
true)
584 void SetIndicesAreGlobal()
586 this->Epetra_CrsGraph::SetIndicesAreGlobal(
true);
595 reinit_matrix (
const Epetra_Map &input_row_map,
596 const Epetra_Map &input_col_map,
598 const bool exchange_data,
599 std_cxx11::shared_ptr<Epetra_Map> &column_space_map,
600 std_cxx11::shared_ptr<Epetra_FECrsMatrix> &matrix,
601 std_cxx11::shared_ptr<Epetra_CrsMatrix> &nonlocal_matrix,
602 std_cxx11::shared_ptr<Epetra_Export> &nonlocal_matrix_exporter)
605 nonlocal_matrix.reset();
606 nonlocal_matrix_exporter.reset();
609 static_cast<size_type
>(n_global_elements(input_row_map)));
611 static_cast<size_type
>(n_global_elements(input_col_map)));
613 column_space_map.reset (
new Epetra_Map (input_col_map));
617 if (relevant_rows.size() == 0)
619 relevant_rows.
set_size(n_global_elements(input_row_map));
620 relevant_rows.add_range(0, n_global_elements(input_row_map));
622 relevant_rows.compress();
623 Assert(relevant_rows.n_elements() >=
static_cast<unsigned int>(input_row_map.NumMyElements()),
624 ExcMessage(
"Locally relevant rows of sparsity pattern must contain " 625 "all locally owned rows"));
630 bool have_ghost_rows =
false;
632 std::vector<::types::global_dof_index> indices;
633 relevant_rows.fill_index_vector(indices);
634 Epetra_Map relevant_map (TrilinosWrappers::types::int_type(-1),
635 TrilinosWrappers::types::int_type(relevant_rows.n_elements()),
636 (indices.empty() ? 0 :
637 reinterpret_cast<TrilinosWrappers::types::int_type *
>(&indices[0])),
638 0, input_row_map.Comm());
639 if (relevant_map.SameAs(input_row_map))
640 have_ghost_rows =
false;
642 have_ghost_rows =
true;
645 const unsigned int n_rows = relevant_rows.n_elements();
646 std::vector<TrilinosWrappers::types::int_type> ghost_rows;
647 std::vector<int> n_entries_per_row(input_row_map.NumMyElements());
648 std::vector<int> n_entries_per_ghost_row;
649 for (
unsigned int i=0, own=0; i<n_rows; ++i)
651 const TrilinosWrappers::types::int_type global_row =
652 relevant_rows.nth_index_in_set(i);
653 if (input_row_map.MyGID(global_row))
654 n_entries_per_row[own++] = sparsity_pattern.
row_length(global_row);
655 else if (sparsity_pattern.
row_length(global_row) > 0)
657 ghost_rows.push_back(global_row);
658 n_entries_per_ghost_row.push_back(sparsity_pattern.
row_length(global_row));
662 Epetra_Map off_processor_map(-1, ghost_rows.size(),
663 (ghost_rows.size()>0)?(&ghost_rows[0]):NULL,
664 0, input_row_map.Comm());
666 std_cxx11::shared_ptr<Epetra_CrsGraph> graph;
667 std_cxx11::shared_ptr<Epetra_CrsGraphMod> nonlocal_graph;
668 if (input_row_map.Comm().NumProc() > 1)
670 graph.reset (
new Epetra_CrsGraph (Copy, input_row_map,
671 (n_entries_per_row.size()>0)?(&n_entries_per_row[0]):NULL,
672 exchange_data ?
false :
true));
673 if (have_ghost_rows ==
true)
674 nonlocal_graph.reset (
new Epetra_CrsGraphMod (off_processor_map,
675 &n_entries_per_ghost_row[0]));
678 graph.reset (
new Epetra_CrsGraph (Copy, input_row_map, input_col_map,
679 (n_entries_per_row.size()>0)?(&n_entries_per_row[0]):NULL,
683 std::vector<TrilinosWrappers::types::int_type> row_indices;
685 for (
unsigned int i=0; i<n_rows; ++i)
687 const TrilinosWrappers::types::int_type global_row =
688 relevant_rows.nth_index_in_set(i);
693 row_indices.resize (row_length, -1);
695 row_indices[col] = sparsity_pattern.
column_number(global_row, col);
697 if (input_row_map.MyGID(global_row))
698 graph->InsertGlobalIndices (global_row, row_length, &row_indices[0]);
701 Assert(nonlocal_graph.get() != 0, ExcInternalError());
702 nonlocal_graph->InsertGlobalIndices (global_row, row_length,
708 if (nonlocal_graph.get() != 0)
714 nonlocal_graph->SetIndicesAreGlobal();
715 Assert(nonlocal_graph->IndicesAreGlobal() ==
true,
717 nonlocal_graph->FillComplete(input_col_map, input_row_map);
718 nonlocal_graph->OptimizeStorage();
723 Epetra_Export exporter(nonlocal_graph->RowMap(), input_row_map);
724 int ierr = graph->Export(*nonlocal_graph, exporter, Add);
726 Assert (ierr==0, ExcTrilinosError(ierr));
729 nonlocal_matrix.reset (
new Epetra_CrsMatrix(Copy, *nonlocal_graph));
732 graph->FillComplete(input_col_map, input_row_map);
733 graph->OptimizeStorage();
736 n_global_cols(*graph)));
738 matrix.reset (
new Epetra_FECrsMatrix(Copy, *graph,
false));
744 template <
typename SparsityPatternType>
748 const Epetra_Map rows (static_cast<TrilinosWrappers::types::int_type>(sparsity_pattern.n_rows()),
751 const Epetra_Map columns (static_cast<TrilinosWrappers::types::int_type>(sparsity_pattern.n_cols()),
755 reinit_matrix (rows, columns, sparsity_pattern,
false,
762 template <
typename SparsityPatternType>
765 const SparsityPatternType &sparsity_pattern,
766 const bool exchange_data)
768 reinit_matrix (input_map, input_map, sparsity_pattern, exchange_data,
778 template <
typename SparsityPatternType>
781 const IndexSet &col_parallel_partitioning,
782 const SparsityPatternType &sparsity_pattern,
783 const MPI_Comm &communicator,
784 const bool exchange_data)
790 reinit_matrix (row_map, col_map, sparsity_pattern, exchange_data,
802 template <
typename SparsityPatternType>
805 const Epetra_Map &col_map,
806 const SparsityPatternType &sparsity_pattern,
807 const bool exchange_data)
809 reinit_matrix (row_map, col_map, sparsity_pattern, exchange_data,
831 matrix.reset (
new Epetra_FECrsMatrix
848 if (
this == &sparse_matrix)
854 matrix.reset (
new Epetra_FECrsMatrix
869 template <
typename number>
872 const IndexSet &col_parallel_partitioning,
873 const ::SparseMatrix<number> &dealii_sparse_matrix,
874 const MPI_Comm &communicator,
875 const double drop_tolerance,
876 const bool copy_values,
877 const ::SparsityPattern *use_this_sparsity)
879 if (copy_values ==
false)
883 if (use_this_sparsity == 0)
884 reinit (row_parallel_partitioning, col_parallel_partitioning,
885 dealii_sparse_matrix.get_sparsity_pattern(),
886 communicator,
false);
888 reinit (row_parallel_partitioning, col_parallel_partitioning,
889 *use_this_sparsity, communicator,
false);
893 const size_type n_rows = dealii_sparse_matrix.m();
898 const ::SparsityPattern &sparsity_pattern =
899 (use_this_sparsity!=0)? *use_this_sparsity :
900 dealii_sparse_matrix.get_sparsity_pattern();
902 if (matrix.get() == 0 ||
906 reinit (row_parallel_partitioning, col_parallel_partitioning,
907 sparsity_pattern, communicator,
false);
915 size_type maximum_row_length = matrix->MaxNumEntries();
916 std::vector<size_type> row_indices (maximum_row_length);
917 std::vector<TrilinosScalar> values (maximum_row_length);
921 if (row_parallel_partitioning.
is_element(row) ==
true)
924 sparsity_pattern.begin(row);
925 typename ::SparseMatrix<number>::const_iterator it =
926 dealii_sparse_matrix.begin(row);
928 if (sparsity_pattern.n_rows() == sparsity_pattern.n_cols())
932 if (std::fabs(it->value()) > drop_tolerance)
934 values[col] = it->value();
935 row_indices[col++] = it->column();
941 while (it != dealii_sparse_matrix.end(row) &&
942 select_index != sparsity_pattern.end(row))
944 while (select_index->
column() < it->column() &&
945 select_index != sparsity_pattern.end(row))
947 while (it->column() < select_index->
column() &&
948 it != dealii_sparse_matrix.end(row))
951 if (it == dealii_sparse_matrix.end(row))
953 if (std::fabs(it->value()) > drop_tolerance)
955 values[col] = it->value();
956 row_indices[col++] = it->column();
961 set (row, col,
reinterpret_cast<size_type *
>(&row_indices[0]),
969 template <
typename number>
972 const double drop_tolerance,
973 const bool copy_values,
974 const ::SparsityPattern *use_this_sparsity)
976 reinit (complete_index_set(dealii_sparse_matrix.m()),
977 complete_index_set(dealii_sparse_matrix.n()),
978 dealii_sparse_matrix, MPI_COMM_SELF, drop_tolerance,
979 copy_values, use_this_sparsity);
984 template <
typename number>
987 const ::SparseMatrix<number> &dealii_sparse_matrix,
988 const double drop_tolerance,
989 const bool copy_values,
990 const ::SparsityPattern *use_this_sparsity)
993 MPI_COMM_SELF, drop_tolerance, copy_values, use_this_sparsity);
998 template <
typename number>
1001 const Epetra_Map &input_col_map,
1002 const ::SparseMatrix<number> &dealii_sparse_matrix,
1003 const double drop_tolerance,
1004 const bool copy_values,
1005 const ::SparsityPattern *use_this_sparsity)
1008 dealii_sparse_matrix, MPI_COMM_SELF,
1009 drop_tolerance, copy_values, use_this_sparsity);
1016 const bool copy_values)
1018 Assert (input_matrix.Filled()==
true,
1019 ExcMessage(
"Input CrsMatrix has not called FillComplete()!"));
1023 const Epetra_CrsGraph *graph = &input_matrix.Graph();
1028 matrix.reset (
new Epetra_FECrsMatrix(Copy, *graph,
false));
1032 if (copy_values ==
true)
1036 const TrilinosScalar *in_values = input_matrix[0];
1037 TrilinosScalar *values = (*matrix)[0];
1038 const size_type my_nonzeros = input_matrix.NumMyNonzeros();
1039 std::memcpy (&values[0], &in_values[0],
1040 my_nonzeros*
sizeof (TrilinosScalar));
1056 if ((operation==::VectorOperation::add) ||
1057 (operation==::VectorOperation::unknown))
1059 else if (operation==::VectorOperation::insert)
1065 ((
last_action == Add) && (operation!=::VectorOperation::insert))
1067 ((
last_action == Insert) && (operation!=::VectorOperation::add)),
1068 ExcMessage(
"Operation and argument to compress() do not match"));
1092 ierr = matrix->OptimizeStorage ();
1114 matrix->FillComplete();
1123 const TrilinosScalar new_diag_value)
1125 Assert (matrix->Filled()==
true, ExcMatrixNotCompressed());
1130 matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(row));
1133 TrilinosScalar *values;
1136 const int ierr = matrix->ExtractMyRowView(local_row, num_entries,
1137 values, col_indices);
1141 ExcTrilinosError(ierr));
1143 int *diag_find = std::find(col_indices,col_indices+num_entries,local_row);
1144 int diag_index = (int)(diag_find - col_indices);
1146 for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
1147 if (diag_index != j || new_diag_value == 0)
1150 if (diag_find && std::fabs(values[diag_index]) == 0.0 &&
1151 new_diag_value != 0.0)
1152 values[diag_index] = new_diag_value;
1160 const TrilinosScalar new_diag_value)
1162 for (
size_type row=0; row<rows.size(); ++row)
1174 int trilinos_i = matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(i)),
1175 trilinos_j = matrix->LCID(static_cast<TrilinosWrappers::types::int_type>(j));
1176 TrilinosScalar value = 0.;
1184 if (trilinos_i == -1)
1194 Assert (matrix->Filled(), ExcMatrixNotCompressed());
1198 int nnz_present = matrix->NumMyEntries(trilinos_i);
1201 TrilinosScalar *values;
1207 int ierr = matrix->ExtractMyRowView(trilinos_i, nnz_extracted,
1208 values, col_indices);
1210 Assert (ierr==0, ExcTrilinosError(ierr));
1212 Assert (nnz_present == nnz_extracted,
1213 ExcDimensionMismatch(nnz_present, nnz_extracted));
1219 int *el_find = std::find(col_indices, col_indices + nnz_present, trilinos_j);
1221 int local_col_index = (int)(el_find - col_indices);
1231 if (local_col_index == nnz_present)
1233 Assert (
false, ExcInvalidIndex (i,j));
1236 value = values[local_col_index];
1250 int trilinos_i = matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(i)),
1251 trilinos_j = matrix->LCID(static_cast<TrilinosWrappers::types::int_type>(j));
1252 TrilinosScalar value = 0.;
1261 if ((trilinos_i == -1 ) || (trilinos_j == -1))
1268 Assert (matrix->Filled(), ExcMatrixNotCompressed());
1272 int nnz_present = matrix->NumMyEntries(trilinos_i);
1275 TrilinosScalar *values;
1280 int ierr = matrix->ExtractMyRowView(trilinos_i, nnz_extracted,
1281 values, col_indices);
1283 Assert (ierr==0, ExcTrilinosError(ierr));
1285 Assert (nnz_present == nnz_extracted,
1286 ExcDimensionMismatch(nnz_present, nnz_extracted));
1291 int *el_find = std::find(col_indices, col_indices + nnz_present, trilinos_j);
1293 int local_col_index = (int)(el_find - col_indices);
1303 if (local_col_index == nnz_present)
1306 value = values[local_col_index];
1317 Assert (
m() ==
n(), ExcNotQuadratic());
1338 Assert (row <
m(), ExcInternalError());
1343 int local_row = matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(row));
1350 int ierr = matrix->NumMyRowEntries (local_row, ncols);
1361 const std::vector<size_type> &col_indices,
1363 const bool elide_zero_values)
1365 Assert (row_indices.size() == values.
m(),
1366 ExcDimensionMismatch(row_indices.size(), values.
m()));
1367 Assert (col_indices.size() == values.
n(),
1368 ExcDimensionMismatch(col_indices.size(), values.
n()));
1370 for (
size_type i=0; i<row_indices.size(); ++i)
1371 set (row_indices[i], col_indices.size(), &col_indices[0], &values(i,0),
1379 const std::vector<size_type> &col_indices,
1380 const std::vector<TrilinosScalar> &values,
1381 const bool elide_zero_values)
1383 Assert (col_indices.size() == values.size(),
1384 ExcDimensionMismatch(col_indices.size(), values.size()));
1386 set (row, col_indices.size(), &col_indices[0], &values[0],
1396 const TrilinosScalar *values,
1397 const bool elide_zero_values)
1407 Assert (ierr == 0, ExcTrilinosError(ierr));
1412 TrilinosWrappers::types::int_type *col_index_ptr;
1413 TrilinosScalar *col_value_ptr;
1414 TrilinosWrappers::types::int_type n_columns;
1416 TrilinosScalar short_val_array[100];
1417 TrilinosWrappers::types::int_type short_index_array[100];
1418 std::vector<TrilinosScalar> long_val_array;
1419 std::vector<TrilinosWrappers::types::int_type> long_index_array;
1425 if (elide_zero_values ==
false)
1427 col_index_ptr = (TrilinosWrappers::types::int_type *)col_indices;
1428 col_value_ptr =
const_cast<TrilinosScalar *
>(values);
1437 long_val_array.resize(n_cols);
1438 long_index_array.resize(n_cols);
1439 col_index_ptr = &long_index_array[0];
1440 col_value_ptr = &long_val_array[0];
1444 col_index_ptr = &short_index_array[0];
1445 col_value_ptr = &short_val_array[0];
1451 const double value = values[j];
1455 col_index_ptr[n_columns] = col_indices[j];
1456 col_value_ptr[n_columns] = value;
1461 Assert(n_columns <= (TrilinosWrappers::types::int_type)n_cols, ExcInternalError());
1473 if (matrix->RowMap().MyGID(static_cast<TrilinosWrappers::types::int_type>(row)) ==
true)
1475 if (matrix->Filled() ==
false)
1477 ierr = matrix->Epetra_CrsMatrix::InsertGlobalValues(
1478 static_cast<TrilinosWrappers::types::int_type>(row),
1479 static_cast<int>(n_columns),const_cast<double *>(col_value_ptr),
1489 ierr = matrix->Epetra_CrsMatrix::ReplaceGlobalValues(row, n_columns,
1503 if (matrix->Filled() ==
false)
1505 ierr = matrix->InsertGlobalValues (1,
1506 (TrilinosWrappers::types::int_type *)&row,
1507 n_columns, col_index_ptr,
1509 Epetra_FECrsMatrix::ROW_MAJOR);
1514 ierr = matrix->ReplaceGlobalValues (1,
1515 (TrilinosWrappers::types::int_type *)&row,
1516 n_columns, col_index_ptr,
1518 Epetra_FECrsMatrix::ROW_MAJOR);
1526 Assert (ierr <= 0, ExcAccessToNonPresentElement(row, col_index_ptr[0]));
1535 const bool elide_zero_values)
1537 Assert (indices.size() == values.
m(),
1538 ExcDimensionMismatch(indices.size(), values.
m()));
1539 Assert (values.
m() == values.
n(), ExcNotQuadratic());
1541 for (
size_type i=0; i<indices.size(); ++i)
1542 add (indices[i], indices.size(), &indices[0], &values(i,0),
1550 const std::vector<size_type> &col_indices,
1552 const bool elide_zero_values)
1554 Assert (row_indices.size() == values.
m(),
1555 ExcDimensionMismatch(row_indices.size(), values.
m()));
1556 Assert (col_indices.size() == values.
n(),
1557 ExcDimensionMismatch(col_indices.size(), values.
n()));
1559 for (
size_type i=0; i<row_indices.size(); ++i)
1560 add (row_indices[i], col_indices.size(), &col_indices[0],
1561 &values(i,0), elide_zero_values);
1568 const std::vector<size_type> &col_indices,
1569 const std::vector<TrilinosScalar> &values,
1570 const bool elide_zero_values)
1572 Assert (col_indices.size() == values.size(),
1573 ExcDimensionMismatch(col_indices.size(), values.size()));
1575 add (row, col_indices.size(), &col_indices[0], &values[0],
1585 const TrilinosScalar *values,
1586 const bool elide_zero_values,
1596 matrix->RowMap(),
false);
1603 TrilinosWrappers::types::int_type *col_index_ptr;
1604 TrilinosScalar *col_value_ptr;
1605 TrilinosWrappers::types::int_type n_columns;
1607 double short_val_array[100];
1608 TrilinosWrappers::types::int_type short_index_array[100];
1609 std::vector<TrilinosScalar> long_val_array;
1610 std::vector<TrilinosWrappers::types::int_type> long_index_array;
1615 if (elide_zero_values ==
false)
1617 col_index_ptr = (TrilinosWrappers::types::int_type *)col_indices;
1618 col_value_ptr =
const_cast<TrilinosScalar *
>(values);
1631 long_val_array.resize(n_cols);
1632 long_index_array.resize(n_cols);
1633 col_index_ptr = &long_index_array[0];
1634 col_value_ptr = &long_val_array[0];
1638 col_index_ptr = &short_index_array[0];
1639 col_value_ptr = &short_val_array[0];
1645 const double value = values[j];
1650 col_index_ptr[n_columns] = col_indices[j];
1651 col_value_ptr[n_columns] = value;
1656 Assert(n_columns <= (TrilinosWrappers::types::int_type)n_cols, ExcInternalError());
1663 if (matrix->RowMap().MyGID(static_cast<TrilinosWrappers::types::int_type>(row)) ==
true)
1665 ierr = matrix->Epetra_CrsMatrix::SumIntoGlobalValues(row, n_columns,
1676 ExcMessage(
"Attempted to write into off-processor matrix row " 1677 "that has not be specified as being writable upon " 1693 ierr = matrix->SumIntoGlobalValues (1,
1694 (TrilinosWrappers::types::int_type *)&row, n_columns,
1697 Epetra_FECrsMatrix::ROW_MAJOR);
1703 std::cout <<
"------------------------------------------" 1705 std::cout <<
"Got error " << ierr <<
" in row " << row
1706 <<
" of proc " << matrix->RowMap().Comm().MyPID()
1707 <<
" when trying to add the columns:" << std::endl;
1708 for (TrilinosWrappers::types::int_type i=0; i<n_columns; ++i)
1709 std::cout << col_index_ptr[i] <<
" ";
1710 std::cout << std::endl << std::endl;
1711 std::cout <<
"Matrix row " 1712 << (matrix->RowMap().MyGID(static_cast<TrilinosWrappers::types::int_type>(row)) ==
false ?
"(nonlocal part)" :
"")
1713 <<
" has the following indices:" << std::endl;
1714 std::vector<TrilinosWrappers::types::int_type> indices;
1715 const Epetra_CrsGraph *graph =
1717 matrix->RowMap().MyGID(static_cast<TrilinosWrappers::types::int_type>(row)) ==
false) ?
1720 indices.resize(graph->NumGlobalIndices(static_cast<TrilinosWrappers::types::int_type>(row)));
1722 graph->ExtractGlobalRowCopy(static_cast<TrilinosWrappers::types::int_type>(row),
1723 indices.size(), n_indices, &indices[0]);
1724 AssertDimension(static_cast<unsigned int>(n_indices), indices.size());
1726 for (TrilinosWrappers::types::int_type i=0; i<n_indices; ++i)
1727 std::cout << indices[i] <<
" ";
1728 std::cout << std::endl << std::endl;
1730 ExcAccessToNonPresentElement(row, col_index_ptr[0]));
1733 Assert (ierr >= 0, ExcTrilinosError(ierr));
1741 Assert (d==0, ExcScalarAssignmentOnlyForZeroValue());
1742 compress (::VectorOperation::unknown);
1744 const int ierr = matrix->PutScalar(d);
1763 ExcMessage(
"Can only add matrices with same distribution of rows"));
1765 ExcMessage(
"Addition of matrices only allowed if matrices are " 1766 "filled, i.e., compress() has been called"));
1768 const std::pair<size_type, size_type>
1770 const bool same_col_map = matrix->ColMap().SameAs(rhs.
matrix->ColMap());
1773 for (
size_type row=local_range.first; row < local_range.second; ++row)
1775 const int row_local =
1776 matrix->RowMap().LID(static_cast<TrilinosWrappers::types::int_type>(row));
1782 int n_entries, rhs_n_entries;
1783 TrilinosScalar *value_ptr, *rhs_value_ptr;
1784 int *index_ptr, *rhs_index_ptr;
1785 ierr = rhs.
matrix->ExtractMyRowView (row_local, rhs_n_entries,
1786 rhs_value_ptr, rhs_index_ptr);
1788 Assert (ierr == 0, ExcTrilinosError(ierr));
1790 ierr = matrix->ExtractMyRowView (row_local, n_entries, value_ptr,
1792 Assert (ierr == 0, ExcTrilinosError(ierr));
1793 bool expensive_checks = (n_entries != rhs_n_entries || !same_col_map);
1794 if (!expensive_checks)
1798 expensive_checks = std::memcmp(static_cast<void *>(index_ptr),
1799 static_cast<void *>(rhs_index_ptr),
1800 sizeof(
int)*n_entries) != 0;
1801 if (!expensive_checks)
1802 for (
int i=0; i<n_entries; ++i)
1803 value_ptr[i] += rhs_value_ptr[i] * factor;
1809 if (expensive_checks)
1811 for (
int i=0; i<rhs_n_entries; ++i)
1813 if (rhs_value_ptr[i] == 0.)
1815 const TrilinosWrappers::types::int_type rhs_global_col =
1816 global_column_index(*rhs.
matrix, rhs_index_ptr[i]);
1817 int local_col = matrix->ColMap().LID(rhs_global_col);
1819 index_ptr+n_entries,
1821 Assert(local_index != index_ptr + n_entries &&
1822 *local_index == local_col,
1823 ExcMessage(
"Adding the entries from the other matrix " 1824 "failed, because the sparsity pattern " 1825 "of that matrix includes more elements than the " 1826 "calling matrix, which is not allowed."));
1827 value_ptr[local_index-index_ptr] += factor * rhs_value_ptr[i];
1845 if (!matrix->UseTranspose())
1847 ierr = matrix->SetUseTranspose (
true);
1852 ierr = matrix->SetUseTranspose (
false);
1862 const int ierr = matrix->Scale (a);
1863 Assert (ierr == 0, ExcTrilinosError(ierr));
1874 Assert (a !=0, ExcDivideByZero());
1876 const TrilinosScalar factor = 1./a;
1878 const int ierr = matrix->Scale (factor);
1879 Assert (ierr == 0, ExcTrilinosError(ierr));
1890 Assert (matrix->Filled(), ExcMatrixNotCompressed());
1891 return matrix->NormOne();
1899 Assert (matrix->Filled(), ExcMatrixNotCompressed());
1900 return matrix->NormInf();
1908 Assert (matrix->Filled(), ExcMatrixNotCompressed());
1909 return matrix->NormFrobenius();
1918 template <
typename VectorType>
1920 void check_vector_map_equality(
const Epetra_CrsMatrix &,
1927 void check_vector_map_equality(
const Epetra_CrsMatrix &
m,
1932 ExcMessage (
"Column map of matrix does not fit with vector map!"));
1934 ExcMessage (
"Row map of matrix does not fit with vector map!"));
1943 template <
typename VectorType>
1946 const VectorType &src)
const 1948 Assert (&src != &dst, ExcSourceEqualsDestination());
1949 Assert (matrix->Filled(), ExcMatrixNotCompressed());
1953 internal::SparseMatrix::check_vector_map_equality(*matrix, src, dst);
1954 const size_type dst_local_size = dst.end() - dst.begin();
1955 AssertDimension (dst_local_size, static_cast<size_type>(matrix->RangeMap().NumMyPoints()));
1956 const size_type src_local_size = src.end() - src.begin();
1957 AssertDimension (src_local_size, static_cast<size_type>(matrix->DomainMap().NumMyPoints()));
1959 Epetra_MultiVector tril_dst (View, matrix->RangeMap(), dst.begin(),
1961 Epetra_MultiVector tril_src (View, matrix->DomainMap(),
1962 const_cast<TrilinosScalar *
>(src.begin()),
1965 const int ierr = matrix->Multiply (
false, tril_src, tril_dst);
1966 Assert (ierr == 0, ExcTrilinosError(ierr));
1972 template <
typename VectorType>
1975 const VectorType &src)
const 1977 Assert (&src != &dst, ExcSourceEqualsDestination());
1978 Assert (matrix->Filled(), ExcMatrixNotCompressed());
1980 internal::SparseMatrix::check_vector_map_equality(*matrix, dst, src);
1981 const size_type dst_local_size = dst.end() - dst.begin();
1982 AssertDimension (dst_local_size, static_cast<size_type>(matrix->DomainMap().NumMyPoints()));
1983 const size_type src_local_size = src.end() - src.begin();
1984 AssertDimension (src_local_size, static_cast<size_type>(matrix->RangeMap().NumMyPoints()));
1986 Epetra_MultiVector tril_dst (View, matrix->DomainMap(), dst.begin(),
1988 Epetra_MultiVector tril_src (View, matrix->RangeMap(),
1989 const_cast<double *
>(src.begin()),
1992 const int ierr = matrix->Multiply (
true, tril_src, tril_dst);
1993 Assert (ierr == 0, ExcTrilinosError(ierr));
1999 template <
typename VectorType>
2002 const VectorType &src)
const 2004 Assert (&src != &dst, ExcSourceEqualsDestination());
2012 temp_vector.
reinit(dst.end()-dst.begin(),
true);
2015 vmult (temp_vector,
static_cast<const ::Vector<TrilinosScalar>&
>(src_view));
2016 if (dst_view.size() > 0)
2017 dst_view += temp_vector;
2022 template <
typename VectorType>
2025 const VectorType &src)
const 2027 Assert (&src != &dst, ExcSourceEqualsDestination());
2035 temp_vector.
reinit(dst.end()-dst.begin(),
true);
2038 Tvmult (temp_vector,
static_cast<const ::Vector<TrilinosScalar>&
>(src_view));
2039 if (dst_view.size() > 0)
2040 dst_view += temp_vector;
2048 Assert (matrix->RowMap().SameAs(matrix->DomainMap()),
2052 temp_vector.
reinit(v,
true);
2054 vmult (temp_vector, v);
2055 return temp_vector*v;
2064 Assert (matrix->RowMap().SameAs(matrix->DomainMap()),
2068 temp_vector.
reinit(v,
true);
2070 vmult (temp_vector, v);
2071 return u*temp_vector;
2092 typedef ::types::global_dof_index
size_type;
2098 const bool transpose_left)
2100 #ifdef DEAL_II_WITH_64BIT_INDICES 2101 Assert(
false,ExcNotImplemented())
2103 const bool use_vector = (V.
size() == inputright.
m() ?
true :
false);
2104 if (transpose_left ==
false)
2106 Assert (inputleft.
n() == inputright.
m(),
2107 ExcDimensionMismatch(inputleft.
n(), inputright.
m()));
2109 ExcMessage (
"Parallel partitioning of A and B does not fit."));
2113 Assert (inputleft.
m() == inputright.
m(),
2114 ExcDimensionMismatch(inputleft.
m(), inputright.
m()));
2116 ExcMessage (
"Parallel partitioning of A and B does not fit."));
2127 Teuchos::RCP<Epetra_CrsMatrix> mod_B;
2128 if (use_vector ==
false)
2130 mod_B = Teuchos::rcp(const_cast<Epetra_CrsMatrix *>
2136 mod_B = Teuchos::rcp(
new Epetra_CrsMatrix
2142 ExcMessage (
"Parallel distribution of matrix B and vector V " 2143 "does not match."));
2146 for (
int i=0; i<local_N; ++i)
2149 double *new_data, *B_data;
2150 mod_B->ExtractMyRowView (i, N_entries, new_data);
2153 for (TrilinosWrappers::types::int_type j=0; j<N_entries; ++j)
2154 new_data[j] = value * B_data[j];
2164 ML_Comm_Create(&comm);
2166 const Epetra_MpiComm *epcomm =
dynamic_cast<const Epetra_MpiComm *
>(&(inputleft.
trilinos_matrix().Comm()));
2168 if (epcomm) ML_Comm_Set_UsrComm(comm,epcomm->Comm());
2170 ML_Operator *A_ = ML_Operator_Create(comm);
2171 ML_Operator *B_ = ML_Operator_Create(comm);
2172 ML_Operator *C_ = ML_Operator_Create(comm);
2175 if (transpose_left ==
false)
2176 ML_Operator_WrapEpetraCrsMatrix
2185 ExcMessage(
"Matrix must be partitioned contiguously between procs."));
2186 for (
unsigned int i=0; i<inputleft.
local_size(); ++i)
2188 int num_entries, * indices;
2191 Assert (num_entries >= 0, ExcInternalError());
2192 #ifndef DEAL_II_WITH_64BIT_INDICES 2194 for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
2195 sparsity_transposed.add (inputleft.
trilinos_matrix().ColMap().GID(indices[j]),
2199 for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
2200 sparsity_transposed.add (inputleft.
trilinos_matrix().ColMap().GID64(indices[j]),
2205 sparsity_transposed.compress();
2206 transposed_mat.
reinit (sparsity_transposed);
2207 for (
unsigned int i=0; i<inputleft.
local_size(); ++i)
2209 int num_entries, * indices;
2213 Assert (num_entries >= 0, ExcInternalError());
2214 #ifndef DEAL_II_WITH_64BIT_INDICES 2216 for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
2221 for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
2226 transposed_mat.
compress(VectorOperation::insert);
2227 ML_Operator_WrapEpetraCrsMatrix
2231 ML_Operator_WrapEpetraCrsMatrix(mod_B.get(),B_,
false);
2241 ML_Operator *Btmp, *Ctmp, *Ctmp2, *tptr;
2242 ML_CommInfoOP *getrow_comm;
2244 TrilinosWrappers::types::int_type N_input_vector = B_->invec_leng;
2245 getrow_comm = B_->getrow->pre_comm;
2246 if ( getrow_comm != NULL)
2247 for (TrilinosWrappers::types::int_type i = 0; i < getrow_comm->N_neighbors; i++)
2248 for (TrilinosWrappers::types::int_type j = 0; j < getrow_comm->neighbors[i].N_send; j++)
2249 AssertThrow (getrow_comm->neighbors[i].send_list[j] < N_input_vector,
2250 ExcInternalError());
2252 ML_create_unique_col_id(N_input_vector, &(B_->getrow->loc_glob_map),
2253 getrow_comm, &max_per_proc, B_->comm);
2254 B_->getrow->use_loc_glob_map = ML_YES;
2255 if (A_->getrow->pre_comm != NULL)
2256 ML_exchange_rows( B_, &Btmp, A_->getrow->pre_comm);
2260 ML_matmat_mult(A_, Btmp , &Ctmp);
2264 ML_free(B_->getrow->loc_glob_map);
2265 B_->getrow->loc_glob_map = NULL;
2266 B_->getrow->use_loc_glob_map = ML_NO;
2267 if (A_->getrow->pre_comm != NULL)
2270 while ( (tptr!= NULL) && (tptr->sub_matrix != B_))
2271 tptr = tptr->sub_matrix;
2272 if (tptr != NULL) tptr->sub_matrix = NULL;
2273 ML_RECUR_CSR_MSRdata_Destroy(Btmp);
2274 ML_Operator_Destroy(&Btmp);
2278 if (A_->getrow->post_comm != NULL)
2279 ML_exchange_rows(Ctmp, &Ctmp2, A_->getrow->post_comm);
2283 ML_back_to_csrlocal(Ctmp2, C_, max_per_proc);
2285 ML_RECUR_CSR_MSRdata_Destroy (Ctmp);
2286 ML_Operator_Destroy (&Ctmp);
2288 if (A_->getrow->post_comm != NULL)
2290 ML_RECUR_CSR_MSRdata_Destroy(Ctmp2);
2291 ML_Operator_Destroy (&Ctmp2);
2296 Epetra_CrsMatrix *C_mat;
2297 ML_Operator2EpetraCrsMatrix(C_, C_mat);
2298 C_mat->FillComplete();
2299 C_mat->OptimizeStorage();
2304 ML_Operator_Destroy (&A_);
2305 ML_Operator_Destroy (&B_);
2306 ML_Operator_Destroy (&C_);
2307 ML_Comm_Destroy (&comm);
2317 #ifdef DEAL_II_WITH_64BIT_INDICES 2318 Assert(
false,ExcNotImplemented())
2320 internals::perform_mmult (*
this, B, C, V,
false);
2330 #ifdef DEAL_II_WITH_64BIT_INDICES 2331 Assert(
false,ExcNotImplemented())
2333 internals::perform_mmult (*
this, B, C, V,
true);
2341 Assert (
false, ExcNotImplemented());
2351 const bool print_detailed_trilinos_information)
const 2353 if (print_detailed_trilinos_information ==
true)
2361 for (
int i=0; i<matrix->NumMyRows(); ++i)
2363 matrix->ExtractMyRowView (i, num_entries, values, indices);
2364 for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
2365 out <<
"(" << global_row_index(*matrix,i) <<
"," 2366 << global_column_index(*matrix,indices[j]) <<
") " 2367 << values[j] << std::endl;
2379 size_type static_memory =
sizeof(
this) +
sizeof (*matrix)
2380 +
sizeof(*matrix->Graph().DataPtr());
2381 return ((
sizeof(TrilinosScalar)+
sizeof(TrilinosWrappers::types::int_type))*
2382 matrix->NumMyNonzeros() +
sizeof(int)*
local_size() + static_memory);
2390 return matrix->DomainMap();
2398 return matrix->RangeMap();
2406 return matrix->RowMap();
2414 return matrix->ColMap();
2422 #ifdef DEAL_II_WITH_MPI 2424 const Epetra_MpiComm *mpi_comm
2425 =
dynamic_cast<const Epetra_MpiComm *
>(&matrix->RangeMap().Comm());
2426 return mpi_comm->Comm();
2429 return MPI_COMM_SELF;
2439 #include "trilinos_sparse_matrix.inst" 2453 const ::SparsityPattern &,
2462 const ::SparsityPattern &,
2472 const ::SparsityPattern &,
2493 const ::Vector<double> &)
const;
2496 const ::parallel::distributed::Vector<double> &)
const;
2508 const ::Vector<double> &)
const;
2511 const ::parallel::distributed::Vector<double> &)
const;
2523 const ::Vector<double> &)
const;
2526 const ::parallel::distributed::Vector<double> &)
const;
2538 const ::Vector<double> &)
const;
2541 const ::parallel::distributed::Vector<double> &)
const;
2544 DEAL_II_NAMESPACE_CLOSE
2546 #endif // DEAL_II_WITH_TRILINOS std_cxx11::shared_ptr< std::vector< TrilinosScalar > > value_cache
TrilinosScalar diag_element(const size_type i) const
Iterator lower_bound(Iterator first, Iterator last, const T &val)
TrilinosScalar el(const size_type i, const size_type j) const
TrilinosScalar frobenius_norm() const
std_cxx11::shared_ptr< std::vector< size_type > > colnum_cache
#define AssertDimension(dim1, dim2)
void clear_row(const size_type row, const TrilinosScalar new_diag_value=0)
size_type column_number(const size_type row, const size_type index) const
void vmult_add(VectorType &dst, const VectorType &src) const
Epetra_CombineMode last_action
const Epetra_FECrsGraph & trilinos_sparsity_pattern() const
::ExceptionBase & ExcMessage(std::string arg1)
::types::global_dof_index size_type
#define AssertIndexRange(index, range)
TrilinosScalar residual(VectorBase &dst, const VectorBase &x, const VectorBase &b) const
std::pair< size_type, size_type > local_range() const
std_cxx11::shared_ptr< Epetra_CrsGraph > nonlocal_graph
TrilinosScalar operator()(const size_type i, const size_type j) const
void reinit(const VectorBase &v, const bool omit_zeroing_entries=false)
#define AssertThrow(cond, exc)
const Epetra_Comm & comm_self()
void clear_rows(const std::vector< size_type > &rows, const TrilinosScalar new_diag_value=0)
void Tvmult(VectorType &dst, const VectorType &src) const
SparseMatrix & operator=(const double d)
std_cxx11::shared_ptr< Epetra_FECrsMatrix > matrix
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
SparseMatrix & operator*=(const TrilinosScalar factor)
const Epetra_Map & col_partitioner() const DEAL_II_DEPRECATED
std_cxx11::shared_ptr< Epetra_CrsMatrix > nonlocal_matrix
void add(const size_type i, const size_type j, const TrilinosScalar value)
#define Assert(cond, exc)
std_cxx11::shared_ptr< Epetra_Map > column_space_map
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const VectorBase &V=VectorBase()) const
const Epetra_Map & vector_partitioner() const
size_type memory_consumption() const
const IndexSet & row_index_set() const
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
void mmult(SparseMatrix &C, const SparseMatrix &B, const VectorBase &V=VectorBase()) const
const Epetra_CrsMatrix & trilinos_matrix() const
void vmult(VectorType &dst, const VectorType &src) const
size_type n_nonzero_elements() const
MPI_Comm get_mpi_communicator() const
const Epetra_Map & domain_partitioner() const DEAL_II_DEPRECATED
unsigned int row_length(const size_type row) const
void set_size(const size_type size)
void Tvmult_add(VectorType &dst, const VectorType &src) const
const Epetra_Map & row_partitioner() const DEAL_II_DEPRECATED
void reinit(const SparsityPatternType &sparsity_pattern)
TrilinosScalar matrix_norm_square(const VectorBase &v) const
size_type row_length(const size_type row) const
TrilinosScalar linfty_norm() const
void copy_from(const SparseMatrix &source)
SparseMatrix & operator/=(const TrilinosScalar factor)
std::pair< size_type, size_type > local_range() const
const Epetra_MultiVector & trilinos_vector() const
TrilinosScalar matrix_scalar_product(const VectorBase &u, const VectorBase &v) const
const Epetra_CrsGraph & trilinos_sparsity_pattern() const
const Epetra_Map & domain_partitioner() const DEAL_II_DEPRECATED
bool is_element(const size_type index) const
void reinit(const size_type n, const bool omit_zeroing_entries=false)
void compress(::VectorOperation::values operation)
void reinit(const size_type m, const size_type n, const size_type n_entries_per_row=0)
void set(const size_type i, const size_type j, const TrilinosScalar value)
#define AssertIsFinite(number)
const Epetra_Map & range_partitioner() const DEAL_II_DEPRECATED
unsigned int local_size() const
real_type l2_norm() const
TrilinosScalar l1_norm() const
std_cxx11::shared_ptr< Epetra_Export > nonlocal_matrix_exporter