10 #ifndef STOKHOS_SPARSE_3_TENSOR_UTILITIES_HPP
11 #define STOKHOS_SPARSE_3_TENSOR_UTILITIES_HPP
19 #include "EpetraExt_RowMatrixOut.h"
31 template <
typename ordinal_type,
typename value_type>
52 for (
typename Cijk_type::k_iterator k_it=Cijk.
k_begin();
53 k_it!=Cijk.
k_end(); ++k_it) {
54 for (
typename Cijk_type::kj_iterator j_it = Cijk.
j_begin(k_it);
55 j_it != Cijk.
j_end(k_it); ++j_it) {
57 for (
typename Cijk_type::kji_iterator i_it = Cijk.
i_begin(j_it);
58 i_it != Cijk.
i_end(j_it); ++i_it) {
79 template <
typename ordinal_type,
typename value_type>
93 for (
typename Cijk_type::k_iterator k_it=Cijk.
k_begin();
94 k_it!=Cijk.
k_end(); ++k_it) {
95 for (
typename Cijk_type::kj_iterator j_it = Cijk.
j_begin(k_it);
96 j_it != Cijk.
j_end(k_it); ++j_it) {
98 for (
typename Cijk_type::kji_iterator i_it = Cijk.
i_begin(j_it);
99 i_it != Cijk.
i_end(j_it); ++i_it) {
112 template <
typename ordinal_type,
typename value_type>
118 const std::string& file)
125 EpetraExt::RowMatrixToMatrixMarketFile(file.c_str(), mat);
128 template <
typename ordinal_type,
typename value_type>
133 const std::string& file)
140 EpetraExt::RowMatrixToMatrixMarketFile(file.c_str(), mat);
145 #endif // SPARSE_3_TENSOR_UTILITIES_HPP
k_iterator k_begin() const
Iterator pointing to first k entry.
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
void sparse3Tensor2MatrixMarket(const Stokhos::OrthogPolyBasis< ordinal_type, value_type > &basis, const Stokhos::Sparse3Tensor< ordinal_type, value_type > &Cijk, const Epetra_Comm &comm, const std::string &file)
kj_iterator j_begin(const k_iterator &k) const
Iterator pointing to first j entry for given k.
int InsertGlobalIndices(int_type GlobalRow, int NumIndices, int_type *Indices)
int FillComplete(bool OptimizeDataStorage=true)
kj_iterator j_end(const k_iterator &k) const
Iterator pointing to last j entry for given k.
int PutScalar(double ScalarConstant)
Abstract base class for multivariate orthogonal polynomials.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
k_iterator k_end() const
Iterator pointing to last k entry.
Stokhos::Sparse3Tensor< int, double > Cijk_type
Teuchos::RCP< Epetra_CrsGraph > sparse3Tensor2CrsGraph(const Stokhos::OrthogPolyBasis< ordinal_type, value_type > &basis, const Stokhos::Sparse3Tensor< ordinal_type, value_type > &Cijk, const Epetra_Comm &comm)
Build an Epetra_CrsGraph from a sparse 3 tensor.
kji_iterator i_begin(const kj_iterator &j) const
Iterator pointing to first i entry for given j and k.
kji_iterator i_end(const kj_iterator &j) const
Iterator pointing to last i entry for given j and k.
virtual ordinal_type size() const =0
Return total size of basis.