16 #ifndef _ZOLTAN2_MODELHELPERS_HPP_
17 #define _ZOLTAN2_MODELHELPERS_HPP_
23 template <
typename User>
24 RCP<Tpetra::CrsMatrix<int,
29 const RCP<
const Comm<int> > comm,
37 typedef int nonzero_t;
38 typedef Tpetra::CrsMatrix<nonzero_t,lno_t,gno_t,node_t> sparse_matrix_type;
39 typedef Tpetra::Map<lno_t, gno_t, node_t> map_type;
41 const GST dummy = Teuchos::OrdinalTraits<GST>::invalid ();
44 if (ia->availAdjs(sourcetarget, through)) {
52 const offset_t *offsets=NULL;
53 const gno_t *adjacencyIds=NULL;
54 ia->getAdjsView(sourcetarget, through, offsets, adjacencyIds);
56 const gno_t *Ids=NULL;
57 ia->getIDsViewOf(sourcetarget, Ids);
59 const gno_t *throughIds=NULL;
60 ia->getIDsViewOf(through, throughIds);
62 size_t LocalNumIDs = ia->getLocalNumOf(sourcetarget);
68 RCP<const map_type> sourcetargetMapG;
69 RCP<const map_type> throughMapG;
72 size_t LocalNumOfThrough = ia->getLocalNumOf(through);
77 min[0] = std::numeric_limits<gno_t>::max();
78 for (
size_t i = 0; i < LocalNumIDs; ++i) {
79 if (Ids[i] < min[0]) {
82 size_t ncols = offsets[i+1] - offsets[i];
83 if (ncols > maxcols) maxcols = ncols;
87 min[1] = std::numeric_limits<gno_t>::max();
88 for (
size_t i = 0; i < LocalNumOfThrough; ++i) {
89 if (throughIds[i] < min[1]) {
90 min[1] = throughIds[i];
95 Teuchos::reduceAll<int, gno_t>(*comm, Teuchos::REDUCE_MIN, 2, min, gmin);
98 ArrayView<const gno_t> sourceTargetGIDs(Ids, LocalNumIDs);
99 sourcetargetMapG = rcp(
new map_type(dummy,
100 sourceTargetGIDs, gmin[0], comm));
114 ArrayView<const gno_t> throughGIDs(throughIds, LocalNumOfThrough);
115 throughMapG = rcp (
new map_type(dummy,
116 throughGIDs, gmin[1], comm));
119 RCP<const map_type> oneToOneTMap =
120 Tpetra::createOneToOne<lno_t, gno_t, node_t>(throughMapG);
126 RCP<sparse_matrix_type> adjsMatrix;
129 Array<size_t> rowlens(LocalNumIDs);
130 for (
size_t localElement = 0; localElement < LocalNumIDs; localElement++)
131 rowlens[localElement] = offsets[localElement+1] - offsets[localElement];
132 adjsMatrix = rcp (
new sparse_matrix_type (sourcetargetMapG,
135 Array<nonzero_t> justOneA(maxcols, 1);
136 ArrayView<const gno_t> adjacencyIdsAV(adjacencyIds, offsets[LocalNumIDs]);
138 for (
size_t localElement=0; localElement<LocalNumIDs; ++localElement){
140 size_t ncols = offsets[localElement+1] - offsets[localElement];
141 adjsMatrix->insertGlobalValues(Ids[localElement],
142 adjacencyIdsAV(offsets[localElement], ncols),
147 adjsMatrix->fillComplete (oneToOneTMap,
148 adjsMatrix->getRowMap());
151 RCP<sparse_matrix_type> secondAdjs =
152 rcp (
new sparse_matrix_type(adjsMatrix->getRowMap(),0));
153 Tpetra::MatrixMatrix::Multiply(*adjsMatrix,
false,*adjsMatrix,
157 return RCP<sparse_matrix_type>();
160 template <
typename User>
163 const RCP<
const Comm<int> > comm,
174 typedef int nonzero_t;
175 typedef Tpetra::CrsMatrix<nonzero_t,lno_t,gno_t,node_t> sparse_matrix_type;
180 size_t LocalNumIDs = ia->getLocalNumOf(sourcetarget);
181 lno_t *start =
new lno_t [LocalNumIDs+1];
183 if (secondAdjs!=RCP<sparse_matrix_type>()) {
187 gno_t
const *Ids=NULL;
188 ia->getIDsViewOf(sourcetarget, Ids);
190 std::vector<gno_t> adj;
192 for (
size_t localElement=0; localElement<LocalNumIDs; ++localElement){
193 start[localElement] = nadj;
194 const gno_t globalRow = Ids[localElement];
195 size_t NumEntries = secondAdjs->getNumEntriesInGlobalRow (globalRow);
196 typename sparse_matrix_type::nonconst_global_inds_host_view_type Indices(
"Indices", NumEntries);
197 typename sparse_matrix_type::nonconst_values_host_view_type Values(
"Values", NumEntries);
198 secondAdjs->getGlobalRowCopy (globalRow,Indices,Values,NumEntries);
200 for (
size_t j = 0; j < NumEntries; ++j) {
201 if(globalRow != Indices[j]) {
202 adj.push_back(Indices[j]);
209 start[LocalNumIDs] = nadj;
211 gno_t *adj_ =
new gno_t [nadj];
213 for (
size_t i=0; i < nadj; i++) {
217 offsets = arcp<offset_t>(start, 0, LocalNumIDs+1,
true);
218 adjacencyIds = arcp<gno_t>(adj_, 0, nadj,
true);
222 for (
size_t i = 0; i <= LocalNumIDs; i++)
225 offsets = arcp<offset_t>(start, 0, LocalNumIDs+1,
true);
226 adjacencyIds = Teuchos::null;
RCP< Tpetra::CrsMatrix< int, typename MeshAdapter< User >::lno_t, typename MeshAdapter< User >::gno_t, typename MeshAdapter< User >::node_t > > get2ndAdjsMatFromAdjs(const Teuchos::RCP< const MeshAdapter< User > > &ia, const RCP< const Comm< int > > comm, Zoltan2::MeshEntityType sourcetarget, Zoltan2::MeshEntityType through)
void get2ndAdjsViewFromAdjs(const Teuchos::RCP< const MeshAdapter< User > > &ia, const RCP< const Comm< int > > comm, Zoltan2::MeshEntityType sourcetarget, Zoltan2::MeshEntityType through, Teuchos::ArrayRCP< typename MeshAdapter< User >::offset_t > &offsets, Teuchos::ArrayRCP< typename MeshAdapter< User >::gno_t > &adjacencyIds)
Defines the MeshAdapter interface.
MeshAdapter defines the interface for mesh input.
map_t::global_ordinal_type gno_t
typename InputTraits< User >::node_t node_t
typename InputTraits< User >::gno_t gno_t
map_t::local_ordinal_type lno_t
typename InputTraits< User >::offset_t offset_t
MeshEntityType
Enumerate entity types for meshes: Regions, Faces, Edges, or Vertices.
Tpetra::global_size_t global_size_t
typename InputTraits< User >::lno_t lno_t