Zoltan2
Loading...
Searching...
No Matches
Zoltan2_ModelHelpers.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Zoltan2: A package of combinatorial algorithms for scientific computing
4//
5// Copyright 2012 NTESS and the Zoltan2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
15
16#ifndef _ZOLTAN2_MODELHELPERS_HPP_
17#define _ZOLTAN2_MODELHELPERS_HPP_
18
19namespace Zoltan2 {
20
21//GFD this declaration is really messy is there a better way? I couldn't typedef outside since
22// there is no user type until the function.
23template <typename User>
24RCP<Tpetra::CrsMatrix<int,
28get2ndAdjsMatFromAdjs(const Teuchos::RCP<const MeshAdapter<User> > &ia,
29 const RCP<const Comm<int> > comm,
30 Zoltan2::MeshEntityType sourcetarget,
32 typedef typename MeshAdapter<User>::gno_t gno_t;
33 typedef typename MeshAdapter<User>::lno_t lno_t;
34 typedef typename MeshAdapter<User>::node_t node_t;
35 typedef typename MeshAdapter<User>::offset_t offset_t;
36
37 typedef int nonzero_t; // adjacency matrix doesn't need scalar_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;
40 typedef Tpetra::global_size_t GST;
41 const GST dummy = Teuchos::OrdinalTraits<GST>::invalid ();
42
43/* Find the adjacency for a nodal based decomposition */
44 if (ia->availAdjs(sourcetarget, through)) {
45 using Teuchos::Array;
46 using Teuchos::as;
47 using Teuchos::RCP;
48 using Teuchos::rcp;
49
50 // Get node-element connectivity
51
52 const offset_t *offsets=NULL;
53 const gno_t *adjacencyIds=NULL;
54 ia->getAdjsView(sourcetarget, through, offsets, adjacencyIds);
55
56 const gno_t *Ids=NULL;
57 ia->getIDsViewOf(sourcetarget, Ids);
58
59 const gno_t *throughIds=NULL;
60 ia->getIDsViewOf(through, throughIds);
61
62 size_t LocalNumIDs = ia->getLocalNumOf(sourcetarget);
63
64 /***********************************************************************/
65 /************************* BUILD MAPS FOR ADJS *************************/
66 /***********************************************************************/
67
68 RCP<const map_type> sourcetargetMapG;
69 RCP<const map_type> throughMapG;
70
71 // count owned nodes
72 size_t LocalNumOfThrough = ia->getLocalNumOf(through);
73
74 // Build a list of the global sourcetarget ids...
75 gno_t min[2];
76 size_t maxcols = 0;
77 min[0] = std::numeric_limits<gno_t>::max();
78 for (size_t i = 0; i < LocalNumIDs; ++i) {
79 if (Ids[i] < min[0]) {
80 min[0] = Ids[i];
81 }
82 size_t ncols = offsets[i+1] - offsets[i];
83 if (ncols > maxcols) maxcols = ncols;
84 }
85
86 // min(throughIds[i])
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];
91 }
92 }
93
94 gno_t gmin[2];
95 Teuchos::reduceAll<int, gno_t>(*comm, Teuchos::REDUCE_MIN, 2, min, gmin);
96
97 //Generate Map for sourcetarget.
98 ArrayView<const gno_t> sourceTargetGIDs(Ids, LocalNumIDs);
99 sourcetargetMapG = rcp(new map_type(dummy,
100 sourceTargetGIDs, gmin[0], comm));
101
102 //Create a new map with IDs uniquely assigned to ranks (oneToOneSTMap)
103 /*RCP<const map_type> oneToOneSTMap =
104 Tpetra::createOneToOne<lno_t, gno_t, node_t>(sourcetargetMapG);*/
105
106 //Generate Map for through.
107// TODO
108// TODO Could check for max through id as well, and if all through ids are
109// TODO in gmin to gmax, then default constructors works below.
110// TODO Otherwise, may need a constructor that is not one-to-one containing
111// TODO all through entities on processor, followed by call to createOneToOne
112// TODO
113
114 ArrayView<const gno_t> throughGIDs(throughIds, LocalNumOfThrough);
115 throughMapG = rcp (new map_type(dummy,
116 throughGIDs, gmin[1], comm));
117
118 //Create a new map with IDs uniquely assigned to ranks (oneToOneTMap)
119 RCP<const map_type> oneToOneTMap =
120 Tpetra::createOneToOne<lno_t, gno_t, node_t>(throughMapG);
121
122 /***********************************************************************/
123 /************************* BUILD GRAPH FOR ADJS ************************/
124 /***********************************************************************/
125
126 RCP<sparse_matrix_type> adjsMatrix;
127
128 // Construct Tpetra::CrsGraph objects.
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,//oneToOneSTMap,
133 rowlens()));
134
135 Array<nonzero_t> justOneA(maxcols, 1);
136 ArrayView<const gno_t> adjacencyIdsAV(adjacencyIds, offsets[LocalNumIDs]);
137
138 for (size_t localElement=0; localElement<LocalNumIDs; ++localElement){
139 // Insert all columns for global row Ids[localElement]
140 size_t ncols = offsets[localElement+1] - offsets[localElement];
141 adjsMatrix->insertGlobalValues(Ids[localElement],
142 adjacencyIdsAV(offsets[localElement], ncols),
143 justOneA(0, ncols));
144 }// *** source loop ***
145
146 //Fill-complete adjs Graph
147 adjsMatrix->fillComplete (oneToOneTMap, //throughMapG,
148 adjsMatrix->getRowMap());
149
150 // Form 2ndAdjs
151 RCP<sparse_matrix_type> secondAdjs =
152 rcp (new sparse_matrix_type(adjsMatrix->getRowMap(),0));
153 Tpetra::MatrixMatrix::Multiply(*adjsMatrix,false,*adjsMatrix,
154 true,*secondAdjs);
155 return secondAdjs;
156 }
157 return RCP<sparse_matrix_type>();
158}
159
160template <typename User>
162 const Teuchos::RCP<const MeshAdapter<User> > &ia,
163 const RCP<const Comm<int> > comm,
164 Zoltan2::MeshEntityType sourcetarget,
166 Teuchos::ArrayRCP<typename MeshAdapter<User>::offset_t> &offsets,
167 Teuchos::ArrayRCP<typename MeshAdapter<User>::gno_t> &adjacencyIds)
168{
169 typedef typename MeshAdapter<User>::gno_t gno_t;
170 typedef typename MeshAdapter<User>::lno_t lno_t;
171 typedef typename MeshAdapter<User>::offset_t offset_t;
172 typedef typename MeshAdapter<User>::node_t node_t;
173
174 typedef int nonzero_t; // adjacency matrix doesn't need scalar_t
175 typedef Tpetra::CrsMatrix<nonzero_t,lno_t,gno_t,node_t> sparse_matrix_type;
176
177 RCP<sparse_matrix_type> secondAdjs = get2ndAdjsMatFromAdjs(ia,comm,sourcetarget,through);
178
179 /* Allocate memory necessary for the adjacency */
180 size_t LocalNumIDs = ia->getLocalNumOf(sourcetarget);
181 lno_t *start = new lno_t [LocalNumIDs+1];
182
183 if (secondAdjs!=RCP<sparse_matrix_type>()) {
184
185 size_t nadj = 0;
186
187 gno_t const *Ids=NULL;
188 ia->getIDsViewOf(sourcetarget, Ids);
189
190 std::vector<gno_t> adj;
191
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);
199
200 for (size_t j = 0; j < NumEntries; ++j) {
201 if(globalRow != Indices[j]) {
202 adj.push_back(Indices[j]);
203 nadj++;;
204 }
205 }
206 }
207
208 Ids = NULL;
209 start[LocalNumIDs] = nadj;
210
211 gno_t *adj_ = new gno_t [nadj];
212
213 for (size_t i=0; i < nadj; i++) {
214 adj_[i] = adj[i];
215 }
216
217 offsets = arcp<offset_t>(start, 0, LocalNumIDs+1, true);
218 adjacencyIds = arcp<gno_t>(adj_, 0, nadj, true);
219 }
220 else {
221 // No adjacencies could be computed; return no edges and valid offsets array
222 for (size_t i = 0; i <= LocalNumIDs; i++)
223 start[i] = 0;
224
225 offsets = arcp<offset_t>(start, 0, LocalNumIDs+1, true);
226 adjacencyIds = Teuchos::null;
227 }
228
229 //return nadj;
230}
231
232}
233
234#endif
typename Zoltan2::InputTraits< ztcrsmatrix_t >::node_t node_t
Defines the MeshAdapter interface.
typename InputTraits< User >::node_t node_t
typename InputTraits< User >::lno_t lno_t
typename InputTraits< User >::gno_t gno_t
typename InputTraits< User >::offset_t offset_t
MeshAdapter defines the interface for mesh input.
map_t::local_ordinal_type lno_t
map_t::global_ordinal_type gno_t
Created by mbenlioglu on Aug 31, 2020.
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)
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)
MeshEntityType
Enumerate entity types for meshes: Regions, Faces, Edges, or Vertices.