16 #include <deal.II/lac/trilinos_index_access.h> 17 #include <deal.II/lac/trilinos_sparse_matrix.h> 19 #ifdef DEAL_II_WITH_TRILINOS 21 # include <deal.II/base/utilities.h> 22 # include <deal.II/lac/sparse_matrix.h> 23 # include <deal.II/lac/trilinos_sparsity_pattern.h> 24 # include <deal.II/lac/sparsity_pattern.h> 25 # include <deal.II/lac/dynamic_sparsity_pattern.h> 26 # include <deal.II/lac/sparsity_tools.h> 27 # include <deal.II/lac/la_parallel_vector.h> 28 # include <deal.II/lac/trilinos_index_access.h> 29 # include <deal.II/lac/trilinos_precondition.h> 31 # include <Epetra_Export.h> 32 # include <ml_epetra_utils.h> 33 # include <ml_struct.h> 34 # include <Teuchos_RCP.hpp> 36 #include <boost/container/small_vector.hpp> 39 DEAL_II_NAMESPACE_OPEN
45 template <
typename VectorType>
46 double *begin(VectorType &V)
51 template <
typename VectorType>
52 const double *begin(
const VectorType &V)
57 template <
typename VectorType>
58 double *end(VectorType &V)
63 template <
typename VectorType>
64 const double *end(
const VectorType &V)
69 #ifdef DEAL_II_WITH_MPI 120 TrilinosWrappers::types::int_type colnums =
matrix->
n();
133 ExtractGlobalRowCopy((TrilinosWrappers::types::int_type)this->
a_row,
136 reinterpret_cast<TrilinosWrappers::types::int_type *>(&((*
colnum_cache)[0])));
164 column_space_map (new Epetra_Map (0, 0,
166 matrix (new Epetra_FECrsMatrix(View, *column_space_map,
167 *column_space_map, 0)),
179 column_space_map (new Epetra_Map (input_map)),
180 matrix (new Epetra_FECrsMatrix(Copy, *column_space_map,
189 const std::vector<unsigned int> &n_entries_per_row)
191 column_space_map (new Epetra_Map (input_map)),
192 matrix (new Epetra_FECrsMatrix
193 (Copy, *column_space_map,
194 (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
203 const Epetra_Map &input_col_map,
206 column_space_map (new Epetra_Map (input_col_map)),
207 matrix (new Epetra_FECrsMatrix(Copy, input_row_map,
216 const Epetra_Map &input_col_map,
217 const std::vector<unsigned int> &n_entries_per_row)
219 column_space_map (new Epetra_Map (input_col_map)),
220 matrix (new Epetra_FECrsMatrix(Copy, input_row_map,
221 (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
231 const unsigned int n_max_entries_per_row)
244 matrix (new Epetra_FECrsMatrix(Copy,
248 n_max_entries_per_row,
258 const std::vector<unsigned int> &n_entries_per_row)
262 matrix (new Epetra_FECrsMatrix(Copy,
266 (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
275 const MPI_Comm &communicator,
276 const unsigned int n_max_entries_per_row)
278 column_space_map (new Epetra_Map(parallel_partitioning.
279 make_trilinos_map(communicator, false))),
280 matrix (new Epetra_FECrsMatrix(Copy,
282 n_max_entries_per_row,
291 const MPI_Comm &communicator,
292 const std::vector<unsigned int> &n_entries_per_row)
294 column_space_map (new Epetra_Map(parallel_partitioning.
295 make_trilinos_map(communicator, false))),
296 matrix (new Epetra_FECrsMatrix(Copy,
298 (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
307 const IndexSet &col_parallel_partitioning,
308 const MPI_Comm &communicator,
311 column_space_map (new Epetra_Map(col_parallel_partitioning.
312 make_trilinos_map(communicator, false))),
313 matrix (new Epetra_FECrsMatrix(Copy,
314 row_parallel_partitioning.
315 make_trilinos_map(communicator, false),
316 n_max_entries_per_row,
325 const IndexSet &col_parallel_partitioning,
326 const MPI_Comm &communicator,
327 const std::vector<unsigned int> &n_entries_per_row)
329 column_space_map (new Epetra_Map(col_parallel_partitioning.
330 make_trilinos_map(communicator, false))),
331 matrix (new Epetra_FECrsMatrix(Copy,
332 row_parallel_partitioning.
333 make_trilinos_map(communicator, false),
334 (int *)const_cast<unsigned int *>(&(n_entries_per_row[0])),
344 column_space_map (new Epetra_Map (sparsity_pattern.domain_partitioner())),
345 matrix (new Epetra_FECrsMatrix(Copy,
346 sparsity_pattern.trilinos_sparsity_pattern(),
352 ExcMessage(
"The Trilinos sparsity pattern has not been compressed."));
360 column_space_map(std::move(other.column_space_map)),
361 matrix(std::move(other.matrix)),
362 nonlocal_matrix(std::move(other.nonlocal_matrix)),
363 nonlocal_matrix_exporter(std::move(other.nonlocal_matrix_exporter)),
364 last_action(other.last_action),
365 compressed(other.compressed)
367 other.last_action = Zero;
368 other.compressed =
false;
385 bool needs_deep_copy =
390 if (!needs_deep_copy)
392 const std::pair<size_type, size_type>
400 const int row_local =
401 matrix->RowMap().LID(static_cast<TrilinosWrappers::types::int_type>(row));
403 int n_entries, rhs_n_entries;
404 TrilinosScalar *value_ptr, *rhs_value_ptr;
405 int *index_ptr, *rhs_index_ptr;
406 int ierr = rhs.
matrix->ExtractMyRowView (row_local, rhs_n_entries,
407 rhs_value_ptr, rhs_index_ptr);
411 ierr =
matrix->ExtractMyRowView (row_local, n_entries, value_ptr,
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];
434 matrix = std_cxx14::make_unique<Epetra_FECrsMatrix>(*rhs.
matrix);
449 template <
typename SparsityPatternType>
451 reinit_matrix (
const Epetra_Map &input_row_map,
452 const Epetra_Map &input_col_map,
453 const SparsityPatternType &sparsity_pattern,
454 const bool exchange_data,
455 std::unique_ptr<Epetra_Map> &column_space_map,
456 std::unique_ptr<Epetra_FECrsMatrix> &matrix,
457 std::unique_ptr<Epetra_CrsMatrix> &nonlocal_matrix,
458 std::unique_ptr<Epetra_Export> &nonlocal_matrix_exporter)
462 nonlocal_matrix.reset();
463 nonlocal_matrix_exporter.reset();
465 if (input_row_map.Comm().MyPID() == 0)
473 column_space_map = std_cxx14::make_unique<Epetra_Map>(input_col_map);
482 trilinos_sparsity.
reinit (input_row_map, input_col_map,
483 sparsity_pattern, exchange_data);
484 matrix = std_cxx14::make_unique<Epetra_FECrsMatrix>
485 (Copy, trilinos_sparsity.trilinos_sparsity_pattern(),
false);
492 std::vector<int> n_entries_per_row(last_row-first_row);
494 for (size_type row=first_row; row<last_row; ++row)
495 n_entries_per_row[row-first_row] = sparsity_pattern.row_length(row);
511 std::unique_ptr<Epetra_CrsGraph> graph;
512 if (input_row_map.Comm().NumProc() > 1)
513 graph = std_cxx14::make_unique<Epetra_CrsGraph>(Copy, input_row_map,
514 n_entries_per_row.data(),
true);
516 graph = std_cxx14::make_unique<Epetra_CrsGraph>(Copy, input_row_map, input_col_map,
517 n_entries_per_row.data(),
true);
524 std::vector<TrilinosWrappers::types::int_type> row_indices;
526 for (size_type row=first_row; row<last_row; ++row)
528 const int row_length = sparsity_pattern.row_length(row);
532 row_indices.resize (row_length, -1);
534 typename SparsityPatternType::iterator p = sparsity_pattern.begin(row);
535 for (size_type col=0; p != sparsity_pattern.end(row); ++p, ++col)
536 row_indices[col] = p->column();
538 graph->Epetra_CrsGraph::InsertGlobalIndices (row, row_length,
546 graph->FillComplete(input_col_map, input_row_map);
547 graph->OptimizeStorage();
555 matrix = std_cxx14::make_unique<Epetra_FECrsMatrix>(Copy, *graph,
false);
567 class Epetra_CrsGraphMod :
public Epetra_CrsGraph
570 Epetra_CrsGraphMod (
const Epetra_Map &row_map,
571 const int *n_entries_per_row)
573 Epetra_CrsGraph(Copy, row_map, n_entries_per_row, true)
576 void SetIndicesAreGlobal()
578 this->Epetra_CrsGraph::SetIndicesAreGlobal(
true);
587 reinit_matrix (
const Epetra_Map &input_row_map,
588 const Epetra_Map &input_col_map,
590 const bool exchange_data,
591 std::unique_ptr<Epetra_Map> &column_space_map,
592 std::unique_ptr<Epetra_FECrsMatrix> &matrix,
593 std::unique_ptr<Epetra_CrsMatrix> &nonlocal_matrix,
594 std::unique_ptr<Epetra_Export> &nonlocal_matrix_exporter)
597 nonlocal_matrix.reset();
598 nonlocal_matrix_exporter.reset();
605 column_space_map = std_cxx14::make_unique<Epetra_Map>(input_col_map);
609 if (relevant_rows.size() == 0)
614 relevant_rows.compress();
615 Assert(relevant_rows.n_elements() >=
static_cast<unsigned int>(input_row_map.NumMyElements()),
616 ExcMessage(
"Locally relevant rows of sparsity pattern must contain " 617 "all locally owned rows"));
622 bool have_ghost_rows =
false;
624 std::vector<::types::global_dof_index> indices;
625 relevant_rows.fill_index_vector(indices);
626 Epetra_Map relevant_map (TrilinosWrappers::types::int_type(-1),
627 TrilinosWrappers::types::int_type(relevant_rows.n_elements()),
628 (indices.empty() ? nullptr :
629 reinterpret_cast<TrilinosWrappers::types::int_type *
>(indices.data())),
630 0, input_row_map.Comm());
631 if (relevant_map.SameAs(input_row_map))
632 have_ghost_rows =
false;
634 have_ghost_rows =
true;
637 const unsigned int n_rows = relevant_rows.n_elements();
638 std::vector<TrilinosWrappers::types::int_type> ghost_rows;
639 std::vector<int> n_entries_per_row(input_row_map.NumMyElements());
640 std::vector<int> n_entries_per_ghost_row;
641 for (
unsigned int i=0, own=0; i<n_rows; ++i)
643 const TrilinosWrappers::types::int_type global_row =
644 relevant_rows.nth_index_in_set(i);
645 if (input_row_map.MyGID(global_row))
646 n_entries_per_row[own++] = sparsity_pattern.
row_length(global_row);
647 else if (sparsity_pattern.
row_length(global_row) > 0)
649 ghost_rows.push_back(global_row);
650 n_entries_per_ghost_row.push_back(sparsity_pattern.
row_length(global_row));
654 Epetra_Map off_processor_map(-1, ghost_rows.size(),
655 (ghost_rows.size()>0)?(ghost_rows.data()):
nullptr,
656 0, input_row_map.Comm());
658 std::unique_ptr<Epetra_CrsGraph> graph;
659 std::unique_ptr<Epetra_CrsGraphMod> nonlocal_graph;
660 if (input_row_map.Comm().NumProc() > 1)
662 graph = std_cxx14::make_unique<Epetra_CrsGraph>(Copy, input_row_map,
663 (n_entries_per_row.size()>0)?(n_entries_per_row.data()):
nullptr,
664 exchange_data ?
false :
true);
665 if (have_ghost_rows ==
true)
666 nonlocal_graph = std_cxx14::make_unique<Epetra_CrsGraphMod>(off_processor_map,
667 n_entries_per_ghost_row.data());
670 graph = std_cxx14::make_unique<Epetra_CrsGraph>(Copy, input_row_map, input_col_map,
671 (n_entries_per_row.size()>0)?(n_entries_per_row.data()):
nullptr,
675 std::vector<TrilinosWrappers::types::int_type> row_indices;
677 for (
unsigned int i=0; i<n_rows; ++i)
679 const TrilinosWrappers::types::int_type global_row =
680 relevant_rows.nth_index_in_set(i);
681 const int row_length = sparsity_pattern.
row_length(global_row);
685 row_indices.resize (row_length, -1);
686 for (
int col=0; col < row_length; ++col)
687 row_indices[col] = sparsity_pattern.
column_number(global_row, col);
689 if (input_row_map.MyGID(global_row))
690 graph->InsertGlobalIndices (global_row, row_length, row_indices.data());
694 nonlocal_graph->InsertGlobalIndices (global_row, row_length,
700 if (nonlocal_graph.get() !=
nullptr)
706 nonlocal_graph->SetIndicesAreGlobal();
707 Assert(nonlocal_graph->IndicesAreGlobal() ==
true,
709 nonlocal_graph->FillComplete(input_col_map, input_row_map);
710 nonlocal_graph->OptimizeStorage();
715 Epetra_Export exporter(nonlocal_graph->RowMap(), input_row_map);
716 int ierr = graph->Export(*nonlocal_graph, exporter, Add);
721 nonlocal_matrix = std_cxx14::make_unique<Epetra_CrsMatrix>(Copy, *nonlocal_graph);
724 graph->FillComplete(input_col_map, input_row_map);
725 graph->OptimizeStorage();
730 matrix = std_cxx14::make_unique<Epetra_FECrsMatrix>(Copy, *graph,
false);
736 template <
typename SparsityPatternType>
740 const Epetra_Map rows (static_cast<TrilinosWrappers::types::int_type>(sparsity_pattern.n_rows()),
743 const Epetra_Map columns (static_cast<TrilinosWrappers::types::int_type>(sparsity_pattern.n_cols()),
747 reinit_matrix (rows, columns, sparsity_pattern,
false,
754 template <
typename SparsityPatternType>
757 const SparsityPatternType &sparsity_pattern,
758 const bool exchange_data)
760 reinit_matrix (input_map, input_map, sparsity_pattern, exchange_data,
770 template <
typename SparsityPatternType>
773 const IndexSet &col_parallel_partitioning,
774 const SparsityPatternType &sparsity_pattern,
775 const MPI_Comm &communicator,
776 const bool exchange_data)
782 reinit_matrix (row_map, col_map, sparsity_pattern, exchange_data,
794 template <
typename SparsityPatternType>
797 const Epetra_Map &col_map,
798 const SparsityPatternType &sparsity_pattern,
799 const bool exchange_data)
801 reinit_matrix (row_map, col_map, sparsity_pattern, exchange_data,
823 matrix = std_cxx14::make_unique<Epetra_FECrsMatrix>
840 if (
this == &sparse_matrix)
846 matrix = std_cxx14::make_unique<Epetra_FECrsMatrix>
861 template <
typename number>
864 const IndexSet &col_parallel_partitioning,
865 const ::SparseMatrix<number> &dealii_sparse_matrix,
866 const MPI_Comm &communicator,
867 const double drop_tolerance,
868 const bool copy_values,
869 const ::SparsityPattern *use_this_sparsity)
871 if (copy_values ==
false)
875 if (use_this_sparsity ==
nullptr)
876 reinit (row_parallel_partitioning, col_parallel_partitioning,
877 dealii_sparse_matrix.get_sparsity_pattern(),
878 communicator,
false);
880 reinit (row_parallel_partitioning, col_parallel_partitioning,
881 *use_this_sparsity, communicator,
false);
885 const size_type n_rows = dealii_sparse_matrix.m();
890 const ::SparsityPattern &sparsity_pattern =
891 (use_this_sparsity!=
nullptr)? *use_this_sparsity :
892 dealii_sparse_matrix.get_sparsity_pattern();
894 if (
matrix.get() ==
nullptr ||
898 reinit (row_parallel_partitioning, col_parallel_partitioning,
899 sparsity_pattern, communicator,
false);
908 std::vector<size_type> row_indices (maximum_row_length);
909 std::vector<TrilinosScalar> values (maximum_row_length);
913 if (row_parallel_partitioning.
is_element(row) ==
true)
916 sparsity_pattern.begin(row);
917 typename ::SparseMatrix<number>::const_iterator it =
918 dealii_sparse_matrix.begin(row);
920 if (sparsity_pattern.n_rows() == sparsity_pattern.n_cols())
924 if (std::fabs(it->value()) > drop_tolerance)
926 values[col] = it->value();
927 row_indices[col++] = it->column();
933 while (it != dealii_sparse_matrix.end(row) &&
934 select_index != sparsity_pattern.end(row))
936 while (select_index->
column() < it->column() &&
937 select_index != sparsity_pattern.end(row))
939 while (it->column() < select_index->
column() &&
940 it != dealii_sparse_matrix.end(row))
943 if (it == dealii_sparse_matrix.end(row))
945 if (std::fabs(it->value()) > drop_tolerance)
947 values[col] = it->value();
948 row_indices[col++] = it->column();
953 set (row, col,
reinterpret_cast<size_type *
>(row_indices.data()),
954 values.data(),
false);
961 template <
typename number>
964 const double drop_tolerance,
965 const bool copy_values,
966 const ::SparsityPattern *use_this_sparsity)
968 reinit (complete_index_set(dealii_sparse_matrix.m()),
969 complete_index_set(dealii_sparse_matrix.n()),
970 dealii_sparse_matrix, MPI_COMM_SELF, drop_tolerance,
971 copy_values, use_this_sparsity);
976 template <
typename number>
979 const ::SparseMatrix<number> &dealii_sparse_matrix,
980 const double drop_tolerance,
981 const bool copy_values,
982 const ::SparsityPattern *use_this_sparsity)
985 MPI_COMM_SELF, drop_tolerance, copy_values, use_this_sparsity);
990 template <
typename number>
993 const Epetra_Map &input_col_map,
994 const ::SparseMatrix<number> &dealii_sparse_matrix,
995 const double drop_tolerance,
996 const bool copy_values,
997 const ::SparsityPattern *use_this_sparsity)
1000 dealii_sparse_matrix, MPI_COMM_SELF,
1001 drop_tolerance, copy_values, use_this_sparsity);
1008 const bool copy_values)
1010 Assert (input_matrix.Filled()==
true,
1011 ExcMessage(
"Input CrsMatrix has not called FillComplete()!"));
1013 column_space_map = std_cxx14::make_unique<Epetra_Map>(input_matrix.DomainMap());
1015 const Epetra_CrsGraph *graph = &input_matrix.Graph();
1020 matrix = std_cxx14::make_unique<Epetra_FECrsMatrix>(Copy, *graph,
false);
1024 if (copy_values ==
true)
1028 const TrilinosScalar *in_values = input_matrix[0];
1029 TrilinosScalar *values = (*matrix)[0];
1030 const size_type my_nonzeros = input_matrix.NumMyNonzeros();
1031 std::memcpy (values, in_values,
1032 my_nonzeros*
sizeof (TrilinosScalar));
1060 ExcMessage(
"Operation and argument to compress() do not match"));
1084 ierr =
matrix->OptimizeStorage ();
1115 const TrilinosScalar new_diag_value)
1122 matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(row));
1125 TrilinosScalar *values;
1128 const int ierr =
matrix->ExtractMyRowView(local_row, num_entries,
1129 values, col_indices);
1135 int *diag_find = std::find(col_indices,col_indices+num_entries,local_row);
1136 int diag_index = (int)(diag_find - col_indices);
1138 for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
1139 if (diag_index != j || new_diag_value == 0)
1142 if (diag_find && std::fabs(values[diag_index]) == 0.0 &&
1143 new_diag_value != 0.0)
1144 values[diag_index] = new_diag_value;
1152 const TrilinosScalar new_diag_value)
1154 for (
size_type row=0; row<rows.size(); ++row)
1166 int trilinos_i =
matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(i)),
1167 trilinos_j =
matrix->LCID(static_cast<TrilinosWrappers::types::int_type>(j));
1168 TrilinosScalar value = 0.;
1176 if (trilinos_i == -1)
1190 int nnz_present =
matrix->NumMyEntries(trilinos_i);
1193 TrilinosScalar *values;
1199 int ierr =
matrix->ExtractMyRowView(trilinos_i, nnz_extracted,
1200 values, col_indices);
1204 Assert (nnz_present == nnz_extracted,
1211 int *el_find = std::find(col_indices, col_indices + nnz_present, trilinos_j);
1213 int local_col_index = (int)(el_find - col_indices);
1223 if (local_col_index == nnz_present)
1228 value = values[local_col_index];
1242 int trilinos_i =
matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(i)),
1243 trilinos_j =
matrix->LCID(static_cast<TrilinosWrappers::types::int_type>(j));
1244 TrilinosScalar value = 0.;
1253 if ((trilinos_i == -1 ) || (trilinos_j == -1))
1264 int nnz_present =
matrix->NumMyEntries(trilinos_i);
1267 TrilinosScalar *values;
1272 int ierr =
matrix->ExtractMyRowView(trilinos_i, nnz_extracted,
1273 values, col_indices);
1277 Assert (nnz_present == nnz_extracted,
1283 int *el_find = std::find(col_indices, col_indices + nnz_present, trilinos_j);
1285 int local_col_index = (int)(el_find - col_indices);
1295 if (local_col_index == nnz_present)
1298 value = values[local_col_index];
1335 int local_row =
matrix->LRID(static_cast<TrilinosWrappers::types::int_type>(row));
1342 int ierr =
matrix->NumMyRowEntries (local_row, ncols);
1353 const std::vector<size_type> &col_indices,
1355 const bool elide_zero_values)
1357 Assert (row_indices.size() == values.
m(),
1359 Assert (col_indices.size() == values.
n(),
1362 for (
size_type i=0; i<row_indices.size(); ++i)
1363 set (row_indices[i], col_indices.size(), col_indices.data(), &values(i,0),
1371 const std::vector<size_type> &col_indices,
1372 const std::vector<TrilinosScalar> &values,
1373 const bool elide_zero_values)
1375 Assert (col_indices.size() == values.size(),
1378 set (row, col_indices.size(), col_indices.data(), values.data(),
1388 const TrilinosScalar *values,
1389 const bool elide_zero_values)
1404 TrilinosWrappers::types::int_type *col_index_ptr;
1405 TrilinosScalar *col_value_ptr;
1406 TrilinosWrappers::types::int_type n_columns;
1408 boost::container::small_vector<TrilinosScalar, 200> local_value_array(n_cols);
1409 boost::container::small_vector<TrilinosWrappers::types::int_type, 200> local_index_array(n_cols);
1414 if (elide_zero_values ==
false)
1416 col_index_ptr = (TrilinosWrappers::types::int_type *)col_indices;
1417 col_value_ptr =
const_cast<TrilinosScalar *
>(values);
1424 col_index_ptr = local_index_array.data();
1425 col_value_ptr = local_value_array.data();
1430 const double value = values[j];
1434 col_index_ptr[n_columns] = col_indices[j];
1435 col_value_ptr[n_columns] = value;
1452 if (
matrix->RowMap().MyGID(static_cast<TrilinosWrappers::types::int_type>(row)) ==
true)
1454 if (
matrix->Filled() ==
false)
1456 ierr =
matrix->Epetra_CrsMatrix::InsertGlobalValues(
1457 static_cast<TrilinosWrappers::types::int_type>(row),
1458 static_cast<int>(n_columns),const_cast<double *>(col_value_ptr),
1468 ierr =
matrix->Epetra_CrsMatrix::ReplaceGlobalValues(row, n_columns,
1482 if (
matrix->Filled() ==
false)
1484 ierr =
matrix->InsertGlobalValues (1,
1485 (TrilinosWrappers::types::int_type *)&row,
1486 n_columns, col_index_ptr,
1488 Epetra_FECrsMatrix::ROW_MAJOR);
1493 ierr =
matrix->ReplaceGlobalValues (1,
1494 (TrilinosWrappers::types::int_type *)&row,
1495 n_columns, col_index_ptr,
1497 Epetra_FECrsMatrix::ROW_MAJOR);
1514 const bool elide_zero_values)
1516 Assert (indices.size() == values.
m(),
1520 for (
size_type i=0; i<indices.size(); ++i)
1521 add (indices[i], indices.size(), indices.data(), &values(i,0),
1529 const std::vector<size_type> &col_indices,
1531 const bool elide_zero_values)
1533 Assert (row_indices.size() == values.
m(),
1535 Assert (col_indices.size() == values.
n(),
1538 for (
size_type i=0; i<row_indices.size(); ++i)
1539 add (row_indices[i], col_indices.size(), col_indices.data(),
1540 &values(i,0), elide_zero_values);
1547 const std::vector<size_type> &col_indices,
1548 const std::vector<TrilinosScalar> &values,
1549 const bool elide_zero_values)
1551 Assert (col_indices.size() == values.size(),
1554 add (row, col_indices.size(), col_indices.data(), values.data(),
1564 const TrilinosScalar *values,
1565 const bool elide_zero_values,
1575 matrix->RowMap(),
false);
1582 TrilinosWrappers::types::int_type *col_index_ptr;
1583 TrilinosScalar *col_value_ptr;
1584 TrilinosWrappers::types::int_type n_columns;
1586 boost::container::small_vector<TrilinosScalar, 100> local_value_array(n_cols);
1587 boost::container::small_vector<TrilinosWrappers::types::int_type, 100> local_index_array(n_cols);
1592 if (elide_zero_values ==
false)
1594 col_index_ptr = (TrilinosWrappers::types::int_type *)col_indices;
1595 col_value_ptr =
const_cast<TrilinosScalar *
>(values);
1606 col_index_ptr = local_index_array.data();
1607 col_value_ptr = local_value_array.data();
1612 const double value = values[j];
1617 col_index_ptr[n_columns] = col_indices[j];
1618 col_value_ptr[n_columns] = value;
1635 if (
matrix->RowMap().MyGID(static_cast<TrilinosWrappers::types::int_type>(row)) ==
true)
1637 ierr =
matrix->Epetra_CrsMatrix::SumIntoGlobalValues(row, n_columns,
1649 ExcMessage(
"Attempted to write into off-processor matrix row " 1650 "that has not be specified as being writable upon " 1667 ierr =
matrix->SumIntoGlobalValues (1,
1668 (TrilinosWrappers::types::int_type *)&row, n_columns,
1671 Epetra_FECrsMatrix::ROW_MAJOR);
1678 std::cout <<
"------------------------------------------" 1680 std::cout <<
"Got error " << ierr <<
" in row " << row
1681 <<
" of proc " <<
matrix->RowMap().Comm().MyPID()
1682 <<
" when trying to add the columns:" << std::endl;
1683 for (TrilinosWrappers::types::int_type i=0; i<n_columns; ++i)
1684 std::cout << col_index_ptr[i] <<
" ";
1685 std::cout << std::endl << std::endl;
1686 std::cout <<
"Matrix row " 1687 << (
matrix->RowMap().MyGID(static_cast<TrilinosWrappers::types::int_type>(row)) ==
false ?
"(nonlocal part)" :
"")
1688 <<
" has the following indices:" << std::endl;
1689 std::vector<TrilinosWrappers::types::int_type> indices;
1690 const Epetra_CrsGraph *graph =
1692 matrix->RowMap().MyGID(static_cast<TrilinosWrappers::types::int_type>(row)) ==
false) ?
1695 indices.resize(graph->NumGlobalIndices(static_cast<TrilinosWrappers::types::int_type>(row)));
1697 graph->ExtractGlobalRowCopy(static_cast<TrilinosWrappers::types::int_type>(row),
1698 indices.size(), n_indices, indices.data());
1699 AssertDimension(static_cast<unsigned int>(n_indices), indices.size());
1701 for (TrilinosWrappers::types::int_type i=0; i<n_indices; ++i)
1702 std::cout << indices[i] <<
" ";
1703 std::cout << std::endl << std::endl;
1719 const int ierr =
matrix->PutScalar(d);
1738 ExcMessage(
"Can only add matrices with same distribution of rows"));
1740 ExcMessage(
"Addition of matrices only allowed if matrices are " 1741 "filled, i.e., compress() has been called"));
1743 const std::pair<size_type, size_type>
1745 const bool same_col_map =
matrix->ColMap().SameAs(rhs.
matrix->ColMap());
1749 const int row_local =
1750 matrix->RowMap().LID(static_cast<TrilinosWrappers::types::int_type>(row));
1756 int n_entries, rhs_n_entries;
1757 TrilinosScalar *value_ptr, *rhs_value_ptr;
1758 int *index_ptr, *rhs_index_ptr;
1759 int ierr = rhs.
matrix->ExtractMyRowView (row_local, rhs_n_entries,
1760 rhs_value_ptr, rhs_index_ptr);
1764 ierr =
matrix->ExtractMyRowView (row_local, n_entries, value_ptr,
1767 bool expensive_checks = (n_entries != rhs_n_entries || !same_col_map);
1768 if (!expensive_checks)
1772 expensive_checks = std::memcmp(static_cast<void *>(index_ptr),
1773 static_cast<void *>(rhs_index_ptr),
1774 sizeof(
int)*n_entries) != 0;
1775 if (!expensive_checks)
1776 for (
int i=0; i<n_entries; ++i)
1777 value_ptr[i] += rhs_value_ptr[i] * factor;
1783 if (expensive_checks)
1785 for (
int i=0; i<rhs_n_entries; ++i)
1787 if (rhs_value_ptr[i] == 0.)
1789 const TrilinosWrappers::types::int_type rhs_global_col =
1791 int local_col =
matrix->ColMap().LID(rhs_global_col);
1793 index_ptr+n_entries,
1795 Assert(local_index != index_ptr + n_entries &&
1796 *local_index == local_col,
1797 ExcMessage(
"Adding the entries from the other matrix " 1798 "failed, because the sparsity pattern " 1799 "of that matrix includes more elements than the " 1800 "calling matrix, which is not allowed."));
1801 value_ptr[local_index-index_ptr] += factor * rhs_value_ptr[i];
1819 if (!
matrix->UseTranspose())
1821 ierr =
matrix->SetUseTranspose (
true);
1826 ierr =
matrix->SetUseTranspose (
false);
1836 const int ierr =
matrix->Scale (a);
1850 const TrilinosScalar factor = 1./a;
1852 const int ierr =
matrix->Scale (factor);
1865 return matrix->NormOne();
1874 return matrix->NormInf();
1883 return matrix->NormFrobenius();
1890 namespace SparseMatrixImplementation
1892 template <
typename VectorType>
1894 void check_vector_map_equality(
const Epetra_CrsMatrix &,
1901 void check_vector_map_equality(
const Epetra_CrsMatrix &m,
1906 ExcMessage (
"Column map of matrix does not fit with vector map!"));
1908 ExcMessage (
"Row map of matrix does not fit with vector map!"));
1917 template <
typename VectorType>
1920 const VectorType &src)
const 1927 internal::SparseMatrixImplementation::check_vector_map_equality(*
matrix, src, dst);
1928 const size_type dst_local_size = internal::end(dst) - internal::begin(dst);
1930 const size_type src_local_size = internal::end(src) - internal::begin(src);
1933 Epetra_MultiVector tril_dst (View,
matrix->RangeMap(), internal::begin(dst),
1935 Epetra_MultiVector tril_src (View,
matrix->DomainMap(),
1936 const_cast<TrilinosScalar *
>(internal::begin(src)),
1939 const int ierr =
matrix->Multiply (
false, tril_src, tril_dst);
1946 template <
typename VectorType>
1949 const VectorType &src)
const 1954 internal::SparseMatrixImplementation::check_vector_map_equality(*
matrix, dst, src);
1955 const size_type dst_local_size = internal::end(dst) - internal::begin(dst);
1957 const size_type src_local_size = internal::end(src) - internal::begin(src);
1960 Epetra_MultiVector tril_dst (View,
matrix->DomainMap(), internal::begin(dst),
1962 Epetra_MultiVector tril_src (View,
matrix->RangeMap(),
1963 const_cast<double *
>(internal::begin(src)),
1966 const int ierr =
matrix->Multiply (
true, tril_src, tril_dst);
1973 template <
typename VectorType>
1976 const VectorType &src)
const 1982 VectorType tmp_vector;
1983 tmp_vector.reinit(dst,
true);
1984 vmult (tmp_vector, src);
1990 template <
typename VectorType>
1993 const VectorType &src)
const 1999 VectorType tmp_vector;
2000 tmp_vector.reinit(dst,
true);
2001 Tvmult (tmp_vector, src);
2014 temp_vector.
reinit(v,
true);
2016 vmult (temp_vector, v);
2017 return temp_vector*v;
2030 temp_vector.
reinit(v,
true);
2032 vmult (temp_vector, v);
2033 return u*temp_vector;
2054 typedef ::types::global_dof_index size_type;
2060 const bool transpose_left)
2062 #ifdef DEAL_II_WITH_64BIT_INDICES 2065 const bool use_vector = (V.
size() == inputright.
m() ? true :
false);
2066 if (transpose_left ==
false)
2068 Assert (inputleft.
n() == inputright.
m(),
2071 ExcMessage (
"Parallel partitioning of A and B does not fit."));
2075 Assert (inputleft.
m() == inputright.
m(),
2078 ExcMessage (
"Parallel partitioning of A and B does not fit."));
2089 Teuchos::RCP<Epetra_CrsMatrix> mod_B;
2090 if (use_vector ==
false)
2092 mod_B = Teuchos::rcp(const_cast<Epetra_CrsMatrix *>
2098 mod_B = Teuchos::rcp(
new Epetra_CrsMatrix
2104 ExcMessage (
"Parallel distribution of matrix B and vector V " 2105 "does not match."));
2108 for (
int i=0; i<local_N; ++i)
2111 double *new_data, *B_data;
2112 mod_B->ExtractMyRowView (i, N_entries, new_data);
2115 for (TrilinosWrappers::types::int_type j=0; j<N_entries; ++j)
2116 new_data[j] = value * B_data[j];
2126 ML_Comm_Create(&comm);
2128 const Epetra_MpiComm *epcomm =
dynamic_cast<const Epetra_MpiComm *
>(&(inputleft.
trilinos_matrix().Comm()));
2130 if (epcomm) ML_Comm_Set_UsrComm(comm,epcomm->Comm());
2132 ML_Operator *A_ = ML_Operator_Create(comm);
2133 ML_Operator *B_ = ML_Operator_Create(comm);
2134 ML_Operator *C_ = ML_Operator_Create(comm);
2137 if (transpose_left ==
false)
2138 ML_Operator_WrapEpetraCrsMatrix
2147 ExcMessage(
"Matrix must be partitioned contiguously between procs."));
2150 int num_entries, * indices;
2157 for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
2158 sparsity_transposed.add
2164 sparsity_transposed.compress();
2165 transposed_mat.
reinit (sparsity_transposed);
2168 int num_entries, * indices;
2176 for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
2182 ML_Operator_WrapEpetraCrsMatrix
2183 (const_cast<Epetra_CrsMatrix *>(&transposed_mat.trilinos_matrix()),
2186 ML_Operator_WrapEpetraCrsMatrix(mod_B.get(),B_,
false);
2196 ML_Operator *Btmp, *Ctmp, *Ctmp2, *tptr;
2197 ML_CommInfoOP *getrow_comm;
2199 TrilinosWrappers::types::int_type N_input_vector = B_->invec_leng;
2200 getrow_comm = B_->getrow->pre_comm;
2201 if ( getrow_comm !=
nullptr)
2202 for (TrilinosWrappers::types::int_type i = 0; i < getrow_comm->N_neighbors; i++)
2203 for (TrilinosWrappers::types::int_type j = 0; j < getrow_comm->neighbors[i].N_send; j++)
2204 AssertThrow (getrow_comm->neighbors[i].send_list[j] < N_input_vector,
2207 ML_create_unique_col_id(N_input_vector, &(B_->getrow->loc_glob_map),
2208 getrow_comm, &max_per_proc, B_->comm);
2209 B_->getrow->use_loc_glob_map = ML_YES;
2210 if (A_->getrow->pre_comm !=
nullptr)
2211 ML_exchange_rows( B_, &Btmp, A_->getrow->pre_comm);
2215 ML_matmat_mult(A_, Btmp, &Ctmp);
2219 ML_free(B_->getrow->loc_glob_map);
2220 B_->getrow->loc_glob_map =
nullptr;
2221 B_->getrow->use_loc_glob_map = ML_NO;
2222 if (A_->getrow->pre_comm !=
nullptr)
2225 while ( (tptr!=
nullptr) && (tptr->sub_matrix != B_))
2226 tptr = tptr->sub_matrix;
2227 if (tptr !=
nullptr) tptr->sub_matrix =
nullptr;
2228 ML_RECUR_CSR_MSRdata_Destroy(Btmp);
2229 ML_Operator_Destroy(&Btmp);
2233 if (A_->getrow->post_comm !=
nullptr)
2234 ML_exchange_rows(Ctmp, &Ctmp2, A_->getrow->post_comm);
2238 ML_back_to_csrlocal(Ctmp2, C_, max_per_proc);
2240 ML_RECUR_CSR_MSRdata_Destroy (Ctmp);
2241 ML_Operator_Destroy (&Ctmp);
2243 if (A_->getrow->post_comm !=
nullptr)
2245 ML_RECUR_CSR_MSRdata_Destroy(Ctmp2);
2246 ML_Operator_Destroy (&Ctmp2);
2251 Epetra_CrsMatrix *C_mat;
2252 ML_Operator2EpetraCrsMatrix(C_, C_mat);
2253 C_mat->FillComplete(mod_B->DomainMap(),
2257 C_mat->OptimizeStorage();
2262 ML_Operator_Destroy (&A_);
2263 ML_Operator_Destroy (&B_);
2264 ML_Operator_Destroy (&C_);
2265 ML_Comm_Destroy (&comm);
2275 #ifdef DEAL_II_WITH_64BIT_INDICES 2278 internals::perform_mmult (*
this, B, C, V,
false);
2288 #ifdef DEAL_II_WITH_64BIT_INDICES 2291 internals::perform_mmult (*
this, B, C, V,
true);
2309 const bool print_detailed_trilinos_information)
const 2311 if (print_detailed_trilinos_information ==
true)
2319 for (
int i=0; i<
matrix->NumMyRows(); ++i)
2321 matrix->ExtractMyRowView (i, num_entries, values, indices);
2322 for (TrilinosWrappers::types::int_type j=0; j<num_entries; ++j)
2325 << values[j] << std::endl;
2337 size_type static_memory =
sizeof(*this) +
sizeof (*matrix)
2338 +
sizeof(*
matrix->Graph().DataPtr());
2339 return ((
sizeof(TrilinosScalar)+
sizeof(TrilinosWrappers::types::int_type))*
2348 return matrix->DomainMap();
2356 return matrix->RangeMap();
2380 #ifdef DEAL_II_WITH_MPI 2382 const Epetra_MpiComm *mpi_comm
2383 =
dynamic_cast<const Epetra_MpiComm *
>(&
matrix->RangeMap().Comm());
2385 return mpi_comm->Comm();
2388 return MPI_COMM_SELF;
2402 #ifndef DEAL_II_WITH_MPI 2404 make_serial_Epetra_map (
const IndexSet &serial_partitioning)
2407 return Epetra_Map(TrilinosWrappers::types::int_type(serial_partitioning.
size()),
2408 TrilinosWrappers::types::int_type(serial_partitioning.
n_elements()),
2409 0, Epetra_SerialComm());
2414 namespace LinearOperatorImplementation
2418 : use_transpose (false),
2419 #ifdef DEAL_II_WITH_MPI
2420 communicator (MPI_COMM_SELF),
2421 domain_map (
IndexSet().make_trilinos_map(communicator.Comm())),
2422 range_map (
IndexSet().make_trilinos_map(communicator.Comm()))
2431 ExcMessage(
"Uninitialized TrilinosPayload::vmult called" 2432 "(Default constructor)"));
2438 ExcMessage(
"Uninitialized TrilinosPayload::Tvmult called" 2439 "(Default constructor)"));
2445 ExcMessage(
"Uninitialized TrilinosPayload::inv_vmult called" 2446 "(Default constructor)"));
2452 ExcMessage(
"Uninitialized TrilinosPayload::inv_Tvmult called" 2453 "(Default constructor)"));
2461 : use_transpose (matrix_exemplar.trilinos_matrix().UseTranspose()),
2462 #ifdef DEAL_II_WITH_MPI
2463 communicator (matrix_exemplar.get_mpi_communicator()),
2464 domain_map (matrix_exemplar.locally_owned_domain_indices().make_trilinos_map(communicator.Comm())),
2465 range_map (matrix_exemplar.locally_owned_range_indices().make_trilinos_map(communicator.Comm()))
2467 domain_map (
internal::make_serial_Epetra_map(matrix_exemplar.locally_owned_domain_indices())),
2468 range_map (
internal::make_serial_Epetra_map(matrix_exemplar.locally_owned_range_indices()))
2471 vmult = [&matrix_exemplar,&matrix](
Range &tril_dst,
const Domain &tril_src)
2476 internal::check_vector_map_equality(matrix_exemplar.
trilinos_matrix(),
2479 internal::check_vector_map_equality(matrix.trilinos_matrix(),
2481 matrix.trilinos_matrix().UseTranspose());
2483 const int ierr = matrix.trilinos_matrix().Apply (tril_src, tril_dst);
2492 internal::check_vector_map_equality(matrix_exemplar.
trilinos_matrix(),
2495 internal::check_vector_map_equality(matrix.trilinos_matrix(),
2497 !matrix.trilinos_matrix().UseTranspose());
2499 Epetra_CrsMatrix &tril_mtrx_non_const =
const_cast<Epetra_CrsMatrix &
>(matrix.trilinos_matrix());
2500 tril_mtrx_non_const.SetUseTranspose(!matrix.trilinos_matrix().UseTranspose());
2501 const int ierr = matrix.trilinos_matrix().Apply (tril_src, tril_dst);
2503 tril_mtrx_non_const.SetUseTranspose(!matrix.trilinos_matrix().UseTranspose());
2509 ExcMessage(
"Uninitialized TrilinosPayload::inv_vmult called" 2510 "(Matrix constructor with matrix exemplar)"));
2516 ExcMessage(
"Uninitialized TrilinosPayload::inv_Tvmult called" 2517 "(Matrix constructor with matrix exemplar)"));
2525 : use_transpose (matrix_exemplar.trilinos_matrix().UseTranspose()),
2526 #ifdef DEAL_II_WITH_MPI
2527 communicator (matrix_exemplar.get_mpi_communicator()),
2528 domain_map (matrix_exemplar.locally_owned_domain_indices().make_trilinos_map(communicator.Comm())),
2529 range_map (matrix_exemplar.locally_owned_range_indices().make_trilinos_map(communicator.Comm()))
2531 domain_map (
internal::make_serial_Epetra_map(matrix_exemplar.locally_owned_domain_indices())),
2532 range_map (
internal::make_serial_Epetra_map(matrix_exemplar.locally_owned_range_indices()))
2535 vmult = [&matrix_exemplar,&preconditioner](
Range &tril_dst,
const Domain &tril_src)
2540 internal::check_vector_map_equality(matrix_exemplar.
trilinos_matrix(),
2547 const int ierr = preconditioner.
trilinos_operator().ApplyInverse (tril_src,tril_dst);
2551 Tvmult = [&matrix_exemplar,&preconditioner](
Domain &tril_dst,
const Range &tril_src)
2556 internal::check_vector_map_equality(matrix_exemplar.
trilinos_matrix(),
2564 const int ierr = preconditioner.
trilinos_operator().ApplyInverse (tril_src,tril_dst);
2574 internal::check_vector_map_equality(matrix_exemplar.
trilinos_matrix(),
2581 const int ierr = preconditioner.
trilinos_operator().ApplyInverse (tril_src,tril_dst);
2590 internal::check_vector_map_equality(matrix_exemplar.
trilinos_matrix(),
2598 const int ierr = preconditioner.
trilinos_operator().ApplyInverse (tril_src,tril_dst);
2608 : use_transpose (preconditioner_exemplar.trilinos_operator().UseTranspose()),
2609 #ifdef DEAL_II_WITH_MPI
2610 communicator (preconditioner_exemplar.get_mpi_communicator()),
2611 domain_map (preconditioner_exemplar.locally_owned_domain_indices().make_trilinos_map(communicator.Comm())),
2612 range_map (preconditioner_exemplar.locally_owned_range_indices().make_trilinos_map(communicator.Comm()))
2614 domain_map (
internal::make_serial_Epetra_map(preconditioner_exemplar.locally_owned_domain_indices())),
2615 range_map (
internal::make_serial_Epetra_map(preconditioner_exemplar.locally_owned_range_indices()))
2618 vmult = [&preconditioner_exemplar,&preconditioner](
Range &tril_dst,
const Domain &tril_src)
2623 internal::check_vector_map_equality(preconditioner_exemplar.
trilinos_operator(),
2634 Tvmult = [&preconditioner_exemplar,&preconditioner](
Domain &tril_dst,
const Range &tril_src)
2639 internal::check_vector_map_equality(preconditioner_exemplar.
trilinos_operator(),
2652 inv_vmult = [&preconditioner_exemplar,&preconditioner](
Domain &tril_dst,
const Range &tril_src)
2657 internal::check_vector_map_equality(preconditioner_exemplar.
trilinos_operator(),
2664 const int ierr = preconditioner.
trilinos_operator().ApplyInverse (tril_src,tril_dst);
2673 internal::check_vector_map_equality(preconditioner_exemplar.
trilinos_operator(),
2681 const int ierr = preconditioner.
trilinos_operator().ApplyInverse (tril_src,tril_dst);
2690 : vmult (payload.vmult),
2691 Tvmult (payload.Tvmult),
2692 inv_vmult (payload.inv_vmult),
2693 inv_Tvmult (payload.inv_Tvmult),
2694 use_transpose (payload.use_transpose),
2695 communicator (payload.communicator),
2696 domain_map (payload.domain_map),
2697 range_map (payload.range_map)
2706 : use_transpose (false),
2707 communicator (first_op.communicator),
2708 domain_map (second_op.domain_map),
2709 range_map (first_op.range_map)
2720 const Range &tril_src)
2722 tril_dst = tril_src;
2726 const Range &tril_src)
2728 tril_dst = tril_src;
2732 const Range &tril_src)
2734 tril_dst = tril_src;
2738 const Range &tril_src)
2740 tril_dst = tril_src;
2756 const int ierr = tril_dst.PutScalar(0.0);
2764 const int ierr = tril_dst.PutScalar(0.0);
2774 const int ierr = tril_dst.PutScalar(0.0);
2784 const int ierr = tril_dst.PutScalar(0.0);
2823 #ifdef DEAL_II_WITH_MPI 2826 return MPI_COMM_SELF;
2890 return "TrilinosPayload";
2945 ExcMessage(
"Operators are set to work on incompatible IndexSets."));
2947 ExcMessage(
"Operators are set to work on incompatible IndexSets."));
2952 return_op.vmult = [first_op, second_op](Range &tril_dst,
2953 const Domain &tril_src)
2962 i->reinit(
IndexSet(first_op_init_map),
2967 const size_type i_local_size = i->end() - i->begin();
2968 AssertDimension (i_local_size, static_cast<size_type>(first_op_init_map.NumMyPoints()));
2970 AssertDimension (i_local_size, static_cast<size_type>(second_op_init_map.NumMyPoints()));
2971 (void)second_op_init_map;
2972 Intermediate tril_int (View, first_op_init_map,
2973 const_cast<TrilinosScalar *>(i->begin()),
2978 second_op.
Apply(tril_src, tril_int);
2979 first_op.
Apply(tril_src, tril_dst);
2980 const int ierr = tril_dst.Update (1.0, tril_int, 1.0);
2984 return_op.Tvmult = [first_op, second_op](Domain &tril_dst,
2985 const Range &tril_src)
2996 const_cast<TrilinosPayload &
>(first_op).
transpose();
2997 const_cast<TrilinosPayload &
>(second_op).
transpose();
3001 i->reinit(
IndexSet(first_op_init_map),
3006 const size_type i_local_size = i->end() - i->begin();
3007 AssertDimension (i_local_size, static_cast<size_type>(first_op_init_map.NumMyPoints()));
3009 AssertDimension (i_local_size, static_cast<size_type>(second_op_init_map.NumMyPoints()));
3010 (void)second_op_init_map;
3011 Intermediate tril_int (View, first_op_init_map,
3012 const_cast<TrilinosScalar *>(i->begin()),
3017 second_op.
Apply(tril_src, tril_int);
3018 first_op.
Apply(tril_src, tril_dst);
3019 const int ierr = tril_dst.Update (1.0, tril_int, 1.0);
3023 const_cast<TrilinosPayload &
>(first_op).
transpose();
3024 const_cast<TrilinosPayload &
>(second_op).
transpose();
3027 return_op.inv_vmult = [first_op, second_op](Domain &tril_dst,
3028 const Range &tril_src)
3037 i->reinit(
IndexSet(first_op_init_map),
3042 const size_type i_local_size = i->end() - i->begin();
3043 AssertDimension (i_local_size, static_cast<size_type>(first_op_init_map.NumMyPoints()));
3045 AssertDimension (i_local_size, static_cast<size_type>(second_op_init_map.NumMyPoints()));
3046 (void)second_op_init_map;
3047 Intermediate tril_int (View, first_op_init_map,
3048 const_cast<TrilinosScalar *>(i->begin()),
3055 const int ierr = tril_dst.Update (1.0, tril_int, 1.0);
3059 return_op.inv_Tvmult = [first_op, second_op](Range &tril_dst,
3060 const Domain &tril_src)
3071 const_cast<TrilinosPayload &
>(first_op).
transpose();
3072 const_cast<TrilinosPayload &
>(second_op).
transpose();
3076 i->reinit(
IndexSet(first_op_init_map),
3081 const size_type i_local_size = i->end() - i->begin();
3082 AssertDimension (i_local_size, static_cast<size_type>(first_op_init_map.NumMyPoints()));
3084 AssertDimension (i_local_size, static_cast<size_type>(second_op_init_map.NumMyPoints()));
3085 (void)second_op_init_map;
3086 Intermediate tril_int (View, first_op_init_map,
3087 const_cast<TrilinosScalar *>(i->begin()),
3094 const int ierr = tril_dst.Update (1.0, tril_int, 1.0);
3098 const_cast<TrilinosPayload &
>(first_op).
transpose();
3099 const_cast<TrilinosPayload &
>(second_op).
transpose();
3107 TrilinosPayload
operator*(
const TrilinosPayload &first_op,
3108 const TrilinosPayload &second_op)
3115 AssertThrow(first_op.locally_owned_domain_indices() == second_op.locally_owned_range_indices(),
3116 ExcMessage(
"Operators are set to work on incompatible IndexSets."));
3118 TrilinosPayload return_op (first_op, second_op);
3121 return_op.vmult = [first_op, second_op](Range &tril_dst,
3122 const Domain &tril_src)
3130 const Epetra_Map &first_op_init_map = first_op.OperatorDomainMap();
3131 i->reinit(
IndexSet(first_op_init_map),
3132 first_op.get_mpi_communicator(),
3136 const size_type i_local_size = i->end() - i->begin();
3137 AssertDimension (i_local_size, static_cast<size_type>(first_op_init_map.NumMyPoints()));
3138 const Epetra_Map &second_op_init_map = second_op.OperatorRangeMap();
3139 AssertDimension (i_local_size, static_cast<size_type>(second_op_init_map.NumMyPoints()));
3140 (void)second_op_init_map;
3141 Intermediate tril_int (View, first_op_init_map,
3142 const_cast<TrilinosScalar *>(i->begin()),
3147 second_op.Apply(tril_src, tril_int);
3148 first_op.Apply(tril_int, tril_dst);
3151 return_op.Tvmult = [first_op, second_op](Domain &tril_dst,
3152 const Range &tril_src)
3163 const_cast<TrilinosPayload &
>(first_op).
transpose();
3164 const_cast<TrilinosPayload &
>(second_op).
transpose();
3167 const Epetra_Map &first_op_init_map = first_op.OperatorRangeMap();
3168 i->reinit(
IndexSet(first_op_init_map),
3169 first_op.get_mpi_communicator(),
3173 const size_type i_local_size = i->end() - i->begin();
3174 AssertDimension (i_local_size, static_cast<size_type>(first_op_init_map.NumMyPoints()));
3175 const Epetra_Map &second_op_init_map = second_op.OperatorDomainMap();
3176 AssertDimension (i_local_size, static_cast<size_type>(second_op_init_map.NumMyPoints()));
3177 (void)second_op_init_map;
3178 Intermediate tril_int (View, first_op_init_map,
3179 const_cast<TrilinosScalar *>(i->begin()),
3183 first_op.Apply(tril_src, tril_int);
3184 second_op.Apply(tril_int, tril_dst);
3187 const_cast<TrilinosPayload &
>(first_op).
transpose();
3188 const_cast<TrilinosPayload &
>(second_op).
transpose();
3191 return_op.inv_vmult = [first_op, second_op](Domain &tril_dst,
3192 const Range &tril_src)
3200 const Epetra_Map &first_op_init_map = first_op.OperatorRangeMap();
3201 i->reinit(
IndexSet(first_op_init_map),
3202 first_op.get_mpi_communicator(),
3206 const size_type i_local_size = i->end() - i->begin();
3207 AssertDimension (i_local_size, static_cast<size_type>(first_op_init_map.NumMyPoints()));
3208 const Epetra_Map &second_op_init_map = second_op.OperatorDomainMap();
3209 AssertDimension (i_local_size, static_cast<size_type>(second_op_init_map.NumMyPoints()));
3210 (void)second_op_init_map;
3211 Intermediate tril_int (View, first_op_init_map,
3212 const_cast<TrilinosScalar *>(i->begin()),
3217 first_op.ApplyInverse(tril_src, tril_int);
3218 second_op.ApplyInverse(tril_int, tril_dst);
3221 return_op.inv_Tvmult = [first_op, second_op](Range &tril_dst,
3222 const Domain &tril_src)
3233 const_cast<TrilinosPayload &
>(first_op).
transpose();
3234 const_cast<TrilinosPayload &
>(second_op).
transpose();
3237 const Epetra_Map &first_op_init_map = first_op.OperatorDomainMap();
3238 i->reinit(
IndexSet(first_op_init_map),
3239 first_op.get_mpi_communicator(),
3243 const size_type i_local_size = i->end() - i->begin();
3244 AssertDimension (i_local_size, static_cast<size_type>(first_op_init_map.NumMyPoints()));
3245 const Epetra_Map &second_op_init_map = second_op.OperatorRangeMap();
3246 AssertDimension (i_local_size, static_cast<size_type>(second_op_init_map.NumMyPoints()));
3247 (void)second_op_init_map;
3248 Intermediate tril_int (View, first_op_init_map,
3249 const_cast<TrilinosScalar *>(i->begin()),
3256 second_op.ApplyInverse(tril_src, tril_int);
3257 first_op.ApplyInverse(tril_int, tril_dst);
3260 const_cast<TrilinosPayload &
>(first_op).
transpose();
3261 const_cast<TrilinosPayload &
>(second_op).
transpose();
3274 #include "trilinos_sparse_matrix.inst" 3288 const ::SparsityPattern &,
3297 const ::SparsityPattern &,
3307 const ::SparsityPattern &,
3320 const MPI::Vector &)
const;
3323 const ::Vector<double> &)
const;
3326 const ::LinearAlgebra::distributed::Vector<double> &)
const;
3327 #ifdef DEAL_II_WITH_MPI 3330 const ::LinearAlgebra::EpetraWrappers::Vector &)
const;
3334 const MPI::Vector &)
const;
3337 const ::Vector<double> &)
const;
3340 const ::LinearAlgebra::distributed::Vector<double> &)
const;
3341 #ifdef DEAL_II_WITH_MPI 3344 const ::LinearAlgebra::EpetraWrappers::Vector &)
const;
3348 const MPI::Vector &)
const;
3351 const ::Vector<double> &)
const;
3354 const ::LinearAlgebra::distributed::Vector<double> &)
const;
3355 #ifdef DEAL_II_WITH_MPI 3358 const ::LinearAlgebra::EpetraWrappers::Vector &)
const;
3362 const MPI::Vector &)
const;
3365 const ::Vector<double> &)
const;
3368 const ::LinearAlgebra::distributed::Vector<double> &)
const;
3369 #ifdef DEAL_II_WITH_MPI 3372 const ::LinearAlgebra::EpetraWrappers::Vector &)
const;
3376 DEAL_II_NAMESPACE_CLOSE
3378 #endif // DEAL_II_WITH_TRILINOS TrilinosScalar diag_element(const size_type i) const
Iterator lower_bound(Iterator first, Iterator last, const T &val)
static ::ExceptionBase & ExcTrilinosError(int arg1)
TrilinosScalar el(const size_type i, const size_type j) const
types::global_dof_index size_type
virtual bool UseTranspose() const
const Epetra_Map & col_partitioner() const
const Epetra_Map & domain_partitioner() const
TrilinosScalar frobenius_norm() const
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
real_type l2_norm() const
#define AssertDimension(dim1, dim2)
TrilinosPayload null_payload() const
Contents is actually a matrix.
static ::ExceptionBase & ExcSourceEqualsDestination()
static ::ExceptionBase & ExcIO()
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
const Epetra_Map & range_partitioner() const
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
void vmult_add(VectorType &dst, const VectorType &src) const
std::function< void(VectorType &, const VectorType &)> inv_Tvmult
virtual double NormInf() const
Epetra_CombineMode last_action
const Epetra_FECrsGraph & trilinos_sparsity_pattern() const
void reinit(const Vector &v, const bool omit_zeroing_entries=false, const bool allow_different_maps=false)
void set(const size_type i, const size_type j, const number value)
::types::global_dof_index size_type
#define AssertIndexRange(index, range)
TrilinosWrappers::types::int_type min_my_gid(const Epetra_BlockMap &map)
std::pair< size_type, size_type > local_range() const
TrilinosScalar operator()(const size_type i, const size_type j) const
std::function< void(VectorType &, const VectorType &)> Tvmult
virtual int ApplyInverse(const VectorType &Y, VectorType &X) const
static ::ExceptionBase & ExcInvalidIndex(size_type arg1, size_type arg2)
const Epetra_Map & row_partitioner() const
#define AssertThrow(cond, exc)
static ::ExceptionBase & ExcAccessToNonPresentElement(size_type arg1, size_type arg2)
SymmetricTensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > operator+(const SymmetricTensor< rank_, dim, Number > &left, const SymmetricTensor< rank_, dim, OtherNumber > &right)
TrilinosScalar residual(MPI::Vector &dst, const MPI::Vector &x, const MPI::Vector &b) const
const Epetra_Comm & comm_self()
static ::ExceptionBase & ExcTrilinosError(int arg1)
SymmetricTensor< rank_, dim, Number > operator*(const SymmetricTensor< rank_, dim, Number > &t, const Number &factor)
TrilinosScalar matrix_scalar_product(const MPI::Vector &u, const MPI::Vector &v) const
virtual void reinit(const SparsityPattern &sparsity)
std::unique_ptr< Epetra_CrsGraph > nonlocal_graph
static ::ExceptionBase & ExcDivideByZero()
void clear_rows(const std::vector< size_type > &rows, const TrilinosScalar new_diag_value=0)
void Tvmult(VectorType &dst, const VectorType &src) const
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
std::shared_ptr< std::vector< size_type > > colnum_cache
SparseMatrix & operator*=(const TrilinosScalar factor)
std::unique_ptr< Epetra_FECrsMatrix > matrix
SparseMatrix & operator=(const SparseMatrix &)=delete
TrilinosPayload transpose_payload() const
bool in_local_range(const size_type index) const
static ::ExceptionBase & ExcMessage(std::string arg1)
TrilinosPayload identity_payload() const
unsigned int global_dof_index
void add(const size_type i, const size_type j, const TrilinosScalar value)
virtual bool HasNormInf() const
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::unique_ptr< Epetra_Export > nonlocal_matrix_exporter
TrilinosWrappers::types::int_type n_global_cols(const Epetra_CrsGraph &graph)
size_type memory_consumption() const
TrilinosWrappers::types::int_type global_row_index(const Epetra_CrsMatrix &matrix, const ::types::global_dof_index i)
const IndexSet & row_index_set() const
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
Epetra_MpiComm communicator
std::function< void(VectorType &, const VectorType &)> inv_vmult
TrilinosScalar matrix_norm_square(const MPI::Vector &v) const
const Epetra_FEVector & trilinos_vector() const
IndexSet locally_owned_range_indices() const
const Epetra_CrsMatrix & trilinos_matrix() const
SymmetricTensor< rank_, dim, Number > transpose(const SymmetricTensor< rank_, dim, Number > &t)
const Epetra_Map & vector_partitioner() const
TrilinosWrappers::types::int_type global_column_index(const Epetra_CrsMatrix &matrix, const ::types::global_dof_index i)
Epetra_Operator & trilinos_operator() const
virtual int SetUseTranspose(bool UseTranspose)
void vmult(VectorType &dst, const VectorType &src) const
Epetra_MultiVector VectorType
size_type n_nonzero_elements() const
MPI_Comm get_mpi_communicator() const
void reinit(const size_type m, const size_type n, const unsigned int max_per_row)
static ::ExceptionBase & ExcNotQuadratic()
TrilinosWrappers::types::int_type n_global_elements(const Epetra_BlockMap &map)
void swap(Vector< Number > &u, Vector< Number > &v)
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
virtual int Apply(const VectorType &X, VectorType &Y) const
virtual const Epetra_Map & OperatorDomainMap() const
void reinit(const SparsityPatternType &sparsity_pattern)
size_type row_length(const size_type row) const
TrilinosScalar linfty_norm() const
void copy_from(const SparseMatrix &source)
const Epetra_Map & domain_partitioner() const
std::function< void(VectorType &, const VectorType &)> vmult
SparseMatrix & operator/=(const TrilinosScalar factor)
TrilinosWrappers::types::int_type max_my_gid(const Epetra_BlockMap &map)
std::pair< size_type, size_type > local_range() const
static ::ExceptionBase & ExcNotImplemented()
IndexSet locally_owned_domain_indices() const
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
virtual const Epetra_Comm & Comm() const
static ::ExceptionBase & ExcMatrixNotCompressed()
const Epetra_CrsGraph & trilinos_sparsity_pattern() const
bool is_element(const size_type index) const
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
std::unique_ptr< Epetra_Map > column_space_map
void compress(::VectorOperation::values operation)
const Epetra_MultiVector & trilinos_vector() const
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
std::unique_ptr< Epetra_CrsMatrix > nonlocal_matrix
MPI_Comm get_mpi_communicator() const
void set(const size_type i, const size_type j, const TrilinosScalar value)
size_type n_elements() const
#define AssertIsFinite(number)
unsigned int local_size() const
virtual const char * Label() const
TrilinosScalar l1_norm() const
void compress(::VectorOperation::values)
static ::ExceptionBase & ExcInternalError()
virtual const Epetra_Map & OperatorRangeMap() const
std::shared_ptr< std::vector< TrilinosScalar > > value_cache