42 #ifndef BELOS_PCE_TPETRA_ADAPTER_HPP 43 #define BELOS_PCE_TPETRA_ADAPTER_HPP 51 #include <Tpetra_MultiVector.hpp> 52 #include <Tpetra_Operator.hpp> 59 #include <BelosConfigDefs.hpp> 60 #include <BelosTypes.hpp> 61 #include <BelosMultiVecTraits.hpp> 62 #include <BelosOperatorTraits.hpp> 64 #ifdef HAVE_BELOS_TSQR 65 # include <Tpetra_TsqrAdaptor.hpp> 66 #endif // HAVE_BELOS_TSQR 85 template<
class BaseScalar,
class Storage,
class LO,
class GO,
class Node>
86 class MultiVecTraits<BaseScalar,
Tpetra::MultiVector< Sacado::PCE::OrthogPoly<BaseScalar,Storage>,LO,GO,Node> >
89 #ifdef HAVE_BELOS_TPETRA_TIMERS 97 return Teuchos::rcp(
new Tpetra::MultiVector<Scalar,LO,GO,Node>(mv.getMap(),numvecs));
102 return Teuchos::rcp(
new Tpetra::MultiVector<Scalar,LO,GO,Node>( mv ) );
108 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::CloneCopy(mv,index): numvecs must be greater than zero.");
109 #ifdef HAVE_TPETRA_DEBUG 111 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::CloneCopy(mv,index): indices must be >= zero.");
113 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::CloneCopy(mv,index): indices must be < mv.getNumVectors().");
115 for (
typename std::vector<int>::size_type
j=1;
j<index.size(); ++
j) {
116 if (index[
j] != index[
j-1]+1) {
119 return mv.subCopy(stinds);
127 CloneCopy (
const Tpetra::MultiVector<Scalar,LO,GO,Node>& mv,
130 const bool validRange = index.
size() > 0 &&
132 index.
ubound() < GetNumberVecs(mv);
135 std::ostringstream os;
136 os <<
"Belos::MultiVecTraits<Scalar, Tpetra::MultiVector<...> >::" 137 "CloneCopy(mv,index=[" << index.
lbound() <<
", " << index.
ubound()
140 os.str() <<
"Empty index range is not allowed.");
142 os.str() <<
"Index range includes negative " 143 "index/ices, which is not allowed.");
146 std::invalid_argument,
147 os.str() <<
"Index range exceeds number of vectors " 148 << mv.getNumVectors() <<
" in the input multivector.");
150 os.str() <<
"Should never get here!");
152 return mv.subCopy (index);
159 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::CloneView(mv,index): numvecs must be greater than zero.");
160 #ifdef HAVE_TPETRA_DEBUG 162 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::CloneView(mv,index): indices must be >= zero.");
163 TEUCHOS_TEST_FOR_EXCEPTION( (
size_t)*std::max_element(index.begin(),index.end()) >= mv.getNumVectors(), std::invalid_argument,
164 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::CloneView(mv,index): indices must be < mv.getNumVectors().");
166 for (
typename std::vector<int>::size_type
j=1;
j<index.size(); ++
j) {
167 if (index[
j] != index[
j-1]+1) {
170 return mv.subViewNonConst(stinds);
185 const int numCols =
static_cast<int> (mv.getNumVectors());
186 const bool validRange = index.
size() > 0 &&
190 std::ostringstream os;
191 os <<
"Belos::MultiVecTraits<Scalar, Tpetra::MultiVector<...> >::" 192 "CloneViewNonConst(mv,index=[" << index.
lbound() <<
", " 193 << index.
ubound() <<
"]): ";
195 os.str() <<
"Empty index range is not allowed.");
197 os.str() <<
"Index range includes negative " 198 "index/ices, which is not allowed.");
200 os.str() <<
"Index range exceeds number of " 201 "vectors " << numCols <<
" in the input " 204 os.str() <<
"Should never get here!");
206 return mv.subViewNonConst (index);
213 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::CloneView(mv,index): numvecs must be greater than zero.");
214 #ifdef HAVE_TPETRA_DEBUG 216 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::CloneView(mv,index): indices must be >= zero.");
217 TEUCHOS_TEST_FOR_EXCEPTION( (
size_t)*std::max_element(index.begin(),index.end()) >= mv.getNumVectors(), std::invalid_argument,
218 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::CloneView(mv,index): indices must be < mv.getNumVectors().");
220 for (
typename std::vector<int>::size_type
j=1;
j<index.size(); ++
j) {
221 if (index[
j] != index[
j-1]+1) {
224 return mv.subView(stinds);
232 CloneView (
const Tpetra::MultiVector<Scalar,LO,GO,Node>& mv,
238 const int numCols =
static_cast<int> (mv.getNumVectors());
239 const bool validRange = index.
size() > 0 &&
243 std::ostringstream os;
244 os <<
"Belos::MultiVecTraits<Scalar, Tpetra::MultiVector<...> >::" 245 "CloneView(mv, index=[" << index.
lbound() <<
", " 246 << index.
ubound() <<
"]): ";
248 os.str() <<
"Empty index range is not allowed.");
250 os.str() <<
"Index range includes negative " 251 "index/ices, which is not allowed.");
253 os.str() <<
"Index range exceeds number of " 254 "vectors " << numCols <<
" in the input " 257 os.str() <<
"Should never get here!");
259 return mv.subView (index);
263 {
return static_cast<ptrdiff_t
>(mv.getGlobalLength()); }
265 static int GetNumberVecs(
const Tpetra::MultiVector<Scalar,LO,GO,Node>& mv )
266 {
return mv.getNumVectors(); }
269 {
return mv.isConstantStride(); }
273 Scalar beta, Tpetra::MultiVector<Scalar,LO,GO,Node>& mv )
276 for (
int i=0; i<
B.numRows(); i++)
277 for (
int j=0;
j<
B.numCols();
j++)
279 MvTimesMatAddMv(alpha,
A, B_pce, beta, mv);
283 Scalar beta, Tpetra::MultiVector<Scalar,LO,GO,Node>& mv )
285 #ifdef HAVE_BELOS_TPETRA_TIMERS 290 Tpetra::Map<LO,GO,Node> LocalMap(
B.numRows(), 0, Teuchos::rcpFromRef< const Teuchos::Comm<int> >(scomm), Tpetra::LocallyReplicated);
294 Tpetra::MultiVector<Scalar,LO,GO,Node> B_mv(Teuchos::rcpFromRef(LocalMap),Bvalues,
B.stride(),
B.numCols());
299 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 )
304 static void MvScale ( Tpetra::MultiVector<Scalar,LO,GO,Node>& mv,
Scalar alpha )
307 static void MvScale ( Tpetra::MultiVector<Scalar,LO,GO,Node>& mv,
const std::vector<BaseScalar>& alphas )
309 std::vector<Scalar> alphas_pce(alphas.size());
310 for (
int i=0; i<alphas.size(); i++) alphas_pce[i] = alphas[i];
311 mv.scale(alphas_pce);
313 static void MvScale ( Tpetra::MultiVector<Scalar,LO,GO,Node>& mv,
const std::vector<Scalar>& alphas )
314 { mv.scale(alphas); }
316 static void MvTransMv(
Scalar alpha,
const Tpetra::MultiVector<Scalar,LO,GO,Node>&
A,
const Tpetra::MultiVector<Scalar,LO,GO,Node>&
B,
Teuchos::SerialDenseMatrix<int,BaseScalar>&
C)
319 MvTransMv(alpha,
A,
B, C_pce);
320 for (
int i=0; i<
C.numRows(); i++)
321 for (
int j=0;
j<
C.numCols();
j++)
322 C(i,
j) = C_pce(i,
j).coeff(0);
324 static void MvTransMv(
Scalar alpha,
const Tpetra::MultiVector<Scalar,LO,GO,Node>&
A,
const Tpetra::MultiVector<Scalar,LO,GO,Node>&
B,
Teuchos::SerialDenseMatrix<int,Scalar>&
C)
326 #ifdef HAVE_BELOS_TPETRA_TIMERS 334 const int numRowsC =
C.numRows(),
335 numColsC =
C.numCols(),
336 strideC =
C.stride();
339 Tpetra::Map<LO,GO,Node> LocalMap(numRowsC, 0, Teuchos::rcpFromRef<
const Teuchos::Comm<int> >(scomm), Tpetra::LocallyReplicated);
341 const bool INIT_TO_ZERO =
true;
342 Tpetra::MultiVector<Scalar,LO,GO,Node> C_mv(Teuchos::rcpFromRef(LocalMap),numColsC, INIT_TO_ZERO);
349 if (pcomm->getSize() == 1) {
352 C_mv.get1dCopy(C_view,strideC);
357 if (strideC == numRowsC) {
359 Teuchos::reduceAll<int,Scalar>(*pcomm,Teuchos::REDUCE_SUM,numColsC*numRowsC,C_mv_view.
getRawPtr(),C_view.getRawPtr());
364 Teuchos::reduceAll<int,Scalar>(*pcomm,Teuchos::REDUCE_SUM,numColsC*numRowsC,C_mv_view.
getRawPtr(),destBuff.
getRawPtr());
365 for (
int j=0;
j < numColsC; ++
j) {
366 for (
int i=0; i < numRowsC; ++i) {
367 C_view[strideC*
j+i] = destBuff[numRowsC*
j+i];
374 static void MvDot(
const Tpetra::MultiVector<Scalar,LO,GO,Node>&
A,
const Tpetra::MultiVector<Scalar,LO,GO,Node>&
B, std::vector<BaseScalar> &dots)
377 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::MvDot(A,B,dots): A and B must have the same number of vectors.");
378 #ifdef HAVE_TPETRA_DEBUG 380 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::MvDot(A,B,dots): dots must have room for all dot products.");
386 for (
unsigned int i=0; i<
A.getNumVectors(); i++)
387 dots[i] = pce_dots[i].coeff(0);
392 #ifdef HAVE_TPETRA_DEBUG 394 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::MvNorm(mv,normvec): normvec must have room for all norms.");
400 mv.norm1(av(0,mv.getNumVectors()));
406 mv.norm2(av(0,mv.getNumVectors()));
412 mv.normInf(av(0,mv.getNumVectors()));
421 static void SetBlock(
const Tpetra::MultiVector<Scalar,LO,GO,Node>&
A,
const std::vector<int>& index, Tpetra::MultiVector<Scalar,LO,GO,Node>& mv )
423 #ifdef HAVE_TPETRA_DEBUG 425 "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::SetBlock(A,index,mv): index must be the same size as A.");
428 if ((
typename std::vector<int>::size_type)
A.getNumVectors() > index.size()) {
435 mvsub = Teuchos::null;
439 SetBlock (
const Tpetra::MultiVector<Scalar,LO,GO,Node>&
A,
441 Tpetra::MultiVector<Scalar,LO,GO,Node>& mv)
451 const bool overflow = maxInt <
A.getNumVectors() && maxInt < mv.getNumVectors();
454 std::ostringstream os;
455 os <<
"Belos::MultiVecTraits<Scalar, Tpetra::MultiVector<Scalar, ..." 456 "> >::SetBlock(A, index=[" << index.
lbound() <<
", " 457 << index.
ubound() <<
"], mv): ";
459 os.str() <<
"Number of columns in the input multi" 460 "vector 'A' (a size_t) overflows int.");
462 os.str() <<
"Number of columns in the output multi" 463 "vector 'mv' (a size_t) overflows int.");
466 const int numColsA =
static_cast<int> (
A.getNumVectors());
467 const int numColsMv =
static_cast<int> (mv.getNumVectors());
469 const bool validIndex = index.
lbound() >= 0 && index.
ubound() < numColsMv;
471 const bool validSource = index.
size() <= numColsA;
473 if (! validIndex || ! validSource)
475 std::ostringstream os;
476 os <<
"Belos::MultiVecTraits<Scalar, Tpetra::MultiVector<Scalar, ..." 477 "> >::SetBlock(A, index=[" << index.
lbound() <<
", " 478 << index.
ubound() <<
"], mv): ";
480 os.str() <<
"Range lower bound must be nonnegative.");
482 os.str() <<
"Range upper bound must be less than " 483 "the number of columns " << numColsA <<
" in the " 484 "'mv' output argument.");
486 os.str() <<
"Range must have no more elements than" 487 " the number of columns " << numColsA <<
" in the " 488 "'A' input argument.");
498 if (index.
lbound() == 0 && index.
ubound()+1 == numColsMv)
499 mv_view = Teuchos::rcpFromRef (mv);
501 mv_view = CloneViewNonConst (mv, index);
507 if (index.
size() == numColsA)
508 A_view = Teuchos::rcpFromRef (
A);
521 Assign (
const Tpetra::MultiVector<Scalar,LO,GO,Node>&
A,
522 Tpetra::MultiVector<Scalar,LO,GO,Node>& mv)
532 const bool overflow = maxInt <
A.getNumVectors() && maxInt < mv.getNumVectors();
535 std::ostringstream os;
536 os <<
"Belos::MultiVecTraits<Scalar, Tpetra::MultiVector<Scalar, ..." 537 "> >::Assign(A, mv): ";
539 os.str() <<
"Number of columns in the input multi" 540 "vector 'A' (a size_t) overflows int.");
542 os.str() <<
"Number of columns in the output multi" 543 "vector 'mv' (a size_t) overflows int.");
547 const int numColsA =
static_cast<int> (
A.getNumVectors());
548 const int numColsMv =
static_cast<int> (mv.getNumVectors());
549 if (numColsA > numColsMv)
551 std::ostringstream os;
552 os <<
"Belos::MultiVecTraits<Scalar, Tpetra::MultiVector<Scalar, ..." 553 "> >::Assign(A, mv): ";
555 os.str() <<
"Input multivector 'A' has " 556 << numColsA <<
" columns, but output multivector " 557 "'mv' has only " << numColsMv <<
" columns.");
565 if (numColsA == numColsMv)
576 static void MvRandom( Tpetra::MultiVector<Scalar,LO,GO,Node>& mv )
582 { mv.putScalar(alpha); }
584 static void MvPrint(
const Tpetra::MultiVector<Scalar,LO,GO,Node>& mv, std::ostream& os )
590 #ifdef HAVE_BELOS_TSQR 591 typedef Tpetra::TsqrAdaptor< Tpetra::MultiVector< Scalar, LO, GO, Node > > tsqr_adaptor_type;
595 #endif // HAVE_BELOS_TSQR 605 template <
class BaseScalar,
class Storage,
class LO,
class GO,
class Node>
606 class OperatorTraits <BaseScalar,
Tpetra::MultiVector<Sacado::PCE::OrthogPoly<BaseScalar,Storage>,LO,GO,Node>, Tpetra::Operator<Sacado::PCE::OrthogPoly<BaseScalar,Storage>,LO,GO,Node> >
611 Apply (
const Tpetra::Operator<Scalar,LO,GO,Node>& Op,
612 const Tpetra::MultiVector<Scalar,LO,GO,Node>& X,
613 Tpetra::MultiVector<Scalar,LO,GO,Node>& Y,
614 ETrans trans=NOTRANS)
631 const std::string otName =
"Belos::OperatorTraits<" + scalarName
632 +
"," + loName +
"," + goName +
"," + nodeName +
">";
634 "get here; fell through a switch statement. " 635 "Please report this bug to the Belos developers.");
642 return Op.hasTransposeApply ();
static void MvPrint(const Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, std::ostream &os)
Belos::OperatorTraits< BaseScalar, Tpetra::MultiVector< Sacado::PCE::OrthogPoly< BaseScalar, Storage >, LO, GO, Node >, Tpetra::Operator< Sacado::PCE::OrthogPoly< BaseScalar, Storage >, LO, GO, Node > >::HasApplyTranspose static bool HasApplyTranspose(const Tpetra::Operator< Scalar, LO, GO, Node > &Op)
static void MvInit(Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, Scalar alpha=Teuchos::ScalarTraits< Scalar >::zero())
static Teuchos::RCP< Tpetra::MultiVector< Scalar, LO, GO, Node > > CloneCopy(const Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const std::vector< int > &index)
static void Assign(const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, Tpetra::MultiVector< Scalar, LO, GO, Node > &mv)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
static void MvTransMv(Scalar alpha, const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, const Tpetra::MultiVector< Scalar, LO, GO, Node > &B, Teuchos::SerialDenseMatrix< int, BaseScalar > &C)
static void MvRandom(Tpetra::MultiVector< Scalar, LO, GO, Node > &mv)
Sacado::PCE::OrthogPoly< BaseScalar, Storage > Scalar
static void MvTimesMatAddMv(Scalar alpha, const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, const Teuchos::SerialDenseMatrix< int, Scalar > &B, Scalar beta, Tpetra::MultiVector< Scalar, LO, GO, Node > &mv)
static void SetBlock(const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, const std::vector< int > &index, Tpetra::MultiVector< Scalar, LO, GO, Node > &mv)
static void MvScale(Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, Scalar alpha)
static ptrdiff_t GetGlobalLength(const Tpetra::MultiVector< Scalar, LO, GO, Node > &mv)
static Teuchos::RCP< const Tpetra::MultiVector< Scalar, LO, GO, Node > > CloneView(const Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const Teuchos::Range1D &index)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
static Teuchos::RCP< Tpetra::MultiVector< Scalar, LO, GO, Node > > CloneCopy(const Tpetra::MultiVector< Scalar, LO, GO, Node > &mv)
static Teuchos::RCP< Tpetra::MultiVector< Scalar, LO, GO, Node > > Clone(const Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const int numvecs)
static void MvScale(Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const std::vector< BaseScalar > &alphas)
static void MvTransMv(Scalar alpha, const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, const Tpetra::MultiVector< Scalar, LO, GO, Node > &B, Teuchos::SerialDenseMatrix< int, Scalar > &C)
static void MvTimesMatAddMv(Scalar alpha, const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, const Teuchos::SerialDenseMatrix< int, BaseScalar > &B, Scalar beta, Tpetra::MultiVector< Scalar, LO, GO, Node > &mv)
static void SetBlock(const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, const Teuchos::Range1D &index, Tpetra::MultiVector< Scalar, LO, GO, Node > &mv)
Belos::OperatorTraits< BaseScalar, Tpetra::MultiVector< Sacado::PCE::OrthogPoly< BaseScalar, Storage >, LO, GO, Node >, Tpetra::Operator< Sacado::PCE::OrthogPoly< BaseScalar, Storage >, LO, GO, Node > >::Apply 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 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)
static int GetNumberVecs(const Tpetra::MultiVector< Scalar, LO, GO, Node > &mv)
static Teuchos::RCP< Tpetra::MultiVector< Scalar, LO, GO, Node > > CloneViewNonConst(Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const std::vector< int > &index)
static Teuchos::RCP< Tpetra::MultiVector< Scalar, LO, GO, Node > > CloneCopy(const Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const Teuchos::Range1D &index)
static void MvNorm(const Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, std::vector< typename Teuchos::ScalarTraits< BaseScalar >::magnitudeType > &normvec, NormType type=TwoNorm)
static bool HasConstantStride(const Tpetra::MultiVector< Scalar, LO, GO, Node > &mv)
static void MvDot(const Tpetra::MultiVector< Scalar, LO, GO, Node > &A, const Tpetra::MultiVector< Scalar, LO, GO, Node > &B, std::vector< BaseScalar > &dots)
static Teuchos::RCP< const Tpetra::MultiVector< Scalar, LO, GO, Node > > CloneView(const Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const std::vector< int > &index)
static Teuchos::RCP< Tpetra::MultiVector< Scalar, LO, GO, Node > > CloneViewNonConst(Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const Teuchos::Range1D &index)
Belos::OperatorTraits< BaseScalar, Tpetra::MultiVector< Sacado::PCE::OrthogPoly< BaseScalar, Storage >, LO, GO, Node >, Tpetra::Operator< Sacado::PCE::OrthogPoly< BaseScalar, Storage >, LO, GO, Node > >::Scalar Sacado::PCE::OrthogPoly< BaseScalar, Storage > Scalar
static void MvScale(Tpetra::MultiVector< Scalar, LO, GO, Node > &mv, const std::vector< Scalar > &alphas)
static std::string name()