Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Classes | Static Public Member Functions | List of all members
Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag > Class Template Reference

#include <Ifpack2_BlockTriDiContainer_decl.hpp>

Inheritance diagram for Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >:
Inheritance graph
[legend]

Classes

struct  ApplyParameters
 Input arguments to applyInverseJacobi More...
 

Public Member Functions

Constructor and destructor
 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. More...
 
 BlockTriDiContainer (const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< local_ordinal_type > > &partitions, bool overlapCommAndComp=false, bool useSequentialMethod=false)
 Constructor for bare construction. More...
 
 ~BlockTriDiContainer () override
 Destructor (declared virtual for memory safety of derived classes). More...
 
Get and set methods
void setParameters (const Teuchos::ParameterList &List) override
 Set all necessary parameters. More...
 
void clearBlocks () override
 
Mathematical functions
void initialize () override
 Do all set-up operations that only require matrix structure. More...
 
void compute () override
 Extract the local tridiagonal block and prepare the solver. More...
 
void applyInverseJacobi (const mv_type &X, mv_type &Y, scalar_type dampingFactor, bool zeroStartingSolution=false, int numSweeps=1) const override
 
ComputeParameters createDefaultComputeParameters () const
 Create a ComputeParameters struct with default values. More...
 
void compute (const ComputeParameters &input)
 Extract the local tridiagonal block and prepare the solver. More...
 
ApplyParameters createDefaultApplyParameters () const
 Create an ApplyParameters struct with default values. More...
 
int applyInverseJacobi (const mv_type &X, mv_type &Y, const ApplyParameters &input) const
 Compute Y := D^{-1} (X - R*Y). More...
 
const magnitude_type getNorms0 () const
 If a norm-based method was used, return a L2 norm of all rhs at the first iteration; otherwise return a minus one indicating norm is not requested. More...
 
const magnitude_type getNormsFinal () const
 If a norm-based method was used, return a L2 norm of all rhs; otherwise return zero. More...
 
void apply (host_view_type X, host_view_type Y, int blockIndex, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const override
 
void weightedApply (host_view_type X, host_view_type Y, host_view_type W, int blockIndex, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const override
 
Miscellaneous methods
std::ostream & print (std::ostream &os) const override
 Print information about this object to the given output stream. More...
 
Implementation of Teuchos::Describable
std::string description () const override
 A one-line description of this object. More...
 
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. More...
 
- Public Member Functions inherited from Ifpack2::Container< MatrixType >
 Container (const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< LO > > &partitions, bool pointIndexed)
 Constructor. More...
 
virtual ~Container ()
 Destructor. More...
 
Teuchos::ArrayView< const LO > getBlockRows (int blockIndex) const
 Local indices of the rows of the input matrix that belong to this block. More...
 
void setBlockSizes (const Teuchos::Array< Teuchos::Array< LO > > &partitions)
 Initialize arrays with information about block sizes. More...
 
bool isInitialized () const
 Whether the container has been successfully initialized. More...
 
bool isComputed () const
 Whether the container has been successfully computed. More...
 
virtual void apply (HostView X, HostView Y, int blockIndex, Teuchos::ETransp mode=Teuchos::NO_TRANS, SC alpha=Teuchos::ScalarTraits< SC >::one(), SC beta=Teuchos::ScalarTraits< SC >::zero()) const =0
 Compute Y := alpha * M^{-1} X + beta*Y. More...
 
virtual void weightedApply (HostView X, HostView Y, HostView D, int blockIndex, Teuchos::ETransp mode=Teuchos::NO_TRANS, SC alpha=Teuchos::ScalarTraits< SC >::one(), SC beta=Teuchos::ScalarTraits< SC >::zero()) const =0
 Compute Y := alpha * diag(D) * M^{-1} (diag(D) * X) + beta*Y. More...
 
virtual void applyInverseJacobi (const mv_type &, mv_type &, SC dampingFactor, bool, int) const =0
 Compute Y := (1 - a) Y + a D^{-1} (X - R*Y). More...
 
virtual void applyMV (mv_type &X, mv_type &Y) const
 Wrapper for apply with MultiVector. More...
 
virtual void weightedApplyMV (mv_type &X, mv_type &Y, vector_type &W) const
 Wrapper for weightedApply with MultiVector. More...
 

Static Public Member Functions

static std::string getName ()
 Get the name of this container type for Details::constructContainer() More...
 
- Static Public Member Functions inherited from Ifpack2::Container< MatrixType >
static std::string getName ()
 

Additional Inherited Members

- Protected Types inherited from Ifpack2::Container< MatrixType >
using ISC = typename Kokkos::Details::ArithTraits< SC >::val_type
 Internal representation of Scalar in Kokkos::View. More...
 
using HostView = typename mv_type::dual_view_type::t_host
 
- Protected Member Functions inherited from Ifpack2::Container< MatrixType >
virtual void DoGSBlock (HostView X, HostView Y, HostView Y2, HostView Resid, SC dampingFactor, LO i) const
 Do one step of Gauss-Seidel on block i (used by DoGaussSeidel and DoSGS) More...
 
- Protected Attributes inherited from Ifpack2::Container< MatrixType >
Teuchos::RCP< const
row_matrix_type > 
inputMatrix_
 The input matrix to the constructor. More...
 
Teuchos::RCP< const
crs_matrix_type > 
inputCrsMatrix_
 The input matrix, dynamic cast to CrsMatrix. May be null. More...
 
Teuchos::RCP< const
block_crs_matrix_type > 
inputBlockMatrix_
 The input matrix, dynamic cast to BlockCrsMatrix. May be null. More...
 
int numBlocks_
 The number of blocks (partitions) in the container. More...
 
Teuchos::Array< LO > blockRows_
 Local indices of the rows of the input matrix that belong to this block. More...
 
Teuchos::Array< LO > blockSizes_
 Number of rows in each block. More...
 
Teuchos::Array< LO > blockOffsets_
 Starting index in blockRows_ of local row indices for each block. More...
 
Teuchos::RCP< vector_type > Diag_
 Diagonal elements. More...
 
bool IsParallel_
 Whether the problem is distributed across multiple MPI processes. More...
 
LO NumLocalRows_
 Number of local rows in input matrix. More...
 
GO NumGlobalRows_
 Number of global rows in input matrix. More...
 
GO NumGlobalNonzeros_
 Number of nonzeros in input matrix. More...
 
bool hasBlockCrs_
 Whether the input matrix is a BlockCRS matrix. More...
 
int bcrsBlockSize_
 If hasBlockCrs_, the number of DOFs per vertex. Otherwise 1. More...
 
bool pointIndexed_
 (If hasBlockCrs_) Whether the blocks are described using sub-block row indices instead of full block rows. More...
 
LO scalarsPerRow_
 
bool IsInitialized_
 If true, the container has been successfully initialized. More...
 
bool IsComputed_
 If true, the container has been successfully computed. More...
 

Detailed Description

template<typename MatrixType>
class Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >

Partial specialization with SIMD<ScalarType> internally used. The code does not support Sacado UQ types. As UQ type also uses SIMD like structure, it needs an additional specialization.

Constructor & Destructor Documentation

template<typename MatrixType >
Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::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.

matrix [in] The original input matrix. This Container will construct a local diagonal block from the rows given by localRows.

Parameters
partitions[in] The set of (local) rows assigned to this container. partitions[i] == j, where i (from 0 to getNumRows() - 1) indicates the SparseContainer's row, and j indicates the local row in the calling process. partitions.size() gives the number of rows in the local matrix on each process. This may be different on different processes. If partitions is empty, then a local block-Jacobi preconditioner is formed; this is equivalent to every part of partitions having size 1.
template<typename MatrixType >
Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::BlockTriDiContainer ( const Teuchos::RCP< const row_matrix_type > &  matrix,
const Teuchos::Array< Teuchos::Array< local_ordinal_type > > &  partitions,
bool  overlapCommAndComp = false,
bool  useSequentialMethod = false 
)

Constructor for bare construction.

Constructor for users who want to use BlockTriDiContainer directly. See main constructor for documentation. This constructor removes the Container-general arguments that are not used in BlockTriDiContainer.

Parameters
overlapCommAndComp[in] Overlap communication and computation. This is not always better; it depends on (at least) the MPI implementation and the machine architecture. Defaults to false. It has to be specified at construction because the core data structures depend on whether overlapCommAndComp is true or false.
useSequentialMethod[in] This is for development and testing purposes only.
template<typename MatrixType >
Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::~BlockTriDiContainer ( )
override

Destructor (declared virtual for memory safety of derived classes).

Member Function Documentation

template<typename MatrixType >
void Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::setParameters ( const Teuchos::ParameterList List)
overridevirtual

Set all necessary parameters.

Implements Ifpack2::Container< MatrixType >.

template<typename MatrixType >
void Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::initialize ( )
overridevirtual

Do all set-up operations that only require matrix structure.

Implements Ifpack2::Container< MatrixType >.

template<typename MatrixType >
void Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::compute ( )
overridevirtual

Extract the local tridiagonal block and prepare the solver.

Implements Ifpack2::Container< MatrixType >.

template<typename MatrixType >
BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::ComputeParameters Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::createDefaultComputeParameters ( ) const

Create a ComputeParameters struct with default values.

template<typename MatrixType >
void Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::compute ( const ComputeParameters &  input)

Extract the local tridiagonal block and prepare the solver.

This version of applyInverseJacobi is meant to be called by direct users of this class, rather than by BlockRelaxation.

Parameters
addRadiallyToDiagonal[in] Add a constant to each diagonal entry of the matrix. It is added according to the (possibly complex) sign of the digonal entry at the time that the entry is the pivot in the LU factorization, such that the entry is moved radially outward in the complex plane. N.B. that this constant modifies the matrix in the linear equation, not simply the diagonal preconditioner.
template<typename MatrixType >
BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::ApplyParameters Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::createDefaultApplyParameters ( ) const

Create an ApplyParameters struct with default values.

template<typename MatrixType >
int Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::applyInverseJacobi ( const mv_type &  X,
mv_type &  Y,
const ApplyParameters &  input 
) const

Compute Y := D^{-1} (X - R*Y).

This version of applyInverseJacobi is meant to be called by direct users of this class, rather than by BlockRelaxation. It supports 2-norm-based termination. It returns the number of sweeps performed.

template<typename MatrixType >
const BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::magnitude_type Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::getNorms0 ( ) const

If a norm-based method was used, return a L2 norm of all rhs at the first iteration; otherwise return a minus one indicating norm is not requested.

template<typename MatrixType >
const BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::magnitude_type Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::getNormsFinal ( ) const

If a norm-based method was used, return a L2 norm of all rhs; otherwise return zero.

template<typename MatrixType >
void Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::apply ( host_view_type  X,
host_view_type  Y,
int  blockIndex,
Teuchos::ETransp  mode = Teuchos::NO_TRANS,
scalar_type  alpha = Teuchos::ScalarTraits<scalar_type>::one(),
scalar_type  beta = Teuchos::ScalarTraits<scalar_type>::zero() 
) const
override

Compute Y := (1 - a) Y + a D^{-1} (X - R*Y). Not supported. Call applyInverseJacobi instead.

template<typename MatrixType >
void Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::weightedApply ( host_view_type  X,
host_view_type  Y,
host_view_type  W,
int  blockIndex,
Teuchos::ETransp  mode = Teuchos::NO_TRANS,
scalar_type  alpha = Teuchos::ScalarTraits<scalar_type>::one(),
scalar_type  beta = Teuchos::ScalarTraits<scalar_type>::zero() 
) const
override

Compute Y := alpha * diag(D) * M^{-1} (diag(D) * X) + beta*Y. Not supported.

template<typename MatrixType >
std::ostream & Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::print ( std::ostream &  os) const
overridevirtual

Print information about this object to the given output stream.

operator<< uses this method.

Implements Ifpack2::Container< MatrixType >.

template<typename MatrixType >
std::string Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::description ( ) const
overridevirtual

A one-line description of this object.

Reimplemented from Teuchos::Describable.

template<typename MatrixType >
void Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::describe ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
overridevirtual

Print the object with some verbosity level to the given FancyOStream.

Reimplemented from Teuchos::Describable.

template<typename MatrixType >
std::string Ifpack2::BlockTriDiContainer< MatrixType, BlockTriDiContainerDetails::ImplSimdTag >::getName ( )
static

Get the name of this container type for Details::constructContainer()


The documentation for this class was generated from the following files: