14 #ifndef _ZOLTAN2_IDENTIFIERMODEL_HPP_
15 #define _ZOLTAN2_IDENTIFIERMODEL_HPP_
33 template <
typename Adapter>
38 #ifndef DOXYGEN_SHOULD_SKIP_THIS
39 typedef typename Adapter::scalar_t scalar_t;
54 const RCP<const Environment> &env,
55 const RCP<
const Comm<int> > &comm,
modelFlag_t &modelFlags);
64 return numGlobalIdentifiers_;
76 ArrayView<input_t> &wgts)
const
78 Ids = ArrayView<const gno_t>();
79 wgts = weights_.view(0, nUserWeights_);
82 Ids = ArrayView<const gno_t>(
83 reinterpret_cast<const gno_t*
>(gids_.getRawPtr()), n);
89 Kokkos::View<const gno_t *, typename node_t::device_type> &Ids,
90 Kokkos::View<scalar_t **, typename node_t::device_type> &wgts)
const {
92 adapter_->getIDsKokkosView(Ids);
93 adapter_->getWeightsKokkosView(wgts);
109 gno_t numGlobalIdentifiers_;
110 const RCP<const Environment> env_;
111 const RCP<const Comm<int> > comm_;
112 ArrayRCP<const gno_t> gids_;
113 const RCP<const Adapter> adapter_;
115 ArrayRCP<input_t> weights_;
119 template <
typename Adapter>
121 const RCP<const Adapter> &ia,
122 const RCP<const Environment> &env,
123 const RCP<
const Comm<int> > &comm,
125 numGlobalIdentifiers_(), env_(env), comm_(comm),
126 gids_(), adapter_(ia), nUserWeights_(0), weights_()
129 size_t nLocalIds = ia->getLocalNumIDs();
130 gno_t lsum = nLocalIds;
131 reduceAll<int, gno_t>(*comm_, Teuchos::REDUCE_SUM, 1, &lsum,
132 &numGlobalIdentifiers_);
135 int tmp = ia->getNumWeightsPerID();
137 Teuchos::reduceAll<int, int>(*comm, Teuchos::REDUCE_MAX, 1,
138 &tmp, &nUserWeights_);
143 Array<const scalar_t *> wgts(nUserWeights_, (
const scalar_t *)NULL);
144 Array<int> wgtStrides(nUserWeights_, 0);
146 if (nUserWeights_ > 0){
147 input_t *w =
new input_t [nUserWeights_];
148 weights_ = arcp<input_t>(w, 0, nUserWeights_);
151 const gno_t *gids=NULL;
155 ia->getIDsView(gids);
156 for (
int idx=0; idx < nUserWeights_; idx++)
157 ia->getWeightsView(wgts[idx], wgtStrides[idx], idx);
162 gids_ = arcp(gids, 0, nLocalIds,
false);
164 if (nUserWeights_ > 0){
165 for (
int idx=0; idx < nUserWeights_; idx++){
166 ArrayRCP<const scalar_t> wgtArray(wgts[idx], 0,
167 nLocalIds*wgtStrides[idx],
false);
168 weights_[idx] = input_t(wgtArray, wgtStrides[idx]);
173 env_->memory(
"After construction of identifier model");
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
Defines the Model interface.
std::bitset< NUM_MODEL_FLAGS > modelFlag_t
map_t::global_ordinal_type gno_t
The StridedData class manages lists of weights or coordinates.
map_t::local_ordinal_type lno_t
size_t getLocalNumObjects() const
Return the local number of objects.
size_t getIdentifierListKokkos(Kokkos::View< const gno_t *, typename node_t::device_type > &Ids, Kokkos::View< scalar_t **, typename node_t::device_type > &wgts) const
The base class for all model classes.
size_t getLocalNumIdentifiers() const
IdentifierModel defines the interface for all identifier models.
Tpetra::global_size_t global_size_t
size_t getGlobalNumObjects() const
Return the global number of objects.
global_size_t getGlobalNumIdentifiers() const
size_t getIdentifierList(ArrayView< const gno_t > &Ids, ArrayView< input_t > &wgts) const
IdentifierModel(const RCP< const Adapter > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &modelFlags)
Constructor.
This file defines the StridedData class.