Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Types | Static Public Member Functions | List of all members
Tpetra::MatrixMarket::Reader< SparseMatrixType > Class Template Reference

Matrix Market file reader for CrsMatrix and MultiVector. More...

#include <MatrixMarket_Tpetra.hpp>

Public Types

typedef SparseMatrixType sparse_matrix_type
 This class' template parameter; a specialization of CrsMatrix. More...
 
typedef
SparseMatrixType::scalar_type 
scalar_type
 
typedef
SparseMatrixType::local_ordinal_type 
local_ordinal_type
 
typedef
SparseMatrixType::global_ordinal_type 
global_ordinal_type
 
typedef SparseMatrixType::node_type node_type
 The fourth template parameter of CrsMatrix and MultiVector. More...
 
typedef CrsGraph
< local_ordinal_type,
global_ordinal_type, node_type
sparse_graph_type
 The CrsGraph specialization associated with SparseMatrixType. More...
 
typedef MultiVector
< scalar_type,
local_ordinal_type,
global_ordinal_type, node_type
multivector_type
 The MultiVector specialization associated with SparseMatrixType. More...
 
typedef Vector< scalar_type,
local_ordinal_type,
global_ordinal_type, node_type
vector_type
 The Vector specialization associated with SparseMatrixType. More...
 

Static Public Member Functions

static Teuchos::RCP
< sparse_graph_type
readSparseGraphFile (const std::string &filename, const Teuchos::RCP< const Teuchos::Comm< int > > &comm, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
 Read sparse graph from the given Matrix Market file. More...
 
static Teuchos::RCP
< sparse_graph_type
readSparseGraphFile (const std::string &filename, const Teuchos::RCP< const Teuchos::Comm< int > > &pComm, const Teuchos::RCP< Teuchos::ParameterList > &constructorParams, const Teuchos::RCP< Teuchos::ParameterList > &fillCompleteParams, const bool tolerant=false, const bool debug=false)
 Read sparse graph from the given Matrix Market file. More...
 
static Teuchos::RCP
< sparse_graph_type
readSparseGraphFile (const std::string &filename, const Teuchos::RCP< const map_type > &rowMap, Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
 Read sparse graph from the given Matrix Market file, with provided Maps. More...
 
static Teuchos::RCP
< sparse_graph_type
readSparseGraph (std::istream &in, const Teuchos::RCP< const Teuchos::Comm< int > > &pComm, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
 Read sparse graph from the given Matrix Market input stream. More...
 
static Teuchos::RCP
< sparse_graph_type
readSparseGraph (std::istream &in, const Teuchos::RCP< const Teuchos::Comm< int > > &pComm, const Teuchos::RCP< Teuchos::ParameterList > &constructorParams, const Teuchos::RCP< Teuchos::ParameterList > &fillCompleteParams, const bool tolerant=false, const bool debug=false)
 Read sparse graph from the given Matrix Market input stream. More...
 
static Teuchos::RCP
< sparse_graph_type
readSparseGraph (std::istream &in, const Teuchos::RCP< const map_type > &rowMap, Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
 Read sparse matrix from the given Matrix Market input stream, with provided Maps. More...
 
static Teuchos::RCP
< sparse_matrix_type
readSparseFile (const std::string &filename, const Teuchos::RCP< const Teuchos::Comm< int > > &pComm, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
 Read sparse matrix from the given Matrix Market file. More...
 
static Teuchos::RCP
< sparse_matrix_type
readSparseFile (const std::string &filename, const Teuchos::RCP< const Teuchos::Comm< int > > &pComm, const Teuchos::RCP< Teuchos::ParameterList > &constructorParams, const Teuchos::RCP< Teuchos::ParameterList > &fillCompleteParams, const bool tolerant=false, const bool debug=false)
 Read sparse matrix from the given Matrix Market file. More...
 
static Teuchos::RCP
< sparse_matrix_type
readSparseFile (const std::string &filename, const Teuchos::RCP< const map_type > &rowMap, Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
 Read sparse matrix from the given Matrix Market file, with provided Maps. More...
 
static Teuchos::RCP
< sparse_matrix_type
readSparse (std::istream &in, const Teuchos::RCP< const Teuchos::Comm< int > > &pComm, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
 Read sparse matrix from the given Matrix Market input stream. More...
 
static Teuchos::RCP
< sparse_matrix_type
readSparse (std::istream &in, const Teuchos::RCP< const Teuchos::Comm< int > > &pComm, const Teuchos::RCP< Teuchos::ParameterList > &constructorParams, const Teuchos::RCP< Teuchos::ParameterList > &fillCompleteParams, const bool tolerant=false, const bool debug=false)
 Read sparse matrix from the given Matrix Market input stream. More...
 
static Teuchos::RCP
< sparse_matrix_type
readSparse (std::istream &in, const Teuchos::RCP< const map_type > &rowMap, Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const bool callFillComplete=true, const bool tolerant=false, const bool debug=false)
 Read sparse matrix from the given Matrix Market input stream, with provided Maps. More...
 
static Teuchos::RCP
< multivector_type
readDenseFile (const std::string &filename, const Teuchos::RCP< const comm_type > &comm, Teuchos::RCP< const map_type > &map, const bool tolerant=false, const bool debug=false)
 Read dense matrix (as a MultiVector) from the given Matrix Market file. More...
 
static Teuchos::RCP< vector_typereadVectorFile (const std::string &filename, const Teuchos::RCP< const comm_type > &comm, Teuchos::RCP< const map_type > &map, const bool tolerant=false, const bool debug=false)
 Read a Vector from the given Matrix Market file. More...
 
static Teuchos::RCP
< multivector_type
readDense (std::istream &in, const Teuchos::RCP< const comm_type > &comm, Teuchos::RCP< const map_type > &map, const bool tolerant=false, const bool debug=false)
 Read dense matrix (as a MultiVector) from the given Matrix Market input stream. More...
 
static Teuchos::RCP< vector_typereadVector (std::istream &in, const Teuchos::RCP< const comm_type > &comm, Teuchos::RCP< const map_type > &map, const bool tolerant=false, const bool debug=false)
 Read Vector from the given Matrix Market input stream. More...
 
static Teuchos::RCP< const
map_type
readMapFile (const std::string &filename, const Teuchos::RCP< const comm_type > &comm, const bool tolerant=false, const bool debug=false)
 Read Map (as a MultiVector) from the given Matrix Market file. More...
 
static Teuchos::RCP< const
map_type
readMap (std::istream &in, const Teuchos::RCP< const comm_type > &comm, const bool tolerant=false, const bool debug=false)
 Read Map (as a MultiVector) from the given input stream. More...
 
static Teuchos::RCP< const
map_type
readMap (std::istream &in, const Teuchos::RCP< const comm_type > &comm, const Teuchos::RCP< Teuchos::FancyOStream > &err, const bool tolerant=false, const bool debug=false)
 Read Map (as a MultiVector) from the given input stream, with optional debugging output stream. More...
 

Detailed Description

template<class SparseMatrixType>
class Tpetra::MatrixMarket::Reader< SparseMatrixType >

Matrix Market file reader for CrsMatrix and MultiVector.

Author
Mark Hoemmen

The Matrix Market (see its web site for details) defines a human-readable text file format for storing sparse and dense matrices. This class defines methods for reading sparse and dense matrices from a Matrix Market file or input stream. It represents sparse matrices as CrsMatrix and dense vectors and matrices as MultiVector. Reader can also read a Map (in the format produced by Writer) from a file or input stream.

All methods of this class only open files or read from input streams on MPI Process 0, with respect to the MPI communicator over which the given CrsMatrix or MultiVector is to be distributed.

We define the MultiVector type returned by readDense() and readDenseFile() using the scalar_type, local_ordinal_type, global_ordinal_type, and node_type typedefs in SparseMatrixType. This ensures that the multivectors returned by those methods have a type compatible with the CrsMatrix sparse matrices returned by readSparse() and readSparseFile(). We do this because the typical use case of Matrix Market files in Trilinos is to test sparse matrix methods, which usually involves reading a sparse matrix A and perhaps also a dense right-hand side b.

Template Parameters
SparseMatrixTypeA specialization of CrsMatrix.

Templating on the specialization of CrsMatrix means that the Reader expects matrix data of a type compatible with the CrsMatrix's scalar_type. In general, Matrix Market files may contain data of integer, real, or complex type. However, the reader methods have to return a CrsMatrix of a specific type, so we require that you declare a Reader with the CrsMatrix type that you want and that you expect the file(s) to contain.

We didn't find any of the alternatives to this approach acceptable. One possibility would have been to have the reader methods return a "container that can hold anything," like a boost::any. However, then you would have to know all five template arguments of the CrsMatrix in order to get the actual CrsMatrix object out. C++ doesn't have algebraic data types (see the Wikipedia entry for a good definition) that are disjoint unions of different types. Thus, we couldn't have had the readers return a CrsMatrix with scalar_type = "int or double or complex<double>." While you can implement such a type in C++ (see e.g., boost::variant), it would not be interchangeable for its component types. This is because it may not have the same memory layout (e.g., copying an array of boost::variant<int, double, complex<double> > bitwise into an array of int may not work).

Definition at line 165 of file MatrixMarket_Tpetra.hpp.

Member Typedef Documentation

template<class SparseMatrixType >
typedef SparseMatrixType Tpetra::MatrixMarket::Reader< SparseMatrixType >::sparse_matrix_type

This class' template parameter; a specialization of CrsMatrix.

Definition at line 168 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
typedef SparseMatrixType::scalar_type Tpetra::MatrixMarket::Reader< SparseMatrixType >::scalar_type

Type of the entries of the sparse matrix. The first template parameter of CrsMatrix and MultiVector.

Definition at line 173 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
typedef SparseMatrixType::local_ordinal_type Tpetra::MatrixMarket::Reader< SparseMatrixType >::local_ordinal_type

Type of the local indices of the sparse matrix. The second template parameter of CrsMatrix and MultiVector.

Definition at line 176 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
typedef SparseMatrixType::global_ordinal_type Tpetra::MatrixMarket::Reader< SparseMatrixType >::global_ordinal_type

Type of the global indices of the sparse matrix.

The third template parameter of CrsMatrix and MultiVector. This is also the type of indices as read from the Matrix Market file. Indices of the sparse matrix are read in as global ordinals, since Matrix Market files represent the whole matrix and don't have a notion of distribution.

Definition at line 185 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
typedef SparseMatrixType::node_type Tpetra::MatrixMarket::Reader< SparseMatrixType >::node_type

The fourth template parameter of CrsMatrix and MultiVector.

Definition at line 187 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
typedef CrsGraph<local_ordinal_type, global_ordinal_type, node_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::sparse_graph_type

The CrsGraph specialization associated with SparseMatrixType.

Definition at line 192 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
typedef MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::multivector_type

The MultiVector specialization associated with SparseMatrixType.

Definition at line 198 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
typedef Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::vector_type

The Vector specialization associated with SparseMatrixType.

Definition at line 204 of file MatrixMarket_Tpetra.hpp.

Member Function Documentation

template<class SparseMatrixType >
static Teuchos::RCP<sparse_graph_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readSparseGraphFile ( const std::string &  filename,
const Teuchos::RCP< const Teuchos::Comm< int > > &  comm,
const bool  callFillComplete = true,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read sparse graph from the given Matrix Market file.

Open the given file on MPI Rank 0 (with respect to the given communicator). The file should contain Matrix Market "coordinate" format sparse matrix data. Read that data on Rank 0, and distribute it to all processors. Return the resulting distributed CrsMatrix.

Note
This is a collective operation. Only Rank 0 opens the file and reads data from it, but all ranks participate and wait for the final result.
Parameters
filename[in] Name of the Matrix Market file.
pComm[in] Communicator containing all processor(s) over which the sparse matrix will be distributed.
callFillComplete[in] Whether to call fillComplete() on the Tpetra::CrsMatrix, after adding all the entries read in from the input stream.
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 1642 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<sparse_graph_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readSparseGraphFile ( const std::string &  filename,
const Teuchos::RCP< const Teuchos::Comm< int > > &  pComm,
const Teuchos::RCP< Teuchos::ParameterList > &  constructorParams,
const Teuchos::RCP< Teuchos::ParameterList > &  fillCompleteParams,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read sparse graph from the given Matrix Market file.

This is a variant of readSparseGraph() that lets you pass parameters to the CrsGraph's constructor and to its fillComplete() method.

Open the given file on Process 0 in the given communicator. The file must contain Matrix Market "coordinate" format sparse matrix data. Read that data on Process 0, and distribute it to all processes. Return the resulting distributed CrsMatrix.

Note
This is a collective operation. Only Process 0 opens the file and reads data from it, but all processes participate and wait for the final result.
Parameters
filename[in] Name of the Matrix Market file.
pComm[in] Communicator containing all process(es) over which the sparse matrix will be distributed.
constructorParams[in/out] Parameters for the CrsMatrix constructor.
fillCompleteParams[in/out] Parameters for CrsMatrix's fillComplete call.
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 1708 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<sparse_graph_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readSparseGraphFile ( const std::string &  filename,
const Teuchos::RCP< const map_type > &  rowMap,
Teuchos::RCP< const map_type > &  colMap,
const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const bool  callFillComplete = true,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read sparse graph from the given Matrix Market file, with provided Maps.

This version of readSparseGraph() requires you to provide a row Map, domain Map, and range Map. You may, if you wish, provide a column Map as well, but this is not required.

Open the given file on Process 0 (with respect to the given Maps' communicator). The file should contain Matrix Market "coordinate" format sparse matrix data. Read that data on Process 0, and distribute it to all processors. Return the resulting distributed CrsMatrix.

Note
This is a collective operation. Only Process 0 opens the file and reads data from it, but all ranks participate and wait for the final result.
Parameters
filename[in] Name of the Matrix Market file.
rowMap[in] The Map over which to distribute rows of the sparse matrix. This must be nonnull.
colMap[in/out] If nonnull: the Map over which to distribute columns of the sparse matrix. If null and if callFillComplete is true, we create this for you.
domainMap[in] The sparse matrix's domain Map. This must be nonnull. It may equal (pointer equality) the row Map, if that would be appropriate for this matrix.
rangeMap[in] The sparse matrix's range Map. This must be nonnull. It may equal (pointer equality) the row Map, if that would be appropriate for this matrix.
callFillComplete[in] Whether to call fillComplete() on the Tpetra::CrsGraph, after adding all the entries read in from the input stream.
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 1787 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<sparse_graph_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readSparseGraph ( std::istream &  in,
const Teuchos::RCP< const Teuchos::Comm< int > > &  pComm,
const bool  callFillComplete = true,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read sparse graph from the given Matrix Market input stream.

The given input stream need only be readable by MPI Rank 0 (with respect to the given communicator). The input stream should contain Matrix Market "coordinate" format sparse matrix data. Read that data on Rank 0, and distribute it to all processors. Return the resulting distributed CrsMatrix.

Note
This is a collective operation. Only Rank 0 reads data from the input stream, but all ranks participate and wait for the final result.
Parameters
in[in] The input stream from which to read.
pComm[in] Communicator containing all processor(s) over which the sparse matrix will be distributed.
callFillComplete[in] Whether to call fillComplete() on the Tpetra::CrsMatrix, after adding all the entries read in from the input stream. (Not calling fillComplete() may be useful if you want to change the matrix after reading it from a file.)
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 1860 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<sparse_graph_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readSparseGraph ( std::istream &  in,
const Teuchos::RCP< const Teuchos::Comm< int > > &  pComm,
const Teuchos::RCP< Teuchos::ParameterList > &  constructorParams,
const Teuchos::RCP< Teuchos::ParameterList > &  fillCompleteParams,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read sparse graph from the given Matrix Market input stream.

This is a variant of readSparse() that lets you pass parameters to the CrsMatrix's constructor and to its fillComplete() method.

The given input stream need only be readable by Process 0 in the given communicator. The input stream must contain Matrix Market "coordinate" format sparse matrix data. Read that data on Process 0, and distribute it to all Processes. Return the resulting distributed CrsMatrix.

Note
This is a collective operation. Only Proces 0 reads data from the input stream, but all processes participate and wait for the final result.
Parameters
in[in] The input stream from which to read.
pComm[in] Communicator containing all process(es) over which the sparse matrix will be distributed.
constructorParams[in/out] Parameters for the CrsMatrix constructor.
fillCompleteParams[in/out] Parameters for CrsMatrix's fillComplete call.
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 1911 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<sparse_graph_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readSparseGraph ( std::istream &  in,
const Teuchos::RCP< const map_type > &  rowMap,
Teuchos::RCP< const map_type > &  colMap,
const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const bool  callFillComplete = true,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read sparse matrix from the given Matrix Market input stream, with provided Maps.

This version of readSparse() requires you to provide a row Map, domain Map, and range Map. You may, if you wish, provide a column Map as well, but this is not required.

The given input stream need only be readable by Process 0 (with respect to the given Maps' communicator). The input stream must contain Matrix Market "coordinate" format sparse matrix data. Read that data on Process 0, and distribute it to all processors. Return the resulting distributed CrsMatrix.

Note
This is a collective operation. Only Process 0 reads data from the input stream, but all processes participate and wait for the final result.
Parameters
in[in/out] The input stream from which to read.
rowMap[in] The Map over which to distribute rows of the sparse matrix. This must be nonnull.
colMap[in/out] If nonnull: the Map over which to distribute columns of the sparse matrix. If null and if callFillComplete is true, we create this for you.
domainMap[in] The sparse matrix's domain Map. This must be nonnull. It may equal (pointer equality) the row Map, if that would be appropriate for this matrix.
rangeMap[in] The sparse matrix's range Map. This must be nonnull. It may equal (pointer equality) the row Map, if that would be appropriate for this matrix.
callFillComplete[in] Whether to call fillComplete() on the Tpetra::CrsMatrix, after adding all the entries read in from the input stream. (Not calling fillComplete() may be useful if you want to change the matrix after reading it from a file.)
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 1970 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<sparse_matrix_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readSparseFile ( const std::string &  filename,
const Teuchos::RCP< const Teuchos::Comm< int > > &  pComm,
const bool  callFillComplete = true,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read sparse matrix from the given Matrix Market file.

Open the given file on MPI Rank 0 (with respect to the given communicator). The file should contain Matrix Market "coordinate" format sparse matrix data. Read that data on Rank 0, and distribute it to all processors. Return the resulting distributed CrsMatrix.

Note
This is a collective operation. Only Rank 0 opens the file and reads data from it, but all ranks participate and wait for the final result.
Parameters
filename[in] Name of the Matrix Market file.
pComm[in] Communicator containing all processor(s) over which the sparse matrix will be distributed.
callFillComplete[in] Whether to call fillComplete() on the Tpetra::CrsMatrix, after adding all the entries read in from the input stream.
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 2013 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<sparse_matrix_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readSparseFile ( const std::string &  filename,
const Teuchos::RCP< const Teuchos::Comm< int > > &  pComm,
const Teuchos::RCP< Teuchos::ParameterList > &  constructorParams,
const Teuchos::RCP< Teuchos::ParameterList > &  fillCompleteParams,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read sparse matrix from the given Matrix Market file.

This is a variant of readSparseFile() that lets you pass parameters to the CrsMatrix's constructor and to its fillComplete() method.

Open the given file on Process 0 in the given communicator. The file must contain Matrix Market "coordinate" format sparse matrix data. Read that data on Process 0, and distribute it to all processes. Return the resulting distributed CrsMatrix.

Note
This is a collective operation. Only Process 0 opens the file and reads data from it, but all processes participate and wait for the final result.
Parameters
filename[in] Name of the Matrix Market file.
pComm[in] Communicator containing all process(es) over which the sparse matrix will be distributed.
constructorParams[in/out] Parameters for the CrsMatrix constructor.
fillCompleteParams[in/out] Parameters for CrsMatrix's fillComplete call.
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 2068 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<sparse_matrix_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readSparseFile ( const std::string &  filename,
const Teuchos::RCP< const map_type > &  rowMap,
Teuchos::RCP< const map_type > &  colMap,
const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const bool  callFillComplete = true,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read sparse matrix from the given Matrix Market file, with provided Maps.

This version of readSparseFile() requires you to provide a row Map, domain Map, and range Map. You may, if you wish, provide a column Map as well, but this is not required.

Open the given file on Process 0 (with respect to the given Maps' communicator). The file should contain Matrix Market "coordinate" format sparse matrix data. Read that data on Process 0, and distribute it to all processors. Return the resulting distributed CrsMatrix.

Note
This is a collective operation. Only Process 0 opens the file and reads data from it, but all ranks participate and wait for the final result.
Parameters
filename[in] Name of the Matrix Market file.
rowMap[in] The Map over which to distribute rows of the sparse matrix. This must be nonnull.
colMap[in/out] If nonnull: the Map over which to distribute columns of the sparse matrix. If null and if callFillComplete is true, we create this for you.
domainMap[in] The sparse matrix's domain Map. This must be nonnull. It may equal (pointer equality) the row Map, if that would be appropriate for this matrix.
rangeMap[in] The sparse matrix's range Map. This must be nonnull. It may equal (pointer equality) the row Map, if that would be appropriate for this matrix.
callFillComplete[in] Whether to call fillComplete() on the Tpetra::CrsMatrix, after adding all the entries read in from the input stream.
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 2122 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<sparse_matrix_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readSparse ( std::istream &  in,
const Teuchos::RCP< const Teuchos::Comm< int > > &  pComm,
const bool  callFillComplete = true,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read sparse matrix from the given Matrix Market input stream.

The given input stream need only be readable by MPI Rank 0 (with respect to the given communicator). The input stream should contain Matrix Market "coordinate" format sparse matrix data. Read that data on Rank 0, and distribute it to all processors. Return the resulting distributed CrsMatrix.

Note
This is a collective operation. Only Rank 0 reads data from the input stream, but all ranks participate and wait for the final result.
Parameters
in[in] The input stream from which to read.
pComm[in] Communicator containing all processor(s) over which the sparse matrix will be distributed.
callFillComplete[in] Whether to call fillComplete() on the Tpetra::CrsMatrix, after adding all the entries read in from the input stream. (Not calling fillComplete() may be useful if you want to change the matrix after reading it from a file.)
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 2193 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<sparse_matrix_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readSparse ( std::istream &  in,
const Teuchos::RCP< const Teuchos::Comm< int > > &  pComm,
const Teuchos::RCP< Teuchos::ParameterList > &  constructorParams,
const Teuchos::RCP< Teuchos::ParameterList > &  fillCompleteParams,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read sparse matrix from the given Matrix Market input stream.

This is a variant of readSparse() that lets you pass parameters to the CrsMatrix's constructor and to its fillComplete() method.

The given input stream need only be readable by Process 0 in the given communicator. The input stream must contain Matrix Market "coordinate" format sparse matrix data. Read that data on Process 0, and distribute it to all Processes. Return the resulting distributed CrsMatrix.

Note
This is a collective operation. Only Proces 0 reads data from the input stream, but all processes participate and wait for the final result.
Parameters
in[in] The input stream from which to read.
pComm[in] Communicator containing all process(es) over which the sparse matrix will be distributed.
constructorParams[in/out] Parameters for the CrsMatrix constructor.
fillCompleteParams[in/out] Parameters for CrsMatrix's fillComplete call.
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 2745 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<sparse_matrix_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readSparse ( std::istream &  in,
const Teuchos::RCP< const map_type > &  rowMap,
Teuchos::RCP< const map_type > &  colMap,
const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const bool  callFillComplete = true,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read sparse matrix from the given Matrix Market input stream, with provided Maps.

This version of readSparse() requires you to provide a row Map, domain Map, and range Map. You may, if you wish, provide a column Map as well, but this is not required.

The given input stream need only be readable by Process 0 (with respect to the given Maps' communicator). The input stream must contain Matrix Market "coordinate" format sparse matrix data. Read that data on Process 0, and distribute it to all processors. Return the resulting distributed CrsMatrix.

Note
This is a collective operation. Only Process 0 reads data from the input stream, but all processes participate and wait for the final result.
Parameters
in[in/out] The input stream from which to read.
rowMap[in] The Map over which to distribute rows of the sparse matrix. This must be nonnull.
colMap[in/out] If nonnull: the Map over which to distribute columns of the sparse matrix. If null and if callFillComplete is true, we create this for you.
domainMap[in] The sparse matrix's domain Map. This must be nonnull. It may equal (pointer equality) the row Map, if that would be appropriate for this matrix.
rangeMap[in] The sparse matrix's range Map. This must be nonnull. It may equal (pointer equality) the row Map, if that would be appropriate for this matrix.
callFillComplete[in] Whether to call fillComplete() on the Tpetra::CrsMatrix, after adding all the entries read in from the input stream. (Not calling fillComplete() may be useful if you want to change the matrix after reading it from a file.)
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 3301 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<multivector_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readDenseFile ( const std::string &  filename,
const Teuchos::RCP< const comm_type > &  comm,
Teuchos::RCP< const map_type > &  map,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read dense matrix (as a MultiVector) from the given Matrix Market file.

Open the given file on MPI Process 0 (with respect to the given communicator). The file should contain Matrix Market "array" format dense matrix data. Read that data on Process 0, and distribute it to all processors. Return the resulting distributed MultiVector.

See documentation of readDense() for details.

Parameters
filename[in] Name of the Matrix Market file from which to read. Both the filename and the file itself are only accessed on Rank 0 of the given communicator.
comm[in] Communicator containing all process(es) over which the dense matrix will be distributed.
map[in/out] On input: if nonnull, the map describing how to distribute the vector (not modified). In this case, the map's communicator and node must equal comm resp. node. If null on input, then on output, a sensible (contiguous and uniformly distributed over the given communicator) map describing the distribution of the output multivector.
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 3922 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<vector_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readVectorFile ( const std::string &  filename,
const Teuchos::RCP< const comm_type > &  comm,
Teuchos::RCP< const map_type > &  map,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read a Vector from the given Matrix Market file.

Same as readDenseFile() but will work with 1D data only.

Open the given file on MPI Process 0 (with respect to the given communicator). The file should contain Matrix Market "array" format dense matrix data with number of columns==1. Read that data on Process 0, and distribute it to all processors. Return the resulting distributed Vector.

See documentation of readVector() for details.

Parameters
filename[in] Name of the Matrix Market file from which to read. Both the filename and the file itself are only accessed on Rank 0 of the given communicator.
comm[in] Communicator containing all process(es) over which the dense matrix will be distributed.
map[in/out] On input: if nonnull, the map describing how to distribute the vector (not modified). In this case, the map's communicator and node must equal comm resp. node. If null on input, then on output, a sensible (contiguous and uniformly distributed over the given communicator) map describing the distribution of the output multivector.
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 3966 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<multivector_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readDense ( std::istream &  in,
const Teuchos::RCP< const comm_type > &  comm,
Teuchos::RCP< const map_type > &  map,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read dense matrix (as a MultiVector) from the given Matrix Market input stream.

The given input stream need only be readable by MPI Rank 0 (with respect to the given communicator). The input stream should contain Matrix Market "array" format dense matrix data. Read that data on Process 0, and distribute it to all processes in the communicator. Return the resulting distributed MultiVector.

Unlike readSparse(), this method allows callers to supply a Map over which to distribute the resulting MultiVector. The Map argument is optional; if null, we construct our own reasonable Map. We let users supply their own Map, because a common case in Tpetra is to read in or construct a sparse matrix first, and then create dense (multi)vectors distributed with the sparse matrix's domain or range Map.

Note
This is a collective operation. Only Process 0 in the communicator opens the file and reads data from it, but all processes in the communicator participate and wait for the final result.
"Tolerant" parsing mode means something different for dense matrices than it does for sparse matrices. Since Matrix Market dense matrix files don't store indices with each value, unlike sparse matrices, we can't determine the matrix dimensions from the data alone. Thus, we require the metadata to include a valid number of rows and columns. "Tolerant" in the dense case refers to the data; in tolerant mode, the number of stored matrix entries may be more or less than the number reported by the metadata (number of rows times number of columns). If more, the extra data are ignored; if less, the remainder is filled in with zeros.
On Map compatibility: Suppose that you write a multivector X to a file. Then, you read it back in as a different multivector Y distributed over the same communicator, but with a Map constructed by the input routine (i.e., a null Map on input to readDenseFile() or readDense()). In that case, the only properties shared by the maps of X and Y are that they have the same communicator and the same number of GIDs. The Maps need not necessarily be compatible (in the sense of Map::isCompatible()), and they certainly need not necessarily be the same Map (in the sense of Map::isSameAs()).
Parameters
in[in] The input stream from which to read. The stream is only accessed on Process 0 of the given communicator.
comm[in] Communicator containing all process(es) over which the dense matrix will be distributed.
map[in/out] On input: if nonnull, the map describing how to distribute the vector (not modified). In this case, the map's communicator and node must equal comm resp. node. If null on input, then on output, a sensible (contiguous and uniformly distributed over the given communicator) map describing the distribution of the output multivector.
tolerant[in] Whether to read the data tolerantly from the stream.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 4047 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<vector_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readVector ( std::istream &  in,
const Teuchos::RCP< const comm_type > &  comm,
Teuchos::RCP< const map_type > &  map,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read Vector from the given Matrix Market input stream.

Definition at line 4061 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<const map_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readMapFile ( const std::string &  filename,
const Teuchos::RCP< const comm_type > &  comm,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read Map (as a MultiVector) from the given Matrix Market file.

Open the given file on MPI Process 0 (with respect to the given communicator). The file should contain Matrix Market "array" format dense matrix data with two columns, as generated by Writer::writeMap() or Writer::writeMapFile(). Read that data on Process 0, and distribute it to all processes. Return the resulting Map.

Parameters
filename[in] Name of the Matrix Market file from which to read. Both the filename and the file itself are only accessed on Process 0 of the given communicator.
comm[in] Communicator containing all process(es) over which the Map will be distributed.
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 4094 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<const map_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readMap ( std::istream &  in,
const Teuchos::RCP< const comm_type > &  comm,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read Map (as a MultiVector) from the given input stream.

Read the given input stream on MPI Process 0 (with respect to the given communicator). The stream should contain Matrix Market "array" format dense matrix data with two columns, as generated by Writer::writeMap() or Writer::writeMapFile(). Distribute the data from Process 0 to all processes. Return the resulting Map.

Parameters
in[in/out] Input stream of Matrix Market data from which to read. This is only accessed on Process 0 of the given communicator.
comm[in] Communicator containing all process(es) over which the Map will be distributed.
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] Whether to produce copious status output useful for Tpetra developers, but probably not useful for anyone else.

Definition at line 5160 of file MatrixMarket_Tpetra.hpp.

template<class SparseMatrixType >
static Teuchos::RCP<const map_type> Tpetra::MatrixMarket::Reader< SparseMatrixType >::readMap ( std::istream &  in,
const Teuchos::RCP< const comm_type > &  comm,
const Teuchos::RCP< Teuchos::FancyOStream > &  err,
const bool  tolerant = false,
const bool  debug = false 
)
inlinestatic

Read Map (as a MultiVector) from the given input stream, with optional debugging output stream.

Warning
We make no promises about backwards compatibility for this method. It may disappear or its interface may change at any time.

Read the given input stream on MPI Process 0 (with respect to the given communicator). The stream should contain Matrix Market "array" format dense matrix data with two columns, as generated by Writer::writeMap() or Writer::writeMapFile(). Distribute the data from Process 0 to all processes. Return the resulting Map.

Parameters
in[in/out] Input stream of Matrix Market data from which to read. This is only accessed on Process 0 of the given communicator.
comm[in] Communicator containing all process(es) over which the Map will be distributed.
err[in] Optional output stream for debugging output. This is only referenced if debug is true.
tolerant[in] Whether to read the data tolerantly from the file.
debug[in] If true, write copious debugging output to err on all processes in comm.

Definition at line 5197 of file MatrixMarket_Tpetra.hpp.


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