Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_Sparse3TensorUtilities.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef STOKHOS_SPARSE_3_TENSOR_UTILITIES_HPP
11 #define STOKHOS_SPARSE_3_TENSOR_UTILITIES_HPP
12 
14 #include "Epetra_Comm.h"
15 #include "Epetra_CrsGraph.h"
16 #include "Epetra_BlockMap.h"
17 #include "Epetra_LocalMap.h"
18 #include "Epetra_CrsMatrix.h"
19 #include "EpetraExt_RowMatrixOut.h"
20 
21 namespace Stokhos {
22 
24 
31  template <typename ordinal_type, typename value_type>
36  const Epetra_Comm& comm)
37  {
39 
40  // Number of stochastic rows
41  ordinal_type num_rows = basis.size();
42 
43  // Replicated local map
44  Epetra_LocalMap map(num_rows, 0, comm);
45 
46  // Graph to be created
48  Teuchos::rcp(new Epetra_CrsGraph(Copy, map, 0));
49 
50  // Loop over Cijk entries including a non-zero in the graph at
51  // indices (i,j) if there is any k for which Cijk is non-zero
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) {
56  ordinal_type j = index(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) {
59  ordinal_type i = index(i_it);
60  graph->InsertGlobalIndices(i, 1, &j);
61  }
62  }
63  }
64 
65  // Sort, remove redundencies, transform to local, ...
66  graph->FillComplete();
67 
68  return graph;
69  }
70 
72 
79  template <typename ordinal_type, typename value_type>
83  const Epetra_BlockMap& map)
84  {
86 
87  // Graph to be created
89  Teuchos::rcp(new Epetra_CrsGraph(Copy, map, 0));
90 
91  // Loop over Cijk entries including a non-zero in the graph at
92  // indices (i,j) if there is any k for which Cijk is non-zero
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) {
97  ordinal_type j = index(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) {
100  ordinal_type i = index(i_it);
101  graph->InsertGlobalIndices(i, 1, &j);
102  }
103  }
104  }
105 
106  // Sort, remove redundencies, transform to local, ...
107  graph->FillComplete();
108 
109  return graph;
110  }
111 
112  template <typename ordinal_type, typename value_type>
113  void
117  const Epetra_Comm& comm,
118  const std::string& file)
119  {
121  Stokhos::sparse3Tensor2CrsGraph(basis, Cijk, comm);
122  Epetra_CrsMatrix mat(Copy, *graph);
123  mat.FillComplete();
124  mat.PutScalar(1.0);
125  EpetraExt::RowMatrixToMatrixMarketFile(file.c_str(), mat);
126  }
127 
128  template <typename ordinal_type, typename value_type>
129  void
132  const Epetra_BlockMap& map,
133  const std::string& file)
134  {
137  Epetra_CrsMatrix mat(Copy, *graph);
138  mat.FillComplete();
139  mat.PutScalar(1.0);
140  EpetraExt::RowMatrixToMatrixMarketFile(file.c_str(), mat);
141  }
142 
143 } // namespace Stokhos
144 
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.