46 #ifndef IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP 47 #define IFPACK2_EXPERIMENTALCRSRBILUK_DECL_HPP 49 #include <Tpetra_Experimental_BlockCrsMatrix.hpp> 51 #include <Ifpack2_RILUK.hpp> 57 template <
class scalar_type,
class impl_scalar_type>
58 struct BlockMatrixOperations
61 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType magnitude_type;
63 typedef Teuchos::ScalarTraits<magnitude_type> STM;
65 void square_matrix_matrix_multiply(
const impl_scalar_type * a,
const impl_scalar_type * b, impl_scalar_type * c,
66 const int nrows,
const impl_scalar_type alpha = STM::one(),
const impl_scalar_type beta = STM::zero() )
const 68 for (
int i = 0; i < nrows*nrows; ++i)
71 for (
int i = 0; i < nrows; ++i)
74 const int ioffset = i*nrows;
75 for (
int k = 0; k < nrows; ++k)
77 const int koffset = k*nrows;
78 const impl_scalar_type val = alpha*a[ioffset+k];
79 for (
int j = 0; j < nrows; ++j)
81 c[ioffset+j] += val*b[koffset+j];
161 template<
class MatrixType>
163 typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type> >
169 typedef typename MatrixType::scalar_type
scalar_type;
173 typedef typename MatrixType::scalar_type impl_scalar_type;
185 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType
magnitude_type;
188 typedef Tpetra::RowMatrix<scalar_type,
194 typedef Tpetra::CrsMatrix<scalar_type,
199 typedef Tpetra::Experimental::BlockCrsMatrix<scalar_type,
202 node_type> block_crs_matrix_type;
204 template <
class NewMatrixType>
friend class RBILUK;
209 RBILUK (
const Teuchos::RCP<const row_matrix_type>& A_in);
217 RBILUK (
const Teuchos::RCP<const block_crs_matrix_type>& A_in);
253 using RILUK<Tpetra::RowMatrix<
typename MatrixType::scalar_type,
254 typename MatrixType::local_ordinal_type,
255 typename MatrixType::global_ordinal_type,
256 typename MatrixType::node_type> >::setMatrix;
281 setMatrix (
const Teuchos::RCP<const block_crs_matrix_type>& A);
288 std::string description ()
const;
324 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
325 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
326 Teuchos::ETransp mode = Teuchos::NO_TRANS,
327 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one (),
328 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero ())
const;
334 Teuchos::RCP<const block_crs_matrix_type> getBlockMatrix ()
const;
337 const block_crs_matrix_type& getLBlock ()
const;
340 const block_crs_matrix_type& getDBlock ()
const;
343 const block_crs_matrix_type& getUBlock ()
const;
346 typedef Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> MV;
347 typedef Teuchos::ScalarTraits<impl_scalar_type> STS;
348 typedef Teuchos::ScalarTraits<magnitude_type> STM;
349 typedef typename block_crs_matrix_type::little_block_type little_block_type;
350 typedef typename block_crs_matrix_type::little_vec_type little_vec_type;
352 void allocate_L_and_U_blocks();
353 void initAllValues (
const block_crs_matrix_type& A);
356 Teuchos::RCP<const row_matrix_type> A_;
359 Teuchos::RCP<const block_crs_matrix_type> A_block_;
362 local_ordinal_type blockSize_;
365 Teuchos::RCP<block_crs_matrix_type> L_block_;
367 Teuchos::RCP<block_crs_matrix_type> U_block_;
369 Teuchos::RCP<block_crs_matrix_type> D_block_;
372 Teuchos::RCP<block_crs_matrix_type> D_block_inverse_;
374 BlockMatrixOperations<scalar_type,impl_scalar_type> blockMatOpts;
376 void square_matrix_matrix_multiply(
const impl_scalar_type * a,
const impl_scalar_type * b, impl_scalar_type * c,
377 const int nrows,
const impl_scalar_type alpha = STM::one(),
const impl_scalar_type beta = STM::zero() )
const 379 for (
int i = 0; i < nrows*nrows; ++i)
382 for (
int i = 0; i < nrows; ++i)
385 const int ioffset = i*nrows;
386 for (
int k = 0; k < nrows; ++k)
388 const int koffset = k*nrows;
389 const impl_scalar_type val = alpha*a[ioffset+k];
390 for (
int j = 0; j < nrows; ++j)
392 c[ioffset+j] += val*b[koffset+j];
Ifpack2 features that are experimental. Use at your own risk.
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:179
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:170
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:185
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition: Ifpack2_RILUK_decl.hpp:243
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:182
Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > crs_matrix_type
Tpetra::CrsMatrix specialization used by this class for representing L and U.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:197
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:176
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization used by this class.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:191
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:72
ILU(k) factorization of a given Tpetra::Experimental::BlockCrsMatrix.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:162