Zoltan2
Loading...
Searching...
No Matches
Zoltan2_Adapter.hpp
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
13
14#ifndef _ZOLTAN2_ADAPTER_HPP_
15#define _ZOLTAN2_ADAPTER_HPP_
16
17#include "Kokkos_StaticCrsGraph.hpp"
18#include <Kokkos_Core.hpp>
19#include <Zoltan2_Standards.hpp>
22#include <fstream>
23
24namespace Zoltan2 {
25
26//template<typename Adapter>
27//class PartitioningSolution;
28
40
49
51public:
52 virtual ~BaseAdapterRoot() = default; // required virtual declaration
53
59 virtual size_t getLocalNumIDs() const = 0;
60
66 virtual int getNumWeightsPerID() const { return 0; };
67};
68
69template <typename User>
71
72public:
79 using host_t = typename Kokkos::HostSpace::memory_space;
80 using device_t = typename node_t::device_type;
81
82 using ConstIdsDeviceView = Kokkos::View<const gno_t *, device_t>;
83 using ConstIdsHostView = typename ConstIdsDeviceView::HostMirror;
84
85 using IdsDeviceView = Kokkos::View<gno_t *, device_t>;
86 using IdsHostView = typename IdsDeviceView::HostMirror;
87
88 using ConstOffsetsDeviceView = Kokkos::View<const offset_t *, device_t>;
89 using ConstOffsetsHostView = typename ConstOffsetsDeviceView::HostMirror;
90
91 using OffsetsDeviceView = Kokkos::View<offset_t *, device_t>;
92 using OffsetsHostView = typename OffsetsDeviceView::HostMirror;
93
94 using ConstScalarsDeviceView = Kokkos::View<const scalar_t *, device_t>;
95 using ConstScalarsHostView = typename ConstScalarsDeviceView::HostMirror;
96
97 using ScalarsDeviceView = Kokkos::View<scalar_t *, device_t>;
98 using ScalarsHostView = typename ScalarsDeviceView::HostMirror;
99
100 using ConstWeightsDeviceView1D = Kokkos::View<const scalar_t *, device_t>;
101 using ConstWeightsHostView1D = typename ConstWeightsDeviceView1D::HostMirror;
102
103 using WeightsDeviceView1D = Kokkos::View<scalar_t *, device_t>;
104 using WeightsHostView1D = typename WeightsDeviceView1D::HostMirror;
105
106 using ConstWeightsDeviceView = Kokkos::View<const scalar_t **, device_t>;
107 using ConstWeightsHostView = typename ConstWeightsDeviceView::HostMirror;
108
109 using WeightsDeviceView = Kokkos::View<scalar_t **, device_t>;
110 using WeightsHostView = typename WeightsDeviceView::HostMirror;
111
114 virtual enum BaseAdapterType adapterType() const = 0;
115
121 virtual void getIDsView(const gno_t *&ids) const {
122 // If adapter does not define getIDsView, getIDsKokkosView is called.
123 // If adapter does not define getIDsKokkosView, getIDsView is called.
124 // Allows forward and backwards compatibility.
125 ConstIdsDeviceView kokkosIds;
126 getIDsKokkosView(kokkosIds);
127 //deep_copy(?)
128 ids = kokkosIds.data();
129 }
130
136 virtual void getIDsKokkosView(ConstIdsDeviceView &ids) const {
137 // If adapter does not define getIDsView, getIDsKokkosView is called.
138 // If adapter does not define getIDsKokkosView, getIDsView is called.
139 // Allows forward and backwards compatibility.
140 const gno_t * ptr_ids;
141 getIDsView(ptr_ids);
142
143 auto non_const_ids = IdsDeviceView("ptr_ids", getLocalNumIDs());
144 auto host_ids = Kokkos::create_mirror_view(non_const_ids);
145 for(size_t i = 0; i < this->getLocalNumIDs(); ++i) {
146 host_ids(i) = ptr_ids[i];
147 }
148 Kokkos::deep_copy(non_const_ids, host_ids);
149 ids = non_const_ids;
150 }
151
156 virtual void getIDsHostView(ConstIdsHostView& hostIds) const {
158 }
159
164 virtual void getIDsDeviceView(ConstIdsDeviceView &deviceIds) const {
166 }
167
169 // * \param wgt on return a pointer to the weights for this idx
170 // * \param stride on return, the value such that
171 // * the \t nth weight should be found at <tt> wgt[n*stride] </tt>.
172 // * \param idx the weight index, zero or greater
173 // * This function or getWeightsKokkosView must be implemented in
174 // * derived adapter if getNumWeightsPerID > 0.
175 // * This function should not be called if getNumWeightsPerID is zero.
176 // */
177 virtual void getWeightsView(const scalar_t *&wgt, int &stride,
178 int idx = 0) const {
179 // If adapter does not define getWeightsView, getWeightsKokkosView is called.
180 // If adapter does not define getWeightsKokkosView, getWeightsView is called.
181 // Allows forward and backwards compatibility.
182 Kokkos::View<scalar_t **, device_t> kokkos_wgts_2d;
183 getWeightsKokkosView(kokkos_wgts_2d);
184 Kokkos::View<scalar_t *, device_t> kokkos_wgts;
185 wgt = Kokkos::subview(kokkos_wgts_2d, Kokkos::ALL, idx).data();
186 stride = 1;
187 }
188
190 // * \param wgt on return a Kokkos view of the weights for this idx
191 // * This function or getWeightsView must be implemented in
192 // * derived adapter if getNumWeightsPerID > 0.
193 // * This function should not be called if getNumWeightsPerID is zero.
194 // */
195 virtual void getWeightsKokkosView(Kokkos::View<scalar_t **,
196 device_t> & wgt) const {
197 // If adapter does not define getWeightsKokkosView, getWeightsView is called.
198 // If adapter does not define getWeightsView, getWeightsKokkosView is called.
199 // Allows forward and backwards compatibility.
200 wgt = Kokkos::View<scalar_t **, device_t>(
202 auto host_wgt = Kokkos::create_mirror_view(wgt);
203 for(int j = 0; j < this->getNumWeightsPerID(); ++j) {
204 const scalar_t * ptr_wgts;
205 int stride;
206 getWeightsView(ptr_wgts, stride, j);
207 size_t i = 0;
208 for(size_t n = 0; n < this->getLocalNumIDs() * stride; n += stride) {
209 host_wgt(i++,j) = ptr_wgts[n];
210 }
211 }
212 Kokkos::deep_copy(wgt, host_wgt);
213 }
214
219 virtual void getWeightsHostView(WeightsHostView1D& hostWgts, int idx = 0) const {
221 }
222
226 virtual void getWeightsHostView(WeightsHostView& hostWgts) const {
228 }
229
234 virtual void getWeightsDeviceView(WeightsDeviceView1D& deviceWgts, int idx = 0) const {
236 }
237
241 virtual void getWeightsDeviceView(WeightsDeviceView& deviceWgts) const {
243 }
244
254 virtual void getPartsView(const part_t *&inputPart) const {
255 // Default behavior: return NULL for inputPart array;
256 // assume input part == rank
257 inputPart = NULL;
258 }
259
260 virtual void getPartsHostView(Kokkos::View<part_t*, host_t> &inputPart) const {
262 }
263
264 virtual void getPartsDeviceView(Kokkos::View<part_t*, device_t> &inputPart) const {
266 }
267
285 template <typename Adapter>
286 void applyPartitioningSolution(const User &in, User *&out,
287 const PartitioningSolution<Adapter> &solution) const {
289 }
290
291protected:
292
293 // Write Chaco-formatted graph and assign files echoing adapter input
294 // This routine is serial and may be slow; use it only for debugging
295 // This function does not write edge info to the graph file, as the
296 // BaseAdapter does not know about edge info; it writes
297 // only the Chaco header and vertex weights (if applicable).
298 void generateWeightFileOnly(const char* fileprefix,
299 const Teuchos::Comm<int> &comm) const;
300
301};
302
303template <typename User>
304class AdapterWithCoords : public BaseAdapter<User>
305{
306public:
309 using host_t = typename Kokkos::HostSpace::memory_space;
310
311 // Coordinates in MJ are LayoutLeft since Tpetra Multivector gives LayoutLeft
312 using CoordsDeviceView = Kokkos::View<scalar_t **, Kokkos::LayoutLeft, device_t>;
313 using CoordsHostView = typename CoordsDeviceView::HostMirror;
314
315public:
316 virtual void getCoordinatesView(const scalar_t *&coords, int &stride,
317 int coordDim) const = 0;
318 virtual void getCoordinatesKokkosView(CoordsDeviceView &elements) const = 0;
319
323 virtual void getCoordinatesDeviceView(CoordsDeviceView &elements) const {
325 }
326};
327
328// Forward declare
329template <typename User>
330class VectorAdapter;
331
332template <typename User, typename UserCoord=User>
334{
335public:
336 virtual void setCoordinateInput(VectorAdapter<UserCoord> *coordData) = 0;
338};
339
340template <typename User>
342 const char *fileprefix,
343 const Teuchos::Comm<int> &comm
344) const
345{
346 int np = comm.getSize();
347 int me = comm.getRank();
348
349 size_t nLocalIDs = this->getLocalNumIDs();
350
351 // Write .graph file: header and weights only (no edges)
352 // Adapters with edges have to implement their own generateFiles function
353 // to provide edge info.
354 {
355 // append suffix to filename
356 std::string filenamestr = fileprefix;
357 filenamestr = filenamestr + ".graph";
358 const char *filename = filenamestr.c_str();
359
360 size_t nGlobalIDs;
361 Teuchos::reduceAll(comm, Teuchos::REDUCE_SUM, 1, &nLocalIDs, &nGlobalIDs);
362
363 int nWgts = this->getNumWeightsPerID();
364
365 for (int p = 0; p < np; p++) {
366
367 // Is it this processor's turn to write to files?
368 if (me == p) {
369
370 std::ofstream fp;
371
372 if (me == 0) {
373 // open file for writing
374 fp.open(filename, std::ios::out);
375 // write Chaco header info
376 // this function assumes no edges
377 fp << nGlobalIDs << " " << 0 << " "
378 << (nWgts ? "010" : "000") << " "
379 << (nWgts > 1 ? std::to_string(nWgts) : " ") << std::endl;
380 }
381 else {
382 // open file for appending
383 fp.open(filename, std::ios::app);
384 }
385
386 if (nWgts) {
387
388 // get weight data
389 const scalar_t **wgts = new const scalar_t *[nWgts];
390 int *strides = new int[nWgts];
391 for (int n = 0; n < nWgts; n++)
392 getWeightsView(wgts[n], strides[n], n);
393
394 // write weights to file
395 for (size_t i = 0; i < nLocalIDs; i++) {
396 for (int n = 0; n < nWgts; n++)
397 fp << wgts[n][i*strides[n]] << " ";
398 fp << "\n";
399 }
400
401 delete [] strides;
402 delete [] wgts;
403 }
404
405 fp.close();
406 }
407
408 comm.barrier();
409 }
410 }
411
412 // write assignments file
413 {
414 std::string filenamestr = fileprefix;
415 filenamestr = filenamestr + ".assign";
416 const char *filename = filenamestr.c_str();
417
418 for (int p = 0; p < np; p++) {
419
420 // Is it this processor's turn to write to files?
421 if (me == p) {
422
423 std::ofstream fp;
424
425 if (me == 0) {
426 // open file for writing
427 fp.open(filename, std::ios::out);
428 }
429 else {
430 // open file for appending
431 fp.open(filename, std::ios::app);
432 }
433
434 const part_t *parts;
435 this->getPartsView(parts);
436
437 for (size_t i = 0; i < nLocalIDs; i++) {
438 fp << (parts != NULL ? parts[i] : me) << "\n";
439 }
440 fp.close();
441 }
442
443 comm.barrier();
444 }
445 }
446}
447
448} //namespace Zoltan2
449
450#endif
#define Z2_THROW_NOT_IMPLEMENTED
Traits for application input objects.
Defines the PartitioningSolution class.
Gathering definitions used in software development.
virtual VectorAdapter< UserCoord > * getCoordinateInput() const =0
virtual void setCoordinateInput(VectorAdapter< UserCoord > *coordData)=0
virtual void getCoordinatesView(const scalar_t *&coords, int &stride, int coordDim) const =0
virtual void getCoordinatesDeviceView(CoordsDeviceView &elements) const
typename BaseAdapter< User >::node_t::device_type device_t
typename CoordsDeviceView::HostMirror CoordsHostView
Kokkos::View< scalar_t **, Kokkos::LayoutLeft, device_t > CoordsDeviceView
typename BaseAdapter< User >::scalar_t scalar_t
virtual void getCoordinatesKokkosView(CoordsDeviceView &elements) const =0
typename Kokkos::HostSpace::memory_space host_t
virtual void getCoordinatesHostView(CoordsHostView &) const
BaseAdapter defines methods required by all Adapters.
virtual ~BaseAdapterRoot()=default
virtual int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
virtual size_t getLocalNumIDs() const =0
Returns the number of objects on this process.
virtual void getWeightsHostView(WeightsHostView &hostWgts) const
Provide a Kokkos view (Host side) of the weights.
typename InputTraits< User >::node_t node_t
Kokkos::View< const scalar_t **, device_t > ConstWeightsDeviceView
Kokkos::View< const scalar_t *, device_t > ConstWeightsDeviceView1D
typename InputTraits< User >::lno_t lno_t
virtual void getPartsHostView(Kokkos::View< part_t *, host_t > &inputPart) const
Kokkos::View< const offset_t *, device_t > ConstOffsetsDeviceView
typename WeightsDeviceView1D::HostMirror WeightsHostView1D
void generateWeightFileOnly(const char *fileprefix, const Teuchos::Comm< int > &comm) const
Kokkos::View< const scalar_t *, device_t > ConstScalarsDeviceView
virtual void getIDsView(const gno_t *&ids) const
Provide a pointer to this process' identifiers.
Kokkos::View< const gno_t *, device_t > ConstIdsDeviceView
virtual void getWeightsHostView(WeightsHostView1D &hostWgts, int idx=0) const
Provide a Kokkos view (Host side) of the weights.
typename OffsetsDeviceView::HostMirror OffsetsHostView
typename ScalarsDeviceView::HostMirror ScalarsHostView
typename ConstOffsetsDeviceView::HostMirror ConstOffsetsHostView
virtual void getPartsView(const part_t *&inputPart) const
Provide pointer to an array containing the input part assignment for each ID. The input part informat...
virtual void getIDsKokkosView(ConstIdsDeviceView &ids) const
Provide a Kokkos view to this process' identifiers.
typename InputTraits< User >::scalar_t scalar_t
typename ConstWeightsDeviceView1D::HostMirror ConstWeightsHostView1D
typename WeightsDeviceView::HostMirror WeightsHostView
typename InputTraits< User >::gno_t gno_t
virtual void getPartsDeviceView(Kokkos::View< part_t *, device_t > &inputPart) const
typename ConstScalarsDeviceView::HostMirror ConstScalarsHostView
typename Kokkos::HostSpace::memory_space host_t
typename ConstIdsDeviceView::HostMirror ConstIdsHostView
Kokkos::View< scalar_t *, device_t > ScalarsDeviceView
Kokkos::View< scalar_t *, device_t > WeightsDeviceView1D
virtual enum BaseAdapterType adapterType() const =0
Returns the type of adapter.
Kokkos::View< gno_t *, device_t > IdsDeviceView
typename InputTraits< User >::offset_t offset_t
virtual void getWeightsDeviceView(WeightsDeviceView1D &deviceWgts, int idx=0) const
Provide a Kokkos view (Device side) of the weights.
virtual void getWeightsView(const scalar_t *&wgt, int &stride, int idx=0) const
Provide pointer to a weight array with stride.
typename IdsDeviceView::HostMirror IdsHostView
virtual void getWeightsKokkosView(Kokkos::View< scalar_t **, device_t > &wgt) const
Provide kokkos view of weights.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
Apply a PartitioningSolution to an input.
typename InputTraits< User >::part_t part_t
Kokkos::View< offset_t *, device_t > OffsetsDeviceView
virtual void getIDsHostView(ConstIdsHostView &hostIds) const
Provide a Kokkos view (Host side) to this process' identifiers.
virtual void getWeightsDeviceView(WeightsDeviceView &deviceWgts) const
Provide a Kokkos view (Device side) of the weights.
typename ConstWeightsDeviceView::HostMirror ConstWeightsHostView
Kokkos::View< scalar_t **, device_t > WeightsDeviceView
typename node_t::device_type device_t
virtual void getIDsDeviceView(ConstIdsDeviceView &deviceIds) const
Provide a Kokkos view (Device side) to this process' identifiers.
A PartitioningSolution is a solution to a partitioning problem.
VectorAdapter defines the interface for vector input.
Created by mbenlioglu on Aug 31, 2020.
BaseAdapterType
An enum to identify general types of adapters.
@ VectorAdapterType
vector data
@ InvalidAdapterType
unused value
@ GraphAdapterType
graph data
@ MatrixAdapterType
matrix data
@ MeshAdapterType
mesh data
@ IdentifierAdapterType
identifier data, just a list of IDs
default_offset_t offset_t
The data type to represent offsets.
default_gno_t gno_t
The ordinal type (e.g., int, long, int64_t) that can represent global counts and identifiers.
default_node_t node_t
The Kokkos node type. This is only meaningful for users of Tpetra objects.
default_lno_t lno_t
The ordinal type (e.g., int, long, int64_t) that represents local counts and local indices.
default_part_t part_t
The data type to represent part numbers.
default_scalar_t scalar_t
The data type for weights and coordinates.