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();
145 if(idx<0 || idx >= numWeights_)
147 std::ostringstream emsg;
148 emsg << __FILE__ <<
":" << __LINE__
149 <<
" Invalid weight index " << idx << std::endl;
150 throw std::runtime_error(emsg.str());
154 weights_[idx].getStridedList(length,
weights, stride);
165 template <
typename Adapter>
169 template <
typename Adapter>
175 RCP<const User> invector_;
176 RCP<const x_mvector_t> vector_;
177 RCP<const Xpetra::Map<lno_t, gno_t, node_t> > map_;
180 ArrayRCP<StridedData<lno_t, scalar_t> > weights_;
187 template <
typename User>
189 const RCP<const User> &invector,
190 std::vector<const scalar_t *> &
weights, std::vector<int> &weightStrides):
191 invector_(invector), vector_(), map_(),
197 RCP<x_mvector_t> tmp =
199 vector_ = rcp_const_cast<const x_mvector_t>(tmp);
203 map_ = vector_->getMap();
205 size_t length = vector_->getLocalLength();
207 if (length > 0 && numWeights_ > 0){
209 for (
int w=0; w < numWeights_; w++){
210 if (weightStrides.size())
211 stride = weightStrides[w];
212 ArrayRCP<const scalar_t> wgtV(
weights[w], 0, stride*length,
false);
213 weights_[w] = input_t(wgtV, stride);
220 template <
typename User>
222 const RCP<const User> &invector):
223 invector_(invector), vector_(), map_(),
224 numWeights_(0), weights_()
227 RCP<x_mvector_t> tmp =
229 vector_ = rcp_const_cast<const x_mvector_t>(tmp);
233 map_ = vector_->getMap();
237 template <
typename User>
239 const scalar_t *&elements,
int &stride,
int idx)
const
244 if (map_->lib() == Xpetra::UseTpetra){
245 const xt_mvector_t *tvector =
246 dynamic_cast<const xt_mvector_t *
>(vector_.get());
248 vecsize = tvector->getLocalLength();
250 ArrayRCP<const scalar_t> data = tvector->getData(idx);
251 elements = data.get();
254 else if (map_->lib() == Xpetra::UseEpetra){
255 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
256 typedef Xpetra::EpetraMultiVectorT<gno_t,node_t> xe_mvector_t;
257 const xe_mvector_t *evector =
258 dynamic_cast<const xe_mvector_t *
>(vector_.get());
260 vecsize = evector->getLocalLength();
262 ArrayRCP<const double> data = evector->getData(idx);
266 elements =
reinterpret_cast<const scalar_t *
>(data.get());
269 throw std::logic_error(
"Epetra requested, but Trilinos is not "
270 "built with Epetra");
274 throw std::logic_error(
"invalid underlying lib");
279 template <
typename User>
280 template <
typename Adapter>
282 const User &in, User *&out,
287 ArrayRCP<gno_t> importList;
291 (solution,
this, importList);
297 importList.getRawPtr());
303 template <
typename User>
304 template <
typename Adapter>
306 const User &in, RCP<User> &out,
311 ArrayRCP<gno_t> importList;
315 (solution,
this, importList);
321 importList.getRawPtr());
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
Helper functions for Partitioning Problems.
This file defines the StridedData class.
Defines the VectorAdapter interface.
Traits of Xpetra classes, including migration method.
InputTraits< User >::part_t part_t
InputTraits< User >::scalar_t scalar_t
InputTraits< User >::lno_t lno_t
InputTraits< User >::gno_t gno_t
A PartitioningSolution is a solution to a partitioning problem.
The StridedData class manages lists of weights or coordinates.
VectorAdapter defines the interface for vector input.
An adapter for Xpetra::MultiVector.
void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const
Provide a pointer to the elements of the specified vector.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
void getIDsView(const gno_t *&ids) const
Provide a pointer to this process' identifiers.
XpetraMultiVectorAdapter(const RCP< const User > &invector, std::vector< const scalar_t * > &weights, std::vector< int > &weightStrides)
Constructor.
int getNumEntriesPerID() const
Return the number of vectors.
~XpetraMultiVectorAdapter()
Destructor.
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
size_t getLocalNumIDs() const
Returns the number of objects on this process.
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...
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...
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.