Zoltan2
Loading...
Searching...
No Matches
rcbTest.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
15#include <Zoltan2_config.h>
21
22using Teuchos::RCP;
23using Teuchos::rcp;
24
25typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> tMVector_t;
27// Need to use zgno_t as zzgid_t in basic types, since zzgid_t=zgno_t in MultiVector
28
29
38 const RCP<const Teuchos::Comm<int> > & comm,
39 int nParts,
40 string &filename,
41 bool doRemap
42)
43{
44 int me = comm->getRank();
45 if (me == 0)
46 std::cout << "Parallel partitioning of " << filename << ".mtx: "
47 << nParts << " parts." << std::endl;
48
49 std::string fname(filename);
50 UserInputForTests uinput(testDataFilePath, fname, comm, true);
51
52 RCP<tMVector_t> coords = uinput.getUICoordinates();
53 if (me == 0)
54 std::cout << "Multivector length = " << coords->getGlobalLength()
55 << " Num vectors = " << coords->getNumVectors() << std::endl;
56
57 RCP<const tMVector_t> coordsConst = rcp_const_cast<const tMVector_t>(coords);
58
60 inputAdapter_t ia(coordsConst);
61 if (me == 0)
62 std::cout << "Adapter constructed" << std::endl;
63
64 Teuchos::ParameterList params("test params");
65 params.set("debug_level", "basic_status");
66 params.set("num_global_parts", nParts);
67 params.set("algorithm", "rcb");
68 params.set("imbalance_tolerance", 1.1);
69 if (doRemap) params.set("remap_parts", true); // bool parameter
70
71#ifdef HAVE_ZOLTAN2_MPI
73 MPI_COMM_WORLD);
74#else
76#endif
77 if (me == 0)
78 std::cout << "Problem constructed" << std::endl;
79
80
81 problem.solve();
82 if (me == 0)
83 std::cout << "Problem solved" << std::endl;
84}
85
86void serialTest(int numParts, bool doRemap)
87{
88 int numCoords = 1000;
89 numParts *= 8;
90
91 std::cout << "Serial partitioning: " << numParts << " parts." << std::endl;
92
93 zgno_t *ids = new zgno_t [numCoords];
94 if (!ids)
95 throw std::bad_alloc();
96 for (int i=0; i < numCoords; i++)
97 ids[i] = i;
98 ArrayRCP<zgno_t> globalIds(ids, 0, numCoords, true);
99
100 Array<ArrayRCP<zscalar_t> > randomCoords(3);
101 UserInputForTests::getUIRandomData(555, numCoords, 0, 10,
102 randomCoords.view(0,3));
103
104 typedef Zoltan2::BasicVectorAdapter<myTypes_t> inputAdapter_t;
105
106 inputAdapter_t ia(numCoords, ids,
107 randomCoords[0].getRawPtr(), randomCoords[1].getRawPtr(),
108 randomCoords[2].getRawPtr(), 1,1,1);
109
110 Teuchos::ParameterList params("test params");
111 params.set("debug_level", "basic_status");
112 params.set("num_global_parts", numParts);
113 params.set("algorithm", "rcb");
114 params.set("imbalance_tolerance", 1.1);
115 if (doRemap) params.set("remap_parts", true); // bool parameter
116
117#ifdef HAVE_ZOLTAN2_MPI
119 &ia, &params, MPI_COMM_SELF);
120#else
121 Zoltan2::PartitioningProblem<inputAdapter_t> serialProblem(&ia, &params);
122#endif
123
124 serialProblem.solve();
125}
126
127void meshCoordinatesTest(const RCP<const Teuchos::Comm<int> > & comm)
128{
129 int xdim = 40;
130 int ydim = 60;
131 int zdim = 20;
132 UserInputForTests uinput(xdim, ydim, zdim, string("Laplace3D"), comm, true, true);
133
134 RCP<tMVector_t> coords = uinput.getUICoordinates();
135
136 size_t localCount = coords->getLocalLength();
137
138 zscalar_t *x=NULL, *y=NULL, *z=NULL;
139 x = coords->getDataNonConst(0).getRawPtr();
140 y = coords->getDataNonConst(1).getRawPtr();
141 z = coords->getDataNonConst(2).getRawPtr();
142
143 const zgno_t *globalIds = coords->getMap()->getLocalElementList().getRawPtr();
144 typedef Zoltan2::BasicVectorAdapter<tMVector_t> inputAdapter_t;
145
146 inputAdapter_t ia(localCount, globalIds, x, y, z, 1, 1, 1);
147
148 Teuchos::ParameterList params("test params");
149 params.set("rectilinear", true); // bool parameter
150
151#ifdef HAVE_ZOLTAN2_MPI
152 Zoltan2::PartitioningProblem<inputAdapter_t> problem(&ia, &params, MPI_COMM_WORLD);
153#else
155#endif
156
157 problem.solve();
158}
159
160int main(int narg, char *arg[])
161{
162 Tpetra::ScopeGuard tscope(&narg, &arg);
163 Teuchos::RCP<const Teuchos::Comm<int> > tcomm = Tpetra::getDefaultComm();
164
165 int rank = tcomm->getRank();
166 int nParts = tcomm->getSize();
167 bool doRemap = false;
168 string filename = "USAir97";
169
170 // Read run-time options.
171 Teuchos::CommandLineProcessor cmdp (false, false);
172 cmdp.setOption("file", &filename, "Name of the Matrix Market file to read");
173 cmdp.setOption("nparts", &nParts, "Number of parts.");
174 cmdp.setOption("remap", "no-remap", &doRemap, "Remap part numbers.");
175 cmdp.parse(narg, arg);
176
177 meshCoordinatesTest(tcomm);
178
179 testFromDataFile(tcomm, nParts, filename, doRemap);
180
181 if (rank == 0)
182 serialTest(nParts, doRemap);
183
184 if (rank == 0)
185 std::cout << "PASS" << std::endl;
186}
#define nParts
Defines the BasicVectorAdapter class.
Defines the PartitioningProblem class.
Defines the PartitioningSolution class.
common code used by tests
float zscalar_t
std::string testDataFilePath(".")
Tpetra::Map ::global_ordinal_type zgno_t
Defines the XpetraMultiVectorAdapter.
int main()
static void getUIRandomData(unsigned int seed, zlno_t length, zscalar_t min, zscalar_t max, ArrayView< ArrayRCP< zscalar_t > > data)
Generate lists of random scalars.
RCP< tMVector_t > getUICoordinates()
A simple class that can be the User template argument for an InputAdapter.
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
PartitioningProblem sets up partitioning problems for the user.
void solve(bool updateInputData=true)
Direct the problem to create a solution.
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tMVector_t
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > myTypes_t
Definition rcbTest.cpp:26
void meshCoordinatesTest(const RCP< const Teuchos::Comm< int > > &comm)
Definition rcbTest.cpp:127
void testFromDataFile(const RCP< const Teuchos::Comm< int > > &comm, int nParts, string &filename, bool doRemap)
Definition rcbTest.cpp:37
void serialTest(int numParts, bool doRemap)
Definition rcbTest.cpp:86