Amesos2 - Direct Sparse Solver Interfaces Version of the Day
Amesos2_TpetraMultiVecAdapter_def.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
18
19#ifndef AMESOS2_TPETRA_MULTIVEC_ADAPTER_DEF_HPP
20#define AMESOS2_TPETRA_MULTIVEC_ADAPTER_DEF_HPP
21
22#include <type_traits>
25
26
27namespace Amesos2 {
28
29 using Tpetra::MultiVector;
30
31 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, class Node >
33 MultiVector<Scalar,
34 LocalOrdinal,
35 GlobalOrdinal,
36 Node> >::MultiVecAdapter( const Teuchos::RCP<multivec_t>& m )
37 : mv_(m)
38 {}
39
40 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, class Node >
41 Teuchos::RCP<
42 MultiVector<Scalar,
43 LocalOrdinal,
44 GlobalOrdinal,
45 Node> >
46 MultiVecAdapter<
47 MultiVector<Scalar,
48 LocalOrdinal,
49 GlobalOrdinal,
50 Node> >::clone() const
51 {
52 using MV = MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
53 Teuchos::RCP<MV> Y (new MV (mv_->getMap(), mv_->getNumVectors(), false));
54 Y->setCopyOrView (Teuchos::View);
55 return Y;
56 }
57
58 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, class Node >
59 typename MultiVecAdapter<
60 MultiVector<Scalar,
61 LocalOrdinal,
62 GlobalOrdinal,
63 Node> >::multivec_t::impl_scalar_type *
64 MultiVecAdapter<
65 MultiVector<Scalar,
66 LocalOrdinal,
67 GlobalOrdinal,
68 Node> >::getMVPointer_impl() const
69 {
70 TEUCHOS_TEST_FOR_EXCEPTION( this->getGlobalNumVectors() != 1,
71 std::invalid_argument,
72 "Amesos2_TpetraMultiVectorAdapter: getMVPointer_impl should only be called for case with a single vector and single MPI process" );
73
74 auto contig_local_view_2d = mv_->getLocalViewHost(Tpetra::Access::ReadWrite);
75 auto contig_local_view_1d = Kokkos::subview (contig_local_view_2d, Kokkos::ALL (), 0);
76 return contig_local_view_1d.data();
77 }
78
79 // TODO Proper type handling:
80 // Consider a MultiVectorTraits class
81 // typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> multivector_type
82 // NOTE: In this class, above already has a typedef multivec_t
83 // typedef typename multivector_type::impl_scalar_type return_scalar_type; // this is the POD type the dual_view_type is templated on
84 // Traits class needed to do this generically for the general MultiVectorAdapter interface
85
86
87 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, class Node >
88 void
89 MultiVecAdapter<
90 MultiVector<Scalar,
91 LocalOrdinal,
92 GlobalOrdinal,
93 Node> >::get1dCopy(const Teuchos::ArrayView<scalar_t>& av,
94 size_t lda,
95 Teuchos::Ptr<
96 const Tpetra::Map<LocalOrdinal,
97 GlobalOrdinal,
98 Node> > distribution_map,
99 EDistribution distribution) const
100 {
101 using Teuchos::as;
102 using Teuchos::RCP;
103 typedef Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> map_type;
104 const size_t num_vecs = getGlobalNumVectors ();
105
106 TEUCHOS_TEST_FOR_EXCEPTION(
107 distribution_map.getRawPtr () == NULL, std::invalid_argument,
108 "Amesos2::MultiVecAdapter::get1dCopy: distribution_map argument is null.");
109 TEUCHOS_TEST_FOR_EXCEPTION(
110 mv_.is_null (), std::logic_error,
111 "Amesos2::MultiVecAdapter::get1dCopy: mv_ is null.");
112 // Check mv_ before getMap(), because the latter calls mv_->getMap().
113 TEUCHOS_TEST_FOR_EXCEPTION(
114 this->getMap ().is_null (), std::logic_error,
115 "Amesos2::MultiVecAdapter::get1dCopy: this->getMap() returns null.");
116
117#ifdef HAVE_AMESOS2_DEBUG
118 const size_t requested_vector_length = distribution_map->getLocalNumElements ();
119 TEUCHOS_TEST_FOR_EXCEPTION(
120 lda < requested_vector_length, std::invalid_argument,
121 "Amesos2::MultiVecAdapter::get1dCopy: On process " <<
122 distribution_map->getComm ()->getRank () << " of the distribution Map's "
123 "communicator, the given stride lda = " << lda << " is not large enough "
124 "for the local vector length " << requested_vector_length << ".");
125 TEUCHOS_TEST_FOR_EXCEPTION(
126 as<size_t> (av.size ()) < as<size_t> ((num_vecs - 1) * lda + requested_vector_length),
127 std::invalid_argument, "Amesos2::MultiVector::get1dCopy: MultiVector "
128 "storage not large enough given leading dimension and number of vectors." );
129#endif // HAVE_AMESOS2_DEBUG
130
131 // Special case when number vectors == 1 and single MPI process
132 if ( num_vecs == 1 && this->getComm()->getRank() == 0 && this->getComm()->getSize() == 1 ) {
133 mv_->get1dCopy (av, lda);
134 }
135 else {
136
137 // (Re)compute the Export object if necessary. If not, then we
138 // don't need to clone distribution_map; we can instead just get
139 // the previously cloned target Map from the Export object.
140 RCP<const map_type> distMap;
141 if (exporter_.is_null () ||
142 ! exporter_->getSourceMap ()->isSameAs (* (this->getMap ())) ||
143 ! exporter_->getTargetMap ()->isSameAs (* distribution_map)) {
144 // Since we're caching the Export object, and since the Export
145 // needs to keep the distribution Map, we have to make a copy of
146 // the latter in order to ensure that it will stick around past
147 // the scope of this function call. (Ptr is not reference
148 // counted.)
149 distMap = rcp(new map_type(*distribution_map));
150 // (Re)create the Export object.
151 exporter_ = rcp (new export_type (this->getMap (), distMap));
152 }
153 else {
154 distMap = exporter_->getTargetMap ();
155 }
156
157 multivec_t redist_mv (distMap, num_vecs);
158
159 // Redistribute the input (multi)vector.
160 redist_mv.doExport (*mv_, *exporter_, Tpetra::REPLACE);
161
162 if ( distribution != CONTIGUOUS_AND_ROOTED ) {
163 // Do this if GIDs contiguous - existing functionality
164 // Copy the imported (multi)vector's data into the ArrayView.
165 redist_mv.get1dCopy (av, lda);
166 }
167 else {
168 // Do this if GIDs not contiguous...
169 // sync is needed for example if mv was updated on device, but will be passed through Amesos2 to solver running on host
170 //redist_mv.template sync < host_execution_space > ();
171 auto contig_local_view_2d = redist_mv.getLocalViewHost(Tpetra::Access::ReadOnly);
172 if ( redist_mv.isConstantStride() ) {
173 for ( size_t j = 0; j < num_vecs; ++j) {
174 auto av_j = av(lda*j, lda);
175 for ( size_t i = 0; i < lda; ++i ) {
176 av_j[i] = contig_local_view_2d(i,j); //lda may not be correct if redist_mv is not constant stride...
177 }
178 }
179 }
180 else {
181 // ... lda should come from Teuchos::Array* allocation,
182 // not the MultiVector, since the MultiVector does NOT
183 // have constant stride in this case.
184 // TODO lda comes from X->getGlobalLength() in solve_impl - should this be changed???
185 const size_t lclNumRows = redist_mv.getLocalLength();
186 for (size_t j = 0; j < redist_mv.getNumVectors(); ++j) {
187 auto av_j = av(lda*j, lclNumRows);
188 auto X_lcl_j_2d = redist_mv.getLocalViewHost(Tpetra::Access::ReadOnly);
189 auto X_lcl_j_1d = Kokkos::subview (X_lcl_j_2d, Kokkos::ALL (), j);
190
191 using val_type = typename std::remove_const<typename decltype( X_lcl_j_1d )::value_type>::type;
192 Kokkos::View<val_type*, Kokkos::HostSpace> umavj ( const_cast< val_type* > ( reinterpret_cast<const val_type*> ( av_j.getRawPtr () ) ), av_j.size () );
193 Kokkos::deep_copy (umavj, X_lcl_j_1d);
194 }
195 }
196 }
197 }
198 }
199
200 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, class Node >
201 template <typename KV>
202 bool
203 MultiVecAdapter<
204 MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>
205 >::get1dCopy_kokkos_view(
206 bool bInitialize,
207 KV& kokkos_view,
208 [[maybe_unused]] size_t lda,
209 Teuchos::Ptr<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > distribution_map,
210 EDistribution distribution
211 ) const {
212 using Teuchos::as;
213 using Teuchos::RCP;
214 typedef Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> map_type;
215 const size_t num_vecs = getGlobalNumVectors ();
216
217 TEUCHOS_TEST_FOR_EXCEPTION(
218 distribution_map.getRawPtr () == NULL, std::invalid_argument,
219 "Amesos2::MultiVecAdapter::get1dCopy_kokkos_view: distribution_map argument is null.");
220 TEUCHOS_TEST_FOR_EXCEPTION(
221 mv_.is_null (), std::logic_error,
222 "Amesos2::MultiVecAdapter::get1dCopy_kokkos_view: mv_ is null.");
223 // Check mv_ before getMap(), because the latter calls mv_->getMap().
224 TEUCHOS_TEST_FOR_EXCEPTION(
225 this->getMap ().is_null (), std::logic_error,
226 "Amesos2::MultiVecAdapter::get1dCopy_kokkos_view: this->getMap() returns null.");
227
228#ifdef HAVE_AMESOS2_DEBUG
229 const size_t requested_vector_length = distribution_map->getLocalNumElements ();
230 TEUCHOS_TEST_FOR_EXCEPTION(
231 lda < requested_vector_length, std::invalid_argument,
232 "Amesos2::MultiVecAdapter::get1dCopy_kokkos_view: On process " <<
233 distribution_map->getComm ()->getRank () << " of the distribution Map's "
234 "communicator, the given stride lda = " << lda << " is not large enough "
235 "for the local vector length " << requested_vector_length << ".");
236
237 // Note do not check size since deep_copy_or_assign_view below will allocate
238 // if necessary - but may just assign ptrs.
239#endif // HAVE_AMESOS2_DEBUG
240
241 // Special case when number vectors == 1 and single MPI process
242 if ( num_vecs == 1 && this->getComm()->getRank() == 0 && this->getComm()->getSize() == 1 ) {
243 if(mv_->isConstantStride()) {
244 bool bAssigned;
245 //deep_copy_or_assign_view(bInitialize, kokkos_view, mv_->getLocalViewDevice(Tpetra::Access::ReadOnly), bAssigned);
246 deep_copy_only(bInitialize, kokkos_view, mv_->getLocalViewDevice(Tpetra::Access::ReadOnly), bAssigned);
247 return bAssigned; // if bAssigned is true we are accessing the mv data directly without a copy
248 }
249 else {
250 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Resolve handling for non-constant stride.");
251 }
252 }
253 else {
254
255 // (Re)compute the Export object if necessary. If not, then we
256 // don't need to clone distribution_map; we can instead just get
257 // the previously cloned target Map from the Export object.
258 RCP<const map_type> distMap;
259 if (exporter_.is_null () ||
260 ! exporter_->getSourceMap ()->isSameAs (* (this->getMap ())) ||
261 ! exporter_->getTargetMap ()->isSameAs (* distribution_map)) {
262 // Since we're caching the Export object, and since the Export
263 // needs to keep the distribution Map, we have to make a copy of
264 // the latter in order to ensure that it will stick around past
265 // the scope of this function call. (Ptr is not reference
266 // counted.)
267 distMap = rcp(new map_type(*distribution_map));
268 // (Re)create the Export object.
269 exporter_ = rcp (new export_type (this->getMap (), distMap));
270 }
271 else {
272 distMap = exporter_->getTargetMap ();
273 }
274
275 multivec_t redist_mv (distMap, num_vecs);
276
277 // Redistribute the input (multi)vector.
278 redist_mv.doExport (*mv_, *exporter_, Tpetra::REPLACE);
279
280 if ( distribution != CONTIGUOUS_AND_ROOTED ) {
281 // Do this if GIDs contiguous - existing functionality
282 // Copy the imported (multi)vector's data into the Kokkos View.
283 bool bAssigned;
284 deep_copy_or_assign_view(bInitialize, kokkos_view, redist_mv.getLocalViewDevice(Tpetra::Access::ReadOnly), bAssigned);
285 return false; // do not return bAssigned because redist_mv was already copied so always return false
286 }
287 else {
288 if(redist_mv.isConstantStride()) {
289 bool bAssigned; // deep_copy_or_assign_view sets true if assigned (no deep copy)
290 deep_copy_or_assign_view(bInitialize, kokkos_view, redist_mv.getLocalViewDevice(Tpetra::Access::ReadOnly), bAssigned);
291 return false; // do not return bAssigned because redist_mv was already copied so always return false
292 }
293 else {
294 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Kokkos adapter non-constant stride not imlemented.");
295 }
296 }
297 }
298 }
299
300 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, class Node >
301 Teuchos::ArrayRCP<Scalar>
302 MultiVecAdapter<
303 MultiVector<Scalar,
304 LocalOrdinal,
305 GlobalOrdinal,
306 Node> >::get1dViewNonConst (bool local)
307 {
308 // FIXME (mfh 22 Jan 2014) When I first found this routine, all of
309 // its code was commented out, and it didn't return anything. The
310 // latter is ESPECIALLY dangerous, given that the return value is
311 // an ArrayRCP. Thus, I added the exception throw below.
312 TEUCHOS_TEST_FOR_EXCEPTION(
313 true, std::logic_error, "Amesos2::MultiVecAdapter::get1dViewNonConst: "
314 "Not implemented.");
315
316 // if( local ){
317 // this->localize();
318 // /* Use the global element list returned by
319 // * mv_->getMap()->getLocalElementList() to get a subCopy of mv_ which we
320 // * assign to l_l_mv_, then return get1dViewNonConst() of l_l_mv_
321 // */
322 // if(l_l_mv_.is_null() ){
323 // Teuchos::ArrayView<const GlobalOrdinal> nodeElements_go
324 // = mv_->getMap()->getLocalElementList();
325 // Teuchos::Array<size_t> nodeElements_st(nodeElements_go.size());
326
327 // // Convert the node element to a list of size_t type objects
328 // typename Teuchos::ArrayView<const GlobalOrdinal>::iterator it_go;
329 // Teuchos::Array<size_t>::iterator it_st = nodeElements_st.begin();
330 // for( it_go = nodeElements_go.begin(); it_go != nodeElements_go.end(); ++it_go ){
331 // *(it_st++) = Teuchos::as<size_t>(*it_go);
332 // }
333
334 // // To be consistent with the globalize() function, get a view of the local mv
335 // l_l_mv_ = l_mv_->subViewNonConst(nodeElements_st);
336
337 // return(l_l_mv_->get1dViewNonConst());
338 // } else {
339 // // We need to re-import values to the local, since changes may have been
340 // // made to the global structure that are not reflected in the local
341 // // view.
342 // Teuchos::ArrayView<const GlobalOrdinal> nodeElements_go
343 // = mv_->getMap()->getLocalElementList();
344 // Teuchos::Array<size_t> nodeElements_st(nodeElements_go.size());
345
346 // // Convert the node element to a list of size_t type objects
347 // typename Teuchos::ArrayView<const GlobalOrdinal>::iterator it_go;
348 // Teuchos::Array<size_t>::iterator it_st = nodeElements_st.begin();
349 // for( it_go = nodeElements_go.begin(); it_go != nodeElements_go.end(); ++it_go ){
350 // *(it_st++) = Teuchos::as<size_t>(*it_go);
351 // }
352
353 // return l_l_mv_->get1dViewNonConst();
354 // }
355 // } else {
356 // if( mv_->isDistributed() ){
357 // this->localize();
358
359 // return l_mv_->get1dViewNonConst();
360 // } else { // not distributed, no need to import
361 // return mv_->get1dViewNonConst();
362 // }
363 // }
364 }
365
366
367 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, class Node>
368 void
369 MultiVecAdapter<
370 MultiVector<Scalar,
371 LocalOrdinal,
372 GlobalOrdinal,
373 Node> >::put1dData(const Teuchos::ArrayView<const scalar_t>& new_data,
374 size_t lda,
375 Teuchos::Ptr<
376 const Tpetra::Map<LocalOrdinal,
377 GlobalOrdinal,
378 Node> > source_map,
379 EDistribution distribution )
380 {
381 using Teuchos::RCP;
382 typedef Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> map_type;
383
384 TEUCHOS_TEST_FOR_EXCEPTION(
385 source_map.getRawPtr () == NULL, std::invalid_argument,
386 "Amesos2::MultiVecAdapter::put1dData: source_map argument is null.");
387 TEUCHOS_TEST_FOR_EXCEPTION(
388 mv_.is_null (), std::logic_error,
389 "Amesos2::MultiVecAdapter::put1dData: the internal MultiVector mv_ is null.");
390 // getMap() calls mv_->getMap(), so test first whether mv_ is null.
391 TEUCHOS_TEST_FOR_EXCEPTION(
392 this->getMap ().is_null (), std::logic_error,
393 "Amesos2::MultiVecAdapter::put1dData: this->getMap() returns null.");
394
395 const size_t num_vecs = getGlobalNumVectors ();
396
397 // Special case when number vectors == 1 and single MPI process
398 if ( num_vecs == 1 && this->getComm()->getRank() == 0 && this->getComm()->getSize() == 1 ) {
399 // num_vecs = 1; stride does not matter
400 auto mv_view_to_modify_2d = mv_->getLocalViewHost(Tpetra::Access::OverwriteAll);
401 for ( size_t i = 0; i < lda; ++i ) {
402 mv_view_to_modify_2d(i,0) = new_data[i]; // Only one vector
403 }
404 }
405 else {
406
407 // (Re)compute the Import object if necessary. If not, then we
408 // don't need to clone source_map; we can instead just get the
409 // previously cloned source Map from the Import object.
410 RCP<const map_type> srcMap;
411 if (importer_.is_null () ||
412 ! importer_->getSourceMap ()->isSameAs (* source_map) ||
413 ! importer_->getTargetMap ()->isSameAs (* (this->getMap ()))) {
414
415 // Since we're caching the Import object, and since the Import
416 // needs to keep the source Map, we have to make a copy of the
417 // latter in order to ensure that it will stick around past the
418 // scope of this function call. (Ptr is not reference counted.)
419 srcMap = rcp(new map_type(*source_map));
420 importer_ = rcp (new import_type (srcMap, this->getMap ()));
421 }
422 else {
423 srcMap = importer_->getSourceMap ();
424 }
425
426 if ( distribution != CONTIGUOUS_AND_ROOTED ) {
427 // Do this if GIDs contiguous - existing functionality
428 // Redistribute the output (multi)vector.
429 const multivec_t source_mv (srcMap, new_data, lda, num_vecs);
430 mv_->doImport (source_mv, *importer_, Tpetra::REPLACE);
431 }
432 else {
433 multivec_t redist_mv (srcMap, num_vecs); // unused for ROOTED case
434 if ( redist_mv.isConstantStride() ) {
435 auto contig_local_view_2d = redist_mv.getLocalViewHost(Tpetra::Access::OverwriteAll);
436 for ( size_t j = 0; j < num_vecs; ++j) {
437 auto av_j = new_data(lda*j, lda);
438 for ( size_t i = 0; i < lda; ++i ) {
439 contig_local_view_2d(i,j) = av_j[i];
440 }
441 }
442 }
443 else {
444 // ... lda should come from Teuchos::Array* allocation,
445 // not the MultiVector, since the MultiVector does NOT
446 // have constant stride in this case.
447 // TODO lda comes from X->getGlobalLength() in solve_impl - should this be changed???
448 const size_t lclNumRows = redist_mv.getLocalLength();
449 for (size_t j = 0; j < redist_mv.getNumVectors(); ++j) {
450 auto av_j = new_data(lda*j, lclNumRows);
451 auto X_lcl_j_2d = redist_mv.getLocalViewHost(Tpetra::Access::ReadOnly);
452 auto X_lcl_j_1d = Kokkos::subview (X_lcl_j_2d, Kokkos::ALL (), j);
453
454 using val_type = typename std::remove_const<typename decltype( X_lcl_j_1d )::value_type>::type;
455 Kokkos::View<val_type*, Kokkos::HostSpace> umavj ( const_cast< val_type* > ( reinterpret_cast<const val_type*> ( av_j.getRawPtr () ) ), av_j.size () );
456 Kokkos::deep_copy (umavj, X_lcl_j_1d);
457 }
458 }
459
460 //typedef typename multivec_t::node_type::memory_space memory_space;
461 //redist_mv.template sync <memory_space> ();
462
463 mv_->doImport (redist_mv, *importer_, Tpetra::REPLACE);
464 }
465 }
466
467 }
468
469 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, class Node>
470 template <typename KV>
471 void
472 MultiVecAdapter<
473 MultiVector<Scalar,
474 LocalOrdinal,
475 GlobalOrdinal,
476 Node> >::put1dData_kokkos_view(KV& kokkos_new_data,
477 size_t lda,
478 Teuchos::Ptr<
479 const Tpetra::Map<LocalOrdinal,
480 GlobalOrdinal,
481 Node> > source_map,
482 EDistribution distribution )
483 {
484 using Teuchos::RCP;
485 typedef Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node> map_type;
486
487 TEUCHOS_TEST_FOR_EXCEPTION(
488 source_map.getRawPtr () == NULL, std::invalid_argument,
489 "Amesos2::MultiVecAdapter::put1dData_kokkos_view: source_map argument is null.");
490 TEUCHOS_TEST_FOR_EXCEPTION(
491 mv_.is_null (), std::logic_error,
492 "Amesos2::MultiVecAdapter::put1dData_kokkos_view: the internal MultiVector mv_ is null.");
493 // getMap() calls mv_->getMap(), so test first whether mv_ is null.
494 TEUCHOS_TEST_FOR_EXCEPTION(
495 this->getMap ().is_null (), std::logic_error,
496 "Amesos2::MultiVecAdapter::put1dData_kokkos_view: this->getMap() returns null.");
497
498 const size_t num_vecs = getGlobalNumVectors ();
499
500 // Special case when number vectors == 1 and single MPI process
501 if ( num_vecs == 1 && this->getComm()->getRank() == 0 && this->getComm()->getSize() == 1 ) {
502
503 // num_vecs = 1; stride does not matter
504
505 // If this is the optimized path then kokkos_new_data will be the dst
506 auto mv_view_to_modify_2d = mv_->getLocalViewDevice(Tpetra::Access::OverwriteAll);
507 //deep_copy_or_assign_view(mv_view_to_modify_2d, kokkos_new_data);
508 deep_copy_only(mv_view_to_modify_2d, kokkos_new_data);
509 }
510 else {
511
512 // (Re)compute the Import object if necessary. If not, then we
513 // don't need to clone source_map; we can instead just get the
514 // previously cloned source Map from the Import object.
515 RCP<const map_type> srcMap;
516 if (importer_.is_null () ||
517 ! importer_->getSourceMap ()->isSameAs (* source_map) ||
518 ! importer_->getTargetMap ()->isSameAs (* (this->getMap ()))) {
519
520 // Since we're caching the Import object, and since the Import
521 // needs to keep the source Map, we have to make a copy of the
522 // latter in order to ensure that it will stick around past the
523 // scope of this function call. (Ptr is not reference counted.)
524 srcMap = rcp(new map_type(*source_map));
525 importer_ = rcp (new import_type (srcMap, this->getMap ()));
526 }
527 else {
528 srcMap = importer_->getSourceMap ();
529 }
530
531 if ( distribution != CONTIGUOUS_AND_ROOTED ) {
532 // Use View scalar type, not MV Scalar because we want Kokkos::complex, not
533 // std::complex to avoid a Kokkos::complex<double> to std::complex<float>
534 // conversion which would require a double copy and fail here. Then we'll be
535 // setup to safely reinterpret_cast complex to std if necessary.
536 typedef typename multivec_t::dual_view_type::t_host::value_type tpetra_mv_view_type;
537 Kokkos::View<tpetra_mv_view_type**,typename KV::array_layout,
538 Kokkos::HostSpace> convert_kokkos_new_data;
539 deep_copy_or_assign_view(convert_kokkos_new_data, kokkos_new_data);
540#ifdef HAVE_TEUCHOS_COMPLEX
541 // convert_kokkos_new_data may be Kokkos::complex and Scalar could be std::complex
542 auto pData = reinterpret_cast<Scalar*>(convert_kokkos_new_data.data());
543#else
544 auto pData = convert_kokkos_new_data.data();
545#endif
546
547 const multivec_t source_mv (srcMap, Teuchos::ArrayView<const scalar_t>(
548 pData, kokkos_new_data.size()), lda, num_vecs);
549 mv_->doImport (source_mv, *importer_, Tpetra::REPLACE);
550 }
551 else {
552 multivec_t redist_mv (srcMap, num_vecs); // unused for ROOTED case
553 // Cuda solvers won't currently use this path since they are just serial
554 // right now, so this mirror should be harmless (and not strictly necessary).
555 // Adding it for future possibilities though we may then refactor this
556 // for better efficiency if the kokkos_new_data view is on device.
557 auto host_kokkos_new_data = Kokkos::create_mirror_view(kokkos_new_data);
558 Kokkos::deep_copy(host_kokkos_new_data, kokkos_new_data);
559 if ( redist_mv.isConstantStride() ) {
560 auto contig_local_view_2d = redist_mv.getLocalViewHost(Tpetra::Access::OverwriteAll);
561 for ( size_t j = 0; j < num_vecs; ++j) {
562 auto av_j = Kokkos::subview(host_kokkos_new_data, Kokkos::ALL, j);
563 for ( size_t i = 0; i < lda; ++i ) {
564 contig_local_view_2d(i,j) = av_j(i);
565 }
566 }
567 }
568 else {
569 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Kokkos adapter "
570 "CONTIGUOUS_AND_ROOTED not implemented for put1dData_kokkos_view "
571 "with non constant stride.");
572 }
573
574 //typedef typename multivec_t::node_type::memory_space memory_space;
575 //redist_mv.template sync <memory_space> ();
576
577 mv_->doImport (redist_mv, *importer_, Tpetra::REPLACE);
578 }
579 }
580
581 }
582
583 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, class Node >
584 template<typename KV, typename host_ordinal_type_array>
585 LocalOrdinal
586 MultiVecAdapter<
587 MultiVector<Scalar,
588 LocalOrdinal,
589 GlobalOrdinal,
590 Node> >::gather (KV& kokkos_new_view,
591 host_ordinal_type_array &perm_g2l,
592 host_ordinal_type_array &recvCountRows,
593 host_ordinal_type_array &recvDisplRows,
594 EDistribution distribution) const
595 {
596 auto nCols = this->mv_->getNumVectors();
597 auto nRows = this->mv_->getGlobalLength();
598 auto comm = this->mv_->getMap()->getComm();
599 auto myRank = comm->getRank();
600 if (myRank == 0) {
601 Kokkos::resize(kokkos_new_view, nRows, nCols);
602 if (perm_g2l.extent(0) == nRows) {
603 Kokkos::resize(this->buf_, nRows, 1);
604 } else {
605 Kokkos::resize(this->buf_, 0, 1);
606 }
607 }
608 {
609 auto nRows_l = this->mv_->getLocalLength();
610 auto lclMV_d = this->mv_->getLocalViewDevice(Tpetra::Access::ReadOnly);
611 auto lclMV = Kokkos::create_mirror_view(lclMV_d);
612 Kokkos::deep_copy(lclMV, lclMV_d);
613 for (size_t j=0; j<nCols; j++) {
614 // lclMV with ReadOnly = const sendbuf
615 const Scalar * sendbuf = reinterpret_cast<const Scalar*> (nRows_l > 0 ? &lclMV(0,j) : lclMV.data());
616 Scalar * recvbuf = reinterpret_cast<Scalar*> (this->buf_.extent(0) > 0 || myRank != 0 ? this->buf_.data() : &kokkos_new_view(0,j));
617 Teuchos::gatherv<int, Scalar> (sendbuf, nRows_l,
618 recvbuf, recvCountRows.data(), recvDisplRows.data(),
619 0, *comm);
620 if (myRank == 0 && this->buf_.extent(0) > 0) {
621 for (global_size_t i=0; i<nRows; i++) kokkos_new_view(perm_g2l(i),j) = this->buf_(i,0);
622 }
623 }
624 }
625 return 0;
626 }
627
628 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, class Node >
629 template<typename KV, typename host_ordinal_type_array>
630 LocalOrdinal
631 MultiVecAdapter<
632 MultiVector<Scalar,
633 LocalOrdinal,
634 GlobalOrdinal,
635 Node> >::scatter (KV& kokkos_new_view,
636 host_ordinal_type_array &perm_g2l,
637 host_ordinal_type_array &sendCountRows,
638 host_ordinal_type_array &sendDisplRows,
639 EDistribution distribution) const
640 {
641 auto nCols = this->mv_->getNumVectors();
642 auto nRows = this->mv_->getGlobalLength();
643 auto nRows_l = this->mv_->getLocalLength();
644 auto comm = this->mv_->getMap()->getComm();
645 auto myRank = comm->getRank();
646 {
647 auto lclMV_d = this->mv_->getLocalViewDevice(Tpetra::Access::OverwriteAll);
648 auto lclMV = Kokkos::create_mirror_view(lclMV_d);
649 for (size_t j=0; j<nCols; j++) {
650 if (myRank == 0 && this->buf_.extent(0) > 0) {
651 for (global_size_t i=0; i<nRows; i++) this->buf_(i, 0) = kokkos_new_view(perm_g2l(i),j);
652 }
653 // lclMV with OverwriteAll
654 Scalar * sendbuf = reinterpret_cast<Scalar*> (this->buf_.extent(0) > 0 || myRank != 0 ? this->buf_.data() : &kokkos_new_view(0,j));
655 Scalar * recvbuf = reinterpret_cast<Scalar*> (nRows_l > 0 ? &lclMV(0,j) : lclMV.data());
656 Teuchos::scatterv<int, Scalar> (sendbuf, sendCountRows.data(), sendDisplRows.data(),
657 recvbuf, nRows_l,
658 0, *comm);
659 }
660 Kokkos::deep_copy(lclMV_d, lclMV);
661 }
662 return 0;
663 }
664
665 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, class Node >
666 std::string
667 MultiVecAdapter<
668 MultiVector<Scalar,
669 LocalOrdinal,
670 GlobalOrdinal,
671 Node> >::description() const
672 {
673 std::ostringstream oss;
674 oss << "Amesos2 adapter wrapping: ";
675 oss << mv_->description();
676 return oss.str();
677 }
678
679
680 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, class Node >
681 void
682 MultiVecAdapter<
683 MultiVector<Scalar,
684 LocalOrdinal,
685 GlobalOrdinal,
686 Node> >::describe (Teuchos::FancyOStream& os,
687 const Teuchos::EVerbosityLevel verbLevel) const
688 {
689 mv_->describe (os, verbLevel);
690 }
691
692
693 template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, class Node >
694 const char* MultiVecAdapter<
695 MultiVector<Scalar,
696 LocalOrdinal,
697 GlobalOrdinal,
698 Node> >::name = "Amesos2 adapter for Tpetra::MultiVector";
699
700} // end namespace Amesos2
701
702#endif // AMESOS2_TPETRA_MULTIVEC_ADAPTER_DEF_HPP
Copy or assign views based on memory spaces.
Amesos2::MultiVecAdapter specialization for the Tpetra::MultiVector class.
A templated MultiVector class adapter for Amesos2.
Definition Amesos2_MultiVecAdapter_decl.hpp:142