Zoltan2
Loading...
Searching...
No Matches
BasicIdentifierInput.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// Basic testing of Zoltan2::BasicIdentifierAdapter
12
15
16#include <Teuchos_DefaultComm.hpp>
17#include <Teuchos_RCP.hpp>
18#include <Teuchos_CommHelpers.hpp>
19
20int main(int narg, char *arg[]) {
21 Tpetra::ScopeGuard tscope(&narg, &arg);
22 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
23
24 int rank = comm->getRank();
25 int nprocs = comm->getSize();
26 int fail = 0, gfail = 0;
27
28 // Create global identifiers with weights
29 zlno_t numLocalIds = 10;
30 int nWeights = 2;
31
32 zgno_t *myIds = new zgno_t[numLocalIds];
33 zscalar_t *weights = new zscalar_t [numLocalIds*nWeights];
34 zgno_t myFirstId = rank * numLocalIds * numLocalIds;
35
36 for (zlno_t i=0; i < numLocalIds; i++){
37 myIds[i] = zgno_t(myFirstId+i);
38 weights[i*nWeights] = 1.0;
39 weights[i*nWeights + 1] = (nprocs-rank) / (i+1);
40 }
41
43 std::vector<const zscalar_t *> weightValues;
44 std::vector<int> strides;
45
46 weightValues.push_back(weights);
47 weightValues.push_back(weights + 1);
48 strides.push_back(2);
49 strides.push_back(2);
50
52 weightValues, strides);
53
54 if (!fail && ia.getLocalNumIDs() != size_t(numLocalIds)){
55 fail = 4;
56 }
57 if (!fail && ia.getNumWeightsPerID() != nWeights) {
58 fail = 5;
59 }
60
61 const zgno_t *globalIdsIn;
62 zscalar_t const *weightsIn[2];
63 int weightStridesIn[2];
64
65 ia.getIDsView(globalIdsIn);
66
67 for (int w=0; !fail && w < nWeights; w++) {
68 ia.getWeightsView(weightsIn[w], weightStridesIn[w], w);
69 }
70
71 const zscalar_t *w1 = weightsIn[0];
72 const zscalar_t *w2 = weightsIn[1];
73 int incr1 = weightStridesIn[0];
74 int incr2 = weightStridesIn[1];
75
76 for (zlno_t i=0; !fail && i < numLocalIds; i++){
77 if (globalIdsIn[i] != zgno_t(myFirstId+i)) {
78 fail = 8;
79 }
80 if (!fail && w1[i*incr1] != 1.0) {
81 fail = 9;
82 }
83 if (!fail && w2[i*incr2] != weights[i*nWeights+1]) {
84 fail = 10;
85 }
86 }
87
88 delete [] myIds;
89 delete [] weights;
90
91 gfail = globalFail(*comm, fail);
92 if (gfail) {
93 printFailureCode(*comm, fail); // will exit(1)
94 }
95 if (rank == 0) {
96 std::cout << "PASS" << std::endl;
97 }
98}
99
int globalFail(const Comm< int > &comm, int fail)
void printFailureCode(const Comm< int > &comm, int fail)
Defines the BasicIdentifierAdapter class.
common code used by tests
float zscalar_t
Tpetra::Map ::local_ordinal_type zlno_t
Tpetra::Map ::global_ordinal_type zgno_t
int main()
This class represents a collection of global Identifiers and their associated weights,...
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
void getIDsView(const gno_t *&Ids) const
Provide a pointer to this process' identifiers.
size_t getLocalNumIDs() const
Returns the number of objects on this process.
A simple class that can be the User template argument for an InputAdapter.
static const std::string fail
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > userTypes_t
static ArrayRCP< ArrayRCP< zscalar_t > > weights