Skip to content

Commit

Permalink
Merge remote-tracking branch 'jakob/JD-Zoltan2' into zoltan-testing
Browse files Browse the repository at this point in the history
  • Loading branch information
Matt Bettencourt committed Jul 19, 2021
2 parents 6cee551 + 710c296 commit b5b4845
Show file tree
Hide file tree
Showing 3 changed files with 108 additions and 52 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -62,7 +62,7 @@ namespace Zoltan2 {
* and their associated weights, if any.
*
* The user supplies the identifiers and weights by way of pointers
* to arrays.
* to arrays.
*
* The template parameter (\c User) is a C++ class type which provides the
* actual data types with which the Zoltan2 library will be compiled, through
Expand Down Expand Up @@ -112,40 +112,67 @@ template <typename User>
// The Adapter interface.
////////////////////////////////////////////////////////////////

size_t getLocalNumIDs() const {
size_t getLocalNumIDs() const override {
return idsView_.extent(0);
}

void getIDsView(const gno_t *&ids) const override {
auto kokkosIds = idsView_.view_host();
ids = kokkosIds.data();
}

void getIDsKokkosView(Kokkos::View<const gno_t *,
typename node_t::device_type> &ids) const {
ids = idsView_;
typename node_t::device_type> &ids) const override {
ids = idsView_.template view<typename node_t::device_type>();
}

int getNumWeightsPerID() const {
int getNumWeightsPerID() const override {
return weightsView_.extent(1);
}

void getWeightsView(const scalar_t *&wgt, int &stride,
int idx = 0) const override
{
auto h_wgts_2d = weightsView_.view_host();

wgt = Kokkos::subview(h_wgts_2d, Kokkos::ALL, idx).data();
stride = 1;
}

void getWeightsKokkosView(Kokkos::View<scalar_t **,
typename node_t::device_type> &wgts) const {
wgts = weightsView_;
typename node_t::device_type> &wgts) const override {
wgts = weightsView_.template view<typename node_t::device_type>();
}

private:
Kokkos::View<const gno_t *, typename node_t::device_type> idsView_;
Kokkos::View<scalar_t **, typename node_t::device_type> weightsView_;
Kokkos::DualView<gno_t *> idsView_;
Kokkos::DualView<scalar_t **> weightsView_;
};

////////////////////////////////////////////////////////////////
// Definitions
////////////////////////////////////////////////////////////////

template <typename User>
BasicKokkosIdentifierAdapter<User>::BasicKokkosIdentifierAdapter(
Kokkos::View<gno_t*, typename node_t::device_type> &ids,
Kokkos::View<scalar_t**, typename node_t::device_type> &weights):
idsView_(ids), weightsView_(weights) {
BasicKokkosIdentifierAdapter<User>::BasicKokkosIdentifierAdapter(
Kokkos::View<gno_t *, typename node_t::device_type> &ids,
Kokkos::View<scalar_t **, typename node_t::device_type> &weights)
{
idsView_ = Kokkos::DualView<gno_t *>("idsView_", ids.extent(0));
Kokkos::deep_copy(idsView_.h_view, ids);

weightsView_ = Kokkos::DualView<scalar_t **>("weightsView_", weights.extent(0), weights.extent(1));
Kokkos::deep_copy(weightsView_.h_view, weights);

weightsView_.modify_host();
weightsView_.sync_host();
weightsView_.template sync<typename node_t::device_type>();

idsView_.modify_host();
idsView_.sync_host();
idsView_.template sync<typename node_t::device_type>();
}

} //namespace Zoltan2

#endif
6 changes: 5 additions & 1 deletion packages/zoltan2/example/block/kokkosBlock.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -113,10 +113,14 @@ int main(int narg, char *arg[]) {
const int nWeights = 1;
Kokkos::View<scalar_t **, typename node_t::device_type>
weights("weights", localCount, nWeights);
auto host_weights = Kokkos::create_mirror_view(weights);

for (int index = 0; index < localCount; index++) {
weights(index, 0) = 1; // Error check relies on uniform weights
host_weights(index, 0) = 1; // Error check relies on uniform weights
}

Kokkos::deep_copy(weights, host_weights);

inputAdapter_t ia(globalIds, weights);

/////////////////////////////////////////////////////////////////////////
Expand Down
97 changes: 61 additions & 36 deletions packages/zoltan2/test/core/partition/mj_backwardcompat.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -79,7 +79,7 @@ class OldSchoolVectorAdapterContig : public Zoltan2::VectorAdapter<User>

int getNumWeightsPerID() const { return (weights != NULL); }

void getWeightsView(const scalar_t *&wgt, int &stride, int idx = 0) const
void getWeightsView(const scalar_t *&wgt, int &stride, int idx = 0) const
{ wgt = weights; stride = 1; }

int getNumEntriesPerID() const { return dim; }
Expand Down Expand Up @@ -122,7 +122,7 @@ class OldSchoolVectorAdapterStrided : public Zoltan2::VectorAdapter<User>

int getNumWeightsPerID() const { return (weights != NULL); }

void getWeightsView(const scalar_t *&wgt, int &stride, int idx = 0) const
void getWeightsView(const scalar_t *&wgt, int &stride, int idx = 0) const
{ wgt = weights; stride = 1; }

int getNumEntriesPerID() const { return dim; }
Expand Down Expand Up @@ -163,77 +163,102 @@ class KokkosVectorAdapter : public Zoltan2::VectorAdapter<User>
{
// create kokkos_gids in default memory space
{
typedef Kokkos::View<gno_t *, typename node_t::device_type> view_ids_t;
view_ids_t gids = view_ids_t(
typedef Kokkos::DualView<gno_t *> view_ids_t;
kokkos_gids = view_ids_t(
Kokkos::ViewAllocateWithoutInitializing("gids"), nids);
typename view_ids_t::HostMirror host_gids =
Kokkos::create_mirror_view(gids);

auto host_gids = kokkos_gids.h_view;
for(size_t n = 0; n < nids; ++n) {
host_gids(n) = gids_[n];
}
Kokkos::deep_copy(gids, host_gids);
kokkos_gids = gids; // to const View

kokkos_gids.template modify<typename view_ids_t::host_mirror_space>();
kokkos_gids.sync_host();
kokkos_gids.template sync<typename node_t::device_type>();
}

// create kokkos_weights in default memory space
if(weights_ != NULL)
{
typedef Kokkos::View<scalar_t **, typename node_t::device_type> view_weights_t;
typedef Kokkos::DualView<scalar_t **> view_weights_t;
kokkos_weights = view_weights_t(
Kokkos::ViewAllocateWithoutInitializing("weights"), nids, 0);
typename view_weights_t::HostMirror host_kokkos_weights =
Kokkos::create_mirror_view(kokkos_weights);
auto host_kokkos_weights = kokkos_weights.h_view;
for(size_t n = 0; n < nids; ++n) {
host_kokkos_weights(n,0) = weights_[n];
}
Kokkos::deep_copy(kokkos_weights, host_kokkos_weights);

kokkos_weights.template modify<typename view_weights_t::host_mirror_space>();
kokkos_weights.sync_host();
kokkos_weights.template sync<typename node_t::device_type>();
}

// create kokkos_coords in default memory space
{
typedef Kokkos::View<scalar_t **, Kokkos::LayoutLeft,
typename node_t::device_type> kokkos_coords_t;
typedef Kokkos::DualView<scalar_t **> kokkos_coords_t;
kokkos_coords = kokkos_coords_t(
Kokkos::ViewAllocateWithoutInitializing("coords"), nids, dim);
typename kokkos_coords_t::HostMirror host_kokkos_coords =
Kokkos::create_mirror_view(kokkos_coords);
auto host_kokkos_coords = kokkos_coords.h_view;

for(size_t n = 0; n < nids; ++n) {
for(int idx = 0; idx < dim; ++idx) {
host_kokkos_coords(n,idx) = coords_[n+idx*nids];
}
}
Kokkos::deep_copy(kokkos_coords, host_kokkos_coords);

kokkos_coords.template modify<typename kokkos_coords_t::host_mirror_space>();
kokkos_coords.sync_host();
kokkos_coords.template sync<typename node_t::device_type>();
}
}

size_t getLocalNumIDs() const { return nids; }

void getIDsView(const gno_t *&ids) const override {
auto kokkosIds = kokkos_gids.view_host();
ids = kokkosIds.data();
}

virtual void getIDsKokkosView(Kokkos::View<const gno_t *,
typename node_t::device_type> &ids) const {
ids = kokkos_gids;
ids = kokkos_gids.template view<typename node_t::device_type>();;
}

int getNumWeightsPerID() const { return (kokkos_weights.size() != 0); }
int getNumWeightsPerID() const { return (kokkos_weights.h_view.size() != 0); }

void getWeightsView(const scalar_t *&wgt, int &stride,
int idx = 0) const override
{
auto h_wgts_2d = kokkos_weights.view_host();

wgt = Kokkos::subview(h_wgts_2d, Kokkos::ALL, idx).data();
stride = 1;
}

virtual void getWeightsKokkosView(Kokkos::View<scalar_t **,
typename node_t::device_type> & wgt) const {
wgt = kokkos_weights;
wgt = kokkos_weights.template view<typename node_t::device_type>();;
}

int getNumEntriesPerID() const { return dim; }

void getEntriesView(const scalar_t *&elements,
int &stride, int idx = 0) const override {
elements = kokkos_coords.view_host().data();
stride = 1;
}

virtual void getEntriesKokkosView(Kokkos::View<scalar_t **,
Kokkos::LayoutLeft, typename node_t::device_type> & coo) const {
coo = kokkos_coords;
coo = kokkos_coords.template view<typename node_t::device_type>();
}

private:
const size_t nids;
Kokkos::View<const gno_t *, typename node_t::device_type> kokkos_gids;
Kokkos::DualView<gno_t *> kokkos_gids;
const int dim;
Kokkos::View<scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type>
kokkos_coords;
Kokkos::View<scalar_t **, typename node_t::device_type> kokkos_weights;
Kokkos::DualView<scalar_t **> kokkos_coords;
Kokkos::DualView<scalar_t **> kokkos_weights;
};

//////////////////////////////////////////////
Expand Down Expand Up @@ -281,7 +306,7 @@ int run_test_strided_versus_contig(const std::string & algorithm) {
globalId_t *globalIds = new globalId_t [localCount];
globalId_t offset = rank * localCount;
for (size_t i=0; i < localCount; i++) globalIds[i] = offset++;

///////////////////////////////////////////////////////////////////////
// Create parameters for an MJ problem

Expand All @@ -296,12 +321,12 @@ int run_test_strided_versus_contig(const std::string & algorithm) {
// Test one: No weights

// Partition using strided coords
stridedAdapter_t *ia1 =
stridedAdapter_t *ia1 =
new stridedAdapter_t(localCount,globalIds,dim,cStrided);

Zoltan2::PartitioningProblem<stridedAdapter_t> *problem1 =
new Zoltan2::PartitioningProblem<stridedAdapter_t>(ia1, &params);

problem1->solve();

quality_t *metricObject1 = new quality_t(ia1, &params, comm,
Expand All @@ -325,7 +350,7 @@ int run_test_strided_versus_contig(const std::string & algorithm) {

Zoltan2::PartitioningProblem<contigAdapter_t> *problem2 =
new Zoltan2::PartitioningProblem<contigAdapter_t>(ia2, &params);

problem2->solve();

// Partition using contiguous coords to generate kokkos adapter
Expand All @@ -344,7 +369,7 @@ int run_test_strided_versus_contig(const std::string & algorithm) {
(problem2->getSolution().getPartListView()[i] !=
problem3->getSolution().getPartListView()[i]))
{
std::cout << rank << " Error: differing parts for index " << i
std::cout << rank << " Error: differing parts for index " << i
<< problem1->getSolution().getPartListView()[i] << " "
<< problem2->getSolution().getPartListView()[i] << " "
<< problem3->getSolution().getPartListView()[i] << std::endl;
Expand All @@ -362,11 +387,11 @@ int run_test_strided_versus_contig(const std::string & algorithm) {
delete ia1;
delete ia2;
delete ia3;

///////////////////////////////////////////////////////////////////////
// Test two: weighted
// Create a Zoltan2 input adapter that includes weights.

scalar_t *weights = new scalar_t [localCount];
for (size_t i=0; i < localCount; i++) weights[i] = 1 + scalar_t(rank);

Expand Down Expand Up @@ -397,18 +422,18 @@ int run_test_strided_versus_contig(const std::string & algorithm) {
ia2 = new contigAdapter_t(localCount, globalIds, dim, cContig, weights);

problem2 = new Zoltan2::PartitioningProblem<contigAdapter_t>(ia2, &params);

problem2->solve();

// compare strided vs contiguous
ndiff = 0;
for (size_t i = 0; i < localCount; i++) {
if (problem1->getSolution().getPartListView()[i] !=
if (problem1->getSolution().getPartListView()[i] !=
problem2->getSolution().getPartListView()[i]) {
std::cout << rank << " Error: differing parts for index " << i
std::cout << rank << " Error: differing parts for index " << i
<< problem1->getSolution().getPartListView()[i] << " "
<< problem2->getSolution().getPartListView()[i] << std::endl;

ndiff++;
}
}
Expand Down

0 comments on commit b5b4845

Please sign in to comment.