Amesos2 - Direct Sparse Solver Interfaces Version of the Day
Public Types | Static Public Attributes | Friends | List of all members
Amesos2::CssMKL< Matrix, Vector > Class Template Reference

Amesos2 interface to the CssMKL package. More...

#include <Amesos2_CssMKL_decl.hpp>

Inheritance diagram for Amesos2::CssMKL< Matrix, Vector >:
Inheritance graph
[legend]
Collaboration diagram for Amesos2::CssMKL< Matrix, Vector >:
Collaboration graph
[legend]

Public Types

typedef CssMKL< Matrix, Vector > type
 
typedef SolverCore< Amesos2::CssMKL, Matrix, Vector > super_type
 
typedef super_type::scalar_type scalar_type
 
typedef super_type::local_ordinal_type local_ordinal_type
 
typedef super_type::global_ordinal_type global_ordinal_type
 
typedef super_type::global_size_type global_size_type
 
typedef super_type::node_type node_type
 
typedef Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > map_type
 
typedef TypeMap< Amesos2::PardisoMKL, scalar_type > type_map
 
typedef type_map::type solver_scalar_type
 
typedef type_map::magnitude_type solver_magnitude_type
 
typedef TypeMap< Amesos2::PardisoMKL, global_ordinal_type >::type int_t
 
typedef FunctionMap< Amesos2::CssMKL, int_t > function_map
 
typedef Kokkos::DefaultHostExecutionSpace HostExecSpaceType
 
typedef Kokkos::View< int_t *, HostExecSpaceType > host_size_type_array
 
typedef Kokkos::View< int_t *, HostExecSpaceType > host_ordinal_type_array
 
typedef Kokkos::View< solver_scalar_type *, HostExecSpaceType > host_value_type_array
 
typedef Amesos2::CssMKL< Matrix, Vector > solver_type
 
typedef Matrix matrix_type
 
typedef Vector vector_type
 

Public Member Functions

Mathematical functions
super_typepreOrdering () override
 Pre-orders the matrix A for minimal fill-in.
 
super_typesymbolicFactorization () override
 Performs symbolic factorization on the matrix A.
 
super_typenumericFactorization () override
 Performs numeric factorization on the matrix A.
 
void solve () override
 Solves $ A X = B$ (or $ A^T X = B$ )
 
void solve (const Teuchos::Ptr< Vector > X, const Teuchos::Ptr< const Vector > B) const override
 Solve $ A X = B$ using the given X and B vectors.
 
void solve (Vector *X, const Vector *B) const override
 Solve $ A X = B$ using the given X and B vectors.
 
int solve_ir (const Teuchos::Ptr< Vector > X, const Teuchos::Ptr< const Vector > B, const int maxNumIters, const bool verbose) const
 
int solve_ir (Vector *X, const Vector *B, const int maxNumIters, const bool verbose) const
 
int solve_ir (const int maxNumIters, const bool verbose)
 
bool matrixShapeOK () override
 Returns true if the solver can handle this matrix shape.
 
void setA (const Teuchos::RCP< const Matrix > a, EPhase keep_phase=CLEAN) override
 Sets the matrix A of this solver.
 
void setA (const Matrix *a, EPhase keep_phase=CLEAN) override
 Sets the matrix A of this solver.
 
void setX (const Teuchos::RCP< Vector > x) override
 Sets the LHS vector X.
 
void setX (Vector *x) override
 Sets the LHS vector X using a raw pointer.
 
const Teuchos::RCP< Vector > getX () override
 Returns the vector that is the LHS of the linear system.
 
Vector * getXRaw () override
 Returns a raw pointer to the LHS of the linear system.
 
void setB (const Teuchos::RCP< const Vector > b) override
 Sets the RHS vector B.
 
void setB (const Vector *b) override
 Sets the RHS vector B using a raw pointer.
 
const Teuchos::RCP< const Vector > getB () override
 Returns the vector that is the RHS of the linear system.
 
const Vector * getBRaw () override
 Returns a raw pointer to the RHS of the linear system.
 
Parameter methods
super_typesetParameters (const Teuchos::RCP< Teuchos::ParameterList > &parameterList) override
 Set/update internal variables and solver options.
 
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters () const override
 Return a const parameter list of all of the valid parameters that this->setParameterList(...) will accept.
 
void setParameterList (const Teuchos::RCP< Teuchos::ParameterList > &parameterList)
 Set or update internal variables and solver options.
 
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList ()
 This is a empty stub.
 
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList ()
 This is an empty stub.
 

Static Public Attributes

static const char * name = "CSSMKL"
 The name of this solver interface.
 

Friends

class SolverCore< Amesos2::CssMKL, Matrix, Vector >
 

Constructor/Destructor methods

host_value_type_array nzvals_view_
 Stores the values of the nonzero entries for CssMKL.
 
host_value_type_array nzvals_temp_
 
host_ordinal_type_array colind_view_
 Stores the location in Ai_ and Aval_ that starts row j.
 
host_size_type_array rowptr_view_
 Stores the row indices of the nonzero entries.
 
Teuchos::Array< solver_scalar_type > xvals_
 Persisting, contiguous, 1D store for X.
 
Teuchos::Array< solver_scalar_type > bvals_
 Persisting, contiguous, 1D store for B.
 
void * pt_ [64]
 CssMKL internal data address pointer.
 
int_t mtype_
 The matrix type. We deal only with unsymmetrix matrices.
 
int_t n_
 Number of equations in the sparse linear system.
 
Teuchos::Array< int_t > perm_
 Permutation vector.
 
int_t nrhs_
 number of righthand-side vectors
 
bool css_initialized_
 
bool is_contiguous_
 
int_t msglvl_
 The messaging level. Set to 1 if you wish for Pardiso MKL to print statistical info.
 
int_t iparm_ [64]
 
MPI_Fint CssComm_
 
Teuchos::RCP< const map_type > css_rowmap_
 
Teuchos::RCP< const map_type > css_contig_rowmap_
 
Teuchos::RCP< const map_type > css_contig_colmap_
 
static const int_t maxfct_ = 1
 
static const int_t mnum_ = 1
 
static const bool complex_
 
 CssMKL (Teuchos::RCP< const Matrix > A, Teuchos::RCP< Vector > X, Teuchos::RCP< const Vector > B)
 Initialize from Teuchos::RCP.
 
 ~CssMKL ()
 Destructor.
 
int preOrdering_impl ()
 Performs pre-ordering on the matrix to increase efficiency.
 
int symbolicFactorization_impl ()
 Perform symbolic factorization of the matrix using CssMKL.
 
int numericFactorization_impl ()
 CssMKL specific numeric factorization.
 
int solve_impl (const Teuchos::Ptr< MultiVecAdapter< Vector > > X, const Teuchos::Ptr< const MultiVecAdapter< Vector > > B) const
 CssMKL specific solve.
 
bool matrixShapeOK_impl () const
 Determines whether the shape of the matrix is OK for this solver.
 
void setParameters_impl (const Teuchos::RCP< Teuchos::ParameterList > &parameterList)
 
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters_impl () const
 
bool loadA_impl (EPhase current_phase)
 Reads matrix data into internal structures.
 
void check_css_mkl_error (EPhase phase, int_t error) const
 Throws an appropriate runtime error in the event that error < 0 .
 
void set_css_mkl_matrix_type (int_t mtype=0)
 
void set_css_mkl_default_parameters (void *pt[], int_t iparm[]) const
 

Accessor methods

Teuchos::RCP< const Teuchos::Comm< int > > getComm () const override
 Returns a pointer to the Teuchos::Comm communicator with this operator.
 
StatusgetStatus () const override
 Returns a reference to this solver's internal status object.
 
std::string description () const override
 Returns a short description of this Solver.
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override
 
void printTiming (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
 Prints timing information about the current solver.
 
void getTiming (Teuchos::ParameterList &timingParameterList) const override
 Extracts timing information from the current solver.
 
void loadA (EPhase current_phase)
 Refresh this solver's internal data about A.
 
typedef Kokkos::View< scalar_type *, HostExecSpaceType > host_scalar_type_array
 
Teuchos::RCP< const MatrixAdapter< Matrix > > matrixA_
 The LHS operator.
 
bool matrix_loaded_
 
Teuchos::RCP< Vector > multiVecX_
 The LHS vector/multi-vector.
 
Teuchos::RCP< const Vector > multiVecB_
 The RHS vector/multi-vector.
 
global_size_type globalNumRows_
 Number of global rows in matrixA_.
 
global_size_type globalNumCols_
 Number of global columns in matrixA_.
 
global_size_type globalNumNonZeros_
 Number of global non-zero values in matrixA_.
 
global_size_type rowIndexBase_
 Index base of rowmap of matrixA_.
 
global_size_type columnIndexBase_
 Index base of column map of matrixA_.
 
Status status_
 Holds status information about a solver.
 
Control control_
 Parameters for solving.
 
Timers timers_
 Various timing statistics.
 
int rank_
 The MPI rank of this image.
 
bool root_
 If true, then this is the root processor.
 
int nprocs_
 Number of process images in the matrix communicator.
 
Teuchos::RCP< const map_type > contig_rowmap_
 
Teuchos::RCP< const map_type > contig_colmap_
 
host_ordinal_type_array perm_g2l
 
host_ordinal_type_array recvCountRows
 
host_ordinal_type_array recvDisplRows
 
host_ordinal_type_array recvCounts
 
host_ordinal_type_array recvDispls
 
host_ordinal_type_array transpose_map
 
host_scalar_type_array nzvals_t
 
void setNnzLU (size_t nnz)
 Set the number of non-zero values in the $L$ and $U$ factors.
 

Detailed Description

template<class Matrix, class Vector>
class Amesos2::CssMKL< Matrix, Vector >

Amesos2 interface to the CssMKL package.

This class provides access to the Pardiso (MKL version 10.3 and compatible) sparse direct solver with out-of-core solve support. Access is provided for float and double scalar types, in both real and complex. Access to to Pardiso's 64-bit integer routines is also provided.

Constructor & Destructor Documentation

◆ CssMKL()

template<class Matrix , class Vector >
Amesos2::CssMKL< Matrix, Vector >::CssMKL ( Teuchos::RCP< const Matrix >  A,
Teuchos::RCP< Vector >  X,
Teuchos::RCP< const Vector >  B 
)

Member Function Documentation

◆ preOrdering_impl()

template<class Matrix , class Vector >
int Amesos2::CssMKL< Matrix, Vector >::preOrdering_impl ( )
private

Performs pre-ordering on the matrix to increase efficiency.

CssMKL does reordering internally during symbolic factorization. Please refer to the "IPARM(2)" parameter for some reordering options.

◆ symbolicFactorization_impl()

template<class Matrix , class Vector >
int Amesos2::CssMKL< Matrix, Vector >::symbolicFactorization_impl ( )
private

Perform symbolic factorization of the matrix using CssMKL.

Called the sequence before numericFactorization.

Exceptions
std::runtime_errorCssMKL is not able to factor the matrix.

◆ numericFactorization_impl()

template<class Matrix , class Vector >
int Amesos2::CssMKL< Matrix, Vector >::numericFactorization_impl ( )
private

CssMKL specific numeric factorization.

Exceptions
std::runtime_errorCssMKL is not able to factor the matrix

◆ solve_impl()

template<class Matrix , class Vector >
int Amesos2::CssMKL< Matrix, Vector >::solve_impl ( const Teuchos::Ptr< MultiVecAdapter< Vector > >  X,
const Teuchos::Ptr< const MultiVecAdapter< Vector > >  B 
) const
private

CssMKL specific solve.

Uses the symbolic and numeric factorizations, along with the RHS vector B to solve the sparse system of equations.

The solution of the system is placed in X.

Exceptions
std::runtime_errorCssMKL is not able to solve the system.

◆ matrixShapeOK_impl()

template<class Matrix , class Vector >
bool Amesos2::CssMKL< Matrix, Vector >::matrixShapeOK_impl ( ) const
private

Determines whether the shape of the matrix is OK for this solver.

Pardiso MKL handles square matrices.

◆ setParameters_impl()

template<class Matrix , class Vector >
void Amesos2::CssMKL< Matrix, Vector >::setParameters_impl ( const Teuchos::RCP< Teuchos::ParameterList > &  parameterList)
private

The Pardiso MKL parameters that are currently recognized are:

  • "IPARM(2)"
  • "IPARM(4)"
  • "IPARM(8)"
  • "IPARM(10)"
  • "IPARM(12)"
  • "IPARM(18)"
  • "IPARM(24)"
  • "IPARM(25)"
  • "IPARM(60)"

Please see the Pardiso MKL documentation for a summary of the meaning and valid values for each parameter.

◆ getValidParameters_impl()

template<class Matrix , class Vector >
Teuchos::RCP< const Teuchos::ParameterList > Amesos2::CssMKL< Matrix, Vector >::getValidParameters_impl ( ) const
private
Returns
a const Teuchos::ParameterList of all valid parameters (set to their default values) for this solver.

◆ loadA_impl()

template<class Matrix , class Vector >
bool Amesos2::CssMKL< Matrix, Vector >::loadA_impl ( EPhase  current_phase)
private

Reads matrix data into internal structures.

Parameters
[in]current_phasean indication of which solution phase this load is being performed for.
Returns
true if the matrix was loaded, false if not

References Amesos2::CONTIGUOUS_AND_ROOTED, Amesos2::DISTRIBUTED_NO_OVERLAP, Amesos2::ROOTED, and Amesos2::SORTED_INDICES.

◆ check_css_mkl_error()

template<class Matrix , class Vector >
void Amesos2::CssMKL< Matrix, Vector >::check_css_mkl_error ( EPhase  phase,
int_t  error 
) const
private

Throws an appropriate runtime error in the event that error < 0 .

Parameters
phasethe phase for which this error is being checked. The meaning of a particular error value may depend on which phase was last performed
errorthe error value returned by CssMKL for the given phase.

We broadcast the input value from the rank=0 image to all others before checking the value. Before doing this we convert the error into an int value which allow us to easily broadcast its value to all process images without having to enable Teuchos long long support in the case where the user is making use of pardiso_64. The valid values of error certainly fit within an int.

◆ set_css_mkl_matrix_type()

template<class Matrix , class Vector >
void Amesos2::CssMKL< Matrix, Vector >::set_css_mkl_matrix_type ( int_t  mtype = 0)
private

Sets the internal mtype_ member. Errors are thrown for unacceptable scalar/mtype combinations.

Parameters
mtypethe type of the matrix. This may come as input from the interface user, or may be set to the default value in case mtype == 0 on entry to this function.

Referenced by Amesos2::CssMKL< Matrix, Vector >::CssMKL().

◆ preOrdering()

Solver< Matrix, Vector > & Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::preOrdering ( void  )
overridevirtualinherited

Pre-orders the matrix A for minimal fill-in.

Rearranges the rows and columns of the matrix A to minimize the amount of fill-in of the non-zero entries of the matrix. Pre-ordering may or may not be supported by the underlying solver. If not supported, a call to this method simply does nothing.

Returns
a reference to this .
See also
symbolicFactorization(), numericFactorization(), and solve()

Implements Amesos2::Solver< Matrix, Vector >.

◆ symbolicFactorization()

Solver< Matrix, Vector > & Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::symbolicFactorization ( void  )
overridevirtualinherited

Performs symbolic factorization on the matrix A.

In addition to performing symbolic factorization on the matrix A, the call to symbolicFactorization() implies that no change will be made to the non-zero structure of the underlying matrix without a subsequent call to symbolicFactorization().

Precondition
Postcondition
Returns
a reference to this .
See also
preOrdering(), numericFactorization(), and solve()

Implements Amesos2::Solver< Matrix, Vector >.

◆ numericFactorization()

Solver< Matrix, Vector > & Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::numericFactorization ( void  )
overridevirtualinherited

Performs numeric factorization on the matrix A.

In addition to performing numeric factorization on the matrix A, the call to numericFactorization() implies that no change will be made to the underlying matrix values without a subsequent call to numericFactorization().

Precondition
  • The non-zero structure of the matrix should not have changed since the last call to symbolicFactorization(). Other changes can have arbitrary consequences.
  • The distribution of the matrix should not have changed since the last call to symbolicFactorization().
Postcondition
Numeric factorization will be performed (or marked to be performed) allowing solve() to be performed correctly despite a potential change in the matrix values (though not in the non-zero structure).
Returns
a reference to this
See also
preOrdering(), symbolicFactorization(), and solve()

Implements Amesos2::Solver< Matrix, Vector >.

◆ solve() [1/3]

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::solve ( void  )
overridevirtualinherited

Solves $ A X = B$ (or $ A^T X = B$ )

Precondition
Postcondition
X will be set such that $ A X = B$ (or $ A^T X = B$ ), within the limits of the accuracy of the underlying solver.
Returns
void
See also
preOrdering(), symbolicFactorization(), and numericFactorization()

Implements Amesos2::Solver< Matrix, Vector >.

◆ solve() [2/3]

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::solve ( const Teuchos::Ptr< Vector >  X,
const Teuchos::Ptr< const Vector >  B 
) const
overridevirtualinherited

Solve $ A X = B$ using the given X and B vectors.

This overload of solve uses the given X and B vectors when solving. This X and B are used in place of any X and B that were given upon construction of the Amesos2 solver instance and are used only for this solve.

If a permanent change of X and B are required, see the setX() and setB() methods.

Postcondition
  • The (multi)vector X contains the solution to the system
  • The X and B given at construction time (if any) are unchanged.

Implements Amesos2::Solver< Matrix, Vector >.

◆ solve() [3/3]

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::solve ( Vector *  X,
const Vector *  B 
) const
overridevirtualinherited

Solve $ A X = B$ using the given X and B vectors.

This overload of solve uses the given X and B vectors when solving. This X and B are used in place of any X and B that were given upon construction of the Amesos2 solver instance and are used only for this solve.

If a permanent change of X and B are required, see the setX() and setB() methods.

Postcondition
  • The (multi)vector X contains the solution to the system
  • The X and B given at construction time (if any) are unchanged.

Implements Amesos2::Solver< Matrix, Vector >.

◆ matrixShapeOK()

bool Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::matrixShapeOK ( void  )
overridevirtualinherited

Returns true if the solver can handle this matrix shape.

Returns true if the matrix shape is one that the underlying concrete sparse direct solver can handle. Classes that work only on square matrices should return false for rectangular matrices. Classes that work only on symmetric matrices would return false for non-symmetric matrices. etc.

Implements Amesos2::Solver< Matrix, Vector >.

◆ setA() [1/2]

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::setA ( const Teuchos::RCP< const Matrix >  a,
EPhase  keep_phase = CLEAN 
)
overridevirtualinherited

Sets the matrix A of this solver.

Parameters
[in]aAn RCP to a matrix will will be used for future computation steps
[in]keep_phaseThis parameter tells the solver what state it should keep. For example, you may want to replace the matrix but keep the symbolic factorization because you know the structure of the new matrix is the same as the structure of the old matrix. In this case you would pass Amesos2::SYMBFACT as this parameter.

The default value for the second parameter is Amesos2::CLEAN, which means that the internal state of the solver will be completely reset. It will be as if no previous computational steps were performed.

Implements Amesos2::Solver< Matrix, Vector >.

◆ setA() [2/2]

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::setA ( const Matrix *  a,
EPhase  keep_phase = CLEAN 
)
inlineoverridevirtualinherited

Sets the matrix A of this solver.

Parameters
[in]aAn raw C pointer to a matrix will will be used for future computation steps.
[in]keep_phaseThis parameter tells the solver what state it should keep. For example, you may want to replace the matrix but keep the symbolic factorization because you know the structure of the new matrix is the same as the structure of the old matrix. In this case you would pass Amesos2::SYMBFACT as this parameter.

The default value for the second parameter is Amesos2::CLEAN, which means that the internal state of the solver will be completely reset. It will be as if no previous computational steps were performed.

Implements Amesos2::Solver< Matrix, Vector >.

◆ setX() [1/2]

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::setX ( const Teuchos::RCP< Vector >  x)
inlineoverridevirtualinherited

Sets the LHS vector X.

Implements Amesos2::Solver< Matrix, Vector >.

◆ setX() [2/2]

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::setX ( Vector *  x)
inlineoverridevirtualinherited

Sets the LHS vector X using a raw pointer.

Implements Amesos2::Solver< Matrix, Vector >.

◆ getX()

const Teuchos::RCP< Vector > Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::getX ( void  )
inlineoverridevirtualinherited

Returns the vector that is the LHS of the linear system.

Implements Amesos2::Solver< Matrix, Vector >.

◆ getXRaw()

Vector * Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::getXRaw ( void  )
inlineoverridevirtualinherited

Returns a raw pointer to the LHS of the linear system.

Implements Amesos2::Solver< Matrix, Vector >.

◆ setB() [1/2]

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::setB ( const Teuchos::RCP< const Vector >  b)
inlineoverridevirtualinherited

Sets the RHS vector B.

Implements Amesos2::Solver< Matrix, Vector >.

◆ setB() [2/2]

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::setB ( const Vector *  b)
inlineoverridevirtualinherited

Sets the RHS vector B using a raw pointer.

Implements Amesos2::Solver< Matrix, Vector >.

◆ getB()

const Teuchos::RCP< const Vector > Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::getB ( void  )
inlineoverridevirtualinherited

Returns the vector that is the RHS of the linear system.

Implements Amesos2::Solver< Matrix, Vector >.

◆ getBRaw()

const Vector * Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::getBRaw ( void  )
inlineoverridevirtualinherited

Returns a raw pointer to the RHS of the linear system.

Implements Amesos2::Solver< Matrix, Vector >.

◆ setParameters()

Solver< Matrix, Vector > & Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::setParameters ( const Teuchos::RCP< Teuchos::ParameterList > &  parameterList)
overridevirtualinherited

Set/update internal variables and solver options.

The setParameters method is consistent over all concrete solvers. It accepts general status and control parameters, as well as parameters specific to a given solver. If the solver does not recognize the parameter, then it will simply be ignored

Note
The ParameterList must be named "Amesos2". A list with any other name will be ignored.
Postcondition
  • Internal variables controlling the factorization and solve will be updated and take effect on all subsequent calls to numericFactorization() and solve().
  • All parameters whose value is to differ from the default values must be included in parameterList. Parameters not specified in parameterList revert to their default values.
Returns
a reference to this

Implements Amesos2::Solver< Matrix, Vector >.

◆ getValidParameters()

Teuchos::RCP< const Teuchos::ParameterList > Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::getValidParameters ( void  ) const
overridevirtualinherited

Return a const parameter list of all of the valid parameters that this->setParameterList(...) will accept.

Note
Check the documentation for your concrete solver to see a complete list of the values that each parameter may take. A solver may also recognize multiple data types as arguments for a particular parameters (eg. recognizing "YES" and "NO" as well as true and false ).

Implements Amesos2::Solver< Matrix, Vector >.

◆ setParameterList()

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::setParameterList ( const Teuchos::RCP< Teuchos::ParameterList > &  parameterList)
inlineinherited

Set or update internal variables and solver options.

Redefined from Teuchos::ParameterListAcceptor

Note
Alias for setParameters()
Parameters
[in]parameterList

◆ getNonconstParameterList()

Teuchos::RCP< Teuchos::ParameterList > Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::getNonconstParameterList ( )
inlineinherited

This is a empty stub.

Returns

◆ unsetParameterList()

Teuchos::RCP< Teuchos::ParameterList > Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::unsetParameterList ( )
inlineinherited

This is an empty stub.

Returns

◆ getComm()

Teuchos::RCP< const Teuchos::Comm< int > > Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::getComm ( void  ) const
inlineoverridevirtualinherited

Returns a pointer to the Teuchos::Comm communicator with this operator.

Implements Amesos2::Solver< Matrix, Vector >.

◆ getStatus()

Status & Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::getStatus ( ) const
inlineoverridevirtualinherited

Returns a reference to this solver's internal status object.

Implements Amesos2::Solver< Matrix, Vector >.

◆ description()

std::string Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::description ( void  ) const
overridevirtualinherited

Returns a short description of this Solver.

Implements Amesos2::Solver< Matrix, Vector >.

◆ describe()

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::describe ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
overridevirtualinherited

Prints the status information about the current solver with some level of verbosity

Implements Amesos2::Solver< Matrix, Vector >.

◆ printTiming()

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::printTiming ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel 
) const
overridevirtualinherited

Prints timing information about the current solver.

The Amesos2::SolverCore base class takes care of tracking total time spent in the Amesos2 interface. Concrete solver interface class are responsible for reporting other timing statistics, which include time spent in:

  • Redistribution of matrix objects,
  • Conversion of matrix objects to solver-specific formats,
  • Redistribution of multi-vector objects,
  • Conversion of multi-vector objects to solver formats,
  • TPL symbolic factorizations,
  • TPL numeric factorizations, and
  • TPL solves

Implements Amesos2::Solver< Matrix, Vector >.

◆ getTiming()

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::getTiming ( Teuchos::ParameterList &  timingParameterList) const
overridevirtualinherited

Extracts timing information from the current solver.

Results are placed into the parameter list timingParameterList.

Parameters
[out]timingParameterListAccepts timing information from the current solver

Implements Amesos2::Solver< Matrix, Vector >.

◆ loadA()

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::loadA ( EPhase  current_phase)
privateinherited

Refresh this solver's internal data about A.

Called whenever it would be necessary to refresh a solver's internal storage of the matrix A, which is whenever a phase is called that is equal to or below the current call.

For example, say a user has just previously called solve(), then calls numericFactorization(). Since the solve phase is greater than the numericFactorization phase, this is an indication that the internal store of A needs refreshing, since the user (assuming the user know what she's doing) wouldn't otherwise need to call numericFactorization following a solve.

◆ setNnzLU()

void Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::setNnzLU ( size_t  nnz)
inlineprotectedinherited

Set the number of non-zero values in the $L$ and $U$ factors.

Concrete solver classes may call this method if they wish to (or are able to) report the number of conbined non-zero count for the $L$ and $U$ factors.

Member Data Documentation

◆ iparm_

template<class Matrix , class Vector >
int_t Amesos2::CssMKL< Matrix, Vector >::iparm_[64]
private

CssMKL parameter vector. Note that the documentation uses 1-based indexing, but our interface must use 0-based indexing

Referenced by Amesos2::CssMKL< Matrix, Vector >::CssMKL().

◆ complex_

template<class Matrix , class Vector >
const bool Amesos2::CssMKL< Matrix, Vector >::complex_
staticprivate
Initial value:
= Meta::or_<std::is_same_v<solver_scalar_type, PMKL::_MKL_Complex8>,
std::is_same_v<solver_scalar_type, PMKL::_DOUBLE_COMPLEX_t>>::value

◆ matrix_loaded_

bool Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::matrix_loaded_
protectedinherited

If true indicates that the current matrix A has been loaded into internal solver structures.

◆ multiVecB_

Teuchos::RCP<const Vector> Amesos2::SolverCore< Amesos2::CssMKL , Matrix, Vector >::multiVecB_
protectedinherited

The RHS vector/multi-vector.

We point to a const Vector because Amesos2 should never directly modify B.


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