50 #ifndef _ZOLTAN2_XPETRAMULTIVECTORADAPTER_HPP_ 51 #define _ZOLTAN2_XPETRAMULTIVECTORADAPTER_HPP_ 58 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA) 59 #include <Xpetra_EpetraMultiVector.hpp> 61 #include <Xpetra_TpetraMultiVector.hpp> 82 template <
typename User>
86 #ifndef DOXYGEN_SHOULD_SKIP_THIS 93 typedef User userCoord_t;
95 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
96 typedef Xpetra::TpetraMultiVector<
120 std::vector<const scalar_t *> &
weights, std::vector<int> &weightStrides);
138 ids = map_->getNodeElementList().getRawPtr();
142 Kokkos::View<const gno_t *, typename node_t::device_type> &ids)
const {
143 if (map_->lib() == Xpetra::UseTpetra) {
144 using device_type =
typename node_t::device_type;
145 const xt_mvector_t *tvector =
146 dynamic_cast<const xt_mvector_t *
>(vector_.get());
154 ids = Kokkos::create_mirror_view_and_copy(device_type(),
155 tvector->getTpetra_MultiVector()->getMap()->getMyGlobalIndices());
157 else if (map_->lib() == Xpetra::UseEpetra) {
158 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA) 162 throw std::logic_error(
"Epetra requested, but Trilinos is not " 163 "built with Epetra");
167 throw std::logic_error(
"getIDsKokkosView called but not on Tpetra or Epetra!");
175 if(idx<0 || idx >= numWeights_)
177 std::ostringstream emsg;
178 emsg << __FILE__ <<
":" << __LINE__
179 <<
" Invalid weight index " << idx << std::endl;
180 throw std::runtime_error(emsg.str());
184 weights_[idx].getStridedList(length,
weights, stride);
188 typename node_t::device_type> &wgt)
const {
189 typedef Kokkos::View<scalar_t**, typename node_t::device_type> view_t;
190 wgt = view_t(
"wgts", vector_->getLocalLength(), numWeights_);
191 typename view_t::HostMirror host_wgt = Kokkos::create_mirror_view(wgt);
192 for(
int idx = 0; idx < numWeights_; ++idx) {
196 weights_[idx].getStridedList(length,
weights, stride);
197 size_t fill_index = 0;
198 for(
size_t n = 0; n < length; n += stride) {
199 host_wgt(fill_index++,idx) =
weights[n];
202 Kokkos::deep_copy(wgt, host_wgt);
215 Kokkos::View<
scalar_t **, Kokkos::LayoutLeft,
216 typename node_t::device_type> & elements)
const;
218 template <
typename Adapter>
222 template <
typename Adapter>
228 RCP<const User> invector_;
229 RCP<const x_mvector_t> vector_;
230 RCP<const Xpetra::Map<lno_t, gno_t, node_t> > map_;
233 ArrayRCP<StridedData<lno_t, scalar_t> > weights_;
240 template <
typename User>
242 const RCP<const User> &invector,
243 std::vector<const scalar_t *> &
weights, std::vector<int> &weightStrides):
244 invector_(invector), vector_(), map_(),
250 RCP<x_mvector_t> tmp =
252 vector_ = rcp_const_cast<
const x_mvector_t>(tmp);
256 map_ = vector_->getMap();
258 size_t length = vector_->getLocalLength();
260 if (length > 0 && numWeights_ > 0){
262 for (
int w=0; w < numWeights_; w++){
263 if (weightStrides.size())
264 stride = weightStrides[w];
265 ArrayRCP<const scalar_t> wgtV(
weights[w], 0, stride*length,
false);
266 weights_[w] = input_t(wgtV, stride);
273 template <
typename User>
275 const RCP<const User> &invector):
276 invector_(invector), vector_(), map_(),
277 numWeights_(0), weights_()
280 RCP<x_mvector_t> tmp =
282 vector_ = rcp_const_cast<
const x_mvector_t>(tmp);
286 map_ = vector_->getMap();
290 template <
typename User>
292 const scalar_t *&elements,
int &stride,
int idx)
const 297 if (map_->lib() == Xpetra::UseTpetra){
298 const xt_mvector_t *tvector =
299 dynamic_cast<const xt_mvector_t *
>(vector_.get());
301 vecsize = tvector->getLocalLength();
303 ArrayRCP<const scalar_t> data = tvector->getData(idx);
304 elements = data.get();
307 else if (map_->lib() == Xpetra::UseEpetra){
308 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA) 309 typedef Xpetra::EpetraMultiVectorT<gno_t,node_t> xe_mvector_t;
310 const xe_mvector_t *evector =
311 dynamic_cast<const xe_mvector_t *
>(vector_.get());
313 vecsize = evector->getLocalLength();
315 ArrayRCP<const double> data = evector->getData(idx);
319 elements =
reinterpret_cast<const scalar_t *
>(data.get());
322 throw std::logic_error(
"Epetra requested, but Trilinos is not " 323 "built with Epetra");
327 throw std::logic_error(
"invalid underlying lib");
332 template <
typename User>
335 Kokkos::View<scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type> & elements)
const 337 if (map_->lib() == Xpetra::UseTpetra){
338 const xt_mvector_t *tvector =
339 dynamic_cast<const xt_mvector_t *
>(vector_.get());
341 Kokkos::View<scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type> view2d =
342 tvector->getTpetra_MultiVector()->template getLocalView<node_t>(Tpetra::Access::ReadWrite);
345 else if (map_->lib() == Xpetra::UseEpetra){
346 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA) 347 typedef Xpetra::EpetraMultiVectorT<gno_t,node_t> xe_mvector_t;
348 const xe_mvector_t *evector =
349 dynamic_cast<const xe_mvector_t *
>(vector_.get());
351 Kokkos::View<scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type>
352 (
"elements", evector->getLocalLength(), evector->getNumVectors());
353 if(evector->getLocalLength() > 0) {
354 for(
size_t idx = 0; idx < evector->getNumVectors(); ++idx) {
355 const scalar_t * ptr;
357 getEntriesView(ptr, stride, idx);
358 for(
size_t n = 0; n < evector->getLocalLength(); ++n) {
359 elements(n, idx) = ptr[n];
364 throw std::logic_error(
"Epetra requested, but Trilinos is not " 365 "built with Epetra");
369 throw std::logic_error(
"getEntriesKokkosView called but not using Tpetra or Epetra!");
374 template <
typename User>
375 template <
typename Adapter>
377 const User &in, User *&out,
382 ArrayRCP<gno_t> importList;
386 (solution,
this, importList);
392 importList.getRawPtr());
398 template <
typename User>
399 template <
typename Adapter>
401 const User &in, RCP<User> &out,
406 ArrayRCP<gno_t> importList;
410 (solution,
this, importList);
416 importList.getRawPtr());
InputTraits< User >::scalar_t scalar_t
Created by mbenlioglu on Aug 31, 2020.
Helper functions for Partitioning Problems.
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
InputTraits< User >::gno_t gno_t
virtual void getIDsKokkosView(Kokkos::View< const gno_t *, typename node_t::device_type > &ids) const
Provide a Kokkos view to this process' identifiers.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
map_t::global_ordinal_type gno_t
InputTraits< User >::lno_t lno_t
Defines the VectorAdapter interface.
static RCP< User > doMigration(const User &from, size_t numLocalRows, const gno_t *myNewRows)
Migrate the object Given a user object and a new row distribution, create and return a new user objec...
Traits of Xpetra classes, including migration method.
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.
size_t getImportList(const PartitioningSolution< SolutionAdapter > &solution, const DataAdapter *const data, ArrayRCP< typename DataAdapter::gno_t > &imports)
From a PartitioningSolution, get a list of IDs to be imported. Assumes part numbers in PartitioningSo...
size_t getLocalNumIDs() const
Returns the number of objects on this process.
void getWeightsKokkos2dView(Kokkos::View< scalar_t **, typename node_t::device_type > &wgt) const
void getIDsKokkosView(Kokkos::View< const gno_t *, typename node_t::device_type > &ids) const
A PartitioningSolution is a solution to a partitioning problem.
int getNumEntriesPerID() const
Return the number of vectors.
void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const
void getEntriesKokkosView(Kokkos::View< scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type > &elements) const
VectorAdapter defines the interface for vector input.
The StridedData class manages lists of weights or coordinates.
InputTraits< User >::node_t node_t
~XpetraMultiVectorAdapter()
Destructor.
InputTraits< User >::part_t part_t
An adapter for Xpetra::MultiVector.
void getIDsView(const gno_t *&ids) const
XpetraMultiVectorAdapter(const RCP< const User > &invector, std::vector< const scalar_t *> &weights, std::vector< int > &weightStrides)
Constructor.
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
This file defines the StridedData class.
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t