26 using Teuchos::ArrayRCP;
35 ArrayRCP<ArrayRCP<zzpart_t> > &idList, ArrayRCP<ArrayRCP<zscalar_t> > &sizeList)
37 ArrayRCP<zzpart_t> *idArrays =
new ArrayRCP<zzpart_t> [wdim];
38 ArrayRCP<zscalar_t> *sizeArrays =
new ArrayRCP<zscalar_t> [wdim];
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);
45 idList = arcp(idArrays, 0, wdim,
true);
46 sizeList = arcp(sizeArrays, 0, wdim,
true);
49 int main(
int narg,
char *arg[])
51 Tpetra::ScopeGuard tscope(&narg, &arg);
52 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
54 int nprocs = comm->getSize();
55 int rank = comm->getRank();
62 int numIdsPerProc = 10;
63 int maxNumWeights = 3;
64 int maxNumPartSizes = nprocs;
65 int *lengths =
new int [maxNumWeights];
69 for (
int w=0; w < maxNumWeights; w++){
70 idLists[w] =
new zzpart_t [maxNumPartSizes];
71 sizeLists[w] =
new zscalar_t [maxNumPartSizes];
82 for (
int i=0, x=rank*numIdsPerProc; i < numIdsPerProc; i++)
90 int numGlobalParts = nprocs;
93 ArrayRCP<ArrayRCP<zzpart_t> > ids;
94 ArrayRCP<ArrayRCP<zscalar_t> > sizes;
96 memset(lengths, 0,
sizeof(
int) * maxNumWeights);
100 sizeLists[0][0] = rank%2 + 1.0;
102 makeArrays(1, lengths, idLists, sizeLists, ids, sizes);
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];
113 for (
int i=0; i < numGlobalParts; i++)
114 normalizedPartSizes[i] /= sumSizes;
119 RCP<Zoltan2::PartitioningSolution<idInput_t> > solution;
126 ids.view(0,nWeights),
127 sizes.view(0,nWeights)));
129 catch (std::exception &e){
137 if (solution->getTargetGlobalNumberOfParts() != size_t(numGlobalParts))
140 if (!fail && solution->getLocalNumberOfParts() != 1)
143 if (!fail && !solution->oneToOnePartDistribution())
146 if (!fail && solution->getPartDistribution() != NULL)
149 if (!fail && solution->getProcDistribution() != NULL)
153 ((nprocs>1 && solution->criteriaHasUniformPartSizes(0)) ||
154 (nprocs==1 && !solution->criteriaHasUniformPartSizes(0))) )
158 for (
int partId=0; !fail && partId < numGlobalParts; partId++){
159 zscalar_t psize = solution->getCriteriaPartSize(0, partId);
161 if ( psize < normalizedPartSizes[partId] - epsilon ||
162 psize > normalizedPartSizes[partId] + epsilon )
167 delete [] normalizedPartSizes;
177 for (
int i=0; i < numIdsPerProc; i++){
178 partAssignments[i] = myGids[i] % numGlobalParts;
180 ArrayRCP<zzpart_t> partList = arcp(partAssignments, 0, numIdsPerProc);
183 solution->setParts(partList);
185 catch (std::exception &e){
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))
211 std::cout <<
"PASS" << std::endl;
225 for (
int w=0; w < maxNumWeights; w++){
226 delete [] idLists[w];
227 delete [] sizeLists[w];
void printFailureCode(const Comm< int > &comm, int fail)
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > zzuser_t
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
int main(int narg, char **arg)
Defines the PartitioningSolution class.
idInput_t::part_t zzpart_t
common code used by tests
This class represents a collection of global Identifiers and their associated weights, if any.
void makeArrays(int wdim, int *lens, zzpart_t **ids, zscalar_t **sizes, ArrayRCP< ArrayRCP< zzpart_t > > &idList, ArrayRCP< ArrayRCP< zscalar_t > > &sizeList)
A PartitioningSolution is a solution to a partitioning problem.
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
static const std::string fail
InputTraits< User >::part_t part_t
Defines the BasicIdentifierAdapter class.
int globalFail(const Comm< int > &comm, int fail)
Zoltan2::BasicIdentifierAdapter< zzuser_t > idInput_t
Tpetra::Map::global_ordinal_type zgno_t