16#ifndef dealii_trilinos_sparse_matrix_h
17# define dealii_trilinos_sparse_matrix_h
22# ifdef DEAL_II_WITH_TRILINOS
36# include <Epetra_Comm.h>
37# include <Epetra_CrsGraph.h>
38# include <Epetra_Export.h>
39# include <Epetra_FECrsMatrix.h>
40# include <Epetra_Map.h>
41# include <Epetra_MultiVector.h>
42# include <Epetra_Operator.h>
46# include <type_traits>
48# ifdef DEAL_II_WITH_MPI
49# include <Epetra_MpiComm.h>
52# include <Epetra_SerialComm.h>
59template <
typename MatrixType>
62template <
typename number>
74 template <
bool Constness>
99 <<
"You tried to access row " << arg1
100 <<
" of a distributed sparsity pattern, "
101 <<
" but only rows " << arg2 <<
" through " << arg3
102 <<
" are stored locally and can be accessed.");
200 template <
bool Constess>
239 template <
bool Other>
336 friend class Reference;
352 template <
bool Constness>
376 template <
bool Other>
434 <<
"Attempt to access element " << arg2 <<
" of row "
435 << arg1 <<
" which doesn't have that many elements.");
443 template <
bool Other>
523 <<
"You tried to access row " << arg1
524 <<
" of a non-contiguous locally owned row set."
525 <<
" The row " << arg1
526 <<
" is not stored locally and can't be accessed.");
577 const unsigned int n_max_entries_per_row);
588 const std::vector<unsigned int> &n_entries_per_row);
632 template <
typename SparsityPatternType>
634 reinit(
const SparsityPatternType &sparsity_pattern);
682 template <
typename number>
684 reinit(const ::SparseMatrix<number> &dealii_sparse_matrix,
685 const double drop_tolerance = 1
e-13,
686 const bool copy_values =
true,
687 const ::SparsityPattern * use_this_sparsity =
nullptr);
695 reinit(
const Epetra_CrsMatrix &input_matrix,
const bool copy_values =
true);
715 const MPI_Comm & communicator = MPI_COMM_WORLD,
716 const unsigned int n_max_entries_per_row = 0);
727 const std::vector<unsigned int> &n_entries_per_row);
744 const IndexSet &col_parallel_partitioning,
745 const MPI_Comm &communicator = MPI_COMM_WORLD,
746 const size_type n_max_entries_per_row = 0);
763 const IndexSet & col_parallel_partitioning,
765 const std::vector<unsigned int> &n_entries_per_row);
787 template <
typename SparsityPatternType>
790 const SparsityPatternType &sparsity_pattern,
791 const MPI_Comm & communicator = MPI_COMM_WORLD,
792 const bool exchange_data =
false);
806 template <
typename SparsityPatternType>
807 typename std::enable_if<
808 !std::is_same<SparsityPatternType,
811 const IndexSet & col_parallel_partitioning,
812 const SparsityPatternType &sparsity_pattern,
813 const MPI_Comm & communicator = MPI_COMM_WORLD,
814 const bool exchange_data =
false);
832 template <
typename number>
835 const ::SparseMatrix<number> &dealii_sparse_matrix,
836 const MPI_Comm & communicator = MPI_COMM_WORLD,
837 const double drop_tolerance = 1
e-13,
838 const bool copy_values =
true,
839 const ::SparsityPattern * use_this_sparsity =
nullptr);
854 template <
typename number>
857 const IndexSet & col_parallel_partitioning,
858 const ::SparseMatrix<number> &dealii_sparse_matrix,
859 const MPI_Comm & communicator = MPI_COMM_WORLD,
860 const double drop_tolerance = 1
e-13,
861 const bool copy_values =
true,
862 const ::SparsityPattern * use_this_sparsity =
nullptr);
900 std::pair<size_type, size_type>
1061 set(
const std::vector<size_type> & indices,
1063 const bool elide_zero_values =
false);
1071 set(
const std::vector<size_type> & row_indices,
1072 const std::vector<size_type> & col_indices,
1074 const bool elide_zero_values =
false);
1105 const std::vector<size_type> & col_indices,
1106 const std::vector<TrilinosScalar> &
values,
1107 const bool elide_zero_values =
false);
1136 template <
typename Number>
1142 const bool elide_zero_values =
false);
1175 add(
const std::vector<size_type> & indices,
1177 const bool elide_zero_values =
true);
1185 add(
const std::vector<size_type> & row_indices,
1186 const std::vector<size_type> & col_indices,
1188 const bool elide_zero_values =
true);
1205 const std::vector<size_type> & col_indices,
1206 const std::vector<TrilinosScalar> &
values,
1207 const bool elide_zero_values =
true);
1227 const bool elide_zero_values =
true,
1228 const bool col_indices_are_sorted =
false);
1308 clear_rows(
const std::vector<size_type> &rows,
1404 template <
typename VectorType>
1405 typename std::enable_if<std::is_same<
typename VectorType::value_type,
1407 vmult(VectorType &dst,
const VectorType &src)
const;
1415 template <
typename VectorType>
1416 typename std::enable_if<!std::is_same<
typename VectorType::value_type,
1418 vmult(VectorType &dst,
const VectorType &src)
const;
1434 template <
typename VectorType>
1435 typename std::enable_if<std::is_same<
typename VectorType::value_type,
1437 Tvmult(VectorType &dst,
const VectorType &src)
const;
1445 template <
typename VectorType>
1446 typename std::enable_if<!std::is_same<
typename VectorType::value_type,
1448 Tvmult(VectorType &dst,
const VectorType &src)
const;
1459 template <
typename VectorType>
1461 vmult_add(VectorType &dst,
const VectorType &src)
const;
1473 template <
typename VectorType>
1475 Tvmult_add(VectorType &dst,
const VectorType &src)
const;
1629 const Epetra_CrsMatrix &
1636 const Epetra_CrsGraph &
1785 print(std::ostream &out,
1786 const bool write_extended_trilinos_info =
false)
const;
1798 <<
"An error with error number " << arg1
1799 <<
" occurred while calling a Trilinos function");
1807 <<
"The entry with index <" << arg1 <<
',' << arg2
1808 <<
"> does not exist.");
1814 "You are attempting an operation on two matrices that "
1815 "are the same object, but the operation requires that the "
1816 "two objects are in fact different.");
1831 <<
"You tried to access element (" << arg1 <<
"/" << arg2
1833 <<
" of a distributed matrix, but only rows in range ["
1834 << arg3 <<
"," << arg4
1835 <<
"] are stored locally and can be accessed.");
1843 <<
"You tried to access element (" << arg1 <<
"/" << arg2
1845 <<
" of a sparse matrix, but it appears to not"
1846 <<
" exist in the Trilinos sparsity pattern.");
1934 const Epetra_MultiVector &src,
1935 const Epetra_MultiVector &dst,
1940 Assert(src.Map().SameAs(mtrx.DomainMap()) ==
true,
1942 "Column map of matrix does not fit with vector map!"));
1943 Assert(dst.Map().SameAs(mtrx.RangeMap()) ==
true,
1944 ExcMessage(
"Row map of matrix does not fit with vector map!"));
1948 Assert(src.Map().SameAs(mtrx.RangeMap()) ==
true,
1950 "Column map of matrix does not fit with vector map!"));
1951 Assert(dst.Map().SameAs(mtrx.DomainMap()) ==
true,
1952 ExcMessage(
"Row map of matrix does not fit with vector map!"));
1961 const Epetra_MultiVector &src,
1962 const Epetra_MultiVector &dst,
1967 Assert(src.Map().SameAs(op.OperatorDomainMap()) ==
true,
1969 "Column map of operator does not fit with vector map!"));
1970 Assert(dst.Map().SameAs(op.OperatorRangeMap()) ==
true,
1972 "Row map of operator does not fit with vector map!"));
1976 Assert(src.Map().SameAs(op.OperatorRangeMap()) ==
true,
1978 "Column map of operator does not fit with vector map!"));
1979 Assert(dst.Map().SameAs(op.OperatorDomainMap()) ==
true,
1981 "Row map of operator does not fit with vector map!"));
1988 namespace LinearOperatorImplementation
2115 template <
typename Solver,
typename Preconditioner>
2116 typename std::enable_if<
2117 std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
2119 Preconditioner>::value,
2140 template <
typename Solver,
typename Preconditioner>
2141 typename std::enable_if<
2142 !(std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
2144 Preconditioner>::value),
2305 virtual const char *
2306 Label()
const override;
2315 virtual const Epetra_Comm &
2316 Comm()
const override;
2325 virtual const Epetra_Map &
2336 virtual const Epetra_Map &
2351# ifdef DEAL_II_WITH_MPI
2424 visit_present_row();
2429 AccessorBase::row()
const
2437 AccessorBase::column()
const
2440 return (*colnum_cache)[a_index];
2445 AccessorBase::index()
const
2452 inline Accessor<true>::Accessor(MatrixType *
matrix,
2459 template <
bool Other>
2460 inline Accessor<true>::Accessor(
const Accessor<Other> &other)
2461 : AccessorBase(other)
2466 Accessor<true>::value()
const
2469 return (*value_cache)[a_index];
2473 inline Accessor<false>::Reference::Reference(
const Accessor<false> &acc)
2474 : accessor(const_cast<Accessor<false> &>(acc))
2478 inline Accessor<false>::Reference::operator
TrilinosScalar()
const
2480 return (*accessor.value_cache)[accessor.a_index];
2483 inline const Accessor<false>::Reference &
2484 Accessor<false>::Reference::operator=(
const TrilinosScalar n)
const
2486 (*accessor.value_cache)[accessor.a_index] = n;
2487 accessor.matrix->set(accessor.row(),
2494 inline const Accessor<false>::Reference &
2495 Accessor<false>::Reference::operator+=(
const TrilinosScalar n)
const
2497 (*accessor.value_cache)[accessor.a_index] += n;
2498 accessor.matrix->set(accessor.row(),
2505 inline const Accessor<false>::Reference &
2506 Accessor<false>::Reference::operator-=(
const TrilinosScalar n)
const
2508 (*accessor.value_cache)[accessor.a_index] -= n;
2509 accessor.matrix->set(accessor.row(),
2516 inline const Accessor<false>::Reference &
2517 Accessor<false>::Reference::operator*=(
const TrilinosScalar n)
const
2519 (*accessor.value_cache)[accessor.a_index] *= n;
2520 accessor.matrix->set(accessor.row(),
2527 inline const Accessor<false>::Reference &
2528 Accessor<false>::Reference::operator/=(
const TrilinosScalar n)
const
2530 (*accessor.value_cache)[accessor.a_index] /= n;
2531 accessor.matrix->set(accessor.row(),
2538 inline Accessor<false>::Accessor(MatrixType *
matrix,
2541 : AccessorBase(
matrix, row, index)
2545 inline Accessor<false>::Reference
2546 Accessor<false>::value()
const
2554 template <
bool Constness>
2555 inline Iterator<Constness>::Iterator(MatrixType *
matrix,
2558 : accessor(
matrix, row, index)
2562 template <
bool Constness>
2563 template <
bool Other>
2564 inline Iterator<Constness>::Iterator(
const Iterator<Other> &other)
2565 : accessor(other.accessor)
2569 template <
bool Constness>
2570 inline Iterator<Constness> &
2581 if (accessor.a_index >= accessor.colnum_cache->size())
2583 accessor.a_index = 0;
2586 while ((accessor.a_row < accessor.matrix->m()) &&
2587 ((accessor.matrix->in_local_range(accessor.a_row) ==
false) ||
2588 (accessor.matrix->row_length(accessor.a_row) == 0)))
2591 accessor.visit_present_row();
2597 template <
bool Constness>
2598 inline Iterator<Constness>
2601 const Iterator<Constness> old_state = *
this;
2608 template <
bool Constness>
2616 template <
bool Constness>
2617 inline const Accessor<Constness> *Iterator<Constness>::operator->()
const
2624 template <
bool Constness>
2628 return (accessor.a_row == other.accessor.a_row &&
2629 accessor.a_index == other.accessor.a_index);
2634 template <
bool Constness>
2638 return !(*
this == other);
2643 template <
bool Constness>
2647 return (accessor.row() < other.accessor.row() ||
2648 (accessor.row() == other.accessor.row() &&
2649 accessor.index() < other.accessor.index()));
2653 template <
bool Constness>
2657 return (other < *
this);
2675 return const_iterator(
this, m(), 0);
2684 if (in_local_range(r) && (row_length(r) > 0))
2685 return const_iterator(
this, r, 0);
2701 if (in_local_range(i) && (row_length(i) > 0))
2702 return const_iterator(
this, i, 0);
2722 return iterator(
this, m(), 0);
2731 if (in_local_range(r) && (row_length(r) > 0))
2732 return iterator(
this, r, 0);
2748 if (in_local_range(i) && (row_length(i) > 0))
2749 return iterator(
this, i, 0);
2762# ifndef DEAL_II_WITH_64BIT_INDICES
2767 end =
matrix->RowMap().MaxMyGID64() + 1;
2789 SparseMatrix::set<TrilinosScalar>(
const size_type row,
2793 const bool elide_zero_values);
2797 template <
typename Number>
2803 const bool elide_zero_values)
2805 std::vector<TrilinosScalar> trilinos_values(n_cols);
2808 row, n_cols, col_indices, trilinos_values.data(), elide_zero_values);
2820 set(i, 1, &j, &value,
false);
2828 const bool elide_zero_values)
2834 for (
size_type i = 0; i < indices.size(); ++i)
2859 if (last_action == Insert)
2862 ierr =
matrix->GlobalAssemble(*column_space_map,
2875 add(i, 1, &j, &value,
false);
2885# ifndef DEAL_II_WITH_64BIT_INDICES
2886 return matrix->NumGlobalRows();
2888 return matrix->NumGlobalRows64();
2909 return matrix->NumMyRows();
2914 inline std::pair<SparseMatrix::size_type, SparseMatrix::size_type>
2918# ifndef DEAL_II_WITH_64BIT_INDICES
2923 end =
matrix->RowMap().MaxMyGID64() + 1;
2934# ifndef DEAL_II_WITH_64BIT_INDICES
2935 return matrix->NumGlobalNonzeros();
2937 return matrix->NumGlobalNonzeros64();
2943 template <
typename SparsityPatternType>
2946 const SparsityPatternType &sparsity_pattern,
2948 const bool exchange_data)
2950 reinit(parallel_partitioning,
2951 parallel_partitioning,
2959 template <
typename number>
2962 const ::SparseMatrix<number> &sparse_matrix,
2964 const double drop_tolerance,
2965 const bool copy_values,
2966 const ::SparsityPattern * use_this_sparsity)
2970 reinit(parallel_partitioning,
2971 parallel_partitioning,
2980 inline const Epetra_CrsMatrix &
2983 return static_cast<const Epetra_CrsMatrix &
>(*matrix);
2988 inline const Epetra_CrsGraph &
3029 namespace LinearOperatorImplementation
3031 template <
typename Solver,
typename Preconditioner>
3032 typename std::enable_if<
3033 std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
3035 Preconditioner>::value,
3036 TrilinosPayload>::type
3037 TrilinosPayload::inverse_payload(
3039 const Preconditioner &preconditioner)
const
3041 const auto &payload = *
this;
3043 TrilinosPayload return_op(payload);
3047 return_op.inv_vmult = [payload, &solver, &preconditioner](
3048 TrilinosPayload::Domain & tril_dst,
3049 const TrilinosPayload::Range &tril_src) {
3052 Assert(&tril_src != &tril_dst,
3057 !payload.UseTranspose());
3058 solver.solve(payload, tril_dst, tril_src, preconditioner);
3061 return_op.inv_Tvmult = [payload, &solver, &preconditioner](
3062 TrilinosPayload::Range & tril_dst,
3063 const TrilinosPayload::Domain &tril_src) {
3066 Assert(&tril_src != &tril_dst,
3071 payload.UseTranspose());
3073 const_cast<TrilinosPayload &
>(payload).
transpose();
3074 solver.solve(payload, tril_dst, tril_src, preconditioner);
3075 const_cast<TrilinosPayload &
>(payload).
transpose();
3080 if (return_op.UseTranspose() ==
true)
3081 std::swap(return_op.inv_vmult, return_op.inv_Tvmult);
3086 template <
typename Solver,
typename Preconditioner>
3087 typename std::enable_if<
3088 !(std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
3090 Preconditioner>::value),
3091 TrilinosPayload>::type
3092 TrilinosPayload::inverse_payload(Solver &,
const Preconditioner &)
const
3094 TrilinosPayload return_op(*
this);
3096 return_op.inv_vmult = [](TrilinosPayload::Domain &,
3097 const TrilinosPayload::Range &) {
3099 ExcMessage(
"Payload inv_vmult disabled because of "
3100 "incompatible solver/preconditioner choice."));
3103 return_op.inv_Tvmult = [](TrilinosPayload::Range &,
3104 const TrilinosPayload::Domain &) {
3106 ExcMessage(
"Payload inv_vmult disabled because of "
3107 "incompatible solver/preconditioner choice."));
3117 SparseMatrix::set<TrilinosScalar>(
const size_type row,
3121 const bool elide_zero_values);
bool operator!=(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
bool operator==(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
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
Reference(const Accessor< false > &accessor)
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(const Accessor< Other > &a)
Accessor(MatrixType *matrix, const size_type row, const size_type index)
TrilinosScalar value() const
TrilinosScalar value() const
const Accessor< Constness > & operator*() const
const Accessor< Constness > * operator->() const
bool operator>(const Iterator< Constness > &) const
typename Accessor< Constness >::MatrixType MatrixType
Iterator(const Iterator< Other > &other)
::types::global_dof_index size_type
bool operator==(const Iterator< Constness > &) const
bool operator!=(const Iterator< Constness > &) const
Iterator(MatrixType *matrix, const size_type row, const size_type index)
Accessor< Constness > accessor
Iterator< Constness > & operator++()
bool operator<(const Iterator< Constness > &) const
Iterator< Constness > operator++(int)
void set(const size_type i, const size_type j, const TrilinosScalar value)
std::enable_if<!std::is_same< typenameVectorType::value_type, TrilinosScalar >::value >::type vmult(VectorType &dst, const VectorType &src) const
std::unique_ptr< Epetra_Map > column_space_map
std::enable_if< std::is_same< typenameVectorType::value_type, TrilinosScalar >::value >::type vmult(VectorType &dst, const VectorType &src) const
MPI_Comm get_mpi_communicator() const
SparseMatrix & operator*=(const TrilinosScalar factor)
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
std::unique_ptr< Epetra_Export > nonlocal_matrix_exporter
TrilinosScalar l1_norm() const
std::unique_ptr< Epetra_FECrsMatrix > matrix
std::enable_if<!std::is_same< typenameVectorType::value_type, TrilinosScalar >::value >::type Tvmult(VectorType &dst, const VectorType &src) const
::types::global_dof_index size_type
TrilinosScalar linfty_norm() const
const Epetra_CrsMatrix & trilinos_matrix() const
size_type memory_consumption() const
TrilinosScalar matrix_norm_square(const MPI::Vector &v) const
Epetra_CombineMode last_action
void reinit(const IndexSet ¶llel_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)
IndexSet locally_owned_range_indices() const
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
void clear_row(const size_type row, const TrilinosScalar new_diag_value=0)
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
void compress(::VectorOperation::values operation)
SparseMatrix & operator/=(const TrilinosScalar factor)
void Tvmult_add(VectorType &dst, const VectorType &src) const
void clear_rows(const std::vector< size_type > &rows, const TrilinosScalar new_diag_value=0)
bool is_compressed() const
TrilinosScalar el(const size_type i, const size_type j) const
IndexSet locally_owned_domain_indices() const
TrilinosScalar residual(MPI::Vector &dst, const MPI::Vector &x, const MPI::Vector &b) const
bool in_local_range(const size_type index) const
void copy_from(const SparseMatrix &source)
unsigned int row_length(const size_type row) const
TrilinosScalar frobenius_norm() const
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)
void reinit(const IndexSet ¶llel_partitioning, const SparsityPatternType &sparsity_pattern, const MPI_Comm &communicator=MPI_COMM_WORLD, const bool exchange_data=false)
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)
std::enable_if< std::is_same< typenameVectorType::value_type, TrilinosScalar >::value >::type Tvmult(VectorType &dst, const VectorType &src) const
TrilinosScalar operator()(const size_type i, const size_type j) const
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)
size_type n_nonzero_elements() const
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
std::unique_ptr< Epetra_CrsMatrix > nonlocal_matrix
TrilinosScalar value_type
Epetra_MpiComm communicator
std::function< void(VectorType &, const VectorType &)> inv_Tvmult
virtual int SetUseTranspose(bool UseTranspose) override
virtual bool UseTranspose() const override
virtual const Epetra_Map & OperatorDomainMap() const override
IndexSet locally_owned_range_indices() const
TrilinosPayload transpose_payload() const
Epetra_MultiVector VectorType
MPI_Comm get_mpi_communicator() const
std::enable_if< std::is_base_of< TrilinosWrappers::SolverBase, Solver >::value &&std::is_base_of< TrilinosWrappers::PreconditionBase, Preconditioner >::value, TrilinosPayload >::type inverse_payload(Solver &, const Preconditioner &) const
virtual ~TrilinosPayload() override=default
virtual int ApplyInverse(const VectorType &Y, VectorType &X) const override
std::function< void(VectorType &, const VectorType &)> Tvmult
virtual int Apply(const VectorType &X, VectorType &Y) const override
std::function< void(VectorType &, const VectorType &)> vmult
virtual const Epetra_Map & OperatorRangeMap() const override
std::function< void(VectorType &, const VectorType &)> inv_vmult
virtual const char * Label() const override
std::enable_if<!(std::is_base_of< TrilinosWrappers::SolverBase, Solver >::value &&std::is_base_of< TrilinosWrappers::PreconditionBase, Preconditioner >::value), TrilinosPayload >::type inverse_payload(Solver &, const Preconditioner &) const
IndexSet locally_owned_domain_indices() const
virtual const Epetra_Comm & Comm() const override
virtual bool HasNormInf() const override
TrilinosPayload identity_payload() const
virtual double NormInf() const override
TrilinosPayload null_payload() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
__global__ void set(Number *val, const Number s, const size_type N)
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)
Expression operator>(const Expression &lhs, const Expression &rhs)
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)
@ matrix
Contents is actually a matrix.
SparseMatrix< double > SparseMatrix
types::global_dof_index size_type
void swap(MemorySpaceData< Number, MemorySpace > &, MemorySpaceData< Number, MemorySpace > &)
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
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::int_type n_global_elements(const Epetra_BlockMap &map)
void copy(const T *begin, const T *end, U *dest)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
unsigned int global_dof_index
static const bool zero_addition_can_be_elided
SynchronousIterators< Iterators > operator++(SynchronousIterators< Iterators > &a)
bool operator<(const SynchronousIterators< Iterators > &a, const SynchronousIterators< Iterators > &b)