10#ifndef IFPACK2_BLOCKTRIDICONTAINER_DEF_HPP
11#define IFPACK2_BLOCKTRIDICONTAINER_DEF_HPP
13#include <Teuchos_Details_MpiTypeTraits.hpp>
15#include <Tpetra_Distributor.hpp>
16#include <Tpetra_BlockMultiVector.hpp>
17#include <Tpetra_BlockCrsMatrix_Helpers.hpp>
19#include <Kokkos_ArithTraits.hpp>
20#include <KokkosBatched_Util.hpp>
21#include <KokkosBatched_Vector.hpp>
22#include <KokkosBatched_AddRadial_Decl.hpp>
23#include <KokkosBatched_AddRadial_Impl.hpp>
24#include <KokkosBatched_Gemm_Decl.hpp>
25#include <KokkosBatched_Gemm_Serial_Impl.hpp>
26#include <KokkosBatched_Gemv_Decl.hpp>
27#include <KokkosBatched_Trsm_Decl.hpp>
28#include <KokkosBatched_Trsm_Serial_Impl.hpp>
29#include <KokkosBatched_Trsv_Decl.hpp>
30#include <KokkosBatched_Trsv_Serial_Impl.hpp>
31#include <KokkosBatched_LU_Decl.hpp>
32#include <KokkosBatched_LU_Serial_Impl.hpp>
35#include "Ifpack2_BlockTriDiContainer_impl.hpp"
46 template <
typename MatrixType>
50 const Teuchos::RCP<const import_type>& importer,
51 const bool overlapCommAndComp,
52 const bool useSeqMethod,
54 const bool explicitConversion)
56 IFPACK2_BLOCKHELPER_TIMER_WITH_FENCE(
"BlockTriDiContainer::initInternal", initInternal,
typename BlockHelperDetails::ImplType<MatrixType>::execution_space);
60 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::createImpl", createImpl);
62 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
69 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::setA", setA);
70 if (explicitConversion) {
71 impl_->A = Teuchos::rcp_dynamic_cast<const block_crs_matrix_type>(matrix);
72 if (impl_->A.is_null()) {
73 TEUCHOS_TEST_FOR_EXCEPT_MSG
74 (block_size == -1,
"A pointwise matrix and block_size = -1 were given as inputs.");
76 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::setA::convertToBlockCrsMatrix", convertToBlockCrsMatrix);
77 impl_->A = Tpetra::convertToBlockCrsMatrix(*Teuchos::rcp_dynamic_cast<const crs_matrix_type>(matrix), block_size,
true);
78 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
85 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
88 impl_->tpetra_importer = Teuchos::null;
89 impl_->async_importer = Teuchos::null;
93 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::createBlockCrsTpetraImporter useSeqMethod", useSeqMethod);
94 if (importer.is_null())
97 impl_->tpetra_importer = importer;
98 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
102 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::createBlockCrsTpetraImporter", createBlockCrsTpetraImporter);
106 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
115 impl_->overlap_communication_and_computation = overlapCommAndComp;
118 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::createZ", createZ);
119 impl_->Z =
typename impl_type::tpetra_multivector_type();
120 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
123 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::createW", createW);
124 impl_->W =
typename impl_type::impl_scalar_type_1d_view();
125 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
128 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
131 template <
typename MatrixType>
136 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::clearInternal", clearInternal);
138 using part_interface_type = BlockHelperDetails::PartInterface<MatrixType>;
143 impl_->A = Teuchos::null;
144 impl_->tpetra_importer = Teuchos::null;
145 impl_->async_importer = Teuchos::null;
147 impl_->Z =
typename impl_type::tpetra_multivector_type();
148 impl_->W =
typename impl_type::impl_scalar_type_1d_view();
150 impl_->part_interface = part_interface_type();
151 impl_->block_tridiags = block_tridiags_type();
152 impl_->a_minus_d = amd_type();
153 impl_->work =
typename impl_type::vector_type_1d_view();
154 impl_->norm_manager = norm_manager_type();
156 impl_ = Teuchos::null;
157 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
160 template <
typename MatrixType>
163 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
164 const Teuchos::RCP<const import_type>& importer,
166 :
Container<MatrixType>(matrix, partitions, pointIndexed), partitions_(partitions)
168 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::BlockTriDiContainer",
BlockTriDiContainer);
169 const bool useSeqMethod =
false;
170 const bool overlapCommAndComp =
false;
171 initInternal(matrix, importer, overlapCommAndComp, useSeqMethod);
172 n_subparts_per_part_ = -1;
174 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
177 template <
typename MatrixType>
180 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
181 const int n_subparts_per_part,
182 const bool overlapCommAndComp,
183 const bool useSeqMethod,
184 const int block_size,
185 const bool explicitConversion)
186 :
Container<MatrixType>(matrix, partitions, false), partitions_(partitions)
188 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::BlockTriDiContainer",
BlockTriDiContainer);
189 initInternal(matrix, Teuchos::null, overlapCommAndComp, useSeqMethod, block_size, explicitConversion);
190 n_subparts_per_part_ = n_subparts_per_part;
191 block_size_ = block_size;
192 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
195 template <
typename MatrixType>
201 template <
typename MatrixType>
206 if (List.isType<
int>(
"partitioner: subparts per part"))
207 n_subparts_per_part_ = List.get<
int>(
"partitioner: subparts per part");
208 if (List.isType<
int>(
"partitioner: block size"))
209 block_size_ = List.get<
int>(
"partitioner: block size");
212 template <
typename MatrixType>
217 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::initialize",
initialize);
220 auto bA = Teuchos::rcp_dynamic_cast<const block_crs_matrix_type>(impl_->A);
222 TEUCHOS_TEST_FOR_EXCEPT_MSG
223 (block_size_ == -1,
"A pointwise matrix and block_size = -1 were given as inputs.");
225 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::initialize::getBlockCrsGraph", getBlockCrsGraph);
226 auto A = Teuchos::rcp_dynamic_cast<const crs_matrix_type>(impl_->A);
227 impl_->blockGraph = Tpetra::getBlockCrsGraph(*A, block_size_,
true);
228 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
232 impl_->blockGraph = Teuchos::rcpFromRef(bA->getCrsGraph());
237 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::createPartInterfaceBlockTridiagsNormManager", createPartInterfaceBlockTridiagsNormManager);
241 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
247 TEUCHOS_ASSERT(!impl_->A.is_null());
249 bool useSeqMethod = !impl_->tpetra_importer.is_null();
253 impl_->part_interface,
254 impl_->block_tridiags,
256 impl_->overlap_communication_and_computation,
257 impl_->async_importer, useSeqMethod);
259 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
262 template <
typename MatrixType>
267 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::compute",
compute);
275 impl_->part_interface, impl_->block_tridiags,
276 Kokkos::ArithTraits<magnitude_type>::zero());
279 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
282 template <
typename MatrixType>
287 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::clearBlocks", clearBlocks);
289 this->IsInitialized_ =
false;
290 this->IsComputed_ =
false;
291 Container<MatrixType>::clearBlocks();
292 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
295 template <
typename MatrixType>
297 BlockTriDiContainer<MatrixType, BlockTriDiContainerDetails::ImplSimdTag>
298 ::applyInverseJacobi (
const mv_type& X, mv_type& Y, scalar_type dampingFactor,
299 bool zeroStartingSolution,
int numSweeps)
const
301 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::applyInverseJacobi",
applyInverseJacobi);
302 const magnitude_type tol = Kokkos::ArithTraits<magnitude_type>::zero();
303 const int check_tol_every = 1;
308 impl_->tpetra_importer,
309 impl_->async_importer,
310 impl_->overlap_communication_and_computation,
311 X, Y, impl_->Z, impl_->W,
312 impl_->part_interface, impl_->block_tridiags, impl_->a_minus_d,
316 zeroStartingSolution,
320 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
323 template <
typename MatrixType>
328 return ComputeParameters();
331 template <
typename MatrixType>
336 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::compute",
compute);
344 impl_->part_interface, impl_->block_tridiags,
345 in.addRadiallyToDiagonal);
348 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
351 template <
typename MatrixType>
357 in.
dampingFactor = Teuchos::ScalarTraits<scalar_type>::one();
361 template <
typename MatrixType>
367 IFPACK2_BLOCKHELPER_TIMER(
"BlockTriDiContainer::applyInverseJacobi",
applyInverseJacobi);
373 impl_->tpetra_importer,
374 impl_->async_importer,
375 impl_->overlap_communication_and_computation,
376 X, Y, impl_->Z, impl_->W,
377 impl_->part_interface, impl_->block_tridiags, impl_->a_minus_d,
386 IFPACK2_BLOCKHELPER_TIMER_FENCE(
typename BlockHelperDetails::ImplType<MatrixType>::execution_space)
390 template <
typename MatrixType>
394 return impl_->norm_manager.getNorms0();
397 template <
typename MatrixType>
401 return impl_->norm_manager.getNormsFinal();
404 template <
typename MatrixType>
407 ::apply (ConstHostView , HostView ,
int , Teuchos::ETransp ,
408 scalar_type , scalar_type )
const
410 TEUCHOS_TEST_FOR_EXCEPT_MSG(
true,
"BlockTriDiContainer::apply is not implemented. You may have reached this message "
411 <<
"because you want to use this container's performance-portable Jacobi iteration. In "
412 <<
"that case, set \"relaxation: type\" to \"MT Split Jacobi\" rather than \"Jacobi\".");
415 template <
typename MatrixType>
419 Teuchos::ETransp , scalar_type , scalar_type )
const
421 TEUCHOS_TEST_FOR_EXCEPT_MSG(
true,
"BlockTriDiContainer::weightedApply is not implemented.");
424 template <
typename MatrixType>
429 Teuchos::FancyOStream fos(Teuchos::rcp(&os,
false));
430 fos.setOutputToRootOnly(0);
435 template <
typename MatrixType>
440 std::ostringstream oss;
441 oss << Teuchos::Describable::description();
444 oss <<
"{status = initialized, computed";
447 oss <<
"{status = initialized, not computed";
451 oss <<
"{status = not initialized, not computed";
458 template <
typename MatrixType>
461 describe (Teuchos::FancyOStream& os,
462 const Teuchos::EVerbosityLevel verbLevel)
const
465 if(verbLevel==Teuchos::VERB_NONE)
return;
466 os <<
"================================================================================" << endl
467 <<
"Ifpack2::BlockTriDiContainer" << endl
468 <<
"Number of blocks = " << this->
numBlocks_ << endl
471 <<
"================================================================================" << endl
475 template <
typename MatrixType>
478 ::getName() {
return "Ifpack2::BlockTriDiContainer::ImplSimdTag"; }
480#define IFPACK2_BLOCKTRIDICONTAINER_INSTANT(S,LO,GO,N) \
481 template class Ifpack2::BlockTriDiContainer< Tpetra::RowMatrix<S, LO, GO, N> >;
Ifpack2::BlockTriDiContainer class declaration.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override
Print the object with some verbosity level to the given FancyOStream.
Definition Ifpack2_BlockTriDiContainer_def.hpp:461
BlockTriDiContainer(const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< local_ordinal_type > > &partitions, const Teuchos::RCP< const import_type > &importer, bool pointIndexed)
Constructor.
Definition Ifpack2_BlockTriDiContainer_def.hpp:162
void applyInverseJacobi(const mv_type &X, mv_type &Y, scalar_type dampingFactor, bool zeroStartingSolution=false, int numSweeps=1) const override
Compute Y := (1 - a) Y + a D^{-1} (X - R*Y).
Definition Ifpack2_BlockTriDiContainer_def.hpp:298
void initialize() override
Do all set-up operations that only require matrix structure.
Definition Ifpack2_BlockTriDiContainer_def.hpp:215
void compute() override
Extract the local tridiagonal block and prepare the solver.
Definition Ifpack2_BlockTriDiContainer_def.hpp:265
Store and solve local block tridiagonal linear problems.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:101
int numBlocks_
The number of blocks (partitions) in the container.
Definition Ifpack2_Container_decl.hpp:259
bool isInitialized() const
Whether the container has been successfully initialized.
Definition Ifpack2_Container_def.hpp:132
Container(const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< LO > > &partitions, bool pointIndexed)
Constructor.
Definition Ifpack2_Container_def.hpp:21
bool IsInitialized_
If true, the container has been successfully initialized.
Definition Ifpack2_Container_decl.hpp:287
bool isComputed() const
Whether the container has been successfully computed.
Definition Ifpack2_Container_def.hpp:137
bool IsComputed_
If true, the container has been successfully computed.
Definition Ifpack2_Container_decl.hpp:289
BlockHelperDetails::PartInterface< MatrixType > createPartInterface(const Teuchos::RCP< const typename BlockHelperDetails::ImplType< MatrixType >::tpetra_row_matrix_type > &A, const Teuchos::RCP< const typename BlockHelperDetails::ImplType< MatrixType >::tpetra_crs_graph_type > &G, const Teuchos::Array< Teuchos::Array< typename BlockHelperDetails::ImplType< MatrixType >::local_ordinal_type > > &partitions, const typename BlockHelperDetails::ImplType< MatrixType >::local_ordinal_type n_subparts_per_part_in)
Definition Ifpack2_BlockTriDiContainer_impl.hpp:1043
int applyInverseJacobi(const Teuchos::RCP< const typename BlockHelperDetails::ImplType< MatrixType >::tpetra_row_matrix_type > &A, const Teuchos::RCP< const typename BlockHelperDetails::ImplType< MatrixType >::tpetra_crs_graph_type > &G, const Teuchos::RCP< const typename BlockHelperDetails::ImplType< MatrixType >::tpetra_import_type > &tpetra_importer, const Teuchos::RCP< AsyncableImport< MatrixType > > &async_importer, const bool overlap_communication_and_computation, const typename BlockHelperDetails::ImplType< MatrixType >::tpetra_multivector_type &X, typename BlockHelperDetails::ImplType< MatrixType >::tpetra_multivector_type &Y, typename BlockHelperDetails::ImplType< MatrixType >::tpetra_multivector_type &Z, typename BlockHelperDetails::ImplType< MatrixType >::impl_scalar_type_1d_view &W, const BlockHelperDetails::PartInterface< MatrixType > &interf, const BlockTridiags< MatrixType > &btdm, const BlockHelperDetails::AmD< MatrixType > &amd, typename BlockHelperDetails::ImplType< MatrixType >::vector_type_1d_view &work, BlockHelperDetails::NormManager< MatrixType > &norm_manager, const typename BlockHelperDetails::ImplType< MatrixType >::impl_scalar_type &damping_factor, bool is_y_zero, const int max_num_sweeps, const typename BlockHelperDetails::ImplType< MatrixType >::magnitude_type tol, const int check_tol_every)
Definition Ifpack2_BlockTriDiContainer_impl.hpp:4875
void performSymbolicPhase(const Teuchos::RCP< const typename BlockHelperDetails::ImplType< MatrixType >::tpetra_row_matrix_type > &A, const Teuchos::RCP< const typename BlockHelperDetails::ImplType< MatrixType >::tpetra_crs_graph_type > &g, const BlockHelperDetails::PartInterface< MatrixType > &interf, BlockTridiags< MatrixType > &btdm, BlockHelperDetails::AmD< MatrixType > &amd, const bool overlap_communication_and_computation, const Teuchos::RCP< AsyncableImport< MatrixType > > &async_importer, bool useSeqMethod)
Definition Ifpack2_BlockTriDiContainer_impl.hpp:1862
Teuchos::RCP< const typename BlockHelperDetails::ImplType< MatrixType >::tpetra_import_type > createBlockCrsTpetraImporter(const Teuchos::RCP< const typename BlockHelperDetails::ImplType< MatrixType >::tpetra_row_matrix_type > &A)
Definition Ifpack2_BlockTriDiContainer_impl.hpp:163
Teuchos::RCP< AsyncableImport< MatrixType > > createBlockCrsAsyncImporter(const Teuchos::RCP< const typename BlockHelperDetails::ImplType< MatrixType >::tpetra_row_matrix_type > &A)
Definition Ifpack2_BlockTriDiContainer_impl.hpp:884
void performNumericPhase(const Teuchos::RCP< const typename BlockHelperDetails::ImplType< MatrixType >::tpetra_row_matrix_type > &A, const Teuchos::RCP< const typename BlockHelperDetails::ImplType< MatrixType >::tpetra_crs_graph_type > &G, const BlockHelperDetails::PartInterface< MatrixType > &interf, BlockTridiags< MatrixType > &btdm, const typename BlockHelperDetails::ImplType< MatrixType >::magnitude_type tiny)
Definition Ifpack2_BlockTriDiContainer_impl.hpp:3599
BlockTridiags< MatrixType > createBlockTridiags(const BlockHelperDetails::PartInterface< MatrixType > &interf)
Definition Ifpack2_BlockTriDiContainer_impl.hpp:1620
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41
Definition Ifpack2_BlockComputeResidualVector.hpp:23
Definition Ifpack2_BlockHelper.hpp:249
Definition Ifpack2_BlockHelper.hpp:351
Definition Ifpack2_BlockTriDiContainer_impl.hpp:1565
forward declaration
Definition Ifpack2_BlockTriDiContainer_impl.hpp:5069
Input arguments to applyInverseJacobi.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:222
bool zeroStartingSolution
Definition Ifpack2_BlockTriDiContainer_decl.hpp:225
scalar_type dampingFactor
Damping factor. Defaults to 1.
Definition Ifpack2_BlockTriDiContainer_decl.hpp:227
int checkToleranceEvery
Definition Ifpack2_BlockTriDiContainer_decl.hpp:248
magnitude_type tolerance
Definition Ifpack2_BlockTriDiContainer_decl.hpp:240
int maxNumSweeps
Definition Ifpack2_BlockTriDiContainer_decl.hpp:230