43 #ifndef IFPACK2_BLOCKTRIDICONTAINER_DECL_HPP
44 #define IFPACK2_BLOCKTRIDICONTAINER_DECL_HPP
49 #include "Ifpack2_config.h"
50 #include "Ifpack2_Container.hpp"
51 #include "Tpetra_MultiVector.hpp"
52 #include "Tpetra_Map.hpp"
53 #include "Tpetra_RowMatrix.hpp"
54 #include "Tpetra_BlockCrsMatrix_decl.hpp"
55 #include <type_traits>
111 namespace BlockTriDiContainerDetails {
116 struct ImplSimdTag {};
117 struct ImplSacadoTag {};
119 template<
typename T>
struct ImplTag {
typedef ImplNotAvailTag type; };
120 template<>
struct ImplTag<float> {
typedef ImplSimdTag type; };
121 template<>
struct ImplTag<double> {
typedef ImplSimdTag type; };
122 template<>
struct ImplTag<std::complex<float> > {
typedef ImplSimdTag type; };
123 template<>
struct ImplTag<std::complex<double> > {
typedef ImplSimdTag type; };
132 template <
typename MatrixType,
141 template <
typename MatrixType>
152 typedef MatrixType matrix_type;
155 typedef typename MatrixType::scalar_type scalar_type;
157 typedef typename Kokkos::ArithTraits<scalar_type>::magnitudeType magnitude_type;
159 typedef typename Container<MatrixType>::local_ordinal_type local_ordinal_type;
161 typedef typename Container<MatrixType>::global_ordinal_type global_ordinal_type;
163 typedef typename Container<MatrixType>::node_type node_type;
165 typedef typename Container<MatrixType>::mv_type mv_type;
166 typedef typename Container<MatrixType>::map_type map_type;
167 typedef typename Container<MatrixType>::vector_type vector_type;
168 typedef typename Container<MatrixType>::import_type import_type;
171 typedef host_view_type HostView;
175 typedef Tpetra::BlockCrsMatrix
176 <scalar_type,local_ordinal_type,global_ordinal_type,node_type> block_crs_matrix_type;
184 typedef typename Container<MatrixType>::row_matrix_type row_matrix_type;
186 static_assert (std::is_same<MatrixType, row_matrix_type>::value,
187 "Ifpack2::BlockTriDiContainer: MatrixType must be a Tpetra::RowMatrix specialization.");
229 bool overlapCommAndComp =
false,
bool useSequentialMethod =
false);
234 struct ComputeParameters {
241 magnitude_type addRadiallyToDiagonal = Kokkos::ArithTraits<magnitude_type>::zero();
245 struct ApplyParameters {
248 bool zeroStartingSolution =
false;
250 scalar_type dampingFactor = Kokkos::ArithTraits<scalar_type>::one();
253 int maxNumSweeps = 1;
263 magnitude_type tolerance = Kokkos::ArithTraits<magnitude_type>::zero();
271 int checkToleranceEvery = 1;
281 void clearBlocks()
override;
288 void initialize ()
override;
291 void compute ()
override;
294 void applyInverseJacobi (
const mv_type& X, mv_type& Y,
295 scalar_type dampingFactor,
296 bool zeroStartingSolution =
false,
297 int numSweeps = 1)
const override;
300 ComputeParameters createDefaultComputeParameters ()
const;
313 void compute (
const ComputeParameters& input);
316 ApplyParameters createDefaultApplyParameters ()
const;
324 int applyInverseJacobi (
const mv_type& X, mv_type& Y,
325 const ApplyParameters& input)
const;
330 const magnitude_type getNorms0 ()
const;
334 const magnitude_type getNormsFinal ()
const;
339 apply (host_view_type X,
349 weightedApply (host_view_type X,
364 std::ostream& print (std::ostream& os)
const override;
371 std::string description ()
const override;
382 static std::string getName();
395 const bool overlapCommAndComp,
396 const bool useSeqMethod);
398 void clearInternal();
408 template <
typename MatrixType>
412 typedef typename MatrixType::scalar_type scalar_type;
413 typedef typename Kokkos::ArithTraits<scalar_type>::magnitudeType magnitude_type;
414 typedef typename Container<MatrixType>::local_ordinal_type local_ordinal_type;
415 typedef typename Container<MatrixType>::global_ordinal_type global_ordinal_type;
417 typedef typename Container<MatrixType>::mv_type mv_type;
418 typedef typename Container<MatrixType>::import_type import_type;
421 typedef typename Container<MatrixType>::row_matrix_type row_matrix_type;
423 static_assert (std::is_same<MatrixType, row_matrix_type>::value,
424 "Ifpack2::BlockTriDiContainer: MatrixType must be a Tpetra::RowMatrix specialization.");
436 void clearBlocks()
override {}
440 void applyInverseJacobi (
const mv_type& X, mv_type& Y,
441 scalar_type dampingFactor,
442 bool zeroStartingSolution =
false,
443 int numSweeps = 1)
const override {}
446 apply (host_view_type X,
454 weightedApply (host_view_type X,
462 std::ostream&
print (std::ostream& os)
const override {
463 return os <<
"Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
466 std::string description ()
const override {
467 return "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
474 out <<
"Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
477 static std::string getName() {
478 return "Ifpack2::BlockTriDiContainer::ImplNotAvailTag";
485 #endif // IFPACK2_BLOCKTRIDICONTAINER_DECL_HPP
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:115
void compute() override
Extract the local diagonal blocks and prepare the solver.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:439
Store and solve local block tridiagonal linear problems.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:134
void initialize() override
Do all set-up operations that only require matrix structure.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:438
#define TEUCHOS_TEST_FOR_EXCEPT_MSG(throw_exception_test, msg)
forward declaration
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:126
typename mv_type::dual_view_type::t_host HostView
Definition: Ifpack2_Container_decl.hpp:139
std::ostream & print(std::ostream &os) const override
Print basic information about the container to os.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:462
void setParameters(const Teuchos::ParameterList &List) override
Set parameters, if any.
Definition: Ifpack2_BlockTriDiContainer_decl.hpp:435
static const EVerbosityLevel verbLevel_default
Interface for creating and solving a set of local linear problems.
Definition: Ifpack2_Container_decl.hpp:112