Zoltan2
Loading...
Searching...
No Matches
Zoltan2_IdentifierModel.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
14#ifndef _ZOLTAN2_IDENTIFIERMODEL_HPP_
15#define _ZOLTAN2_IDENTIFIERMODEL_HPP_
16
17#include <Zoltan2_Adapter.hpp>
18#include <Zoltan2_Model.hpp>
20
21namespace Zoltan2 {
22
33template <typename Adapter>
34class IdentifierModel : public Model<Adapter>
35{
36public:
37
38#ifndef DOXYGEN_SHOULD_SKIP_THIS
39 typedef typename Adapter::scalar_t scalar_t;
40 typedef typename Adapter::gno_t gno_t;
41 typedef typename Adapter::lno_t lno_t;
42 typedef typename Adapter::node_t node_t;
43 typedef StridedData<lno_t, scalar_t> input_t;
44#endif
45
53 IdentifierModel(const RCP<const Adapter> &ia,
54 const RCP<const Environment> &env,
55 const RCP<const Comm<int> > &comm, modelFlag_t &modelFlags);
56
59 inline size_t getLocalNumIdentifiers() const { return gids_.size(); }
60
64 return numGlobalIdentifiers_;
65 }
66
75 inline size_t getIdentifierList(ArrayView<const gno_t> &Ids,
76 ArrayView<input_t> &wgts) const
77 {
78 Ids = ArrayView<const gno_t>();
79 wgts = weights_.view(0, nUserWeights_);
80 size_t n = getLocalNumIdentifiers();
81 if (n){
82 Ids = ArrayView<const gno_t>(
83 reinterpret_cast<const gno_t*>(gids_.getRawPtr()), n);
84 }
85 return n;
86 }
87
89 Kokkos::View<const gno_t *, typename node_t::device_type> &Ids,
90 Kokkos::View<scalar_t **, typename node_t::device_type> &wgts) const {
91 try {
92 adapter_->getIDsKokkosView(Ids);
93 adapter_->getWeightsKokkosView(wgts);
94 }
96
98 }
99
101 // The Model interface.
103
104 inline size_t getLocalNumObjects() const {return getLocalNumIdentifiers();}
105
106 inline size_t getGlobalNumObjects() const {return getGlobalNumIdentifiers();}
107
108private:
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_;
114 int nUserWeights_;
115 ArrayRCP<input_t> weights_;
116};
117
119template <typename Adapter>
121 const RCP<const Adapter> &ia,
122 const RCP<const Environment> &env,
123 const RCP<const Comm<int> > &comm,
124 modelFlag_t &/* modelFlags */):
125 numGlobalIdentifiers_(), env_(env), comm_(comm),
126 gids_(), adapter_(ia), nUserWeights_(0), weights_()
127{
128 // Get the local and global problem size
129 size_t nLocalIds = ia->getLocalNumIDs();
130 gno_t lsum = nLocalIds;
131 reduceAll<int, gno_t>(*comm_, Teuchos::REDUCE_SUM, 1, &lsum,
132 &numGlobalIdentifiers_);
133
134 // Get the number of weights
135 int tmp = ia->getNumWeightsPerID();
136 // Use max number of weights over all processes as nUserWeights_
137 Teuchos::reduceAll<int, int>(*comm, Teuchos::REDUCE_MAX, 1,
138 &tmp, &nUserWeights_);
139
140 // Prepare to store views from input adapter
141 // TODO: Do we have to store these views, or can we get them on an
142 // TODO: as-needed basis?
143 Array<const scalar_t *> wgts(nUserWeights_, (const scalar_t *)NULL);
144 Array<int> wgtStrides(nUserWeights_, 0);
145
146 if (nUserWeights_ > 0){
147 input_t *w = new input_t [nUserWeights_];
148 weights_ = arcp<input_t>(w, 0, nUserWeights_);
149 }
150
151 const gno_t *gids=NULL;
152
153 // Get the input adapter's views
154 try{
155 ia->getIDsView(gids);
156 for (int idx=0; idx < nUserWeights_; idx++)
157 ia->getWeightsView(wgts[idx], wgtStrides[idx], idx);
158 }
160
161 if (nLocalIds){
162 gids_ = arcp(gids, 0, nLocalIds, false);
163
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]);
169 }
170 }
171 }
172
173 env_->memory("After construction of identifier model");
174}
175
176} // namespace Zoltan2
177
178#endif
typename Zoltan2::InputTraits< ztcrsmatrix_t >::node_t node_t
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
Defines the Model interface.
This file defines the StridedData class.
IdentifierModel defines the interface for all identifier models.
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
size_t getLocalNumObjects() const
Return the local number of objects.
global_size_t getGlobalNumIdentifiers() const
IdentifierModel(const RCP< const Adapter > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &modelFlags)
Constructor.
size_t getGlobalNumObjects() const
Return the global number of objects.
size_t getIdentifierList(ArrayView< const gno_t > &Ids, ArrayView< input_t > &wgts) const
The base class for all model classes.
The StridedData class manages lists of weights or coordinates.
map_t::local_ordinal_type lno_t
map_t::global_ordinal_type gno_t
Created by mbenlioglu on Aug 31, 2020.
Tpetra::global_size_t global_size_t
std::bitset< NUM_MODEL_FLAGS > modelFlag_t