57 #include <Teuchos_Comm.hpp>
58 #include <Teuchos_CommHelpers.hpp>
59 #include <Teuchos_DefaultComm.hpp>
60 #include <Teuchos_RCP.hpp>
67 using Teuchos::rcp_const_cast;
68 using Teuchos::rcp_dynamic_cast;
77 typename rowAdapter_t::ConstWeightsHostView1D::execution_space;
79 template <
typename offset_t>
81 const offset_t *offsets,
const zgno_t *edgeIds) {
82 int rank = comm->getRank();
83 int nprocs = comm->getSize();
85 for (
int p = 0; p < nprocs; p++) {
87 std::cout << rank <<
":" << std::endl;
88 for (
zlno_t i = 0; i < nvtx; i++) {
89 std::cout <<
" vertex " << vtxIds[i] <<
": ";
90 for (offset_t j = offsets[i]; j < offsets[i + 1]; j++) {
91 std::cout << edgeIds[j] <<
" ";
93 std::cout << std::endl;
102 template <
typename adapter_t,
typename graph_t>
105 using idsHost_t =
typename adapter_t::ConstIdsHostView;
106 using offsetsHost_t =
typename adapter_t::ConstOffsetsHostView;
108 typename adapter_t::user_t::nonconst_local_inds_host_view_type;
110 const auto nvtx = graph.getLocalNumRows();
111 const auto nedges = graph.getLocalNumEntries();
112 const auto maxNumEntries = graph.getLocalMaxNumRowEntries();
114 typename adapter_t::Base::ConstIdsHostView adjIdsHost_(
"adjIdsHost_", nedges);
115 typename adapter_t::Base::ConstOffsetsHostView offsHost_(
"offsHost_",
118 localInds_t nbors(
"nbors", maxNumEntries);
120 for (
size_t v = 0; v < nvtx; v++) {
121 size_t numColInds = 0;
122 graph.getLocalRowCopy(v, nbors, numColInds);
124 offsHost_(v + 1) = offsHost_(v) + numColInds;
125 for (
size_t e = offsHost_(v), i = 0; e < offsHost_(v + 1); e++) {
126 adjIdsHost_(e) = graph.getColMap()->getGlobalElement(nbors(i++));
130 idsHost_t vtxIdsHost;
131 ia.getVertexIDsHostView(vtxIdsHost);
133 const auto graphIDS = graph.getRowMap()->getLocalElementList();
137 idsHost_t adjIdsHost;
138 offsetsHost_t offsetsHost;
139 ia.getEdgesHostView(offsetsHost, adjIdsHost);
145 template <
typename adapter_t,
typename graph_t>
147 using idsDevice_t =
typename adapter_t::ConstIdsDeviceView;
148 using idsHost_t =
typename adapter_t::ConstIdsHostView;
149 using offsetsDevice_t =
typename adapter_t::ConstOffsetsDeviceView;
150 using offsetsHost_t =
typename adapter_t::ConstOffsetsHostView;
151 using weightsDevice_t =
typename adapter_t::WeightsDeviceView1D;
152 using weightsHost_t =
typename adapter_t::WeightsHostView1D;
154 const auto nVtx = ia.getLocalNumIDs();
163 idsDevice_t vtxIdsDevice;
164 ia.getVertexIDsDeviceView(vtxIdsDevice);
165 idsHost_t vtxIdsHost;
166 ia.getVertexIDsHostView(vtxIdsHost);
174 idsDevice_t adjIdsDevice;
175 offsetsDevice_t offsetsDevice;
177 ia.getEdgesDeviceView(offsetsDevice, adjIdsDevice);
179 idsHost_t adjIdsHost;
180 offsetsHost_t offsetsHost;
181 ia.getEdgesHostView(offsetsHost, adjIdsHost);
190 typename adapter_t::ConstWeightsDeviceView1D{}, 50),
193 weightsDevice_t wgts0(
"wgts0", nVtx);
194 Kokkos::parallel_for(
195 nVtx, KOKKOS_LAMBDA(
const int idx) { wgts0(idx) = idx * 2; });
202 weightsDevice_t wgts1(
"wgts1", nVtx);
203 Kokkos::parallel_for(
204 nVtx, KOKKOS_LAMBDA(
const int idx) { wgts1(idx) = idx * 3; });
212 weightsDevice_t weightsDevice;
215 weightsHost_t weightsHost;
223 weightsDevice_t weightsDevice;
226 weightsHost_t weightsHost;
234 weightsDevice_t wgtsDevice;
238 weightsHost_t wgtsHost;
239 Z2_TEST_THROW(ia.getVertexWeightsHostView(wgtsHost, 2), std::runtime_error);
245 int main(
int narg,
char *arg[]) {
252 Tpetra::ScopeGuard tscope(&narg, &arg);
253 const auto comm = Tpetra::getDefaultComm();
256 Teuchos::ParameterList params;
257 params.set(
"input file",
"simple");
258 params.set(
"file type",
"Chaco");
263 const auto crsGraph = uinput->getUITpetraCrsGraph();
264 const auto rowGraph = rcp_dynamic_cast<
ztrowgraph_t>(crsGraph);
266 const auto nvtx = rowGraph->getLocalNumRows();
271 const int nWeights = 2;
279 auto tpetraCrsGraphInput = rcp(
new crsAdapter_t(crsGraph, nWeights));
284 crsPart_t *p =
new crsPart_t[nvtx];
285 memset(p, 0,
sizeof(crsPart_t) * nvtx);
286 ArrayRCP<crsPart_t> solnParts(p, 0, nvtx,
true);
288 crsSoln_t solution(env, comm, nWeights);
289 solution.setParts(solnParts);
290 tpetraCrsGraphInput->applyPartitioningSolution(*crsGraph, mMigrate,
292 const auto newG = rcp(mMigrate);
297 PrintFromRoot(
"Input adapter for Tpetra::RowGraph migrated to proc 0");
308 auto tpetraRowGraphInput = rcp(
new rowAdapter_t(rowGraph, nWeights));
312 rowPart_t *p =
new rowPart_t[nvtx];
313 memset(p, 0,
sizeof(rowPart_t) * nvtx);
314 ArrayRCP<rowPart_t> solnParts(p, 0, nvtx,
true);
316 rowSoln_t solution(env, comm, nWeights);
317 solution.setParts(solnParts);
320 tpetraRowGraphInput->applyPartitioningSolution(*crsGraph, mMigrate,
322 const auto newG = rcp(mMigrate);
327 PrintFromRoot(
"Input adapter for Tpetra::RowGraph migrated to proc 0");
331 }
catch (std::exception &e) {
332 std::cout << e.what() << std::endl;
void PrintFromRoot(const std::string &message)
Defines TpetraRowGraphAdapter class.
Defines TpetraCrsGraphAdapter class.
Provides access for Zoltan2 to Tpetra::CrsMatrix data.
Provides access for Zoltan2 to Tpetra::CrsGraph data.
int main(int narg, char **arg)
common code used by tests
typename InputTraits< User >::part_t part_t
#define Z2_TEST_COMPARE_ARRAYS(val1, val2)
A PartitioningSolution is a solution to a partitioning problem.
#define Z2_TEST_DEVICE_HOST_VIEWS(deviceView, hostView)
#define Z2_TEST_EQUALITY(val1, val2)
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
Tpetra::Map::local_ordinal_type zlno_t
#define Z2_TEST_THROW(code, ExceptType)
Tpetra::Map::global_ordinal_type zgno_t
#define Z2_TEST_NOTHROW(code)
Provides access for Zoltan2 to Tpetra::RowGraph data.