Zoltan2
Loading...
Searching...
No Matches
PartitioningSolution.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// Test the PartitioningSolution class.
12//
13// Normally a Solution is created by a Problem.
14// We create a few Solutions in this unit test.
15
16// TODO test ::convertSolutionToImportList
17
21
25
26using Teuchos::ArrayRCP;
27using Teuchos::Array;
28using Teuchos::RCP;
29using Teuchos::rcp;
30using Teuchos::arcp;
31
32
33
34void makeArrays(int wdim, int *lens, zzpart_t **ids, zscalar_t **sizes,
35 ArrayRCP<ArrayRCP<zzpart_t> > &idList, ArrayRCP<ArrayRCP<zscalar_t> > &sizeList)
36{
37 ArrayRCP<zzpart_t> *idArrays = new ArrayRCP<zzpart_t> [wdim];
38 ArrayRCP<zscalar_t> *sizeArrays = new ArrayRCP<zscalar_t> [wdim];
39
40 for (int w=0; w < wdim; w++){
41 idArrays[w] = arcp(ids[w], 0, lens[w], false);
42 sizeArrays[w] = arcp(sizes[w], 0, lens[w], false);
43 }
44
45 idList = arcp(idArrays, 0, wdim, true);
46 sizeList = arcp(sizeArrays, 0, wdim, true);
47}
48
49int main(int narg, char *arg[])
50{
51 Tpetra::ScopeGuard tscope(&narg, &arg);
52 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
53
54 int nprocs = comm->getSize();
55 int rank = comm->getRank();
56 int fail=0, gfail=0;
57 double epsilon = 10e-6;
58
60 // Arrays to hold part Ids and part Sizes for each weight
61
62 int numIdsPerProc = 10;
63 int maxNumWeights = 3;
64 int maxNumPartSizes = nprocs;
65 int *lengths = new int [maxNumWeights];
66 zzpart_t **idLists = new zzpart_t * [maxNumWeights];
67 zscalar_t **sizeLists = new zscalar_t * [maxNumWeights];
68
69 for (int w=0; w < maxNumWeights; w++){
70 idLists[w] = new zzpart_t [maxNumPartSizes];
71 sizeLists[w] = new zscalar_t [maxNumPartSizes];
72 }
73
75 // A default environment
76 RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
77
79 // A simple identifier map.
80
81 zgno_t *myGids = new zgno_t [numIdsPerProc];
82 for (int i=0, x=rank*numIdsPerProc; i < numIdsPerProc; i++)
83 myGids[i] = x++;
84
86 // TEST:
87 // One weight, one part per proc.
88 // Some part sizes are 2 and some are 1.
89
90 int numGlobalParts = nprocs;
91 int nWeights = 1;
92
93 ArrayRCP<ArrayRCP<zzpart_t> > ids;
94 ArrayRCP<ArrayRCP<zscalar_t> > sizes;
95
96 memset(lengths, 0, sizeof(int) * maxNumWeights);
97
98 lengths[0] = 1; // We give a size for 1 part.
99 idLists[0][0] = rank; // The part is my part.
100 sizeLists[0][0] = rank%2 + 1.0; // The size is 1.0 or 2.0
101
102 makeArrays(1, lengths, idLists, sizeLists, ids, sizes);
103
104 // Normalized part size for every part, for checking later on
105
106 zscalar_t *normalizedPartSizes = new zscalar_t [numGlobalParts];
107 zscalar_t sumSizes=0;
108 for (int i=0; i < numGlobalParts; i++){
109 normalizedPartSizes[i] = 1.0;
110 if (i % 2) normalizedPartSizes[i] = 2.0;
111 sumSizes += normalizedPartSizes[i];
112 }
113 for (int i=0; i < numGlobalParts; i++)
114 normalizedPartSizes[i] /= sumSizes;
115
117 // Create a solution object with part size information, and check it.
118
119 RCP<Zoltan2::PartitioningSolution<idInput_t> > solution;
120
121 try{
123 env, // application environment info
124 comm, // problem communicator
125 nWeights, // number of weights
126 ids.view(0,nWeights), // part ids
127 sizes.view(0,nWeights))); // part sizes
128 }
129 catch (std::exception &e){
130 fail=1;
131 }
132
133 TEST_FAIL_AND_EXIT(*comm, fail==0, "constructor call 1", 1);
134
135 // Test the Solution queries that are used by algorithms
136
137 if (solution->getTargetGlobalNumberOfParts() != size_t(numGlobalParts))
138 fail=2;
139
140 if (!fail && solution->getLocalNumberOfParts() != 1)
141 fail=3;
142
143 if (!fail && !solution->oneToOnePartDistribution())
144 fail=4;
145
146 if (!fail && solution->getPartDistribution() != NULL)
147 fail=5;
148
149 if (!fail && solution->getProcDistribution() != NULL)
150 fail=6;
151
152 if (!fail &&
153 ((nprocs>1 && solution->criteriaHasUniformPartSizes(0)) ||
154 (nprocs==1 && !solution->criteriaHasUniformPartSizes(0))) )
155 fail=8;
156
157 if (!fail){
158 for (int partId=0; !fail && partId < numGlobalParts; partId++){
159 zscalar_t psize = solution->getCriteriaPartSize(0, partId);
160
161 if ( psize < normalizedPartSizes[partId] - epsilon ||
162 psize > normalizedPartSizes[partId] + epsilon )
163 fail=9;
164 }
165 }
166
167 delete [] normalizedPartSizes;
168
169 gfail = globalFail(*comm, fail);
170 if (gfail){
171 printFailureCode(*comm, fail); // exits after printing "FAIL"
172 }
173
174 // Test the Solution set method that is called by algorithms
175
176 zzpart_t *partAssignments = new zzpart_t [numIdsPerProc];
177 for (int i=0; i < numIdsPerProc; i++){
178 partAssignments[i] = myGids[i] % numGlobalParts; // round robin
179 }
180 ArrayRCP<zzpart_t> partList = arcp(partAssignments, 0, numIdsPerProc);
181
182 try{
183 solution->setParts(partList);
184 }
185 catch (std::exception &e){
186 fail=10;
187 }
188
189 gfail = globalFail(*comm, fail);
190 if (gfail){
191 printFailureCode(*comm, fail); // exits after printing "FAIL"
192 }
193
194 // Test the Solution get methods that may be called by users
195 // or migration functions.
196
197 if (!fail){
198 const zzpart_t *parts = solution->getPartListView();
199 for (int i=0; !fail && i < numIdsPerProc; i++){
200 if (parts[i] != zzpart_t(myGids[i] % numGlobalParts))
201 fail = 13;
202 }
203 }
204
205 gfail = globalFail(*comm, fail);
206 if (gfail){
207 printFailureCode(*comm, fail); // exits after printing "FAIL"
208 }
209
210 if (rank==0)
211 std::cout << "PASS" << std::endl;
212
214 // TODO:
216 // Create a solution object without part size information, and check it.
218 // Test multiple weights.
220 // Test multiple parts per process.
222 // Specify a list of parts of size 0. (The rest should be uniform.)
223
224 delete [] lengths;
225 for (int w=0; w < maxNumWeights; w++){
226 delete [] idLists[w];
227 delete [] sizeLists[w];
228 }
229 delete [] idLists;
230 delete [] sizeLists;
231 delete [] myGids;
232}
int globalFail(const Comm< int > &comm, int fail)
void printFailureCode(const Comm< int > &comm, int fail)
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
idInput_t::part_t zzpart_t
void makeArrays(int wdim, int *lens, zzpart_t **ids, zscalar_t **sizes, ArrayRCP< ArrayRCP< zzpart_t > > &idList, ArrayRCP< ArrayRCP< zscalar_t > > &sizeList)
Zoltan2::BasicIdentifierAdapter< zzuser_t > idInput_t
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > zzuser_t
Defines the BasicIdentifierAdapter class.
Defines the PartitioningSolution class.
common code used by tests
float zscalar_t
Tpetra::Map ::global_ordinal_type zgno_t
int main()
This class represents a collection of global Identifiers and their associated weights,...
A simple class that can be the User template argument for an InputAdapter.
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
A PartitioningSolution is a solution to a partitioning problem.
static const std::string fail
#define epsilon
Definition nd.cpp:47