42 #ifndef IFPACK2_TRIDICONTAINER_DECL_HPP
43 #define IFPACK2_TRIDICONTAINER_DECL_HPP
50 #include "Ifpack2_Details_LapackSupportsScalar.hpp"
51 #include "Tpetra_MultiVector.hpp"
52 #include "Tpetra_Map.hpp"
53 #include "Tpetra_RowMatrix.hpp"
56 #include <type_traits>
105 template<
class MatrixType,
class LocalScalarType,
109 template<
class MatrixType,
class LocalScalarType>
120 typedef MatrixType matrix_type;
122 typedef LocalScalarType local_scalar_type;
125 typedef typename Container<MatrixType>::scalar_type scalar_type;
127 typedef typename Container<MatrixType>::local_ordinal_type local_ordinal_type;
129 typedef typename Container<MatrixType>::global_ordinal_type global_ordinal_type;
131 typedef typename Container<MatrixType>::node_type node_type;
133 typedef typename Container<MatrixType>::mv_type mv_type;
134 typedef typename Container<MatrixType>::map_type map_type;
135 typedef typename Container<MatrixType>::vector_type vector_type;
137 typedef typename Container<MatrixType>::import_type import_type;
139 typedef typename Container<MatrixType>::HostView HostView;
140 typedef Tpetra::MultiVector<local_scalar_type, local_ordinal_type, global_ordinal_type, node_type> local_mv_type;
141 typedef typename Kokkos::View<local_scalar_type**, Kokkos::HostSpace> HostViewLocal;
143 static_assert (std::is_same<MatrixType, Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type>>::value,
144 "Ifpack2::TriDiContainer: MatrixType must be a Tpetra::RowMatrix specialization.");
154 typedef typename Container<MatrixType>::row_matrix_type row_matrix_type;
177 scalar_type DampingFactor);
190 virtual bool isInitialized()
const;
193 virtual bool isComputed()
const;
203 virtual void initialize ();
206 virtual void compute ();
220 void applyImpl (HostViewLocal& X,
225 local_scalar_type alpha,
226 local_scalar_type beta)
const;
230 weightedApply (HostView& X,
246 virtual std::ostream& print (std::ostream& os)
const;
253 virtual std::string description ()
const;
264 static std::string getName();
278 std::vector<Teuchos::SerialTriDiMatrix<int, local_scalar_type>> diagBlocks_;
281 mutable std::vector<HostViewLocal> X_local;
284 mutable std::vector<HostViewLocal> Y_local;
296 local_scalar_type* scalars_;
302 template<
class MatrixType,
class LocalScalarType>
313 typedef MatrixType matrix_type;
315 typedef LocalScalarType local_scalar_type;
318 typedef typename Container<MatrixType>::scalar_type scalar_type;
320 typedef typename Container<MatrixType>::local_ordinal_type local_ordinal_type;
322 typedef typename Container<MatrixType>::global_ordinal_type global_ordinal_type;
324 typedef typename Container<MatrixType>::node_type node_type;
326 typedef typename Container<MatrixType>::mv_type mv_type;
327 typedef typename Container<MatrixType>::map_type map_type;
328 typedef typename Container<MatrixType>::vector_type vector_type;
330 typedef typename Container<MatrixType>::import_type import_type;
332 typedef typename Container<MatrixType>::HostView HostView;
333 typedef Tpetra::MultiVector<local_scalar_type, local_ordinal_type, global_ordinal_type, node_type> local_mv_type;
334 typedef typename Kokkos::View<local_scalar_type**, Kokkos::HostSpace> HostViewLocal;
336 static_assert (std::is_same<MatrixType, Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type>>::value,
337 "Ifpack2::TriDiContainer: MatrixType must be a Tpetra::RowMatrix specialization.");
347 typedef typename Container<MatrixType>::row_matrix_type row_matrix_type;
370 scalar_type DampingFactor);
383 virtual bool isInitialized()
const;
386 virtual bool isComputed()
const;
396 virtual void initialize ();
399 virtual void compute ();
413 void applyImpl (HostViewLocal& X,
418 local_scalar_type alpha,
419 local_scalar_type beta)
const;
423 weightedApply (HostView& X,
439 virtual std::ostream& print (std::ostream& os)
const;
446 virtual std::string description ()
const;
457 static std::string getName();
471 std::vector<Teuchos::SerialTriDiMatrix<int, local_scalar_type>> diagBlocks_;
474 mutable std::vector<HostViewLocal> X_local;
477 mutable std::vector<HostViewLocal> Y_local;
489 local_scalar_type* scalars_;
497 #endif // IFPACK2_TRIDICONTAINER_DECL_HPP
Store and solve a local TriDi linear problem.
Definition: Ifpack2_TriDiContainer_decl.hpp:107
Ifpack2::Container class declaration.
Declaration and definition of the Ifpack2::Details::MultiVectorLocalGatherScatter class...
Ifpack2::Partitioner:
Definition: Ifpack2_Partitioner.hpp:179
static const EVerbosityLevel verbLevel_default
Interface for creating and solving a local linear problem.
Definition: Ifpack2_Container.hpp:114
Type traits class that says whether Teuchos::LAPACK has a valid implementation for the given ScalarTy...
Definition: Ifpack2_Details_LapackSupportsScalar.hpp:17