10 #ifndef IFPACK2_DENSECONTAINER_DECL_HPP 
   11 #define IFPACK2_DENSECONTAINER_DECL_HPP 
   16 #include "Ifpack2_Container.hpp" 
   17 #include "Tpetra_MultiVector.hpp" 
   18 #include "Tpetra_Map.hpp" 
   19 #include "Tpetra_RowMatrix.hpp" 
   21 #include "Tpetra_BlockCrsMatrix_decl.hpp" 
   69 template <
class MatrixType, 
class LocalScalarType>
 
   81   using matrix_type = MatrixType;
 
   83   using LSC = LocalScalarType;
 
   86   using typename Container<MatrixType>::SC;
 
   88   using typename Container<MatrixType>::LO;
 
   90   using typename Container<MatrixType>::GO;
 
   92   using typename Container<MatrixType>::NO;
 
   94   using typename Container<MatrixType>::mv_type;
 
   95   using local_mv_type = Tpetra::MultiVector<LSC, LO, GO, NO>;
 
   96   using typename Container<MatrixType>::map_type;
 
   97   using typename Container<MatrixType>::vector_type;
 
   98   using typename Container<MatrixType>::import_type;
 
  101   using typename ContainerImpl<MatrixType, LocalScalarType>::LISC;
 
  103   using typename ContainerImpl<MatrixType, LocalScalarType>::HostViewLocal;
 
  104   using typename ContainerImpl<MatrixType, LocalScalarType>::HostSubviewLocal;
 
  105   using typename ContainerImpl<MatrixType, LocalScalarType>::ConstHostSubviewLocal;
 
  107   static_assert(std::is_same<MatrixType, Tpetra::RowMatrix<SC, LO, GO, NO>>::value,
 
  108                 "Ifpack2::DenseContainer: Please use MatrixType = Tpetra::RowMatrix.");
 
  118   using typename Container<MatrixType>::row_matrix_type;
 
  120   using block_crs_matrix_type = Tpetra::BlockCrsMatrix<SC, LO, GO, NO>;
 
  170   virtual std::ostream& 
print(std::ostream& os) 
const;
 
  208   solveBlock(ConstHostSubviewLocal X,
 
  213              const LSC beta) 
const;
 
  216   std::vector<Teuchos::SerialDenseMatrix<int, LSC>> diagBlocks_;
 
  222   bool hasBlockCrsMatrix_;
 
  233 #endif  // IFPACK2_DENSECONTAINER_DECL_HPP 
typename Kokkos::ArithTraits< SC >::val_type ISC
Internal representation of Scalar in Kokkos::View. 
Definition: Ifpack2_Container_decl.hpp:104
The implementation of the numerical features of Container (Jacobi, Gauss-Seidel, SGS). This class allows a custom scalar type (LocalScalarType) to be used for storing blocks and solving the block systems. Hiding this template parameter from the Container interface simplifies the BlockRelaxation and ContainerFactory classes. 
Definition: Ifpack2_Container_decl.hpp:311
virtual void compute()
Extract the local diagonal block and prepare the solver. 
Definition: Ifpack2_DenseContainer_def.hpp:76
static std::string getName()
Get the name of this container type for Details::constructContainer() 
Definition: Ifpack2_DenseContainer_def.hpp:366
Store and solve a local dense linear problem. 
Definition: Ifpack2_DenseContainer_decl.hpp:70
DenseContainer(const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< LO >> &partitions, const Teuchos::RCP< const import_type > &importer, bool pointIndexed)
Constructor. 
virtual std::ostream & print(std::ostream &os) const 
Print information about this object to the given output stream. 
Definition: Ifpack2_DenseContainer_def.hpp:314
typename mv_type::dual_view_type::t_host HostView
Definition: Ifpack2_Container_decl.hpp:109
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const 
Print the object with some verbosity level to the given FancyOStream. 
Definition: Ifpack2_DenseContainer_def.hpp:343
virtual ~DenseContainer()
Destructor (declared virtual for memory safety of derived classes). 
Definition: Ifpack2_DenseContainer_def.hpp:58
static const EVerbosityLevel verbLevel_default
virtual std::string description() const 
A one-line description of this object. 
Definition: Ifpack2_DenseContainer_def.hpp:324
virtual void initialize()
Do all set-up operations that only require matrix structure. 
Definition: Ifpack2_DenseContainer_def.hpp:62