Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Belos_TpetraAdapter_MP_Vector.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Stokhos Package
5 // Copyright (2009) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Eric T. Phipps (etphipp@sandia.gov).
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef BELOS_TPETRA_ADAPTER_MP_VECTOR_HPP
43 #define BELOS_TPETRA_ADAPTER_MP_VECTOR_HPP
44 
45 #include "BelosTpetraAdapter.hpp"
48 #include "KokkosBlas.hpp"
49 
50 #ifdef HAVE_BELOS_TSQR
52 #endif // HAVE_BELOS_TSQR
53 
54 namespace Belos {
55 
57  //
58  // Implementation of Belos::MultiVecTraits for Tpetra::MultiVector.
59  //
61 
72  template<class Storage, class LO, class GO, class Node>
73  class MultiVecTraits<typename Storage::value_type,
74  Tpetra::MultiVector< Sacado::MP::Vector<Storage>,
75  LO, GO, Node > > {
76  public:
77  typedef typename Storage::ordinal_type s_ordinal;
78  typedef typename Storage::value_type BaseScalar;
80  typedef Tpetra::MultiVector<Scalar, LO, GO, Node> MV;
81  public:
82  typedef typename Tpetra::MultiVector<Scalar,LO,GO,Node>::dot_type dot_type;
83  typedef typename Tpetra::MultiVector<Scalar,LO,GO,Node>::mag_type mag_type;
84 
90  static Teuchos::RCP<MV> Clone (const MV& X, const int numVecs) {
91  Teuchos::RCP<MV> Y (new MV (X.getMap (), numVecs, false));
92  Y->setCopyOrView (Teuchos::View);
93  return Y;
94  }
95 
97  static Teuchos::RCP<MV> CloneCopy (const MV& X)
98  {
99  // Make a deep copy of X. The one-argument copy constructor
100  // does a shallow copy by default; the second argument tells it
101  // to do a deep copy.
102  Teuchos::RCP<MV> X_copy (new MV (X, Teuchos::Copy));
103  // Make Tpetra::MultiVector use the new view semantics. This is
104  // a no-op for the Kokkos refactor version of Tpetra; it only
105  // does something for the "classic" version of Tpetra. This
106  // shouldn't matter because Belos only handles MV through RCP
107  // and through this interface anyway, but it doesn't hurt to set
108  // it and make sure that it works.
109  X_copy->setCopyOrView (Teuchos::View);
110  return X_copy;
111  }
112 
124  static Teuchos::RCP<MV>
125  CloneCopy (const MV& mv, const std::vector<int>& index)
126  {
127 #ifdef HAVE_TPETRA_DEBUG
128  const char fnName[] = "Belos::MultiVecTraits::CloneCopy(mv,index)";
129  const size_t inNumVecs = mv.getNumVectors ();
131  index.size () > 0 && *std::min_element (index.begin (), index.end ()) < 0,
132  std::runtime_error, fnName << ": All indices must be nonnegative.");
134  index.size () > 0 &&
135  static_cast<size_t> (*std::max_element (index.begin (), index.end ())) >= inNumVecs,
136  std::runtime_error,
137  fnName << ": All indices must be strictly less than the number of "
138  "columns " << inNumVecs << " of the input multivector mv.");
139 #endif // HAVE_TPETRA_DEBUG
140 
141  // Tpetra wants an array of size_t, not of int.
142  Teuchos::Array<size_t> columns (index.size ());
143  for (std::vector<int>::size_type j = 0; j < index.size (); ++j) {
144  columns[j] = index[j];
145  }
146  // mfh 14 Aug 2014: Tpetra already detects and optimizes for a
147  // continuous column index range in MultiVector::subCopy, so we
148  // don't have to check here.
149  Teuchos::RCP<MV> X_copy = mv.subCopy (columns ());
150  X_copy->setCopyOrView (Teuchos::View);
151  return X_copy;
152  }
153 
160  static Teuchos::RCP<MV>
161  CloneCopy (const MV& mv, const Teuchos::Range1D& index)
162  {
163  const bool validRange = index.size() > 0 &&
164  index.lbound() >= 0 &&
165  index.ubound() < GetNumberVecs(mv);
166  if (! validRange) { // invalid range; generate error message
167  std::ostringstream os;
168  os << "Belos::MultiVecTraits::CloneCopy(mv,index=["
169  << index.lbound() << "," << index.ubound() << "]): ";
171  index.size() == 0, std::invalid_argument,
172  os.str() << "Empty index range is not allowed.");
174  index.lbound() < 0, std::invalid_argument,
175  os.str() << "Index range includes negative index/ices, which is not "
176  "allowed.");
178  index.ubound() >= GetNumberVecs(mv), std::invalid_argument,
179  os.str() << "Index range exceeds number of vectors "
180  << mv.getNumVectors() << " in the input multivector.");
181  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
182  os.str() << "Should never get here!");
183  }
184  Teuchos::RCP<MV> X_copy = mv.subCopy (index);
185  X_copy->setCopyOrView (Teuchos::View);
186  return X_copy;
187  }
188 
189 
190  static Teuchos::RCP<MV>
191  CloneViewNonConst (MV& mv, const std::vector<int>& index)
192  {
193 #ifdef HAVE_TPETRA_DEBUG
194  const char fnName[] = "Belos::MultiVecTraits::CloneViewNonConst(mv,index)";
195  const size_t numVecs = mv.getNumVectors ();
197  index.size () > 0 && *std::min_element (index.begin (), index.end ()) < 0,
198  std::invalid_argument,
199  fnName << ": All indices must be nonnegative.");
201  index.size () > 0 &&
202  static_cast<size_t> (*std::max_element (index.begin (), index.end ())) >= numVecs,
203  std::invalid_argument,
204  fnName << ": All indices must be strictly less than the number of "
205  "columns " << numVecs << " in the input MultiVector mv.");
206 #endif // HAVE_TPETRA_DEBUG
207 
208  // Tpetra wants an array of size_t, not of int.
209  Teuchos::Array<size_t> columns (index.size ());
210  for (std::vector<int>::size_type j = 0; j < index.size (); ++j) {
211  columns[j] = index[j];
212  }
213  // mfh 14 Aug 2014: Tpetra already detects and optimizes for a
214  // continuous column index range in
215  // MultiVector::subViewNonConst, so we don't have to check here.
216  Teuchos::RCP<MV> X_view = mv.subViewNonConst (columns ());
217  X_view->setCopyOrView (Teuchos::View);
218  return X_view;
219  }
220 
221  static Teuchos::RCP<MV>
223  {
224  // NOTE (mfh 11 Jan 2011) We really should check for possible
225  // overflow of int here. However, the number of columns in a
226  // multivector typically fits in an int.
227  const int numCols = static_cast<int> (mv.getNumVectors());
228  const bool validRange = index.size() > 0 &&
229  index.lbound() >= 0 && index.ubound() < numCols;
230  if (! validRange) {
231  std::ostringstream os;
232  os << "Belos::MultiVecTraits::CloneViewNonConst(mv,index=["
233  << index.lbound() << ", " << index.ubound() << "]): ";
235  index.size() == 0, std::invalid_argument,
236  os.str() << "Empty index range is not allowed.");
238  index.lbound() < 0, std::invalid_argument,
239  os.str() << "Index range includes negative inde{x,ices}, which is "
240  "not allowed.");
242  index.ubound() >= numCols, std::invalid_argument,
243  os.str() << "Index range exceeds number of vectors " << numCols
244  << " in the input multivector.");
246  true, std::logic_error,
247  os.str() << "Should never get here!");
248  }
249  Teuchos::RCP<MV> X_view = mv.subViewNonConst (index);
250  X_view->setCopyOrView (Teuchos::View);
251  return X_view;
252  }
253 
255  CloneView (const MV& mv, const std::vector<int>& index)
256  {
257 #ifdef HAVE_TPETRA_DEBUG
258  const char fnName[] = "Belos::MultiVecTraits<Scalar, "
259  "Tpetra::MultiVector<...> >::CloneView(mv,index)";
260  const size_t numVecs = mv.getNumVectors ();
262  *std::min_element (index.begin (), index.end ()) < 0,
263  std::invalid_argument,
264  fnName << ": All indices must be nonnegative.");
266  static_cast<size_t> (*std::max_element (index.begin (), index.end ())) >= numVecs,
267  std::invalid_argument,
268  fnName << ": All indices must be strictly less than the number of "
269  "columns " << numVecs << " in the input MultiVector mv.");
270 #endif // HAVE_TPETRA_DEBUG
271 
272  // Tpetra wants an array of size_t, not of int.
273  Teuchos::Array<size_t> columns (index.size ());
274  for (std::vector<int>::size_type j = 0; j < index.size (); ++j) {
275  columns[j] = index[j];
276  }
277  // mfh 14 Aug 2014: Tpetra already detects and optimizes for a
278  // continuous column index range in MultiVector::subView, so we
279  // don't have to check here.
280  Teuchos::RCP<const MV> X_view = mv.subView (columns);
281  Teuchos::rcp_const_cast<MV> (X_view)->setCopyOrView (Teuchos::View);
282  return X_view;
283  }
284 
286  CloneView (const MV& mv, const Teuchos::Range1D& index)
287  {
288  // NOTE (mfh 11 Jan 2011) We really should check for possible
289  // overflow of int here. However, the number of columns in a
290  // multivector typically fits in an int.
291  const int numCols = static_cast<int> (mv.getNumVectors());
292  const bool validRange = index.size() > 0 &&
293  index.lbound() >= 0 && index.ubound() < numCols;
294  if (! validRange) {
295  std::ostringstream os;
296  os << "Belos::MultiVecTraits::CloneView(mv, index=["
297  << index.lbound () << ", " << index.ubound() << "]): ";
298  TEUCHOS_TEST_FOR_EXCEPTION(index.size() == 0, std::invalid_argument,
299  os.str() << "Empty index range is not allowed.");
300  TEUCHOS_TEST_FOR_EXCEPTION(index.lbound() < 0, std::invalid_argument,
301  os.str() << "Index range includes negative index/ices, which is not "
302  "allowed.");
304  index.ubound() >= numCols, std::invalid_argument,
305  os.str() << "Index range exceeds number of vectors " << numCols
306  << " in the input multivector.");
307  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
308  os.str() << "Should never get here!");
309  }
310  Teuchos::RCP<const MV> X_view = mv.subView (index);
311  Teuchos::rcp_const_cast<MV> (X_view)->setCopyOrView (Teuchos::View);
312  return X_view;
313  }
314 
315  static ptrdiff_t GetGlobalLength (const MV& mv) {
316  return static_cast<ptrdiff_t> (mv.getGlobalLength ());
317  }
318 
319  static int GetNumberVecs (const MV& mv) {
320  return static_cast<int> (mv.getNumVectors ());
321  }
322 
323  static bool HasConstantStride (const MV& mv) {
324  return mv.isConstantStride ();
325  }
326 
327  static void
328  MvTimesMatAddMv (const dot_type& alpha,
329  const Tpetra::MultiVector<Scalar,LO,GO,Node>& A,
331  const dot_type& beta,
332  Tpetra::MultiVector<Scalar,LO,GO,Node>& C)
333  {
334  using Teuchos::RCP;
335  using Teuchos::rcp;
336 
337  // Check if numRowsB == numColsB == 1, in which case we can call update()
338  const int numRowsB = B.numRows ();
339  const int numColsB = B.numCols ();
340  const int strideB = B.stride ();
341  if (numRowsB == 1 && numColsB == 1) {
342  C.update (alpha*B(0,0), A, beta);
343  return;
344  }
345 
346  // Ensure A and C have constant stride
347  RCP<const MV> Atmp, Ctmp;
348  if (A.isConstantStride() == false) Atmp = rcp (new MV (A, Teuchos::Copy));
349  else Atmp = rcp(&A,false);
350 
351  if (C.isConstantStride() == false) Ctmp = rcp (new MV (C, Teuchos::Copy));
352  else Ctmp = rcp(&C,false);
353 
354  // Create flattened view's
355  typedef Tpetra::MultiVector<dot_type,LO,GO,Node> FMV;
356  typedef typename FMV::dual_view_type::t_dev flat_view_type;
358  flat_view_type flat_A_view = Atmp->template getLocalView<execution_space>();
359  flat_view_type flat_C_view = Ctmp->template getLocalView<execution_space>();
360 
361  // Create a view for B on the host
362  typedef Kokkos::View<dot_type**, Kokkos::LayoutLeft, Kokkos::HostSpace> b_host_view_type;
363  b_host_view_type B_view_host_input( B.values(), strideB, numColsB);
364  auto B_view_host = Kokkos::subview(B_view_host_input,
365  Kokkos::pair<int,int>(0,numRowsB),
366  Kokkos::pair<int,int>(0,numColsB));
367 
368  // Create view for B on the device -- need to be careful to get the
369  // right stride to match B
370  typedef Kokkos::View<dot_type**, Kokkos::LayoutLeft, execution_space> b_view_type;
371  typedef Kokkos::View<dot_type*, Kokkos::LayoutLeft, execution_space> b_1d_view_type;
372  b_1d_view_type B_1d_view_dev(Kokkos::ViewAllocateWithoutInitializing("B"), numRowsB*numColsB);
373  b_view_type B_view_dev( B_1d_view_dev.data(), numRowsB, numColsB);
374  Kokkos::deep_copy(B_view_dev, B_view_host);
375 
376  // Do local multiply
377  {
378  const char ctransA = 'N', ctransB = 'N';
380  &ctransA, &ctransB,
381  alpha, flat_A_view, B_view_dev, beta, flat_C_view);
382  }
383  // Copy back to C if we made a copy
384  if (C.isConstantStride() == false)
385  C.assign(*Ctmp);
386  }
387 
395  static void
396  MvAddMv (Scalar alpha,
397  const Tpetra::MultiVector<Scalar,LO,GO,Node>& A,
398  Scalar beta,
399  const Tpetra::MultiVector<Scalar,LO,GO,Node>& B,
400  Tpetra::MultiVector<Scalar,LO,GO,Node>& mv)
401  {
402  mv.update (alpha, A, beta, B, Teuchos::ScalarTraits<Scalar>::zero ());
403  }
404 
405  static void
406  MvScale (Tpetra::MultiVector<Scalar,LO,GO,Node>& mv,
407  const Scalar& alpha)
408  {
409  mv.scale (alpha);
410  }
411 
412  static void
413  MvScale (Tpetra::MultiVector<Scalar,LO,GO,Node>& mv,
414  const std::vector<BaseScalar>& alphas)
415  {
416  std::vector<Scalar> alphas_mp(alphas.size());
417  const size_t sz = alphas.size();
418  for (size_t i=0; i<sz; ++i)
419  alphas_mp[i] = alphas[i];
420  mv.scale (alphas_mp);
421  }
422 
423  static void
424  MvScale (Tpetra::MultiVector<Scalar,LO,GO,Node>& mv,
425  const std::vector<Scalar>& alphas)
426  {
427  mv.scale (alphas);
428  }
429 
430  static void
432  const Tpetra::MultiVector<Scalar,LO,GO,Node>& A,
433  const Tpetra::MultiVector<Scalar,LO,GO,Node>& B,
435  {
436  using Teuchos::Comm;
437  using Teuchos::RCP;
438  using Teuchos::rcp;
439  using Teuchos::REDUCE_SUM;
440  using Teuchos::reduceAll;
441 
442  // Check if numRowsC == numColsC == 1, in which case we can call dot()
443  const int numRowsC = C.numRows ();
444  const int numColsC = C.numCols ();
445  const int strideC = C.stride ();
446  if (numRowsC == 1 && numColsC == 1) {
447  if (alpha == Teuchos::ScalarTraits<Scalar>::zero ()) {
448  // Short-circuit, as required by BLAS semantics.
449  C(0,0) = alpha;
450  return;
451  }
452  A.dot (B, Teuchos::ArrayView<dot_type> (C.values (), 1));
453  if (alpha != Teuchos::ScalarTraits<Scalar>::one ()) {
454  C(0,0) *= alpha;
455  }
456  return;
457  }
458 
459  // Ensure A and B have constant stride
460  RCP<const MV> Atmp, Btmp;
461  if (A.isConstantStride() == false) Atmp = rcp (new MV (A, Teuchos::Copy));
462  else Atmp = rcp(&A,false);
463 
464  if (B.isConstantStride() == false) Btmp = rcp (new MV (B, Teuchos::Copy));
465  else Btmp = rcp(&B,false);
466 
467  // Create flattened Kokkos::MultiVector's
468  typedef Tpetra::MultiVector<dot_type,LO,GO,Node> FMV;
469  typedef typename FMV::dual_view_type::t_dev flat_view_type;
471  flat_view_type flat_A_view = Atmp->template getLocalView<execution_space>();
472  flat_view_type flat_B_view = Btmp->template getLocalView<execution_space>();
473 
474  // Create a view for C on the host
475  typedef Kokkos::View<dot_type**, Kokkos::LayoutLeft, Kokkos::HostSpace> c_host_view_type;
476  c_host_view_type C_view_host_input( C.values(), strideC, numColsC);
477  auto C_view_host = Kokkos::subview(C_view_host_input,
478  Kokkos::pair<int,int>(0,numRowsC),
479  Kokkos::pair<int,int>(0,numColsC));
480 
481  // Create view for C on the device -- need to be careful to get the
482  // right stride to match C (allow setting to 0 for first-touch)
483  typedef Kokkos::View<dot_type**, Kokkos::LayoutLeft, execution_space> c_view_type;
484  typedef Kokkos::View<dot_type*, Kokkos::LayoutLeft, execution_space> c_1d_view_type;
485  c_1d_view_type C_1d_view_dev("C", numRowsC*numColsC);
486  c_view_type C_view_dev( C_1d_view_dev.data(), numRowsC, numColsC);
487 
488  // Do local multiply
489  {
490  const char ctransA = 'C', ctransB = 'N';
492  &ctransA, &ctransB,
493  alpha, flat_A_view, flat_B_view,
494  Kokkos::Details::ArithTraits<dot_type>::zero(),
495  C_view_dev);
496  }
497 
498  // reduce across processors -- could check for RDMA
499  RCP<const Comm<int> > pcomm = A.getMap()->getComm ();
500  if (pcomm->getSize () == 1)
501  Kokkos::deep_copy(C_view_host, C_view_dev);
502  else {
503  typedef Kokkos::View<dot_type*, Kokkos::LayoutLeft, Kokkos::HostSpace> c_1d_host_view_type;
504  c_1d_host_view_type C_1d_view_tmp(Kokkos::ViewAllocateWithoutInitializing("C_tmp"), strideC*numColsC);
505  c_host_view_type C_view_tmp( C_1d_view_tmp.data(),
506  strideC, numColsC);
507  Kokkos::deep_copy(Kokkos::subview(C_view_tmp,
508  Kokkos::pair<int,int>(0,numRowsC),
509  Kokkos::pair<int,int>(0,numColsC)) ,
510  C_view_dev);
511  reduceAll<int> (*pcomm, REDUCE_SUM, strideC*numColsC,
512  C_view_tmp.data(),
513  C_view_host.data());
514  }
515  }
516 
518  static void
519  MvDot (const Tpetra::MultiVector<Scalar,LO,GO,Node>& A,
520  const Tpetra::MultiVector<Scalar,LO,GO,Node>& B,
521  std::vector<dot_type>& dots)
522  {
523  const size_t numVecs = A.getNumVectors ();
524 
526  numVecs != B.getNumVectors (), std::invalid_argument,
527  "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::MvDot(A,B,dots): "
528  "A and B must have the same number of columns. "
529  "A has " << numVecs << " column(s), "
530  "but B has " << B.getNumVectors () << " column(s).");
531 #ifdef HAVE_TPETRA_DEBUG
533  dots.size() < numVecs, std::invalid_argument,
534  "Belos::MultiVecTraits<Scalar,Tpetra::MultiVector>::MvDot(A,B,dots): "
535  "The output array 'dots' must have room for all dot products. "
536  "A and B each have " << numVecs << " column(s), "
537  "but 'dots' only has " << dots.size() << " entry(/ies).");
538 #endif // HAVE_TPETRA_DEBUG
539 
541  A.dot (B, av (0, numVecs));
542  }
543 
545  static void
546  MvNorm (const Tpetra::MultiVector<Scalar,LO,GO,Node>& mv,
547  std::vector<mag_type> &normvec,
548  NormType type=TwoNorm)
549  {
550 
551 #ifdef HAVE_TPETRA_DEBUG
552  typedef std::vector<int>::size_type size_type;
554  normvec.size () < static_cast<size_type> (mv.getNumVectors ()),
555  std::invalid_argument,
556  "Belos::MultiVecTraits::MvNorm(mv,normvec): The normvec output "
557  "argument must have at least as many entries as the number of vectors "
558  "(columns) in the MultiVector mv. normvec.size() = " << normvec.size ()
559  << " < mv.getNumVectors() = " << mv.getNumVectors () << ".");
560 #endif
561 
562  Teuchos::ArrayView<mag_type> av(normvec);
563  switch (type) {
564  case OneNorm:
565  mv.norm1(av(0,mv.getNumVectors()));
566  break;
567  case TwoNorm:
568  mv.norm2(av(0,mv.getNumVectors()));
569  break;
570  case InfNorm:
571  mv.normInf(av(0,mv.getNumVectors()));
572  break;
573  default:
574  // Throw logic_error rather than invalid_argument, because if
575  // we get here, it's probably the fault of a Belos solver,
576  // rather than a user giving Belos an invalid input.
578  true, std::logic_error,
579  "Belos::MultiVecTraits::MvNorm: Invalid NormType value " << type
580  << ". Valid values are OneNorm=" << OneNorm << ", TwoNorm="
581  << TwoNorm <<", and InfNorm=" << InfNorm << ". If you are a Belos "
582  "user and have not modified Belos in any way, and you get this "
583  "message, then this is probably a bug in the Belos solver you were "
584  "using. Please report this to the Belos developers.");
585  }
586  }
587 
588  static void
589  SetBlock (const MV& A, const std::vector<int>& index, MV& mv)
590  {
591  using Teuchos::Range1D;
592  using Teuchos::RCP;
593  const size_t inNumVecs = A.getNumVectors ();
594 #ifdef HAVE_TPETRA_DEBUG
596  inNumVecs < static_cast<size_t> (index.size ()), std::invalid_argument,
597  "Belos::MultiVecTraits::SetBlock(A,index,mv): 'index' argument must "
598  "have no more entries as the number of columns in the input MultiVector"
599  " A. A.getNumVectors() = " << inNumVecs << " < index.size () = "
600  << index.size () << ".");
601 #endif // HAVE_TPETRA_DEBUG
602  RCP<MV> mvsub = CloneViewNonConst (mv, index);
603  if (inNumVecs > static_cast<size_t> (index.size ())) {
604  RCP<const MV> Asub = A.subView (Range1D (0, index.size () - 1));
605  ::Tpetra::deep_copy (*mvsub, *Asub);
606  } else {
607  ::Tpetra::deep_copy (*mvsub, A);
608  }
609  }
610 
611  static void
612  SetBlock (const MV& A, const Teuchos::Range1D& index, MV& mv)
613  {
614  // Range1D bounds are signed; size_t is unsigned.
615  // Assignment of Tpetra::MultiVector is a deep copy.
616 
617  // Tpetra::MultiVector::getNumVectors() returns size_t. It's
618  // fair to assume that the number of vectors won't overflow int,
619  // since the typical use case of multivectors involves few
620  // columns, but it's friendly to check just in case.
621  const size_t maxInt =
622  static_cast<size_t> (Teuchos::OrdinalTraits<int>::max ());
623  const bool overflow =
624  maxInt < A.getNumVectors () && maxInt < mv.getNumVectors ();
625  if (overflow) {
626  std::ostringstream os;
627  os << "Belos::MultiVecTraits::SetBlock(A, index=[" << index.lbound ()
628  << ", " << index.ubound () << "], mv): ";
630  maxInt < A.getNumVectors (), std::range_error, os.str () << "Number "
631  "of columns (size_t) in the input MultiVector 'A' overflows int.");
633  maxInt < mv.getNumVectors (), std::range_error, os.str () << "Number "
634  "of columns (size_t) in the output MultiVector 'mv' overflows int.");
635  }
636  // We've already validated the static casts above.
637  const int numColsA = static_cast<int> (A.getNumVectors ());
638  const int numColsMv = static_cast<int> (mv.getNumVectors ());
639  // 'index' indexes into mv; it's the index set of the target.
640  const bool validIndex =
641  index.lbound () >= 0 && index.ubound () < numColsMv;
642  // We can't take more columns out of A than A has.
643  const bool validSource = index.size () <= numColsA;
644 
645  if (! validIndex || ! validSource) {
646  std::ostringstream os;
647  os << "Belos::MultiVecTraits::SetBlock(A, index=[" << index.lbound ()
648  << ", " << index.ubound () << "], mv): ";
650  index.lbound() < 0, std::invalid_argument,
651  os.str() << "Range lower bound must be nonnegative.");
653  index.ubound() >= numColsMv, std::invalid_argument,
654  os.str() << "Range upper bound must be less than the number of "
655  "columns " << numColsA << " in the 'mv' output argument.");
657  index.size() > numColsA, std::invalid_argument,
658  os.str() << "Range must have no more elements than the number of "
659  "columns " << numColsA << " in the 'A' input argument.");
661  true, std::logic_error, "Should never get here!");
662  }
663 
664  // View of the relevant column(s) of the target multivector mv.
665  // We avoid view creation overhead by only creating a view if
666  // the index range is different than [0, (# columns in mv) - 1].
667  Teuchos::RCP<MV> mv_view;
668  if (index.lbound () == 0 && index.ubound () + 1 == numColsMv) {
669  mv_view = Teuchos::rcpFromRef (mv); // Non-const, non-owning RCP
670  } else {
671  mv_view = CloneViewNonConst (mv, index);
672  }
673 
674  // View of the relevant column(s) of the source multivector A.
675  // If A has fewer columns than mv_view, then create a view of
676  // the first index.size() columns of A.
677  Teuchos::RCP<const MV> A_view;
678  if (index.size () == numColsA) {
679  A_view = Teuchos::rcpFromRef (A); // Const, non-owning RCP
680  } else {
681  A_view = CloneView (A, Teuchos::Range1D (0, index.size () - 1));
682  }
683 
684  ::Tpetra::deep_copy (*mv_view, *A_view);
685  }
686 
687  static void Assign (const MV& A, MV& mv)
688  {
689  const char errPrefix[] = "Belos::MultiVecTraits::Assign(A, mv): ";
690 
691  // Range1D bounds are signed; size_t is unsigned.
692  // Assignment of Tpetra::MultiVector is a deep copy.
693 
694  // Tpetra::MultiVector::getNumVectors() returns size_t. It's
695  // fair to assume that the number of vectors won't overflow int,
696  // since the typical use case of multivectors involves few
697  // columns, but it's friendly to check just in case.
698  const size_t maxInt =
699  static_cast<size_t> (Teuchos::OrdinalTraits<int>::max ());
700  const bool overflow =
701  maxInt < A.getNumVectors () && maxInt < mv.getNumVectors ();
702  if (overflow) {
704  maxInt < A.getNumVectors(), std::range_error,
705  errPrefix << "Number of columns in the input multivector 'A' "
706  "(a size_t) overflows int.");
708  maxInt < mv.getNumVectors(), std::range_error,
709  errPrefix << "Number of columns in the output multivector 'mv' "
710  "(a size_t) overflows int.");
712  true, std::logic_error, "Should never get here!");
713  }
714  // We've already validated the static casts above.
715  const int numColsA = static_cast<int> (A.getNumVectors ());
716  const int numColsMv = static_cast<int> (mv.getNumVectors ());
717  if (numColsA > numColsMv) {
719  numColsA > numColsMv, std::invalid_argument,
720  errPrefix << "Input multivector 'A' has " << numColsA << " columns, "
721  "but output multivector 'mv' has only " << numColsMv << " columns.");
722  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "Should never get here!");
723  }
724  if (numColsA == numColsMv) {
725  ::Tpetra::deep_copy (mv, A);
726  } else {
727  Teuchos::RCP<MV> mv_view =
728  CloneViewNonConst (mv, Teuchos::Range1D (0, numColsA - 1));
729  ::Tpetra::deep_copy (*mv_view, A);
730  }
731  }
732 
733 
734  static void MvRandom (MV& mv) {
735  mv.randomize ();
736  }
737 
738  static void
739  MvInit (MV& mv, const Scalar alpha = Teuchos::ScalarTraits<Scalar>::zero ())
740  {
741  mv.putScalar (alpha);
742  }
743 
744  static void MvPrint (const MV& mv, std::ostream& os) {
745  Teuchos::FancyOStream fos (Teuchos::rcpFromRef (os));
746  mv.describe (fos, Teuchos::VERB_EXTREME);
747  }
748 
749 #ifdef HAVE_BELOS_TSQR
750  typedef Tpetra::TsqrAdaptor< Tpetra::MultiVector< Scalar, LO, GO, Node > > tsqr_adaptor_type;
754 #endif // HAVE_BELOS_TSQR
755  };
756 
758  //
759  // Implementation of the Belos::OperatorTraits for Tpetra::Operator.
760  //
762 
764  template <class Storage, class LO, class GO, class Node>
765  class OperatorTraits <typename Storage::value_type,
766  Tpetra::MultiVector<Sacado::MP::Vector<Storage>,
767  LO,GO,Node>,
768  Tpetra::Operator<Sacado::MP::Vector<Storage>,
769  LO,GO,Node> >
770  {
771  public:
773  static void
774  Apply (const Tpetra::Operator<Scalar,LO,GO,Node>& Op,
775  const Tpetra::MultiVector<Scalar,LO,GO,Node>& X,
776  Tpetra::MultiVector<Scalar,LO,GO,Node>& Y,
777  ETrans trans=NOTRANS)
778  {
779  switch (trans) {
780  case NOTRANS:
781  Op.apply(X,Y,Teuchos::NO_TRANS);
782  break;
783  case TRANS:
784  Op.apply(X,Y,Teuchos::TRANS);
785  break;
786  case CONJTRANS:
787  Op.apply(X,Y,Teuchos::CONJ_TRANS);
788  break;
789  default:
790  const std::string scalarName = Teuchos::TypeNameTraits<Scalar>::name();
791  const std::string loName = Teuchos::TypeNameTraits<LO>::name();
792  const std::string goName = Teuchos::TypeNameTraits<GO>::name();
793  const std::string nodeName = Teuchos::TypeNameTraits<Node>::name();
794  const std::string otName = "Belos::OperatorTraits<" + scalarName
795  + "," + loName + "," + goName + "," + nodeName + ">";
796  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, otName << ": Should never "
797  "get here; fell through a switch statement. "
798  "Please report this bug to the Belos developers.");
799  }
800  }
801 
802  static bool
803  HasApplyTranspose (const Tpetra::Operator<Scalar,LO,GO,Node>& Op)
804  {
805  return Op.hasTransposeApply ();
806  }
807  };
808 
809 } // end of Belos namespace
810 
811 #endif
812 // end of file BELOS_TPETRA_ADAPTER_HPP
ScalarType * values() const
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)
Ordinal ubound() 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)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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)
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]).
static Teuchos::RCP< MV > Clone(const MV &X, const int numVecs)
Create a new MultiVector with numVecs columns.
Ordinal lbound() const
OrdinalType numCols() const
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.
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)
Ordinal size() const
OrdinalType stride() const
static std::string name()
OrdinalType numRows() const
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 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