42 #ifndef BELOS_TPETRA_ADAPTER_UQ_PCE_HPP
43 #define BELOS_TPETRA_ADAPTER_UQ_PCE_HPP
45 #include "BelosTpetraAdapter.hpp"
47 #include "KokkosBlas.hpp"
49 #ifdef HAVE_BELOS_TSQR
51 #endif // HAVE_BELOS_TSQR
71 template<
class Storage,
class LO,
class GO,
class Node>
73 Tpetra::MultiVector< Sacado::UQ::PCE<Storage>,
79 typedef Tpetra::MultiVector<Scalar, LO, GO, Node>
MV;
81 typedef typename Tpetra::MultiVector<Scalar,LO,GO,Node>::dot_type
dot_type;
82 typedef typename Tpetra::MultiVector<Scalar,LO,GO,Node>::mag_type
mag_type;
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 ();
328 const Tpetra::MultiVector<Scalar,LO,GO,Node>&
A,
331 Tpetra::MultiVector<Scalar,LO,GO,Node>&
C)
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);
354 typedef Tpetra::MultiVector<dot_type,LO,GO,Node> FMV;
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_input( B.
values(), strideB, numColsB);
363 auto B_view_host = Kokkos::subview( B_view_host_input,
364 Kokkos::pair<int,int>(0, numRowsB),
365 Kokkos::pair<int,int>(0, numColsB));
369 typedef Kokkos::View<dot_type**, Kokkos::LayoutLeft, execution_space> b_view_type;
370 typedef Kokkos::View<dot_type*, Kokkos::LayoutLeft, execution_space> b_1d_view_type;
371 b_1d_view_type B_1d_view_dev(Kokkos::ViewAllocateWithoutInitializing(
"B"), numRowsB*numColsB);
372 b_view_type B_view_dev( B_1d_view_dev.data(), numRowsB, numColsB);
377 const char ctransA =
'N', ctransB =
'N';
380 alpha, flat_A_view, B_view_dev, beta, flat_C_view);
383 if (C.isConstantStride() ==
false)
396 const Tpetra::MultiVector<Scalar,LO,GO,Node>&
A,
398 const Tpetra::MultiVector<Scalar,LO,GO,Node>&
B,
399 Tpetra::MultiVector<Scalar,LO,GO,Node>& mv)
405 MvScale (Tpetra::MultiVector<Scalar,LO,GO,Node>& mv,
412 MvScale (Tpetra::MultiVector<Scalar,LO,GO,Node>& mv,
413 const std::vector<BaseScalar>& alphas)
415 std::vector<Scalar> alphas_mp(alphas.size());
416 const size_t sz = alphas.size();
417 for (
size_t i=0; i<sz; ++i)
418 alphas_mp[i] = alphas[i];
419 mv.scale (alphas_mp);
423 MvScale (Tpetra::MultiVector<Scalar,LO,GO,Node>& mv,
424 const std::vector<Scalar>& alphas)
431 const Tpetra::MultiVector<Scalar,LO,GO,Node>&
A,
432 const Tpetra::MultiVector<Scalar,LO,GO,Node>&
B,
439 using Teuchos::reduceAll;
442 const int numRowsC = C.
numRows ();
443 const int numColsC = C.
numCols ();
444 const int strideC = C.
stride ();
445 if (numRowsC == 1 && numColsC == 1) {
459 RCP<const MV> Atmp, Btmp;
461 else Atmp =
rcp(&A,
false);
464 else Btmp =
rcp(&B,
false);
467 typedef Tpetra::MultiVector<dot_type,LO,GO,Node> FMV;
468 typedef typename FMV::dual_view_type::t_dev flat_view_type;
470 flat_view_type flat_A_view = Atmp->template getLocalView<execution_space>();
471 flat_view_type flat_B_view = Btmp->template getLocalView<execution_space>();
474 typedef Kokkos::View<dot_type**, Kokkos::LayoutLeft, Kokkos::HostSpace> c_host_view_type;
475 c_host_view_type C_view_host_input( C.
values(), strideC, numColsC);
476 auto C_view_host = Kokkos::subview(C_view_host_input,
477 Kokkos::pair<int,int>(0, numRowsC),
478 Kokkos::pair<int,int>(0, numColsC));
482 typedef Kokkos::View<dot_type**, Kokkos::LayoutLeft, execution_space> c_view_type;
483 typedef Kokkos::View<dot_type*, Kokkos::LayoutLeft, execution_space> c_1d_view_type;
484 c_1d_view_type C_1d_view_dev(
"C", numRowsC*numColsC);
485 c_view_type C_view_dev( C_1d_view_dev.data(), numRowsC, numColsC);
489 const char ctransA =
'C', ctransB =
'N';
492 alpha, flat_A_view, flat_B_view,
493 Kokkos::Details::ArithTraits<dot_type>::zero(),
497 RCP<const Comm<int> > pcomm = A.getMap()->getComm ();
498 if (pcomm->getSize () == 1)
501 typedef Kokkos::View<dot_type*, Kokkos::LayoutLeft, Kokkos::HostSpace> c_1d_host_view_type;
502 c_1d_host_view_type C_1d_view_tmp(Kokkos::ViewAllocateWithoutInitializing(
"C_tmp"), strideC*numColsC);
503 c_host_view_type C_view_tmp( C_1d_view_tmp.data(),
506 Kokkos::pair<int,int>(0, numRowsC),
507 Kokkos::pair<int,int>(0, numColsC)),
509 reduceAll<int> (*pcomm, REDUCE_SUM, strideC*numColsC,
517 MvDot (
const Tpetra::MultiVector<Scalar,LO,GO,Node>&
A,
518 const Tpetra::MultiVector<Scalar,LO,GO,Node>&
B,
519 std::vector<dot_type>& dots)
521 const size_t numVecs = A.getNumVectors ();
524 numVecs != B.getNumVectors (), std::invalid_argument,
525 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::MvDot(A,B,dots): "
526 "A and B must have the same number of columns. "
527 "A has " << numVecs <<
" column(s), "
528 "but B has " << B.getNumVectors () <<
" column(s).");
529 #ifdef HAVE_TPETRA_DEBUG
531 dots.size() < numVecs, std::invalid_argument,
532 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::MvDot(A,B,dots): "
533 "The output array 'dots' must have room for all dot products. "
534 "A and B each have " << numVecs <<
" column(s), "
535 "but 'dots' only has " << dots.size() <<
" entry(/ies).");
536 #endif // HAVE_TPETRA_DEBUG
539 A.dot (B, av (0, numVecs));
544 MvNorm (
const Tpetra::MultiVector<Scalar,LO,GO,Node>& mv,
545 std::vector<mag_type> &normvec,
546 NormType type=TwoNorm)
549 #ifdef HAVE_TPETRA_DEBUG
550 typedef std::vector<int>::size_type size_type;
552 normvec.size () <
static_cast<size_type
> (mv.getNumVectors ()),
553 std::invalid_argument,
554 "Belos::MultiVecTraits::MvNorm(mv,normvec): The normvec output "
555 "argument must have at least as many entries as the number of vectors "
556 "(columns) in the MultiVector mv. normvec.size() = " << normvec.size ()
557 <<
" < mv.getNumVectors() = " << mv.getNumVectors () <<
".");
563 mv.norm1(av(0,mv.getNumVectors()));
566 mv.norm2(av(0,mv.getNumVectors()));
569 mv.normInf(av(0,mv.getNumVectors()));
576 true, std::logic_error,
577 "Belos::MultiVecTraits::MvNorm: Invalid NormType value " << type
578 <<
". Valid values are OneNorm=" << OneNorm <<
", TwoNorm="
579 << TwoNorm <<
", and InfNorm=" << InfNorm <<
". If you are a Belos "
580 "user and have not modified Belos in any way, and you get this "
581 "message, then this is probably a bug in the Belos solver you were "
582 "using. Please report this to the Belos developers.");
591 const size_t inNumVecs = A.getNumVectors ();
592 #ifdef HAVE_TPETRA_DEBUG
594 inNumVecs < static_cast<size_t> (index.size ()), std::invalid_argument,
595 "Belos::MultiVecTraits::SetBlock(A,index,mv): 'index' argument must "
596 "have no more entries as the number of columns in the input MultiVector"
597 " A. A.getNumVectors() = " << inNumVecs <<
" < index.size () = "
598 << index.size () <<
".");
599 #endif // HAVE_TPETRA_DEBUG
600 RCP<MV> mvsub = CloneViewNonConst (mv, index);
601 if (inNumVecs > static_cast<size_t> (index.size ())) {
602 RCP<const MV> Asub = A.subView (Range1D (0, index.size () - 1));
619 const size_t maxInt =
621 const bool overflow =
622 maxInt < A.getNumVectors () && maxInt < mv.getNumVectors ();
624 std::ostringstream os;
625 os <<
"Belos::MultiVecTraits::SetBlock(A, index=[" << index.
lbound ()
626 <<
", " << index.
ubound () <<
"], mv): ";
628 maxInt < A.getNumVectors (), std::range_error, os.str () <<
"Number "
629 "of columns (size_t) in the input MultiVector 'A' overflows int.");
631 maxInt < mv.getNumVectors (), std::range_error, os.str () <<
"Number "
632 "of columns (size_t) in the output MultiVector 'mv' overflows int.");
635 const int numColsA =
static_cast<int> (A.getNumVectors ());
636 const int numColsMv =
static_cast<int> (mv.getNumVectors ());
638 const bool validIndex =
641 const bool validSource = index.
size () <= numColsA;
643 if (! validIndex || ! validSource) {
644 std::ostringstream os;
645 os <<
"Belos::MultiVecTraits::SetBlock(A, index=[" << index.
lbound ()
646 <<
", " << index.
ubound () <<
"], mv): ";
648 index.
lbound() < 0, std::invalid_argument,
649 os.str() <<
"Range lower bound must be nonnegative.");
651 index.
ubound() >= numColsMv, std::invalid_argument,
652 os.str() <<
"Range upper bound must be less than the number of "
653 "columns " << numColsA <<
" in the 'mv' output argument.");
655 index.
size() > numColsA, std::invalid_argument,
656 os.str() <<
"Range must have no more elements than the number of "
657 "columns " << numColsA <<
" in the 'A' input argument.");
659 true, std::logic_error,
"Should never get here!");
666 if (index.
lbound () == 0 && index.
ubound () + 1 == numColsMv) {
667 mv_view = Teuchos::rcpFromRef (mv);
669 mv_view = CloneViewNonConst (mv, index);
676 if (index.
size () == numColsA) {
677 A_view = Teuchos::rcpFromRef (A);
687 const char errPrefix[] =
"Belos::MultiVecTraits::Assign(A, mv): ";
696 const size_t maxInt =
698 const bool overflow =
699 maxInt < A.getNumVectors () && maxInt < mv.getNumVectors ();
702 maxInt < A.getNumVectors(), std::range_error,
703 errPrefix <<
"Number of columns in the input multivector 'A' "
704 "(a size_t) overflows int.");
706 maxInt < mv.getNumVectors(), std::range_error,
707 errPrefix <<
"Number of columns in the output multivector 'mv' "
708 "(a size_t) overflows int.");
710 true, std::logic_error,
"Should never get here!");
713 const int numColsA =
static_cast<int> (A.getNumVectors ());
714 const int numColsMv =
static_cast<int> (mv.getNumVectors ());
715 if (numColsA > numColsMv) {
717 numColsA > numColsMv, std::invalid_argument,
718 errPrefix <<
"Input multivector 'A' has " << numColsA <<
" columns, "
719 "but output multivector 'mv' has only " << numColsMv <<
" columns.");
722 if (numColsA == numColsMv) {
739 mv.putScalar (alpha);
747 #ifdef HAVE_BELOS_TSQR
748 typedef Tpetra::TsqrAdaptor< Tpetra::MultiVector< Scalar, LO, GO, Node > > tsqr_adaptor_type;
752 #endif // HAVE_BELOS_TSQR
762 template <
class Storage,
class LO,
class GO,
class Node>
764 Tpetra::MultiVector<Sacado::UQ::PCE<Storage>,
766 Tpetra::Operator<Sacado::UQ::PCE<Storage>,
772 Apply (
const Tpetra::Operator<Scalar,LO,GO,Node>& Op,
773 const Tpetra::MultiVector<Scalar,LO,GO,Node>& X,
774 Tpetra::MultiVector<Scalar,LO,GO,Node>& Y,
775 ETrans trans=NOTRANS)
792 const std::string otName =
"Belos::OperatorTraits<" + scalarName
793 +
"," + loName +
"," + goName +
"," + nodeName +
">";
795 "get here; fell through a switch statement. "
796 "Please report this bug to the Belos developers.");
803 return Op.hasTransposeApply ();
ScalarType * values() const
Tpetra::MultiVector< Scalar, LO, GO, Node > MV
static Teuchos::RCP< MV > CloneCopy(const MV &X)
Create and return a deep copy of X.
static void SetBlock(const MV &A, const std::vector< int > &index, MV &mv)
static bool HasConstantStride(const MV &mv)
Sacado::UQ::PCE< Storage > Scalar
Tpetra::MultiVector< Scalar, LO, GO, Node >::mag_type mag_type
Kokkos::DefaultExecutionSpace execution_space
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
static void MvPrint(const MV &mv, std::ostream &os)
static Teuchos::RCP< const MV > CloneView(const MV &mv, const Teuchos::Range1D &index)
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)
Tpetra::MultiVector< Scalar, LO, GO, Node >::dot_type dot_type
static void MvScale(Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const std::vector< Scalar > &alphas)
static bool HasApplyTranspose(const Tpetra::Operator< Scalar, LO, GO, Node > &Op)
Sacado::UQ::PCE< Storage > Scalar
static void SetBlock(const MV &A, const Teuchos::Range1D &index, MV &mv)
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const std::vector< int > &index)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
static void MvScale(Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const Scalar &alpha)
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 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 Teuchos::RCP< MV > CloneCopy(const MV &mv, const Teuchos::Range1D &index)
Create and return a deep copy of the given columns of mv.
void deep_copy(const Stokhos::CrsMatrix< ValueType, DstDevice, Layout > &dst, const Stokhos::CrsMatrix< ValueType, SrcDevice, Layout > &src)
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
static void MvScale(Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const std::vector< BaseScalar > &alphas)
OrdinalType numCols() const
static void MvRandom(MV &mv)
static int GetNumberVecs(const 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.
std::enable_if< Kokkos::is_view_mp_vector< Kokkos::View< DA, PA...> >::value &&Kokkos::is_view_mp_vector< Kokkos::View< DB, PB...> >::value &&Kokkos::is_view_mp_vector< Kokkos::View< DC, PC...> >::value >::type gemm(const char transA[], const char transB[], typename Kokkos::View< DA, PA...>::const_value_type &alpha, const Kokkos::View< DA, PA...> &A, const Kokkos::View< DB, PB...> &B, typename Kokkos::View< DC, PC...>::const_value_type &beta, const Kokkos::View< DC, PC...> &C)
static Teuchos::RCP< const MV > CloneView(const MV &mv, const std::vector< int > &index)
static Teuchos::RCP< MV > Clone(const MV &X, const int numVecs)
Create a new MultiVector with numVecs columns.
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]).
static void Assign(const MV &A, MV &mv)
Storage::value_type BaseScalar
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const Teuchos::Range1D &index)
static ptrdiff_t GetGlobalLength(const MV &mv)
OrdinalType stride() const
Storage::ordinal_type s_ordinal
static std::string name()
OrdinalType numRows() const
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)