Amesos2 - Direct Sparse Solver Interfaces Version of the Day
Amesos2_EpetraMultiVecAdapter_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Amesos2: Templated Direct Sparse Solver Package
4//
5// Copyright 2011 NTESS and the Amesos2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
19#ifndef AMESOS2_EPETRA_MULTIVEC_ADAPTER_DECL_HPP
20#define AMESOS2_EPETRA_MULTIVEC_ADAPTER_DECL_HPP
21
22#include <Teuchos_RCP.hpp>
23#include <Teuchos_Array.hpp>
24#include <Teuchos_as.hpp>
25
26#include <Tpetra_Vector.hpp>
27#include <Tpetra_Map.hpp>
28
29#include <Epetra_MultiVector.h>
30
31#include "Amesos2_MultiVecAdapter_decl.hpp"
33
34namespace Amesos2 {
35
41 template <>
42 class MultiVecAdapter<Epetra_MultiVector>
43 {
44 public:
45
46 // public type definitions
47 typedef double scalar_t;
48 typedef int local_ordinal_t;
49 typedef Tpetra::Map<>::global_ordinal_type global_ordinal_t;
50 typedef size_t global_size_t;
51 typedef Tpetra::Map<>::node_type node_t;
52 typedef Epetra_MultiVector multivec_t;
53
54 friend Teuchos::RCP<MultiVecAdapter<multivec_t> > createMultiVecAdapter<>(Teuchos::RCP<multivec_t>);
55 friend Teuchos::RCP<const MultiVecAdapter<multivec_t> > createConstMultiVecAdapter<>(Teuchos::RCP<const multivec_t>);
56
57
58 static const char* name;
59
60
61 protected:
64
70 MultiVecAdapter( const Teuchos::RCP<multivec_t>& m );
71
72
73 public:
74
75 ~MultiVecAdapter() = default;
76
77
79 bool isLocallyIndexed() const;
80
81 bool isGloballyIndexed() const;
82
83
84 Teuchos::RCP<Epetra_MultiVector> clone() const;
85
92 Teuchos::RCP<const Tpetra::Map<local_ordinal_t, global_ordinal_t, node_t> >
93 getMap() const;
94
96 const Teuchos::RCP<const Teuchos::Comm<int> > getComm() const;
97
98
100 size_t getLocalLength() const;
101
102
104 size_t getLocalNumVectors() const;
105
106
108 global_size_t getGlobalLength() const;
109
110
112 size_t getGlobalNumVectors() const;
113
114
116 size_t getStride() const;
117
118
120 bool isConstantStride() const;
121
122
124 Teuchos::RCP<const Tpetra::Vector<scalar_t,local_ordinal_t,global_ordinal_t,node_t> >
125 getVector( size_t j ) const;
126
127
135 Teuchos::RCP<Tpetra::Vector<scalar_t,local_ordinal_t,global_ordinal_t,node_t> >
136 getVectorNonConst( size_t j );
137
138
140 double * getMVPointer_impl() const;
141
147 void get1dCopy( const Teuchos::ArrayView<scalar_t>& A,
148 size_t lda,
149 Teuchos::Ptr<
150 const Tpetra::Map<local_ordinal_t,
151 global_ordinal_t,
152 node_t> > distribution_map,
153 EDistribution distribution) const;
154
155 template<typename KV>
156 bool get1dCopy_kokkos_view(bool bInitialize, KV & A,
157 size_t lda,
158 Teuchos::Ptr<
159 const Tpetra::Map<local_ordinal_t,
160 global_ordinal_t,
161 node_t> > distribution_map,
162 EDistribution distribution) const {
163 Kokkos::View<double**, Kokkos::LayoutLeft, Kokkos::HostSpace> host_new_data;
164 get1dCopy_kokkos_view_host(host_new_data, lda, distribution_map, distribution);
165 bool bAssigned;
166 deep_copy_or_assign_view(bInitialize, A, host_new_data, bAssigned);
167 return false; // currently Epetra and prior get1dCopy_kokkos_view_host call cannot get direct assignment so always return false
168 }
169
170 void get1dCopy_kokkos_view_host(
171 Kokkos::View<scalar_t**, Kokkos::LayoutLeft, Kokkos::HostSpace> & new_data,
172 size_t lda,
173 Teuchos::Ptr<
174 const Tpetra::Map<local_ordinal_t,
175 global_ordinal_t,
176 node_t> > distribution_map,
177 EDistribution) const;
178
196 Teuchos::ArrayRCP<scalar_t> get1dViewNonConst( bool local = false );
197
198
209 void put1dData( const Teuchos::ArrayView<const scalar_t>& new_data,
210 size_t lda,
211 Teuchos::Ptr<
212 const Tpetra::Map<local_ordinal_t,
213 global_ordinal_t,
214 node_t> > source_map,
215 EDistribution distribution );
216
217 template<typename KV>
218 void put1dData_kokkos_view(
219 KV & new_data,
220 size_t lda,
221 Teuchos::Ptr<
222 const Tpetra::Map<local_ordinal_t,
223 global_ordinal_t,
224 node_t> > source_map,
225 EDistribution distribution ) {
226 Kokkos::View<scalar_t**, Kokkos::LayoutLeft, Kokkos::HostSpace> host_new_data(
227 Kokkos::ViewAllocateWithoutInitializing("host_new_data"),
228 new_data.extent(0), new_data.extent(1));
229 Kokkos::deep_copy(host_new_data, new_data);
230 put1dData_kokkos_view_host(host_new_data, lda, source_map, distribution);
231 }
232
233 void put1dData_kokkos_view_host(
234 Kokkos::View<scalar_t**, Kokkos::LayoutLeft, Kokkos::HostSpace> & new_data,
235 size_t lda,
236 Teuchos::Ptr<
237 const Tpetra::Map<local_ordinal_t,
238 global_ordinal_t,
239 node_t> > source_map,
240 EDistribution distribution );
241
242 template<typename KV, typename host_ordinal_type_array>
243 int
244 gather (KV& kokkos_new_view,
245 host_ordinal_type_array &perm_g2l,
246 host_ordinal_type_array &recvCountRows,
247 host_ordinal_type_array &recvDisplRows,
248 EDistribution distribution) const
249 {
250 auto comm = this->getComm();
251 int myRank = comm->getRank();
252 int nCols = this->mv_->NumVectors();
253 int nRows = this->mv_->GlobalLength();
254 int nRows_l = this->mv_->MyLength();
255 if (myRank == 0) {
256 Kokkos::resize(kokkos_new_view, nRows, nCols);
257 if (int(perm_g2l.extent(0)) == nRows) {
258 Kokkos::resize(this->buf_, nRows, 1);
259 } else {
260 Kokkos::resize(this->buf_, 0, 1);
261 }
262 }
263 {
264 for (int j=0; j<nCols; j++) {
265 scalar_t * recvbuf = reinterpret_cast<scalar_t*> (myRank != 0 || this->buf_.extent(0) > 0 ? this->buf_.data() : &kokkos_new_view(0,j));
266 Teuchos::gatherv<int, scalar_t> (const_cast<scalar_t*> ((*this->mv_)[j]), nRows_l,
267 recvbuf, recvCountRows.data(), recvDisplRows.data(),
268 0, *comm);
269 if (myRank == 0 && this->buf_.extent(0) > 0) {
270 for (int i=0; i<nRows; i++) kokkos_new_view(perm_g2l(i),j) = this->buf_(i,0);
271 }
272 }
273 }
274 return 0;
275 }
276
277 template<typename KV, typename host_ordinal_type_array>
278 int
279 scatter (KV& kokkos_new_view,
280 host_ordinal_type_array &perm_g2l,
281 host_ordinal_type_array &sendCountRows,
282 host_ordinal_type_array &sendDisplRows,
283 EDistribution distribution) const
284 {
285 auto comm = this->getMap()->getComm();
286 int myRank = comm->getRank();
287 int nCols = this->mv_->NumVectors();
288 int nRows = this->mv_->GlobalLength();
289 int nRows_l = this->mv_->MyLength();
290 {
291 for (int j=0; j<nCols; j++) {
292 if (myRank == 0 && this->buf_.extent(0) > 0) {
293 for (int i=0; i<nRows; i++) this->buf_(i, 0) = kokkos_new_view(perm_g2l(i),j);
294 }
295 scalar_t * sendbuf = reinterpret_cast<scalar_t*> (myRank != 0 || this->buf_.extent(0) > 0 ? this->buf_.data() : &kokkos_new_view(0,j));
296 Teuchos::scatterv<int, scalar_t> (sendbuf, sendCountRows.data(), sendDisplRows.data(),
297 reinterpret_cast<scalar_t*> ((*this->mv_)[j]), nRows_l,
298 0, *comm);
299 }
300 }
301 return 0;
302 }
303
305 std::string description() const;
306
307
309 void describe( Teuchos::FancyOStream& os,
310 const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default) const;
311
312
313 private:
314
316 Teuchos::RCP<multivec_t> mv_;
317 mutable Kokkos::View<scalar_t**, Kokkos::LayoutLeft, Kokkos::HostSpace> buf_;
318
319 mutable Teuchos::RCP<Epetra_Import> importer_;
320 mutable Teuchos::RCP<Epetra_Export> exporter_;
321
322 mutable Teuchos::RCP<const Epetra_BlockMap> mv_map_;
323
324 }; // end class MultiVecAdapter<NewMultiVec>
325
326} // end namespace Amesos2
327
328
329#endif // AMESOS2_EPETRA_MULTIVEC_ADAPTER_DECL_HPP
Copy or assign views based on memory spaces.
EDistribution
Definition Amesos2_TypeDecl.hpp:89
Teuchos::RCP< multivec_t > mv_
The multi-vector this adapter wraps.
Definition Amesos2_EpetraMultiVecAdapter_decl.hpp:316
A templated MultiVector class adapter for Amesos2.
Definition Amesos2_MultiVecAdapter_decl.hpp:142