42 #ifndef BELOS_TPETRA_ADAPTER_MP_VECTOR_HPP 43 #define BELOS_TPETRA_ADAPTER_MP_VECTOR_HPP 45 #include "BelosTpetraAdapter.hpp" 47 #include "Kokkos_MV_GEMM.hpp" 49 #ifdef HAVE_BELOS_TSQR 51 #endif // HAVE_BELOS_TSQR 71 template<
class Storage,
class LO,
class GO,
class Node>
126 #ifdef HAVE_TPETRA_DEBUG 127 const char fnName[] =
"Belos::MultiVecTraits::CloneCopy(mv,index)";
128 const size_t inNumVecs = mv.getNumVectors ();
130 index.size () > 0 && *std::min_element (index.begin (), index.end ()) < 0,
131 std::runtime_error, fnName <<
": All indices must be nonnegative.");
134 static_cast<size_t> (*std::max_element (index.begin (), index.end ())) >= inNumVecs,
136 fnName <<
": All indices must be strictly less than the number of " 137 "columns " << inNumVecs <<
" of the input multivector mv.");
138 #endif // HAVE_TPETRA_DEBUG 142 for (std::vector<int>::size_type
j = 0;
j < index.size (); ++
j) {
143 columns[
j] = index[
j];
162 const bool validRange = index.
size() > 0 &&
164 index.
ubound() < GetNumberVecs(mv);
166 std::ostringstream os;
167 os <<
"Belos::MultiVecTraits::CloneCopy(mv,index=[" 170 index.
size() == 0, std::invalid_argument,
171 os.str() <<
"Empty index range is not allowed.");
173 index.
lbound() < 0, std::invalid_argument,
174 os.str() <<
"Index range includes negative index/ices, which is not " 177 index.
ubound() >= GetNumberVecs(mv), std::invalid_argument,
178 os.str() <<
"Index range exceeds number of vectors " 179 << mv.getNumVectors() <<
" in the input multivector.");
181 os.str() <<
"Should never get here!");
192 #ifdef HAVE_TPETRA_DEBUG 193 const char fnName[] =
"Belos::MultiVecTraits::CloneViewNonConst(mv,index)";
194 const size_t numVecs = mv.getNumVectors ();
196 index.size () > 0 && *std::min_element (index.begin (), index.end ()) < 0,
197 std::invalid_argument,
198 fnName <<
": All indices must be nonnegative.");
201 static_cast<size_t> (*std::max_element (index.begin (), index.end ())) >= numVecs,
202 std::invalid_argument,
203 fnName <<
": All indices must be strictly less than the number of " 204 "columns " << numVecs <<
" in the input MultiVector mv.");
205 #endif // HAVE_TPETRA_DEBUG 209 for (std::vector<int>::size_type
j = 0;
j < index.size (); ++
j) {
210 columns[
j] = index[
j];
226 const int numCols =
static_cast<int> (mv.getNumVectors());
227 const bool validRange = index.
size() > 0 &&
230 std::ostringstream os;
231 os <<
"Belos::MultiVecTraits::CloneViewNonConst(mv,index=[" 234 index.
size() == 0, std::invalid_argument,
235 os.str() <<
"Empty index range is not allowed.");
237 index.
lbound() < 0, std::invalid_argument,
238 os.str() <<
"Index range includes negative inde{x,ices}, which is " 241 index.
ubound() >= numCols, std::invalid_argument,
242 os.str() <<
"Index range exceeds number of vectors " << numCols
243 <<
" in the input multivector.");
245 true, std::logic_error,
246 os.str() <<
"Should never get here!");
256 #ifdef HAVE_TPETRA_DEBUG 257 const char fnName[] =
"Belos::MultiVecTraits<Scalar, " 258 "Tpetra::MultiVector<...> >::CloneView(mv,index)";
259 const size_t numVecs = mv.getNumVectors ();
261 *std::min_element (index.begin (), index.end ()) < 0,
262 std::invalid_argument,
263 fnName <<
": All indices must be nonnegative.");
265 static_cast<size_t> (*std::max_element (index.begin (), index.end ())) >= numVecs,
266 std::invalid_argument,
267 fnName <<
": All indices must be strictly less than the number of " 268 "columns " << numVecs <<
" in the input MultiVector mv.");
269 #endif // HAVE_TPETRA_DEBUG 273 for (std::vector<int>::size_type
j = 0;
j < index.size (); ++
j) {
274 columns[
j] = index[
j];
290 const int numCols =
static_cast<int> (mv.getNumVectors());
291 const bool validRange = index.
size() > 0 &&
294 std::ostringstream os;
295 os <<
"Belos::MultiVecTraits::CloneView(mv, index=[" 298 os.str() <<
"Empty index range is not allowed.");
300 os.str() <<
"Index range includes negative index/ices, which is not " 303 index.
ubound() >= numCols, std::invalid_argument,
304 os.str() <<
"Index range exceeds number of vectors " << numCols
305 <<
" in the input multivector.");
307 os.str() <<
"Should never get here!");
315 return static_cast<ptrdiff_t
> (mv.getGlobalLength ());
319 return static_cast<int> (mv.getNumVectors ());
323 return mv.isConstantStride ();
337 const int numRowsB =
B.numRows ();
338 const int numColsB =
B.numCols ();
339 const int strideB =
B.stride ();
340 if (numRowsB == 1 && numColsB == 1) {
341 C.update (alpha*
B(0,0),
A, beta);
346 RCP<const MV> Atmp, Ctmp;
348 else Atmp =
rcp(&
A,
false);
351 else Ctmp =
rcp(&
C,
false);
355 typedef typename FMV::dual_view_type::t_dev flat_view_type;
357 flat_view_type flat_A_view = Atmp->template getLocalView<execution_space>();
358 flat_view_type flat_C_view = Ctmp->template getLocalView<execution_space>();
361 typedef Kokkos::View<dot_type**, Kokkos::LayoutLeft, Kokkos::HostSpace> b_host_view_type;
362 b_host_view_type B_view_host(
B.values(), strideB, numColsB);
366 typedef Kokkos::View<dot_type**, Kokkos::LayoutLeft, execution_space> b_view_type;
367 typedef Kokkos::View<dot_type*, Kokkos::LayoutLeft, execution_space> b_1d_view_type;
368 b_1d_view_type B_1d_view_dev(Kokkos::ViewAllocateWithoutInitializing(
"B"), strideB*numColsB);
369 b_view_type B_view_dev( B_1d_view_dev.ptr_on_device(), strideB, numColsB);
373 Kokkos::DeviceGEMM<dot_type,execution_space>::GEMM(
375 alpha, flat_A_view, B_view_dev, beta, flat_C_view);
378 if (
C.isConstantStride() ==
false)
408 const std::vector<BaseScalar>& alphas)
410 std::vector<Scalar> alphas_mp(alphas.size());
411 const size_t sz = alphas.size();
412 for (
size_t i=0; i<sz; ++i)
413 alphas_mp[i] = alphas[i];
414 mv.scale (alphas_mp);
419 const std::vector<Scalar>& alphas)
433 using Teuchos::REDUCE_SUM;
437 const int numRowsC =
C.numRows ();
438 const int numColsC =
C.numCols ();
439 const int strideC =
C.stride ();
440 if (numRowsC == 1 && numColsC == 1) {
454 RCP<const MV> Atmp, Btmp;
456 else Atmp =
rcp(&
A,
false);
459 else Btmp =
rcp(&
B,
false);
463 typedef typename FMV::dual_view_type::t_dev flat_view_type;
465 flat_view_type flat_A_view = Atmp->template getLocalView<execution_space>();
466 flat_view_type flat_B_view = Btmp->template getLocalView<execution_space>();
469 typedef Kokkos::View<dot_type**, Kokkos::LayoutLeft, Kokkos::HostSpace> c_host_view_type;
470 c_host_view_type C_view_host(
C.values(), strideC, numColsC);
474 typedef Kokkos::View<dot_type**, Kokkos::LayoutLeft, execution_space> c_view_type;
475 typedef Kokkos::View<dot_type*, Kokkos::LayoutLeft, execution_space> c_1d_view_type;
476 c_1d_view_type C_1d_view_dev(
"C", strideC*numColsC);
477 c_view_type C_view_dev( C_1d_view_dev.ptr_on_device(), strideC, numColsC);
480 Kokkos::DeviceGEMM<dot_type,execution_space>::GEMM(
482 alpha, flat_A_view, flat_B_view,
483 Kokkos::Details::ArithTraits<dot_type>::zero(),
487 RCP<const Comm<int> > pcomm =
A.getMap()->getComm ();
488 if (pcomm->getSize () == 1)
491 typedef Kokkos::View<dot_type*, Kokkos::LayoutLeft, Kokkos::HostSpace> c_1d_host_view_type;
492 c_1d_host_view_type C_1d_view_tmp(Kokkos::ViewAllocateWithoutInitializing(
"C_tmp"), strideC*numColsC);
493 c_host_view_type C_view_tmp( C_1d_view_tmp.ptr_on_device(),
496 reduceAll<int> (*pcomm, REDUCE_SUM, strideC*numColsC,
497 C_view_tmp.ptr_on_device(),
498 C_view_host.ptr_on_device());
506 std::vector<dot_type>& dots)
508 const size_t numVecs =
A.getNumVectors ();
511 numVecs !=
B.getNumVectors (), std::invalid_argument,
512 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::MvDot(A,B,dots): " 513 "A and B must have the same number of columns. " 514 "A has " << numVecs <<
" column(s), " 515 "but B has " <<
B.getNumVectors () <<
" column(s).");
516 #ifdef HAVE_TPETRA_DEBUG 518 dots.size() < numVecs, std::invalid_argument,
519 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::MvDot(A,B,dots): " 520 "The output array 'dots' must have room for all dot products. " 521 "A and B each have " << numVecs <<
" column(s), " 522 "but 'dots' only has " << dots.size() <<
" entry(/ies).");
523 #endif // HAVE_TPETRA_DEBUG 526 A.dot (
B, av (0, numVecs));
532 std::vector<mag_type> &normvec,
533 NormType type=TwoNorm)
536 #ifdef HAVE_TPETRA_DEBUG 537 typedef std::vector<int>::size_type size_type;
539 normvec.size () <
static_cast<size_type
> (mv.getNumVectors ()),
540 std::invalid_argument,
541 "Belos::MultiVecTraits::MvNorm(mv,normvec): The normvec output " 542 "argument must have at least as many entries as the number of vectors " 543 "(columns) in the MultiVector mv. normvec.size() = " << normvec.size ()
544 <<
" < mv.getNumVectors() = " << mv.getNumVectors () <<
".");
550 mv.norm1(av(0,mv.getNumVectors()));
553 mv.norm2(av(0,mv.getNumVectors()));
556 mv.normInf(av(0,mv.getNumVectors()));
563 true, std::logic_error,
564 "Belos::MultiVecTraits::MvNorm: Invalid NormType value " << type
565 <<
". Valid values are OneNorm=" << OneNorm <<
", TwoNorm=" 566 << TwoNorm <<
", and InfNorm=" << InfNorm <<
". If you are a Belos " 567 "user and have not modified Belos in any way, and you get this " 568 "message, then this is probably a bug in the Belos solver you were " 569 "using. Please report this to the Belos developers.");
578 const size_t inNumVecs =
A.getNumVectors ();
579 #ifdef HAVE_TPETRA_DEBUG 581 inNumVecs < static_cast<size_t> (index.size ()), std::invalid_argument,
582 "Belos::MultiVecTraits::SetBlock(A,index,mv): 'index' argument must " 583 "have no more entries as the number of columns in the input MultiVector" 584 " A. A.getNumVectors() = " << inNumVecs <<
" < index.size () = " 585 << index.size () <<
".");
586 #endif // HAVE_TPETRA_DEBUG 587 RCP<MV> mvsub = CloneViewNonConst (mv, index);
588 if (inNumVecs > static_cast<size_t> (index.size ())) {
589 RCP<const MV> Asub =
A.subView (Range1D (0, index.size () - 1));
606 const size_t maxInt =
608 const bool overflow =
609 maxInt <
A.getNumVectors () && maxInt < mv.getNumVectors ();
611 std::ostringstream os;
612 os <<
"Belos::MultiVecTraits::SetBlock(A, index=[" << index.
lbound ()
613 <<
", " << index.
ubound () <<
"], mv): ";
615 maxInt <
A.getNumVectors (), std::range_error, os.str () <<
"Number " 616 "of columns (size_t) in the input MultiVector 'A' overflows int.");
618 maxInt < mv.getNumVectors (), std::range_error, os.str () <<
"Number " 619 "of columns (size_t) in the output MultiVector 'mv' overflows int.");
622 const int numColsA =
static_cast<int> (
A.getNumVectors ());
623 const int numColsMv =
static_cast<int> (mv.getNumVectors ());
625 const bool validIndex =
628 const bool validSource = index.
size () <= numColsA;
630 if (! validIndex || ! validSource) {
631 std::ostringstream os;
632 os <<
"Belos::MultiVecTraits::SetBlock(A, index=[" << index.
lbound ()
633 <<
", " << index.
ubound () <<
"], mv): ";
635 index.
lbound() < 0, std::invalid_argument,
636 os.str() <<
"Range lower bound must be nonnegative.");
638 index.
ubound() >= numColsMv, std::invalid_argument,
639 os.str() <<
"Range upper bound must be less than the number of " 640 "columns " << numColsA <<
" in the 'mv' output argument.");
642 index.
size() > numColsA, std::invalid_argument,
643 os.str() <<
"Range must have no more elements than the number of " 644 "columns " << numColsA <<
" in the 'A' input argument.");
646 true, std::logic_error,
"Should never get here!");
653 if (index.
lbound () == 0 && index.
ubound () + 1 == numColsMv) {
654 mv_view = Teuchos::rcpFromRef (mv);
656 mv_view = CloneViewNonConst (mv, index);
663 if (index.
size () == numColsA) {
664 A_view = Teuchos::rcpFromRef (
A);
674 const char errPrefix[] =
"Belos::MultiVecTraits::Assign(A, mv): ";
683 const size_t maxInt =
685 const bool overflow =
686 maxInt <
A.getNumVectors () && maxInt < mv.getNumVectors ();
689 maxInt <
A.getNumVectors(), std::range_error,
690 errPrefix <<
"Number of columns in the input multivector 'A' " 691 "(a size_t) overflows int.");
693 maxInt < mv.getNumVectors(), std::range_error,
694 errPrefix <<
"Number of columns in the output multivector 'mv' " 695 "(a size_t) overflows int.");
697 true, std::logic_error,
"Should never get here!");
700 const int numColsA =
static_cast<int> (
A.getNumVectors ());
701 const int numColsMv =
static_cast<int> (mv.getNumVectors ());
702 if (numColsA > numColsMv) {
704 numColsA > numColsMv, std::invalid_argument,
705 errPrefix <<
"Input multivector 'A' has " << numColsA <<
" columns, " 706 "but output multivector 'mv' has only " << numColsMv <<
" columns.");
709 if (numColsA == numColsMv) {
726 mv.putScalar (alpha);
734 #ifdef HAVE_BELOS_TSQR 735 typedef Tpetra::TsqrAdaptor< Tpetra::MultiVector< Scalar, LO, GO, Node > > tsqr_adaptor_type;
739 #endif // HAVE_BELOS_TSQR 749 template <
class Storage,
class LO,
class GO,
class Node>
753 Tpetra::Operator<Sacado::MP::Vector<Storage>,
759 Apply (
const Tpetra::Operator<Scalar,LO,GO,Node>& Op,
762 ETrans trans=NOTRANS)
779 const std::string otName =
"Belos::OperatorTraits<" + scalarName
780 +
"," + loName +
"," + goName +
"," + nodeName +
">";
782 "get here; fell through a switch statement. " 783 "Please report this bug to the Belos developers.");
790 return Op.hasTransposeApply ();
Tpetra::MultiVector< Scalar, LO, GO, Node >::mag_type mag_type
static int GetNumberVecs(const MV &mv)
static bool HasApplyTranspose(const Tpetra::Operator< Scalar, LO, GO, Node > &Op)
static void Assign(const MV &A, MV &mv)
static Teuchos::RCP< MV > CloneCopy(const MV &mv, const std::vector< int > &index)
Create and return a deep copy of the given columns of mv.
static void MvScale(Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const std::vector< BaseScalar > &alphas)
Kokkos::DefaultExecutionSpace execution_space
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Storage::ordinal_type s_ordinal
Tpetra::MultiVector< Scalar, LO, GO, Node > MV
Sacado::MP::Vector< Storage > Scalar
static bool HasConstantStride(const MV &mv)
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const Teuchos::Range1D &index)
static void MvTimesMatAddMv(const dot_type &alpha, const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, const Teuchos::SerialDenseMatrix< int, dot_type > &B, const dot_type &beta, Tpetra::MultiVector< Scalar, LO, GO, Node > &C)
static void MvScale(Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const Scalar &alpha)
static void SetBlock(const MV &A, const std::vector< int > &index, MV &mv)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
static void MvPrint(const MV &mv, std::ostream &os)
static void MvScale(Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const std::vector< Scalar > &alphas)
static void MvTransMv(dot_type alpha, const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, const Tpetra::MultiVector< Scalar, LO, GO, Node > &B, Teuchos::SerialDenseMatrix< int, dot_type > &C)
static ptrdiff_t GetGlobalLength(const MV &mv)
TEUCHOS_DEPRECATED void reduceAll(const Comm< Ordinal > &comm, const EReductionType reductType, const Packet &send, Packet *globalReduct)
void deep_copy(const Stokhos::CrsMatrix< ValueType, DstDevice, Layout > &dst, const Stokhos::CrsMatrix< ValueType, SrcDevice, Layout > &src)
static void MvNorm(const Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, std::vector< mag_type > &normvec, NormType type=TwoNorm)
For all columns j of mv, set normvec[j] = norm(mv[j]).
Storage::value_type BaseScalar
Tpetra::MultiVector< Scalar, LO, GO, Node >::dot_type dot_type
static Teuchos::RCP< MV > CloneCopy(const MV &X)
Create and return a deep copy of X.
static Teuchos::RCP< MV > Clone(const MV &X, const int numVecs)
Create a new MultiVector with numVecs columns.
static void Apply(const Tpetra::Operator< Scalar, LO, GO, Node > &Op, const Tpetra::MultiVector< Scalar, LO, GO, Node > &X, Tpetra::MultiVector< Scalar, LO, GO, Node > &Y, ETrans trans=NOTRANS)
static Teuchos::RCP< MV > CloneCopy(const MV &mv, const Teuchos::Range1D &index)
Create and return a deep copy of the given columns of mv.
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const std::vector< int > &index)
Sacado::MP::Vector< Storage > Scalar
static void MvRandom(MV &mv)
static void MvInit(MV &mv, const Scalar alpha=Teuchos::ScalarTraits< Scalar >::zero())
static Teuchos::RCP< const MV > CloneView(const MV &mv, const Teuchos::Range1D &index)
static Teuchos::RCP< const MV > CloneView(const MV &mv, const std::vector< int > &index)
static std::string name()
static void MvDot(const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, const Tpetra::MultiVector< Scalar, LO, GO, Node > &B, std::vector< dot_type > &dots)
For all columns j of A, set dots[j] := A[j]^T * B[j].
static void SetBlock(const MV &A, const Teuchos::Range1D &index, MV &mv)
static void MvAddMv(Scalar alpha, const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, Scalar beta, const Tpetra::MultiVector< Scalar, LO, GO, Node > &B, Tpetra::MultiVector< Scalar, LO, GO, Node > &mv)
mv := alpha*A + beta*B