deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
trilinos_sparse_matrix.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2008 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_trilinos_sparse_matrix_h
16# define dealii_trilinos_sparse_matrix_h
17
18
19# include <deal.II/base/config.h>
20
21# ifdef DEAL_II_WITH_TRILINOS
22
26
35
37# include <Epetra_Comm.h>
38# include <Epetra_CrsGraph.h>
39# include <Epetra_Export.h>
40# include <Epetra_FECrsMatrix.h>
41# include <Epetra_Map.h>
42# include <Epetra_MpiComm.h>
43# include <Epetra_MultiVector.h>
44# include <Epetra_Operator.h>
46
47# include <cmath>
48# include <iterator>
49# include <memory>
50# include <type_traits>
51# include <vector>
52
54
55// forward declarations
56# ifndef DOXYGEN
57template <typename MatrixType>
58class BlockMatrixBase;
59
60template <typename number>
61class SparseMatrix;
62class SparsityPattern;
64
65namespace TrilinosWrappers
66{
67 class SparseMatrix;
68 class SparsityPattern;
69
70 namespace SparseMatrixIterators
71 {
72 template <bool Constness>
73 class Iterator;
74 }
75} // namespace TrilinosWrappers
76# endif
77
78namespace TrilinosWrappers
79{
84 {
89
94 std::size_t,
95 std::size_t,
96 std::size_t,
97 << "You tried to access row " << arg1
98 << " of a distributed sparsity pattern, "
99 << " but only rows " << arg2 << " through " << arg3
100 << " are stored locally and can be accessed.");
101
110 {
111 public:
116
121 const size_type row,
122 const size_type index);
123
128 row() const;
129
134 index() const;
135
140 column() const;
141
142 protected:
154
159
165 void
167
180 std::shared_ptr<std::vector<size_type>> colnum_cache;
181
185 std::shared_ptr<std::vector<TrilinosScalar>> value_cache;
186 };
187
198 template <bool Constess>
199 class Accessor : public AccessorBase
200 {
205 value() const;
206
212 };
213
214
215
219 template <>
220 class Accessor<true> : public AccessorBase
221 {
222 public:
228
234
239 template <bool Other>
241
246 value() const;
247
248 private:
249 // Make iterator class a friend.
250 template <bool>
251 friend class Iterator;
252 };
253
257 template <>
258 class Accessor<false> : public AccessorBase
259 {
261 {
262 public:
267
271 operator TrilinosScalar() const;
272
276 const Reference &
277 operator=(const TrilinosScalar n) const;
278
282 const Reference &
284
288 const Reference &
290
294 const Reference &
296
300 const Reference &
302
303 private:
309 };
310
311 public:
317
323
328 value() const;
329
330 private:
331 // Make iterator class a friend.
332 template <bool>
333 friend class Iterator;
334
335 // Make Reference object a friend.
336 friend class Reference;
337 };
338
352 template <bool Constness>
354 {
355 public:
360
366
372
378
383 Iterator(MatrixType *matrix, const size_type row, const size_type index);
384
388 template <bool Other>
390
396
402
406 const Accessor<Constness> &
407 operator*() const;
408
412 const Accessor<Constness> *
413 operator->() const;
414
419 template <bool OtherConstness>
420 bool
422
426 template <bool OtherConstness>
427 bool
429
435 template <bool OtherConstness>
436 bool
437 operator<(const Iterator<OtherConstness> &) const;
438
442 template <bool OtherConstness>
443 bool
445
450 size_type,
451 size_type,
452 << "Attempt to access element " << arg2 << " of row "
453 << arg1 << " which doesn't have that many elements.");
454
455 private:
460
461 template <bool Other>
462 friend class Iterator;
463 };
464
465 } // namespace SparseMatrixIterators
466} // namespace TrilinosWrappers
467
468DEAL_II_NAMESPACE_CLOSE // Do not convert for module purposes
469
470 namespace std
471{
472 template <bool Constness>
475 {
476 using iterator_category = forward_iterator_tag;
478 typename ::TrilinosWrappers::SparseMatrixIterators::Iterator<
479 Constness>::value_type;
481 typename ::TrilinosWrappers::SparseMatrixIterators::Iterator<
482 Constness>::difference_type;
483 };
484} // namespace std
485
486DEAL_II_NAMESPACE_OPEN // Do not convert for module purposes
487
488
489 namespace TrilinosWrappers
490{
552 {
553 public:
558
563 std::size_t,
564 << "You tried to access row " << arg1
565 << " of a non-contiguous locally owned row set."
566 << " The row " << arg1
567 << " is not stored locally and can't be accessed.");
568
576 struct Traits
577 {
582 static const bool zero_addition_can_be_elided = true;
583 };
584
589
594
599
607 SparseMatrix();
608
617 const size_type n,
618 const unsigned int n_max_entries_per_row);
619
628 const size_type n,
629 const std::vector<unsigned int> &n_entries_per_row);
630
634 SparseMatrix(const SparsityPattern &InputSparsityPattern);
635
640 SparseMatrix(SparseMatrix &&other) noexcept;
641
645 SparseMatrix(const SparseMatrix &) = delete;
646
651 operator=(const SparseMatrix &) = delete;
652
656 virtual ~SparseMatrix() override = default;
657
673 template <typename SparsityPatternType>
674 void
675 reinit(const SparsityPatternType &sparsity_pattern);
676
689 void
690 reinit(const SparsityPattern &sparsity_pattern);
691
700 void
701 reinit(const SparseMatrix &sparse_matrix);
702
723 template <typename number>
724 void
725 reinit(const ::SparseMatrix<number> &dealii_sparse_matrix,
726 const double drop_tolerance = 1e-13,
727 const bool copy_values = true,
728 const ::SparsityPattern *use_this_sparsity = nullptr);
729
735 void
736 reinit(const Epetra_CrsMatrix &input_matrix, const bool copy_values = true);
738
755 SparseMatrix(const IndexSet &parallel_partitioning,
756 const MPI_Comm communicator = MPI_COMM_WORLD,
757 const unsigned int n_max_entries_per_row = 0);
758
766 SparseMatrix(const IndexSet &parallel_partitioning,
767 const MPI_Comm communicator,
768 const std::vector<unsigned int> &n_entries_per_row);
769
784 SparseMatrix(const IndexSet &row_parallel_partitioning,
785 const IndexSet &col_parallel_partitioning,
786 const MPI_Comm communicator = MPI_COMM_WORLD,
787 const size_type n_max_entries_per_row = 0);
788
803 SparseMatrix(const IndexSet &row_parallel_partitioning,
804 const IndexSet &col_parallel_partitioning,
805 const MPI_Comm communicator,
806 const std::vector<unsigned int> &n_entries_per_row);
807
828 template <typename SparsityPatternType>
829 void
830 reinit(const IndexSet &parallel_partitioning,
831 const SparsityPatternType &sparsity_pattern,
832 const MPI_Comm communicator = MPI_COMM_WORLD,
833 const bool exchange_data = false);
834
847 template <typename SparsityPatternType>
848 std::enable_if_t<
849 !std::is_same_v<SparsityPatternType, ::SparseMatrix<double>>>
850 reinit(const IndexSet &row_parallel_partitioning,
851 const IndexSet &col_parallel_partitioning,
852 const SparsityPatternType &sparsity_pattern,
853 const MPI_Comm communicator = MPI_COMM_WORLD,
854 const bool exchange_data = false);
855
872 template <typename number>
873 void
874 reinit(const IndexSet &parallel_partitioning,
875 const ::SparseMatrix<number> &dealii_sparse_matrix,
876 const MPI_Comm communicator = MPI_COMM_WORLD,
877 const double drop_tolerance = 1e-13,
878 const bool copy_values = true,
879 const ::SparsityPattern *use_this_sparsity = nullptr);
880
894 template <typename number>
895 void
896 reinit(const IndexSet &row_parallel_partitioning,
897 const IndexSet &col_parallel_partitioning,
898 const ::SparseMatrix<number> &dealii_sparse_matrix,
899 const MPI_Comm communicator = MPI_COMM_WORLD,
900 const double drop_tolerance = 1e-13,
901 const bool copy_values = true,
902 const ::SparsityPattern *use_this_sparsity = nullptr);
908
913 m() const;
914
919 n() const;
920
929 unsigned int
930 local_size() const;
931
940 std::pair<size_type, size_type>
941 local_range() const;
942
947 bool
948 in_local_range(const size_type index) const;
949
954 std::uint64_t
956
960 unsigned int
961 row_length(const size_type row) const;
962
969 bool
971
978 memory_consumption() const;
979
984 get_mpi_communicator() const;
985
991
1001 SparseMatrix &
1002 operator=(const double d);
1003
1011 void
1012 clear();
1013
1041 void
1043
1065 void
1066 set(const size_type i, const size_type j, const TrilinosScalar value);
1067
1100 void
1101 set(const std::vector<size_type> &indices,
1102 const FullMatrix<TrilinosScalar> &full_matrix,
1103 const bool elide_zero_values = false);
1104
1110 void
1111 set(const std::vector<size_type> &row_indices,
1112 const std::vector<size_type> &col_indices,
1113 const FullMatrix<TrilinosScalar> &full_matrix,
1114 const bool elide_zero_values = false);
1115
1143 void
1144 set(const size_type row,
1145 const std::vector<size_type> &col_indices,
1146 const std::vector<TrilinosScalar> &values,
1147 const bool elide_zero_values = false);
1148
1176 template <typename Number>
1177 void
1178 set(const size_type row,
1179 const size_type n_cols,
1180 const size_type *col_indices,
1181 const Number *values,
1182 const bool elide_zero_values = false);
1183
1193 void
1194 add(const size_type i, const size_type j, const TrilinosScalar value);
1195
1214 void
1215 add(const std::vector<size_type> &indices,
1216 const FullMatrix<TrilinosScalar> &full_matrix,
1217 const bool elide_zero_values = true);
1218
1224 void
1225 add(const std::vector<size_type> &row_indices,
1226 const std::vector<size_type> &col_indices,
1227 const FullMatrix<TrilinosScalar> &full_matrix,
1228 const bool elide_zero_values = true);
1229
1243 void
1244 add(const size_type row,
1245 const std::vector<size_type> &col_indices,
1246 const std::vector<TrilinosScalar> &values,
1247 const bool elide_zero_values = true);
1248
1262 void
1263 add(const size_type row,
1264 const size_type n_cols,
1265 const size_type *col_indices,
1266 const TrilinosScalar *values,
1267 const bool elide_zero_values = true,
1268 const bool col_indices_are_sorted = false);
1269
1273 SparseMatrix &
1274 operator*=(const TrilinosScalar factor);
1275
1279 SparseMatrix &
1280 operator/=(const TrilinosScalar factor);
1281
1285 void
1286 copy_from(const SparseMatrix &source);
1287
1295 void
1296 add(const TrilinosScalar factor, const SparseMatrix &matrix);
1297
1324 void
1325 clear_row(const size_type row, const TrilinosScalar new_diag_value = 0);
1326
1347 void
1349 const TrilinosScalar new_diag_value = 0);
1350
1360 void
1361 transpose();
1362
1368
1378 operator()(const size_type i, const size_type j) const;
1379
1397 el(const size_type i, const size_type j) const;
1398
1406 diag_element(const size_type i) const;
1407
1413
1444 template <typename VectorType>
1445 std::enable_if_t<
1446 std::is_same_v<typename VectorType::value_type, TrilinosScalar>>
1447 vmult(VectorType &dst, const VectorType &src) const;
1448
1455 template <typename VectorType>
1456 std::enable_if_t<
1457 !std::is_same_v<typename VectorType::value_type, TrilinosScalar>>
1458 vmult(VectorType &dst, const VectorType &src) const;
1459
1474 template <typename VectorType>
1475 std::enable_if_t<
1476 std::is_same_v<typename VectorType::value_type, TrilinosScalar>>
1477 Tvmult(VectorType &dst, const VectorType &src) const;
1478
1485 template <typename VectorType>
1486 std::enable_if_t<
1487 !std::is_same_v<typename VectorType::value_type, TrilinosScalar>>
1488 Tvmult(VectorType &dst, const VectorType &src) const;
1489
1499 template <typename VectorType>
1500 void
1501 vmult_add(VectorType &dst, const VectorType &src) const;
1502
1513 template <typename VectorType>
1514 void
1515 Tvmult_add(VectorType &dst, const VectorType &src) const;
1516
1539 matrix_norm_square(const MPI::Vector &v) const;
1540
1561 matrix_scalar_product(const MPI::Vector &u, const MPI::Vector &v) const;
1562
1579 template <typename VectorType>
1581 residual(VectorType &dst, const VectorType &x, const VectorType &b) const;
1582
1597 void
1599 const SparseMatrix &B,
1600 const MPI::Vector &V = MPI::Vector()) const;
1601
1602
1619 void
1621 const SparseMatrix &B,
1622 const MPI::Vector &V = MPI::Vector()) const;
1623
1629
1638 l1_norm() const;
1639
1649 linfty_norm() const;
1650
1656 frobenius_norm() const;
1657
1663
1668 const Epetra_CrsMatrix &
1670
1675 const Epetra_CrsGraph &
1677
1679
1684
1689 IndexSet
1691
1697 IndexSet
1699
1701
1706
1726 begin() const;
1727
1731 iterator
1733
1739 end() const;
1740
1744 iterator
1746
1776 begin(const size_type r) const;
1777
1781 iterator
1783
1794 end(const size_type r) const;
1795
1799 iterator
1800 end(const size_type r);
1801
1807
1813 void
1814 write_ascii();
1815
1823 void
1824 print(std::ostream &out,
1825 const bool write_extended_trilinos_info = false) const;
1826
1837 int,
1838 << "An error with error number " << arg1
1839 << " occurred while calling a Trilinos function. "
1840 "\n\n"
1841 "For historical reasons, many Trilinos functions express "
1842 "errors by returning specific integer values to indicate "
1843 "certain errors. Unfortunately, different Trilinos functions "
1844 "often use the same integer values for different kinds of "
1845 "errors, and in most cases it is also not documented what "
1846 "each error code actually means. As a consequence, it is often "
1847 "difficult to say what a particular error (in this case, "
1848 "the error with integer code '"
1849 << arg1
1850 << "') represents and how one should fix a code to avoid it. "
1851 "The best one can often do is to look up the call stack to "
1852 "see which deal.II function generated the error, and which "
1853 "Trilinos function the error code had originated from; "
1854 "then look up the Trilinos source code of that function (for "
1855 "example on github) to see what code path set that error "
1856 "code. Short of going through all of that, the only other "
1857 "option is to guess the cause of the error from "
1858 "the context in which the error appeared.");
1859
1860
1865 size_type,
1866 size_type,
1867 << "The entry with index <" << arg1 << ',' << arg2
1868 << "> does not exist.");
1869
1874 "You are attempting an operation on two vectors that "
1875 "are the same object, but the operation requires that the "
1876 "two objects are in fact different.");
1877
1882
1887 size_type,
1888 size_type,
1889 size_type,
1890 size_type,
1891 << "You tried to access element (" << arg1 << '/' << arg2
1892 << ')'
1893 << " of a distributed matrix, but only rows in range ["
1894 << arg3 << ',' << arg4
1895 << "] are stored locally and can be accessed.");
1896
1901 size_type,
1902 size_type,
1903 << "You tried to access element (" << arg1 << '/' << arg2
1904 << ')' << " of a sparse matrix, but it appears to not"
1905 << " exist in the Trilinos sparsity pattern.");
1907
1908
1909
1910 protected:
1911 private:
1916 std::unique_ptr<Epetra_Map> column_space_map;
1917
1923 std::unique_ptr<Epetra_FECrsMatrix> matrix;
1924
1930 std::unique_ptr<Epetra_CrsMatrix> nonlocal_matrix;
1931
1935 std::unique_ptr<Epetra_Export> nonlocal_matrix_exporter;
1936
1948 Epetra_CombineMode last_action;
1949
1955
1968 void
1970
1978 void
1980
1981 // To allow calling protected prepare_add() and prepare_set().
1982 friend class BlockMatrixBase<SparseMatrix>;
1983 };
1984
1985
1986
1987 // forwards declarations
1988 class SolverBase;
1989 class PreconditionBase;
1990
1991 namespace internal
1992 {
1993 inline void
1994 check_vector_map_equality(const Epetra_CrsMatrix &mtrx,
1995 const Epetra_MultiVector &src,
1996 const Epetra_MultiVector &dst,
1997 const bool transpose)
1998 {
1999 if (transpose == false)
2000 {
2001 Assert(src.Map().SameAs(mtrx.DomainMap()) == true,
2002 ExcMessage(
2003 "Column map of matrix does not fit with vector map!"));
2004 Assert(dst.Map().SameAs(mtrx.RangeMap()) == true,
2005 ExcMessage("Row map of matrix does not fit with vector map!"));
2006 }
2007 else
2008 {
2009 Assert(src.Map().SameAs(mtrx.RangeMap()) == true,
2010 ExcMessage(
2011 "Column map of matrix does not fit with vector map!"));
2012 Assert(dst.Map().SameAs(mtrx.DomainMap()) == true,
2013 ExcMessage("Row map of matrix does not fit with vector map!"));
2014 }
2015 }
2016
2017 inline void
2019 const Epetra_MultiVector &src,
2020 const Epetra_MultiVector &dst,
2021 const bool transpose)
2022 {
2023 if (transpose == false)
2024 {
2025 Assert(src.Map().SameAs(op.OperatorDomainMap()) == true,
2026 ExcMessage(
2027 "Column map of operator does not fit with vector map!"));
2028 Assert(dst.Map().SameAs(op.OperatorRangeMap()) == true,
2029 ExcMessage(
2030 "Row map of operator does not fit with vector map!"));
2031 }
2032 else
2033 {
2034 Assert(src.Map().SameAs(op.OperatorRangeMap()) == true,
2035 ExcMessage(
2036 "Column map of operator does not fit with vector map!"));
2037 Assert(dst.Map().SameAs(op.OperatorDomainMap()) == true,
2038 ExcMessage(
2039 "Row map of operator does not fit with vector map!"));
2040 }
2041 }
2042
2043
2045 {
2066 {
2067 public:
2071 using VectorType = Epetra_MultiVector;
2072
2077
2082
2087
2092
2101
2105 TrilinosPayload(const TrilinosWrappers::SparseMatrix &matrix_exemplar,
2106 const TrilinosWrappers::SparseMatrix &matrix);
2107
2111 TrilinosPayload(const TrilinosPayload &payload_exemplar,
2112 const TrilinosWrappers::SparseMatrix &matrix);
2113
2118 const TrilinosWrappers::SparseMatrix &matrix_exemplar,
2119 const TrilinosWrappers::PreconditionBase &preconditioner);
2120
2125 const TrilinosWrappers::PreconditionBase &preconditioner_exemplar,
2126 const TrilinosWrappers::PreconditionBase &preconditioner);
2127
2132 const TrilinosPayload &payload_exemplar,
2133 const TrilinosWrappers::PreconditionBase &preconditioner);
2134
2138 TrilinosPayload(const TrilinosPayload &payload);
2139
2147 TrilinosPayload(const TrilinosPayload &first_op,
2148 const TrilinosPayload &second_op);
2149
2153 virtual ~TrilinosPayload() override = default;
2154
2159 identity_payload() const;
2160
2165 null_payload() const;
2166
2171 transpose_payload() const;
2172
2189 template <typename Solver, typename Preconditioner>
2190 std::enable_if_t<
2191 std::is_base_of_v<TrilinosWrappers::SolverBase, Solver> &&
2192 std::is_base_of_v<TrilinosWrappers::PreconditionBase,
2193 Preconditioner>,
2195 inverse_payload(Solver &, const Preconditioner &) const;
2196
2214 template <typename Solver, typename Preconditioner>
2215 std::enable_if_t<
2216 !(std::is_base_of_v<TrilinosWrappers::SolverBase, Solver> &&
2217 std::is_base_of_v<TrilinosWrappers::PreconditionBase,
2218 Preconditioner>),
2220 inverse_payload(Solver &, const Preconditioner &) const;
2221
2223
2228
2234 IndexSet
2236
2242 IndexSet
2244
2248 MPI_Comm
2249 get_mpi_communicator() const;
2250
2257 void
2258 transpose();
2259
2267 std::function<void(VectorType &, const VectorType &)> vmult;
2268
2276 std::function<void(VectorType &, const VectorType &)> Tvmult;
2277
2286 std::function<void(VectorType &, const VectorType &)> inv_vmult;
2287
2296 std::function<void(VectorType &, const VectorType &)> inv_Tvmult;
2297
2299
2304
2311 virtual bool
2312 UseTranspose() const override;
2313
2329 virtual int
2330 SetUseTranspose(bool UseTranspose) override;
2331
2343 virtual int
2344 Apply(const VectorType &X, VectorType &Y) const override;
2345
2364 virtual int
2365 ApplyInverse(const VectorType &Y, VectorType &X) const override;
2367
2372
2379 virtual const char *
2380 Label() const override;
2381
2389 virtual const Epetra_Comm &
2390 Comm() const override;
2391
2399 virtual const Epetra_Map &
2400 OperatorDomainMap() const override;
2401
2410 virtual const Epetra_Map &
2411 OperatorRangeMap() const override;
2413
2414 private:
2424 template <typename EpetraOpType>
2425 TrilinosPayload(EpetraOpType &op,
2426 const bool supports_inverse_operations,
2427 const bool use_transpose,
2428 const MPI_Comm mpi_communicator,
2431
2437
2442 Epetra_MpiComm communicator;
2443
2448 Epetra_Map domain_map;
2449
2454 Epetra_Map range_map;
2455
2464 virtual bool
2465 HasNormInf() const override;
2466
2474 virtual double
2475 NormInf() const override;
2476 };
2477
2483 operator+(const TrilinosPayload &first_op,
2484 const TrilinosPayload &second_op);
2485
2491 operator*(const TrilinosPayload &first_op,
2492 const TrilinosPayload &second_op);
2493
2494 } // namespace LinearOperatorImplementation
2495 } /* namespace internal */
2496
2497
2498
2499 // ----------------------- inline and template functions --------------------
2500
2501# ifndef DOXYGEN
2502
2503 namespace SparseMatrixIterators
2504 {
2506 size_type row,
2507 size_type index)
2508 : matrix(matrix)
2509 , a_row(row)
2510 , a_index(index)
2511 {
2512 visit_present_row();
2513 }
2514
2515
2516 inline AccessorBase::size_type
2517 AccessorBase::row() const
2518 {
2519 Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
2520 return a_row;
2521 }
2522
2523
2524 inline AccessorBase::size_type
2525 AccessorBase::column() const
2526 {
2527 Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
2528 return (*colnum_cache)[a_index];
2529 }
2530
2531
2532 inline AccessorBase::size_type
2533 AccessorBase::index() const
2534 {
2535 Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
2536 return a_index;
2537 }
2538
2539
2540 inline Accessor<true>::Accessor(MatrixType *matrix,
2541 const size_type row,
2542 const size_type index)
2543 : AccessorBase(const_cast<SparseMatrix *>(matrix), row, index)
2544 {}
2545
2546
2547 template <bool Other>
2548 inline Accessor<true>::Accessor(const Accessor<Other> &other)
2549 : AccessorBase(other)
2550 {}
2551
2552
2553 inline TrilinosScalar
2554 Accessor<true>::value() const
2555 {
2556 Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
2557 return (*value_cache)[a_index];
2558 }
2559
2560
2561 inline Accessor<false>::Reference::Reference(const Accessor<false> &acc)
2562 : accessor(const_cast<Accessor<false> &>(acc))
2563 {}
2564
2565
2566 inline Accessor<false>::Reference::operator TrilinosScalar() const
2567 {
2568 return (*accessor.value_cache)[accessor.a_index];
2569 }
2570
2571 inline const Accessor<false>::Reference &
2572 Accessor<false>::Reference::operator=(const TrilinosScalar n) const
2573 {
2574 (*accessor.value_cache)[accessor.a_index] = n;
2575 accessor.matrix->set(accessor.row(),
2576 accessor.column(),
2577 static_cast<TrilinosScalar>(*this));
2578 return *this;
2579 }
2580
2581
2582 inline const Accessor<false>::Reference &
2583 Accessor<false>::Reference::operator+=(const TrilinosScalar n) const
2584 {
2585 (*accessor.value_cache)[accessor.a_index] += n;
2586 accessor.matrix->set(accessor.row(),
2587 accessor.column(),
2588 static_cast<TrilinosScalar>(*this));
2589 return *this;
2590 }
2591
2592
2593 inline const Accessor<false>::Reference &
2594 Accessor<false>::Reference::operator-=(const TrilinosScalar n) const
2595 {
2596 (*accessor.value_cache)[accessor.a_index] -= n;
2597 accessor.matrix->set(accessor.row(),
2598 accessor.column(),
2599 static_cast<TrilinosScalar>(*this));
2600 return *this;
2601 }
2602
2603
2604 inline const Accessor<false>::Reference &
2605 Accessor<false>::Reference::operator*=(const TrilinosScalar n) const
2606 {
2607 (*accessor.value_cache)[accessor.a_index] *= n;
2608 accessor.matrix->set(accessor.row(),
2609 accessor.column(),
2610 static_cast<TrilinosScalar>(*this));
2611 return *this;
2612 }
2613
2614
2615 inline const Accessor<false>::Reference &
2616 Accessor<false>::Reference::operator/=(const TrilinosScalar n) const
2617 {
2618 (*accessor.value_cache)[accessor.a_index] /= n;
2619 accessor.matrix->set(accessor.row(),
2620 accessor.column(),
2621 static_cast<TrilinosScalar>(*this));
2622 return *this;
2623 }
2624
2625
2626 inline Accessor<false>::Accessor(MatrixType *matrix,
2627 const size_type row,
2628 const size_type index)
2629 : AccessorBase(matrix, row, index)
2630 {}
2631
2632
2633 inline Accessor<false>::Reference
2634 Accessor<false>::value() const
2635 {
2636 Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
2637 return {*this};
2638 }
2639
2640
2641
2642 template <bool Constness>
2643 inline Iterator<Constness>::Iterator(MatrixType *matrix,
2644 const size_type row,
2645 const size_type index)
2646 : accessor(matrix, row, index)
2647 {}
2648
2649
2650 template <bool Constness>
2651 template <bool Other>
2652 inline Iterator<Constness>::Iterator(const Iterator<Other> &other)
2653 : accessor(other.accessor)
2654 {}
2655
2656
2657 template <bool Constness>
2658 inline Iterator<Constness> &
2659 Iterator<Constness>::operator++()
2660 {
2661 Assert(accessor.a_row < accessor.matrix->m(), ExcIteratorPastEnd());
2662
2663 ++accessor.a_index;
2664
2665 // If at end of line: do one
2666 // step, then cycle until we
2667 // find a row with a nonzero
2668 // number of entries.
2669 if (accessor.a_index >= accessor.colnum_cache->size())
2670 {
2671 accessor.a_index = 0;
2672 ++accessor.a_row;
2673
2674 while ((accessor.a_row < accessor.matrix->m()) &&
2675 ((accessor.matrix->in_local_range(accessor.a_row) == false) ||
2676 (accessor.matrix->row_length(accessor.a_row) == 0)))
2677 ++accessor.a_row;
2678
2679 accessor.visit_present_row();
2680 }
2681 return *this;
2682 }
2683
2684
2685 template <bool Constness>
2686 inline Iterator<Constness>
2687 Iterator<Constness>::operator++(int)
2688 {
2689 const Iterator<Constness> old_state = *this;
2690 ++(*this);
2691 return old_state;
2692 }
2693
2694
2695
2696 template <bool Constness>
2697 inline const Accessor<Constness> &
2698 Iterator<Constness>::operator*() const
2699 {
2700 return accessor;
2701 }
2702
2703
2704
2705 template <bool Constness>
2706 inline const Accessor<Constness> *
2707 Iterator<Constness>::operator->() const
2708 {
2709 return &accessor;
2710 }
2711
2712
2713
2714 template <bool Constness>
2715 template <bool OtherConstness>
2716 inline bool
2717 Iterator<Constness>::operator==(const Iterator<OtherConstness> &other) const
2718 {
2719 return (accessor.a_row == other.accessor.a_row &&
2720 accessor.a_index == other.accessor.a_index);
2721 }
2722
2723
2724
2725 template <bool Constness>
2726 template <bool OtherConstness>
2727 inline bool
2728 Iterator<Constness>::operator!=(const Iterator<OtherConstness> &other) const
2729 {
2730 return !(*this == other);
2731 }
2732
2733
2734
2735 template <bool Constness>
2736 template <bool OtherConstness>
2737 inline bool
2738 Iterator<Constness>::operator<(const Iterator<OtherConstness> &other) const
2739 {
2740 return (accessor.row() < other.accessor.row() ||
2741 (accessor.row() == other.accessor.row() &&
2742 accessor.index() < other.accessor.index()));
2743 }
2744
2745
2746 template <bool Constness>
2747 template <bool OtherConstness>
2748 inline bool
2749 Iterator<Constness>::operator>(const Iterator<OtherConstness> &other) const
2750 {
2751 return (other < *this);
2752 }
2753
2754 } // namespace SparseMatrixIterators
2755
2756
2757
2759 {
2760 return begin(0);
2761 }
2762
2763
2764
2766 {
2767 return const_iterator(this, m(), 0);
2768 }
2769
2770
2771
2772 inline SparseMatrix::const_iterator SparseMatrix::begin(const size_type r)
2773 const
2774 {
2775 AssertIndexRange(r, m());
2776 if (in_local_range(r) && (row_length(r) > 0))
2777 return const_iterator(this, r, 0);
2778 else
2779 return end(r);
2780 }
2781
2782
2783
2784 inline SparseMatrix::const_iterator SparseMatrix::end(const size_type r) const
2785 {
2786 AssertIndexRange(r, m());
2787
2788 // place the iterator on the first entry
2789 // past this line, or at the end of the
2790 // matrix
2791 for (size_type i = r + 1; i < m(); ++i)
2792 if (in_local_range(i) && (row_length(i) > 0))
2793 return const_iterator(this, i, 0);
2794
2795 // if there is no such line, then take the
2796 // end iterator of the matrix
2797 return end();
2798 }
2799
2800
2801
2803 {
2804 return begin(0);
2805 }
2806
2807
2808
2810 {
2811 return iterator(this, m(), 0);
2812 }
2813
2814
2815
2816 inline SparseMatrix::iterator SparseMatrix::begin(const size_type r)
2817 {
2818 AssertIndexRange(r, m());
2819 if (in_local_range(r) && (row_length(r) > 0))
2820 return iterator(this, r, 0);
2821 else
2822 return end(r);
2823 }
2824
2825
2826
2827 inline SparseMatrix::iterator SparseMatrix::end(const size_type r)
2828 {
2829 AssertIndexRange(r, m());
2830
2831 // place the iterator on the first entry
2832 // past this line, or at the end of the
2833 // matrix
2834 for (size_type i = r + 1; i < m(); ++i)
2835 if (in_local_range(i) && (row_length(i) > 0))
2836 return iterator(this, i, 0);
2837
2838 // if there is no such line, then take the
2839 // end iterator of the matrix
2840 return end();
2841 }
2842
2843
2844
2845 inline bool SparseMatrix::in_local_range(const size_type index) const
2846 {
2848# ifndef DEAL_II_WITH_64BIT_INDICES
2849 begin = matrix->RowMap().MinMyGID();
2850 end = matrix->RowMap().MaxMyGID() + 1;
2851# else
2852 begin = matrix->RowMap().MinMyGID64();
2853 end = matrix->RowMap().MaxMyGID64() + 1;
2854# endif
2855
2856 return ((index >= static_cast<size_type>(begin)) &&
2857 (index < static_cast<size_type>(end)));
2858 }
2859
2860
2861
2862 inline bool SparseMatrix::is_compressed() const
2863 {
2864 return compressed;
2865 }
2866
2867
2868
2869 // Inline the set() and add() functions, since they will be called
2870 // frequently, and the compiler can optimize away some unnecessary loops
2871 // when the sizes are given at compile time.
2872 template <>
2873 void SparseMatrix::set<TrilinosScalar>(const size_type row,
2874 const size_type n_cols,
2875 const size_type *col_indices,
2876 const TrilinosScalar *values,
2877 const bool elide_zero_values);
2878
2879
2880
2881 template <typename Number>
2882 void SparseMatrix::set(const size_type row,
2883 const size_type n_cols,
2884 const size_type *col_indices,
2885 const Number *values,
2886 const bool elide_zero_values)
2887 {
2888 std::vector<TrilinosScalar> trilinos_values(n_cols);
2889 std::copy(values, values + n_cols, trilinos_values.begin());
2890 this->set(
2891 row, n_cols, col_indices, trilinos_values.data(), elide_zero_values);
2892 }
2893
2894
2895
2896 inline void SparseMatrix::set(const size_type i,
2897 const size_type j,
2898 const TrilinosScalar value)
2899 {
2901
2902 set(i, 1, &j, &value, false);
2903 }
2904
2905
2906
2907 inline void SparseMatrix::set(const std::vector<size_type> &indices,
2908 const FullMatrix<TrilinosScalar> &values,
2909 const bool elide_zero_values)
2910 {
2911 Assert(indices.size() == values.m(),
2912 ExcDimensionMismatch(indices.size(), values.m()));
2913 Assert(values.m() == values.n(), ExcNotQuadratic());
2914
2915 for (size_type i = 0; i < indices.size(); ++i)
2916 set(indices[i],
2917 indices.size(),
2918 indices.data(),
2919 &values(i, 0),
2920 elide_zero_values);
2921 }
2922
2923
2924
2925 inline void SparseMatrix::add(const size_type i,
2926 const size_type j,
2927 const TrilinosScalar value)
2928 {
2930
2931 if (value == 0)
2932 {
2933 // we have to check after Insert/Add in any case to be consistent
2934 // with the MPI communication model, but we can save some
2935 // work if the addend is zero. However, these actions are done in case
2936 // we pass on to the other function.
2937
2938 // TODO: fix this (do not run compress here, but fail)
2939 if (last_action == Insert)
2940 {
2941 const int ierr = matrix->GlobalAssemble(*column_space_map,
2942 matrix->RowMap(),
2943 false);
2944
2945 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2946 }
2947
2948 last_action = Add;
2949
2950 return;
2951 }
2952 else
2953 add(i, 1, &j, &value, false);
2954 }
2955
2956
2957
2958 // inline "simple" functions that are called frequently and do only involve
2959 // a call to some Trilinos function.
2961 {
2962# ifndef DEAL_II_WITH_64BIT_INDICES
2963 return matrix->NumGlobalRows();
2964# else
2965 return matrix->NumGlobalRows64();
2966# endif
2967 }
2968
2969
2970
2972 {
2973 // If the matrix structure has not been fixed (i.e., we did not have a
2974 // sparsity pattern), it does not know about the number of columns so we
2975 // must always take this from the additional column space map
2976 Assert(column_space_map.get() != nullptr, ExcInternalError());
2977 return n_global_elements(*column_space_map);
2978 }
2979
2980
2981
2982 inline unsigned int SparseMatrix::local_size() const
2983 {
2984 return matrix->NumMyRows();
2985 }
2986
2987
2988
2989 inline std::pair<SparseMatrix::size_type, SparseMatrix::size_type>
2991 {
2993# ifndef DEAL_II_WITH_64BIT_INDICES
2994 begin = matrix->RowMap().MinMyGID();
2995 end = matrix->RowMap().MaxMyGID() + 1;
2996# else
2997 begin = matrix->RowMap().MinMyGID64();
2998 end = matrix->RowMap().MaxMyGID64() + 1;
2999# endif
3000
3001 return std::make_pair(begin, end);
3002 }
3003
3004
3005
3006 inline std::uint64_t SparseMatrix::n_nonzero_elements() const
3007 {
3008 // Trilinos uses 64bit functions internally for attribute access, which
3009 // return `long long`. They also offer 32bit variants that return `int`,
3010 // however those call the 64bit version and convert the values to 32bit.
3011 // There is no necessity in using the 32bit versions at all.
3012 return static_cast<std::uint64_t>(matrix->NumGlobalNonzeros64());
3013 }
3014
3015
3016
3017 template <typename SparsityPatternType>
3018 inline void SparseMatrix::reinit(const IndexSet &parallel_partitioning,
3019 const SparsityPatternType &sparsity_pattern,
3020 const MPI_Comm communicator,
3021 const bool exchange_data)
3022 {
3023 reinit(parallel_partitioning,
3024 parallel_partitioning,
3025 sparsity_pattern,
3026 communicator,
3027 exchange_data);
3028 }
3029
3030
3031
3032 template <typename number>
3033 inline void SparseMatrix::reinit(
3034 const IndexSet &parallel_partitioning,
3035 const ::SparseMatrix<number> &sparse_matrix,
3036 const MPI_Comm communicator,
3037 const double drop_tolerance,
3038 const bool copy_values,
3039 const ::SparsityPattern *use_this_sparsity)
3040 {
3041 Epetra_Map map =
3042 parallel_partitioning.make_trilinos_map(communicator, false);
3043 reinit(parallel_partitioning,
3044 parallel_partitioning,
3045 sparse_matrix,
3046 drop_tolerance,
3047 copy_values,
3048 use_this_sparsity);
3049 }
3050
3051
3052
3053 inline const Epetra_CrsMatrix &SparseMatrix::trilinos_matrix() const
3054 {
3055 return static_cast<const Epetra_CrsMatrix &>(*matrix);
3056 }
3057
3058
3059
3060 inline const Epetra_CrsGraph &SparseMatrix::trilinos_sparsity_pattern() const
3061 {
3062 return matrix->Graph();
3063 }
3064
3065
3066
3067 inline IndexSet SparseMatrix::locally_owned_domain_indices() const
3068 {
3069 return IndexSet(matrix->DomainMap());
3070 }
3071
3072
3073
3074 inline IndexSet SparseMatrix::locally_owned_range_indices() const
3075 {
3076 return IndexSet(matrix->RangeMap());
3077 }
3078
3079
3080
3081 inline void SparseMatrix::prepare_add()
3082 {
3083 // nothing to do here
3084 }
3085
3086
3087
3088 inline void SparseMatrix::prepare_set()
3089 {
3090 // nothing to do here
3091 }
3092
3093
3094
3095 template <typename VectorType>
3096 inline TrilinosScalar SparseMatrix::residual(VectorType & dst,
3097 const VectorType &x,
3098 const VectorType &b) const
3099 {
3100 vmult(dst, x);
3101 dst -= b;
3102 dst *= -1.;
3103
3104 return dst.l2_norm();
3105 }
3106
3107
3108 namespace internal
3109 {
3110 namespace LinearOperatorImplementation
3111 {
3112 template <typename EpetraOpType>
3113 TrilinosPayload::TrilinosPayload(
3114 EpetraOpType &op,
3115 const bool supports_inverse_operations,
3116 const bool use_transpose,
3117 const MPI_Comm mpi_communicator,
3118 const IndexSet &locally_owned_domain_indices,
3119 const IndexSet &locally_owned_range_indices)
3120 : use_transpose(use_transpose)
3121 , communicator(mpi_communicator)
3122 , domain_map(
3123 locally_owned_domain_indices.make_trilinos_map(communicator.Comm()))
3124 , range_map(
3125 locally_owned_range_indices.make_trilinos_map(communicator.Comm()))
3126 {
3127 vmult = [&op](Range &tril_dst, const Domain &tril_src) {
3128 // Duplicated from TrilinosWrappers::PreconditionBase::vmult
3129 // as well as from TrilinosWrappers::SparseMatrix::Tvmult
3130 Assert(&tril_src != &tril_dst,
3133 tril_src,
3134 tril_dst,
3135 op.UseTranspose());
3136
3137 const int ierr = op.Apply(tril_src, tril_dst);
3138 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
3139 };
3140
3141 Tvmult = [&op](Domain &tril_dst, const Range &tril_src) {
3142 // Duplicated from TrilinosWrappers::PreconditionBase::vmult
3143 // as well as from TrilinosWrappers::SparseMatrix::Tvmult
3144 Assert(&tril_src != &tril_dst,
3147 tril_src,
3148 tril_dst,
3149 !op.UseTranspose());
3150
3151 op.SetUseTranspose(!op.UseTranspose());
3152 const int ierr = op.Apply(tril_src, tril_dst);
3153 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
3154 op.SetUseTranspose(!op.UseTranspose());
3155 };
3156
3157 if (supports_inverse_operations)
3158 {
3159 inv_vmult = [&op](Domain &tril_dst, const Range &tril_src) {
3160 // Duplicated from TrilinosWrappers::PreconditionBase::vmult
3161 // as well as from TrilinosWrappers::SparseMatrix::Tvmult
3162 Assert(
3163 &tril_src != &tril_dst,
3166 tril_src,
3167 tril_dst,
3168 !op.UseTranspose());
3169
3170 const int ierr = op.ApplyInverse(tril_src, tril_dst);
3171 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
3172 };
3173
3174 inv_Tvmult = [&op](Range &tril_dst, const Domain &tril_src) {
3175 // Duplicated from TrilinosWrappers::PreconditionBase::vmult
3176 // as well as from TrilinosWrappers::SparseMatrix::Tvmult
3177 Assert(
3178 &tril_src != &tril_dst,
3181 tril_src,
3182 tril_dst,
3183 op.UseTranspose());
3184
3185 op.SetUseTranspose(!op.UseTranspose());
3186 const int ierr = op.ApplyInverse(tril_src, tril_dst);
3187 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
3188 op.SetUseTranspose(!op.UseTranspose());
3189 };
3190 }
3191 else
3192 {
3193 inv_vmult = [](Domain &, const Range &) {
3194 Assert(false,
3195 ExcMessage(
3196 "Uninitialized TrilinosPayload::inv_vmult called. "
3197 "The operator does not support inverse operations."));
3198 };
3199
3200 inv_Tvmult = [](Range &, const Domain &) {
3201 Assert(false,
3202 ExcMessage(
3203 "Uninitialized TrilinosPayload::inv_Tvmult called. "
3204 "The operator does not support inverse operations."));
3205 };
3206 }
3207 }
3208
3209
3210 template <typename Solver, typename Preconditioner>
3211 std::enable_if_t<
3212 std::is_base_of_v<TrilinosWrappers::SolverBase, Solver> &&
3213 std::is_base_of_v<TrilinosWrappers::PreconditionBase, Preconditioner>,
3216 Solver &solver,
3217 const Preconditioner &preconditioner) const
3218 {
3219 const auto &payload = *this;
3220
3221 TrilinosPayload return_op(payload);
3222
3223 // Capture by copy so the payloads are always valid
3224
3225 return_op.inv_vmult = [payload, &solver, &preconditioner](
3226 TrilinosPayload::Domain &tril_dst,
3227 const TrilinosPayload::Range &tril_src) {
3228 // Duplicated from TrilinosWrappers::PreconditionBase::vmult
3229 // as well as from TrilinosWrappers::SparseMatrix::Tvmult
3230 Assert(&tril_src != &tril_dst,
3233 tril_src,
3234 tril_dst,
3235 !payload.UseTranspose());
3236 solver.solve(payload, tril_dst, tril_src, preconditioner);
3237 };
3238
3239 return_op.inv_Tvmult = [payload, &solver, &preconditioner](
3240 TrilinosPayload::Range &tril_dst,
3241 const TrilinosPayload::Domain &tril_src) {
3242 // Duplicated from TrilinosWrappers::PreconditionBase::vmult
3243 // as well as from TrilinosWrappers::SparseMatrix::Tvmult
3244 Assert(&tril_src != &tril_dst,
3247 tril_src,
3248 tril_dst,
3249 payload.UseTranspose());
3250
3251 const_cast<TrilinosPayload &>(payload).transpose();
3252 solver.solve(payload, tril_dst, tril_src, preconditioner);
3253 const_cast<TrilinosPayload &>(payload).transpose();
3254 };
3255
3256 // If the input operator is already setup for transpose operations, then
3257 // we must do similar with its inverse.
3258 if (return_op.UseTranspose() == true)
3259 std::swap(return_op.inv_vmult, return_op.inv_Tvmult);
3260
3261 return return_op;
3262 }
3263
3264 template <typename Solver, typename Preconditioner>
3265 std::enable_if_t<
3266 !(std::is_base_of_v<TrilinosWrappers::SolverBase, Solver> &&
3267 std::is_base_of_v<TrilinosWrappers::PreconditionBase,
3268 Preconditioner>),
3270 TrilinosPayload::inverse_payload(Solver &, const Preconditioner &) const
3271 {
3272 TrilinosPayload return_op(*this);
3273
3274 return_op.inv_vmult = [](TrilinosPayload::Domain &,
3275 const TrilinosPayload::Range &) {
3276 AssertThrow(false,
3277 ExcMessage("Payload inv_vmult disabled because of "
3278 "incompatible solver/preconditioner choice."));
3279 };
3280
3281 return_op.inv_Tvmult = [](TrilinosPayload::Range &,
3282 const TrilinosPayload::Domain &) {
3283 AssertThrow(false,
3284 ExcMessage("Payload inv_vmult disabled because of "
3285 "incompatible solver/preconditioner choice."));
3286 };
3287
3288 return return_op;
3289 }
3290 } // namespace LinearOperatorImplementation
3291 } // namespace internal
3292
3293 template <>
3294 void SparseMatrix::set<TrilinosScalar>(const size_type row,
3295 const size_type n_cols,
3296 const size_type *col_indices,
3297 const TrilinosScalar *values,
3298 const bool elide_zero_values);
3299# endif // DOXYGEN
3300
3301} /* namespace TrilinosWrappers */
3302
3303
3305
3306
3307# else
3308
3309// Make sure the scripts that create the C++20 module input files have
3310// something to latch on if the preprocessor #ifdef above would
3311// otherwise lead to an empty content of the file.
3314
3315# endif // DEAL_II_WITH_TRILINOS
3316
3317
3318/*----------------------- trilinos_sparse_matrix.h --------------------*/
3319
3320#endif
3321/*----------------------- trilinos_sparse_matrix.h --------------------*/
Epetra_Map make_trilinos_map(const MPI_Comm communicator=MPI_COMM_WORLD, const bool overlapping=false) const
AccessorBase(SparseMatrix *matrix, const size_type row, const size_type index)
std::shared_ptr< std::vector< TrilinosScalar > > value_cache
std::shared_ptr< std::vector< size_type > > colnum_cache
const Reference & operator-=(const TrilinosScalar n) const
const Reference & operator*=(const TrilinosScalar n) const
const Reference & operator+=(const TrilinosScalar n) const
const Reference & operator/=(const TrilinosScalar n) const
const Reference & operator=(const TrilinosScalar n) const
Accessor(MatrixType *matrix, const size_type row, const size_type index)
Accessor(MatrixType *matrix, const size_type row, const size_type index)
const Accessor< Constness > & operator*() const
const Accessor< Constness > * operator->() const
typename Accessor< Constness >::MatrixType MatrixType
bool operator==(const Iterator< OtherConstness > &) const
bool operator!=(const Iterator< OtherConstness > &) const
bool operator<(const Iterator< OtherConstness > &) const
Iterator(const Iterator< Other > &other)
bool operator>(const Iterator< OtherConstness > &) const
Iterator(MatrixType *matrix, const size_type row, const size_type index)
void set(const size_type i, const size_type j, const TrilinosScalar value)
std::unique_ptr< Epetra_Map > column_space_map
TrilinosScalar residual(VectorType &dst, const VectorType &x, const VectorType &b) const
SparseMatrix & operator*=(const TrilinosScalar factor)
SparseMatrixIterators::Iterator< false > iterator
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
std::unique_ptr< Epetra_Export > nonlocal_matrix_exporter
void compress(VectorOperation::values operation)
std::unique_ptr< Epetra_FECrsMatrix > matrix
::types::global_dof_index size_type
const Epetra_CrsMatrix & trilinos_matrix() const
TrilinosScalar matrix_norm_square(const MPI::Vector &v) const
IndexSet locally_owned_range_indices() const
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
std::enable_if_t< std::is_same_v< typename VectorType::value_type, TrilinosScalar > > Tvmult(VectorType &dst, const VectorType &src) const
void clear_row(const size_type row, const TrilinosScalar new_diag_value=0)
void clear_rows(const ArrayView< const size_type > &rows, const TrilinosScalar new_diag_value=0)
const_iterator begin() const
void reinit(const SparsityPatternType &sparsity_pattern)
void vmult_add(VectorType &dst, const VectorType &src) const
SparseMatrix & operator=(const SparseMatrix &)=delete
const_iterator end(const size_type r) const
SparseMatrix & operator/=(const TrilinosScalar factor)
void Tvmult_add(VectorType &dst, const VectorType &src) const
std::enable_if_t< !std::is_same_v< SparsityPatternType, ::SparseMatrix< double > > > reinit(const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const SparsityPatternType &sparsity_pattern, const MPI_Comm communicator=MPI_COMM_WORLD, const bool exchange_data=false)
TrilinosScalar el(const size_type i, const size_type j) const
IndexSet locally_owned_domain_indices() const
bool in_local_range(const size_type index) const
void copy_from(const SparseMatrix &source)
unsigned int row_length(const size_type row) const
std::uint64_t n_nonzero_elements() const
std::enable_if_t< std::is_same_v< typename VectorType::value_type, TrilinosScalar > > vmult(VectorType &dst, const VectorType &src) const
const_iterator begin(const size_type r) const
const Epetra_CrsGraph & trilinos_sparsity_pattern() const
SparseMatrix(const SparseMatrix &)=delete
iterator end(const size_type r)
TrilinosScalar diag_element(const size_type i) const
void add(const size_type i, const size_type j, const TrilinosScalar value)
unsigned int local_size() const
iterator begin(const size_type r)
TrilinosScalar operator()(const size_type i, const size_type j) const
void reinit(const IndexSet &parallel_partitioning, const ::SparseMatrix< number > &dealii_sparse_matrix, const MPI_Comm communicator=MPI_COMM_WORLD, const double drop_tolerance=1e-13, const bool copy_values=true, const ::SparsityPattern *use_this_sparsity=nullptr)
virtual ~SparseMatrix() override=default
void set(const size_type row, const size_type n_cols, const size_type *col_indices, const Number *values, const bool elide_zero_values=false)
SparseMatrixIterators::Iterator< true > const_iterator
TrilinosScalar matrix_scalar_product(const MPI::Vector &u, const MPI::Vector &v) const
const_iterator end() const
void set(const std::vector< size_type > &indices, const FullMatrix< TrilinosScalar > &full_matrix, const bool elide_zero_values=false)
std::pair< size_type, size_type > local_range() const
void reinit(const IndexSet &parallel_partitioning, const SparsityPatternType &sparsity_pattern, const MPI_Comm communicator=MPI_COMM_WORLD, const bool exchange_data=false)
std::unique_ptr< Epetra_CrsMatrix > nonlocal_matrix
std::function< void(VectorType &, const VectorType &)> inv_Tvmult
virtual int ApplyInverse(const VectorType &Y, VectorType &X) const override
virtual int Apply(const VectorType &X, VectorType &Y) const override
std::enable_if_t< std::is_base_of_v< TrilinosWrappers::SolverBase, Solver > &&std::is_base_of_v< TrilinosWrappers::PreconditionBase, Preconditioner >, TrilinosPayload > inverse_payload(Solver &, const Preconditioner &) const
TrilinosPayload(EpetraOpType &op, const bool supports_inverse_operations, const bool use_transpose, const MPI_Comm mpi_communicator, const IndexSet &locally_owned_domain_indices, const IndexSet &locally_owned_range_indices)
std::enable_if_t< !(std::is_base_of_v< TrilinosWrappers::SolverBase, Solver > &&std::is_base_of_v< TrilinosWrappers::PreconditionBase, Preconditioner >), TrilinosPayload > inverse_payload(Solver &, const Preconditioner &) const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition config.h:603
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition config.h:647
DerivativeForm< 1, spacedim, dim, Number > transpose(const DerivativeForm< 1, dim, spacedim, Number > &DF)
static ::ExceptionBase & ExcInvalidIndexWithinRow(size_type arg1, size_type arg2)
#define DeclException0(Exception0)
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
static ::ExceptionBase & ExcAccessToNonlocalRow(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcInvalidIndex(size_type arg1, size_type arg2)
#define Assert(cond, exc)
static ::ExceptionBase & ExcIteratorPastEnd()
static ::ExceptionBase & ExcAccessToNonPresentElement(size_type arg1, size_type arg2)
#define AssertIsFinite(number)
static ::ExceptionBase & ExcAccessToNonlocalRow(std::size_t arg1)
static ::ExceptionBase & ExcBeyondEndOfMatrix()
#define DeclException2(Exception2, type1, type2, outsequence)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcSourceEqualsDestination()
#define DeclException3(Exception3, type1, type2, type3, outsequence)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcNotQuadratic()
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcMatrixNotCompressed()
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcTrilinosError(int arg1)
#define AssertThrow(cond, exc)
const bool IsBlockVector< VectorType >::value
@ matrix
Contents is actually a matrix.
SparseMatrix< double > SparseMatrix
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
TrilinosPayload operator+(const TrilinosPayload &first_op, const TrilinosPayload &second_op)
TrilinosPayload operator*(const TrilinosPayload &first_op, const TrilinosPayload &second_op)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
void check_vector_map_equality(const Epetra_CrsMatrix &mtrx, const Epetra_MultiVector &src, const Epetra_MultiVector &dst, const bool transpose)
TrilinosWrappers::types::int64_type n_global_elements(const Epetra_BlockMap &map)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
STL namespace.
unsigned int global_dof_index
Definition types.h:94
typename ::TrilinosWrappers::SparseMatrixIterators::Iterator< Constness >::difference_type difference_type
typename ::TrilinosWrappers::SparseMatrixIterators::Iterator< Constness >::value_type value_type
double TrilinosScalar
Definition types.h:198