Zoltan2
Loading...
Searching...
No Matches
IdentifierModel.cpp
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
10//
11// Testing of IdentifierModel
12//
13// TODO test with BasicIdentifierAdapter with weights
14
19
20#include <set>
21#include <bitset>
22
23#include <Teuchos_Comm.hpp>
24#include <Teuchos_DefaultComm.hpp>
25#include <Teuchos_ArrayView.hpp>
26#include <Teuchos_OrdinalTraits.hpp>
27
28#include <Tpetra_CrsMatrix.hpp>
29
30using Teuchos::RCP;
31using Teuchos::rcp;
32using Teuchos::Comm;
33
34void testIdentifierModel(std::string fname, zgno_t xdim, zgno_t ydim, zgno_t zdim,
35 const RCP<const Comm<int> > &comm)
36{
37 int rank = comm->getRank();
38 int fail = 0, gfail = 0;
39
40 std::bitset<Zoltan2::NUM_MODEL_FLAGS> modelFlags = 0;
41
42 RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
43
45 // Use an Tpetra::CrsMatrix for the user data.
47 typedef Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t> tcrsMatrix_t;
48
49 UserInputForTests *uinput;
50 if (fname.size() > 0)
51 uinput = new UserInputForTests(testDataFilePath, fname, comm, true);
52 else
53 uinput = new UserInputForTests(xdim,ydim,zdim,string(""),comm, true, true);
54
55 RCP<tcrsMatrix_t > M = uinput->getUITpetraCrsMatrix();
56 zlno_t nLocalIds = M->getLocalNumRows();
57 zgno_t nGlobalIds = M->getGlobalNumRows();
58
59 ArrayView<const zgno_t> idList = M->getRowMap()->getLocalElementList();
60 std::set<zgno_t> idSet(idList.begin(), idList.end());
61
63 // Create an IdentifierModel with this input
65
69
70 RCP<const adapter_t> ia = Teuchos::rcp(new adapter_t(M));
71
73 RCP<const base_adapter_t> base_ia =
74 Teuchos::rcp_dynamic_cast<const base_adapter_t>(ia);
75
76 try{
78 base_ia, env, comm, modelFlags);
79 }
80 catch (std::exception &e){
81 std::cerr << rank << ") " << e.what() << std::endl;
82 fail = 1;
83 }
84
85 gfail = globalFail(*comm, fail);
86
87 if (gfail)
88 printFailureCode(*comm, fail);
89
90 // Test the IdentifierModel interface
91
92 if (model->getLocalNumIdentifiers() != size_t(nLocalIds)) {
93 std::cerr << rank << ") getLocalNumIdentifiers "
94 << model->getLocalNumIdentifiers() << " "
95 << nLocalIds << std::endl;
96 fail = 2;
97 }
98
99 if (!fail && model->getGlobalNumIdentifiers() != size_t(nGlobalIds)) {
100 std::cerr << rank << ") getGlobalNumIdentifiers "
101 << model->getGlobalNumIdentifiers() << " "
102 << nGlobalIds << std::endl;
103 fail = 3;
104 }
105
106 gfail = globalFail(*comm, fail);
107
108 if (gfail)
109 printFailureCode(*comm, fail);
110
111 ArrayView<const zgno_t> gids;
112 ArrayView<input_t> wgts;
113
114 model->getIdentifierList(gids, wgts);
115
116 if (!fail && gids.size() != nLocalIds) {
117 std::cerr << rank << ") getIdentifierList IDs "
118 << gids.size() << " "
119 << nLocalIds << std::endl;
120 fail = 5;
121 }
122
123 if (!fail && wgts.size() != 0) {
124 std::cerr << rank << ") getIdentifierList Weights "
125 << wgts.size() << " "
126 << 0 << std::endl;
127 fail = 6;
128 }
129
130 for (zlno_t i=0; !fail && i < nLocalIds; i++){
131 std::set<zgno_t>::iterator next = idSet.find(gids[i]);
132 if (next == idSet.end()) {
133 std::cerr << rank << ") getIdentifierList gid not found "
134 << gids[i] << std::endl;
135 fail = 7;
136 }
137 }
138
139 Kokkos::View<const zgno_t*, typename tcrsMatrix_t::device_type> gidsKokkos;
140 Kokkos::View<zscalar_t **, typename tcrsMatrix_t::device_type> wgtsKokkos;
141
142 model->getIdentifierListKokkos(gidsKokkos, wgtsKokkos);
143
144 auto gidsKokkosHost = Kokkos::create_mirror_view(gidsKokkos);
145 Kokkos::deep_copy(gidsKokkosHost, gidsKokkos);
146 auto wgtsKokkosHost = Kokkos::create_mirror_view(wgtsKokkos);
147 Kokkos::deep_copy(wgtsKokkosHost, wgtsKokkos);
148
149
150 if (!fail && gidsKokkosHost.extent(0) != static_cast<size_t>(nLocalIds)) {
151 std::cerr << rank << ") getIdentifierList IDs "
152 << gidsKokkosHost.extent(0) << " "
153 << nLocalIds << std::endl;
154 fail = 5;
155 }
156
157 if (!fail && wgtsKokkosHost.extent(1) != 0) {
158 std::cerr << rank << ") getIdentifierList Weights "
159 << wgtsKokkosHost.extent(1) << " "
160 << 0 << std::endl;
161 fail = 6;
162 }
163
164 for (zlno_t i=0; !fail && i < nLocalIds; i++){
165 std::set<zgno_t>::iterator next = idSet.find(gidsKokkosHost(i));
166 if (next == idSet.end()) {
167 std::cerr << rank << ") getIdentifierList gid not found "
168 << gidsKokkosHost(i) << std::endl;
169 fail = 7;
170 }
171 }
172
173
174 gfail = globalFail(*comm, fail);
175
176 if (gfail)
177 printFailureCode(*comm, fail);
178
179 delete model;
180 delete uinput;
181}
182
183int main(int narg, char *arg[])
184{
185 Tpetra::ScopeGuard tscope(&narg, &arg);
186 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
187
188 int rank = comm->getRank();
189
190 string fname("simple");
191
192
193 if (rank == 0){
194 std::cout << fname;
195 }
196 testIdentifierModel(fname, 0,0,0,comm);
197
198 if (rank==0) std::cout << "PASS" << std::endl;
199
200 return 0;
201}
int globalFail(const Comm< int > &comm, int fail)
void printFailureCode(const Comm< int > &comm, int fail)
void testIdentifierModel(std::string fname, zgno_t xdim, zgno_t ydim, zgno_t zdim, const RCP< const Comm< int > > &comm)
Defines the BasicIdentifierAdapter class.
Defines the IdentifierModel interface.
common code used by tests
Tpetra::Map ::local_ordinal_type zlno_t
std::string testDataFilePath(".")
Tpetra::Map ::global_ordinal_type zgno_t
Defines the XpetraCrsMatrixAdapter class.
int main()
RCP< tcrsMatrix_t > getUITpetraCrsMatrix()
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
IdentifierModel defines the interface for all identifier models.
size_t getIdentifierListKokkos(Kokkos::View< const gno_t *, typename node_t::device_type > &Ids, Kokkos::View< scalar_t **, typename node_t::device_type > &wgts) const
global_size_t getGlobalNumIdentifiers() const
size_t getIdentifierList(ArrayView< const gno_t > &Ids, ArrayView< input_t > &wgts) const
MatrixAdapter defines the adapter interface for matrices.
The StridedData class manages lists of weights or coordinates.
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
static const std::string fail
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > tcrsMatrix_t
Zoltan2::BaseAdapter< userTypes_t > base_adapter_t