Amesos2 - Direct Sparse Solver Interfaces  Version of the Day
Classes | Functions
Amesos2 Utilities

Classes

struct  Amesos2::Util::get_cxs_helper< Matrix, S, GO, GS, Op >
 A generic base class for the CRS and CCS helpers. More...
 
struct  Amesos2::Util::get_ccs_helper< Matrix, S, GO, GS >
 A generic helper class for getting a CCS representation of a Matrix. More...
 
struct  Amesos2::Util::get_crs_helper< Matrix, S, GO, GS >
 Similar to get_ccs_helper , but used to get a CRS representation of the given matrix. More...
 

Functions

template<typename LO , typename GO , typename GS , typename Node >
const Teuchos::RCP< const
Tpetra::Map< LO, GO, Node > > 
Amesos2::Util::getGatherMap (const Teuchos::RCP< const Tpetra::Map< LO, GO, Node > > &map)
 Gets a Tpetra::Map described by the EDistribution. More...
 
template<typename LO , typename GO , typename GS , typename Node >
const Teuchos::RCP< const
Tpetra::Map< LO, GO, Node > > 
Amesos2::Util::getDistributionMap (EDistribution distribution, GS num_global_elements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm, GO indexBase=0, const Teuchos::RCP< const Tpetra::Map< LO, GO, Node > > &map=Teuchos::null)
 
template<typename LO , typename GO , typename GS , typename Node >
RCP< Tpetra::Map< LO, GO, Node > > Amesos2::Util::epetra_map_to_tpetra_map (const Epetra_BlockMap &map)
 Transform an Epetra_Map object into a Tpetra::Map.
 
template<typename LO , typename GO , typename GS , typename Node >
RCP< Epetra_Map > Amesos2::Util::tpetra_map_to_epetra_map (const Tpetra::Map< LO, GO, Node > &map)
 Transform a Tpetra::Map object into an Epetra_Map.
 
const RCP< const Teuchos::Comm
< int > > 
Amesos2::Util::to_teuchos_comm (RCP< const Epetra_Comm > c)
 Transform an Epetra_Comm object into a Teuchos::Comm object.
 
const RCP< const Epetra_Comm > Amesos2::Util::to_epetra_comm (RCP< const Teuchos::Comm< int > > c)
 Transfrom a Teuchos::Comm object into an Epetra_Comm object.
 
template<typename Scalar , typename GlobalOrdinal , typename GlobalSizeT >
void Amesos2::Util::transpose (ArrayView< Scalar > vals, ArrayView< GlobalOrdinal > indices, ArrayView< GlobalSizeT > ptr, ArrayView< Scalar > trans_vals, ArrayView< GlobalOrdinal > trans_indices, ArrayView< GlobalSizeT > trans_ptr)
 
template<typename Scalar1 , typename Scalar2 >
void Amesos2::Util::scale (ArrayView< Scalar1 > vals, size_t l, size_t ld, ArrayView< Scalar2 > s)
 Scales a 1-D representation of a multivector. More...
 
template<typename Scalar1 , typename Scalar2 , class BinaryOp >
void Amesos2::Util::scale (ArrayView< Scalar1 > vals, size_t l, size_t ld, ArrayView< Scalar2 > s, BinaryOp binary_op)
 Scales a 1-D representation of a multivector. More...
 
void Amesos2::Util::printLine (Teuchos::FancyOStream &out)
 Prints a line of 70 "-"s on std::cout. More...
 
template<typename Scalar , typename GlobalOrdinal , typename GlobalSizeT >
void Amesos2::Util::transpose (Teuchos::ArrayView< Scalar > vals, Teuchos::ArrayView< GlobalOrdinal > indices, Teuchos::ArrayView< GlobalSizeT > ptr, Teuchos::ArrayView< Scalar > trans_vals, Teuchos::ArrayView< GlobalOrdinal > trans_indices, Teuchos::ArrayView< GlobalSizeT > trans_ptr)
 
template<typename Scalar1 , typename Scalar2 >
void Amesos2::Util::scale (Teuchos::ArrayView< Scalar1 > vals, size_t l, size_t ld, Teuchos::ArrayView< Scalar2 > s)
 
template<typename Scalar1 , typename Scalar2 , class BinaryOp >
void Amesos2::Util::scale (Teuchos::ArrayView< Scalar1 > vals, size_t l, size_t ld, Teuchos::ArrayView< Scalar2 > s, BinaryOp binary_op)
 
template<class row_ptr_view_t , class cols_view_t , class per_view_t >
void Amesos2::Util::reorder (row_ptr_view_t &row_ptr, cols_view_t &cols, per_view_t &perm, per_view_t &peri, size_t &nnz)
 
template<class values_view_t , class row_ptr_view_t , class cols_view_t , class per_view_t >
void Amesos2::Util::reorder_values (values_view_t &values, const row_ptr_view_t &orig_row_ptr, const row_ptr_view_t &new_row_ptr, const cols_view_t &orig_cols, const per_view_t &perm, const per_view_t &peri, size_t nnz)
 
template<class array_view_t , class per_view_t >
void Amesos2::Util::apply_reorder_permutation (const array_view_t &array, array_view_t &permuted_array, const per_view_t &permutation)
 

Detailed Description

Function Documentation

template<typename LO , typename GO , typename GS , typename Node >
const Teuchos::RCP< const Tpetra::Map< LO, GO, Node > > Amesos2::Util::getGatherMap ( const Teuchos::RCP< const Tpetra::Map< LO, GO, Node > > &  map)

Gets a Tpetra::Map described by the EDistribution.

Parameters
distributionThe distribution that the returned map will conform to
num_global_elementsA global_size_t value that gives the number of global elements in the map.
commThe communicator to create the map on.
Template Parameters
LOThe local ordinal type
GOThe global ordinal type
GSThe global size type
NodeThe Kokkos node type
template<typename Scalar , typename GlobalOrdinal , typename GlobalSizeT >
void Amesos2::Util::transpose ( ArrayView< Scalar >  vals,
ArrayView< GlobalOrdinal >  indices,
ArrayView< GlobalSizeT >  ptr,
ArrayView< Scalar >  trans_vals,
ArrayView< GlobalOrdinal >  trans_indices,
ArrayView< GlobalSizeT >  trans_ptr 
)

Transposes the compressed sparse matrix representation.

template<typename Scalar1 , typename Scalar2 >
void Amesos2::Util::scale ( ArrayView< Scalar1 >  vals,
size_t  l,
size_t  ld,
ArrayView< Scalar2 >  s 
)

Scales a 1-D representation of a multivector.

Parameters
[in/out]vals The values of the multi-vector. On exit will contain the scaled values.
[in]lThe length of each vector in the multivector
[in]ldThe leading dimension of the multivector
[in]sContains the scaling factors of the diagonal scaling matrix

The first vector will be scaled by s[0] , the second vector by s[1] , etc.

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

template<typename Scalar1 , typename Scalar2 , class BinaryOp >
void Amesos2::Util::scale ( ArrayView< Scalar1 >  vals,
size_t  l,
size_t  ld,
ArrayView< Scalar2 >  s,
BinaryOp  binary_op 
)

Scales a 1-D representation of a multivector.

Parameters
[in/out]vals The values of the multi-vector. On exit will contain the scaled values.
[in]lThe length of each vector in the multivector
[in]ldThe leading dimension of the multivector
[in]sContains the scaling factors of the diagonal scaling matrix

Scales each vector by diag(s), with the scaling multiplication being performed by the `binary_op' parameter. BinaryOp is some class that defines a operator() method as

Scalar1 operator()(Scalar1 x, Scalar2 y){ }
void Amesos2::Util::printLine ( Teuchos::FancyOStream &  out)

Prints a line of 70 "-"s on std::cout.

Prints a line of 80 "-"s on out.

Referenced by Amesos2::SolverCore< ConcreteSolver, Matrix, Vector >::describe(), and Amesos2::SolverCore< ConcreteSolver, Matrix, Vector >::printTiming().